Jpp  master_rocky-40-g5f0272dcd
the software that should make you happy
JHermiteSpline.hh
Go to the documentation of this file.
1 #ifndef __JTOOLS__JHERMITESPLINE__
2 #define __JTOOLS__JHERMITESPLINE__
3 
4 #include <utility>
5 #include <sstream>
6 #include <cmath>
7 
8 #include "JMath/JZero.hh"
9 #include "JLang/JException.hh"
10 #include "JLang/JClass.hh"
12 #include "JTools/JFunctional.hh"
13 #include "JTools/JDistance.hh"
14 #include "JTools/JResult.hh"
16 #include "JTools/JMapCollection.hh"
17 #include "JTools/JElement.hh"
18 
19 
20 /**
21  * \author mdejong
22  */
23 
24 namespace JTOOLS {}
25 namespace JPP { using namespace JTOOLS; }
26 
27 namespace JTOOLS {
28 
29  using JLANG::JNoValue;
33 
34 
35  /**
36  * Template base class spline interpolations.
37  *
38  * This class implements the JFunctional interface.
39  *
40  * Note that the data structure of the elements in the collection should have the additional methods:
41  * <pre>
42  * ordinate_type getU() const;
43  * void setU(ordinate_type u);
44  * </pre>
45  * to get and set the derivatives, respectively.
46  *
47  * Note that -by default- the compilation is for a monotonous interpolation.
48  */
49  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
51  public JCollection_t<JElement_t, JDistance_t>,
52  public virtual JFunctional<>
53  {
54  public:
55 
56  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
57 
58  typedef typename collection_type::abscissa_type abscissa_type;
59  typedef typename collection_type::ordinate_type ordinate_type;
60  typedef typename collection_type::value_type value_type;
61  typedef typename collection_type::distance_type distance_type;
62 
63  typedef typename collection_type::const_iterator const_iterator;
64  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
65  typedef typename collection_type::iterator iterator;
66  typedef typename collection_type::reverse_iterator reverse_iterator;
67 
69 
70 
71  /**
72  * Determination of derivatives.
73  *
74  * \param monotone monotone
75  */
76  void compile(const bool monotone)
77  {
78  using namespace std;
79 
80  if (this->size() >= 2u) {
81 
82  {
83  iterator j = this->begin(), i = j++;
84 
85  i->setU((j->getY() - i->getY()) / this->getDistance(i->getX(), j->getX()));
86  }
87 
88  {
89  reverse_iterator j = this->rbegin(), i = j++;
90 
91  i->setU((j->getY() - i->getY()) / this->getDistance(i->getX(), j->getX()));
92  }
93 
94  for (iterator k = this->begin(), i = k++, j = k++; k != this->end(); ++i, ++j, ++k) {
95  j->setU(0.5 * ((j->getY() - i->getY()) / this->getDistance(i->getX(), j->getX()) +
96  (k->getY() - j->getY()) / this->getDistance(j->getX(), k->getX())));
97  }
98 
99  if (monotone) {
100 
101  for (iterator j = this->begin(), i = j++; j != this->end(); ++i, ++j) {
102  if (i->getY() == j->getY()) {
103  j->setU(JMATH::zero);
104  }
105  }
106 
107  for (iterator j = this->begin(), i = j++; j != this->end(); ++i, ++j) {
108 
109  const ordinate_type u = (j->getY() - i->getY()) / this->getDistance(i->getX(), j->getX());
110  const ordinate_type w = (i->getU()*i->getU() + j->getU()*j->getU());
111 
112  if (w > 9.0*u*u) {
113 
114  const ordinate_type v = 3.0*u/sqrt(w);
115 
116  i->setU(v*i->getU());
117  j->setU(v*j->getU());
118  }
119  }
120  }
121  }
122  }
123 
124 
125  protected:
126 
127  static abscissa_type h00 (abscissa_type t) { return (1.0 + 2*t) * (1.0 - t) * (1.0 - t); }
128  static abscissa_type h10 (abscissa_type t) { return t * (1.0 - t) * (1.0 - t); }
129  static abscissa_type h01 (abscissa_type t) { return t * t * (3.0 - 2*t); }
130  static abscissa_type h11 (abscissa_type t) { return t * t * (t - 1.0); }
131 
132  static abscissa_type h00p(abscissa_type t) { return 6 * t * (t - 1.0); }
133  static abscissa_type h10p(abscissa_type t) { return t * (3*t - 4.0) + 1.0; }
134  static abscissa_type h01p(abscissa_type t) { return 6 * t * (1.0 -t); }
135  static abscissa_type h11p(abscissa_type t) { return t * (3*t - 2.0); }
136 
137  static abscissa_type H00 (abscissa_type t) { return t * (t * t * (0.5*t - 1.0) + 1.0); }
138  static abscissa_type H10 (abscissa_type t) { return t * t * (t * (0.25*t - 2.0/3.0) + 0.5); }
139  static abscissa_type H01 (abscissa_type t) { return t * t * t * (1.0 - 0.5*t); }
140  static abscissa_type H11 (abscissa_type t) { return t * t * t * (0.25*t - 1.0/3.0); }
141 
142 
143  /**
144  * Default constructor.
145  */
147  {}
148 
149 
150  /**
151  * Determination of derivatives.
152  */
153  virtual void do_compile() override
154  {
155  compile(true);
156  }
157  };
158 
159 
160  /**
161  * Template definition for functional collection with spline interpolation.
162  */
163  template<class JElement_t,
164  template<class, class> class JCollection_t,
165  class JResult_t,
166  class JDistance_t>
168 
169 
170  /**
171  * Template specialisation for functional collection with spline interpolation.
172  */
173  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
174  class JHermiteSplineFunction<JElement_t,
175  JCollection_t,
176  typename JResultType<typename JElement_t::ordinate_type>::result_type,
177  JDistance_t> :
178  public JHermiteSplineCollection<JElement_t, JCollection_t, JDistance_t>,
179  public virtual JFunction<typename JElement_t::abscissa_type,
180  typename JResultType<typename JElement_t::ordinate_type>::result_type>
181  {
182  public:
183 
185 
186  typedef typename collection_type::abscissa_type abscissa_type;
187  typedef typename collection_type::ordinate_type ordinate_type;
188  typedef typename collection_type::value_type value_type;
189  typedef typename collection_type::distance_type distance_type;
190 
191  typedef typename collection_type::const_iterator const_iterator;
192  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
193  typedef typename collection_type::iterator iterator;
194  typedef typename collection_type::reverse_iterator reverse_iterator;
195 
198 
202 
203 
204  /**
205  * Default constructor.
206  */
208  {}
209 
210 
211  /**
212  * Recursive interpolation method implementation.
213  *
214  * \param pX pointer to abscissa values
215  * \return function value
216  */
217  virtual result_type evaluate(const argument_type* pX) const override
218  {
219  const argument_type x = *pX;
220 
221  if (this->size() <= 1u) {
222 
223  try {
224  return this->getExceptionHandler().action();
225  }
226  catch (const JException& error) {
227 
228  std::ostringstream os;
229 
230  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
231 
232  throw JFunctionalException(os.str());
233  }
234  }
235 
236  const_iterator p = this->lower_bound(x);
237 
238  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
239  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
240 
241  try {
242  return this->getExceptionHandler().action();
243  }
244  catch (const JException& error) {
245 
246  std::ostringstream os;
247 
248  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
249  << STREAM("?") << x << " <> "
250  << STREAM("?") << this->begin() ->getX() << ' '
251  << STREAM("?") << this->rbegin()->getX();
252 
253  throw JValueOutOfRange(os.str());
254  }
255  }
256 
257  const_iterator q = p--;
258 
259  const double dx = this->getDistance(p->getX(), q->getX());
260  const double t = this->getDistance(p->getX(), x) / dx;
261 
262  return h00(t)*p->getY() + h10(t)*p->getU()*dx + h01(t)*q->getY() + h11(t)*q->getU()*dx;
263  }
264 
265  protected:
266 
267  using collection_type::h00;
268  using collection_type::h10;
269  using collection_type::h01;
270  using collection_type::h11;
271  };
272 
273 
274  /**
275  * Template specialisation for spline interpolation method with returning JResultDerivative data structure.
276  */
277  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
278  class JHermiteSplineFunction<JElement_t,
279  JCollection_t,
280  JResultDerivative<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
281  JDistance_t> :
282  public JHermiteSplineCollection<JElement_t, JCollection_t, JDistance_t>,
283  public virtual JFunction<typename JElement_t::abscissa_type,
284  JResultDerivative<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
285  {
286  public:
287 
289 
290  typedef typename collection_type::abscissa_type abscissa_type;
291  typedef typename collection_type::ordinate_type ordinate_type;
292  typedef typename collection_type::value_type value_type;
293  typedef typename collection_type::distance_type distance_type;
294 
295  typedef typename collection_type::const_iterator const_iterator;
296  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
297  typedef typename collection_type::iterator iterator;
298  typedef typename collection_type::reverse_iterator reverse_iterator;
299 
302 
306 
307 
308  /**
309  * Default constructor.
310  */
312  {}
313 
314 
315  /**
316  * Recursive interpolation method implementation.
317  *
318  * \param pX pointer to abscissa values
319  * \return function value
320  */
321  virtual result_type evaluate(const argument_type* pX) const override
322  {
323  const argument_type x = *pX;
324 
325  if (this->size() <= 1u) {
326 
327  try {
328  return this->getExceptionHandler().action();
329  }
330  catch (const JException& error) {
331 
332  std::ostringstream os;
333 
334  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
335 
336  throw JFunctionalException(os.str());
337  }
338  }
339 
340  const_iterator p = this->lower_bound(x);
341 
342 
343  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
344  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
345 
346  try {
347  return this->getExceptionHandler().action();
348  }
349  catch (const JException& error) {
350 
351  std::ostringstream os;
352 
353  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
354  << STREAM("?") << x << " <> "
355  << STREAM("?") << this->begin() ->getX() << ' '
356  << STREAM("?") << this->rbegin()->getX();
357 
358  throw JValueOutOfRange(os.str());
359  }
360  }
361 
362  const_iterator q = p--;
363 
364  const double dx = this->getDistance(p->getX(), q->getX());
365  const double t = this->getDistance(p->getX(), x) / dx;
366 
367  result.f = h00 (t)*p->getY() + h10 (t)*p->getU()*dx + h01 (t)*q->getY() + h11 (t)*q->getU()*dx;
368  result.fp = h00p(t)*p->getY()/dx + h10p(t)*p->getU() + h01p(t)*q->getY()/dx + h11p(t)*q->getU();
369 
370  return result;
371  }
372 
373 
374  protected:
375 
376  using collection_type::h00;
377  using collection_type::h10;
378  using collection_type::h01;
379  using collection_type::h11;
380 
381  using collection_type::h00p;
382  using collection_type::h10p;
383  using collection_type::h01p;
384  using collection_type::h11p;
385 
386  private:
388  };
389 
390 
391  /**
392  * Template specialisation for spline interpolation method with returning JResultPDF data structure.
393  *
394  * Note that the data structure of the elements in the collection should have the additional methods:
395  * <pre>
396  * ordinate_type getIntegral() const;
397  * void setIntegral(ordinate_type v);
398  * </pre>
399  * to get and set the integral values, respectively.
400  */
401  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
402  class JHermiteSplineFunction<JElement_t,
403  JCollection_t,
404  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
405  JDistance_t> :
406  public JHermiteSplineCollection<JElement_t, JCollection_t, JDistance_t>,
407  public virtual JFunction<typename JElement_t::abscissa_type,
408  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
409  {
410  public:
411 
413 
414  typedef typename collection_type::abscissa_type abscissa_type;
415  typedef typename collection_type::ordinate_type ordinate_type;
416  typedef typename collection_type::value_type value_type;
417  typedef typename collection_type::distance_type distance_type;
418 
419  typedef typename collection_type::const_iterator const_iterator;
420  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
421  typedef typename collection_type::iterator iterator;
422  typedef typename collection_type::reverse_iterator reverse_iterator;
423 
426 
430 
431 
432  /**
433  * Default constructor.
434  */
436  {}
437 
438 
439  /**
440  * Recursive interpolation method implementation.
441  *
442  * \param pX pointer to abscissa values
443  * \return function value
444  */
445  virtual result_type evaluate(const argument_type* pX) const override
446  {
447  const argument_type x = *pX;
448 
449  if (this->size() <= 1u) {
450 
451  try {
452  return this->getExceptionHandler().action();
453  }
454  catch (const JException& error) {
455 
456  std::ostringstream os;
457 
458  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
459 
460  throw JFunctionalException(os.str());
461  }
462  }
463 
464  const_iterator p = this->lower_bound(x);
465 
466  if (p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) {
467 
468  try {
469 
470  result = this->getExceptionHandler().action();
471 
472  // overwrite integral values
473 
474  result.v = 0;
475  result.V = this->rbegin()->getIntegral();
476 
477  } catch(const JException& error) {
478 
479  std::ostringstream os;
480 
481  os << __FILE__ << ':' << __LINE__ << " abscissa out of range " << STREAM("?") << x << " < " << STREAM("?") << this->begin() ->getX();
482 
483  throw JValueOutOfRange(os.str());
484  }
485 
486  return result;
487 
488  } else if (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision) {
489 
490  try {
491 
492  result = this->getExceptionHandler().action();
493 
494  // overwrite integral values
495 
496  result.v = this->rbegin()->getIntegral();
497  result.V = this->rbegin()->getIntegral();
498 
499  } catch(const JException& error) {
500 
501  std::ostringstream os;
502 
503  os << __FILE__ << ':' << __LINE__ << " abscissa out of range " << STREAM("?") << x << " > " << STREAM("?") << this->rbegin()->getX();
504 
505  throw JValueOutOfRange(os.str());
506  }
507 
508  return result;
509  }
510 
511  const_iterator q = p--;
512 
513  const double dx = this->getDistance(p->getX(), q->getX());
514  const double t = this->getDistance(p->getX(), x) / dx;
515 
516  result.f = h00 (t)*p->getY() + h10 (t)*p->getU()*dx + h01 (t)*q->getY() + h11 (t)*q->getU()*dx;
517  result.fp = h00p(t)*p->getY()/dx + h10p(t)*p->getU() + h01p(t)*q->getY()/dx + h11p(t)*q->getU();
518  result.v = (p->getIntegral() +
519  (H00 (t)*p->getY() + H10 (t)*p->getU()*dx + H01 (t)*q->getY() + H11 (t)*q->getU()*dx)*dx);
520  result.V = this->rbegin()->getIntegral();
521 
522  return result;
523  }
524 
525 
526  protected:
527 
528  using collection_type::h00;
529  using collection_type::h10;
530  using collection_type::h01;
531  using collection_type::h11;
532 
533  using collection_type::h00p;
534  using collection_type::h10p;
535  using collection_type::h01p;
536  using collection_type::h11p;
537 
538  using collection_type::H00;
539  using collection_type::H10;
540  using collection_type::H01;
541  using collection_type::H11;
542 
543  /**
544  * Determination of derivatives.
545  */
546  virtual void do_compile() override
547  {
548  if (!this->empty()) {
549 
551 
552  this->begin()->setIntegral(JMATH::zero);
553 
554  for (iterator j = this->begin(), i = j++; j != this->end(); ++i, ++j) {
555 
556  const double dx = this->getDistance(i->getX(), j->getX());
557  const ordinate_type y = i->getY() + j->getY();
558  const ordinate_type z = i->getU() - j->getU();
559 
560  const ordinate_type v = dx * 0.50 * y;
561  const ordinate_type w = dx * 1.00 * z*dx/12.0;
562 
563  j->setIntegral(i->getIntegral() + v + w);
564  }
565  }
566  }
567 
568  private:
570  };
571 
572 
573  /**
574  * Template class for spline interpolation in 1D
575  *
576  * This class implements the JFunction1D interface.
577  */
578  template<class JElement_t,
579  template<class, class> class JCollection_t,
580  class JResult_t = typename JElement_t::ordinate_type,
583  public JHermiteSplineFunction<JElement_t, JCollection_t, JResult_t, JDistance_t>,
584  public virtual JFunction1D<typename JElement_t::abscissa_type, JResult_t>
585  {
586  public:
587 
589 
594 
599 
601 
605 
606 
607  /**
608  * Default contructor.
609  */
611  {}
612  };
613 
614 
615  /**
616  * Functional map with spline interpolation.
617  */
618  template<class JKey_t,
619  class JValue_t,
620  template<class, class, class> class JMap_t,
621  class JResult_t,
622  class JDistance_t = JDistance<JKey_t> >
624  public JMap_t<JKey_t, JValue_t, JDistance_t>,
625  public JFunction<JKey_t, JResult_t>
626  {
627  public:
628 
629  typedef JMap_t<JKey_t, JValue_t, JDistance_t> collection_type;
631 
632  typedef typename collection_type::abscissa_type abscissa_type;
633  typedef typename collection_type::ordinate_type ordinate_type;
634  typedef typename collection_type::value_type value_type;
635  typedef typename collection_type::distance_type distance_type;
636 
637  typedef typename collection_type::const_iterator const_iterator;
638  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
639  typedef typename collection_type::iterator iterator;
640  typedef typename collection_type::reverse_iterator reverse_iterator;
641 
645 
650 
651 
652  /**
653  * Default constructor.
654  */
656  {}
657 
658 
659  /**
660  * Recursive interpolation method implementation.
661  *
662  * \param pX pointer to abscissa values
663  * \return function value
664  */
665  virtual result_type evaluate(const argument_type* pX) const override
666  {
667  const argument_type x = *pX;
668 
669  ++pX; // next argument value
670 
671  const_iterator p = this->begin();
672 
673  for (typename JHermiteSplineFunction1D_t::iterator q = buffer.begin(); q != buffer.end(); ++q, ++p) {
674  q->getY() = JFunction<argument_type, data_type>::getValue(p->getY(), pX);
675  }
676 
677  buffer.compile();
678 
679  return buffer(x);
680  }
681 
682 
683  private:
684  /**
685  * Function compilation.
686  */
687  virtual void do_compile() override
688  {
689  buffer.clear();
690 
691  for (iterator i = this->begin(); i != this->end(); ++i) {
692  buffer.put(i->getX(), data_type());
693  }
694  }
695 
696 
698  };
699 
700 
701  /**
702  * Conversion of data points to integral values.
703  *
704  * The integration includes the use of 2nd derivatives of the data points of the input spline interpolating function.
705  *
706  * \param input collection
707  * \param output mappable collection
708  * \return integral
709  */
710  template<class JElement_t,
711  template<class, class> class JCollection_t,
712  class JResult_t,
713  class JDistance_t>
714  inline typename JElement_t::ordinate_type
716  typename JMappable<JElement_t>::map_type& output)
717  {
718  typedef typename JElement_t::ordinate_type ordinate_type;
720 
722 
723  if (input.getSize() > 1) {
724 
725  output.put(input.begin()->getX(), V);
726 
727  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
728 
729  const double dx = input.getDistance(i->getX(), j->getX());
730  const ordinate_type y = i->getY() + j->getY();
731  const ordinate_type z = i->getU() - j->getU();
732 
733  const ordinate_type v = dx * 0.50 * y;
734  const ordinate_type w = dx * 1.00 * z*dx/12.0;
735 
736  V += v + w;
737 
738  output.put(j->getX(), V);
739  }
740  }
741 
742  return V;
743  }
744 
745 
746  /**
747  * Conversion of data points to integral values.
748  *
749  * The integration directly uses the integral values of the input spline interpolating function.
750  *
751  * \param input collection
752  * \param output mappable collection
753  * \return integral
754  */
755  template<class JElement_t,
756  template<class, class> class JCollection_t,
757  class JDistance_t>
758  inline typename JElement_t::ordinate_type
759  integrate(const JHermiteSplineFunction1D<JElement_t, JCollection_t, JResultPDF<typename JElement_t::ordinate_type>, JDistance_t>& input,
760  typename JMappable<JElement_t>::map_type& output)
761  {
762  typedef typename JElement_t::ordinate_type ordinate_type;
765 
766  if (input.getSize() > 1) {
767 
768  for (const_iterator i = input.begin(); i != input.end(); ++i) {
769  output.put(i->getX(), i->getIntegral());
770  }
771 
772  return input.rbegin()->getIntegral();
773  }
774 
775  return JMATH::zero;
776  }
777 }
778 
779 #endif
The elements in a collection are sorted according to their abscissa values and a given distance opera...
Exceptions.
This include file containes various data structures that can be used as specific return types for the...
Definition of zero value for any class.
Exception for division by zero.
Definition: JException.hh:288
General exception.
Definition: JException.hh:24
Exception for a functional operation.
Definition: JException.hh:144
Exception for missing value.
Definition: JException.hh:216
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:180
Template definition of function object interface.
Definition: JFunctional.hh:72
JArgument_t argument_type
Definition: JFunctional.hh:84
const JExceptionHandler & getExceptionHandler() const
Get exception handler.
Definition: JFunctional.hh:292
JResult_t result_type
Definition: JFunctional.hh:87
static result_type getValue(const JFunctional &function, const argument_type *pX)
Recursive function value evaluation.
Definition: JFunctional.hh:107
Template base class spline interpolations.
collection_type::ordinate_type ordinate_type
static abscissa_type h01(abscissa_type t)
static abscissa_type h00p(abscissa_type t)
collection_type::const_reverse_iterator const_reverse_iterator
collection_type::reverse_iterator reverse_iterator
collection_type::const_iterator const_iterator
static abscissa_type h11p(abscissa_type t)
JHermiteSplineCollection()
Default constructor.
static abscissa_type h10p(abscissa_type t)
static abscissa_type h00(abscissa_type t)
static abscissa_type h01p(abscissa_type t)
static abscissa_type H00(abscissa_type t)
JCollection_t< JElement_t, JDistance_t > collection_type
collection_type::distance_type distance_type
static abscissa_type h11(abscissa_type t)
collection_type::iterator iterator
virtual void do_compile() override
Determination of derivatives.
static abscissa_type H01(abscissa_type t)
collection_type::value_type value_type
void compile(const bool monotone)
Determination of derivatives.
static abscissa_type H10(abscissa_type t)
static abscissa_type H11(abscissa_type t)
collection_type::abscissa_type abscissa_type
static abscissa_type h10(abscissa_type t)
Template class for spline interpolation in 1D.
collection_type::abscissa_type abscissa_type
JFunction1D< abscissa_type, JResult_t > function_type
collection_type::const_reverse_iterator const_reverse_iterator
collection_type::iterator iterator
JHermiteSplineCollection< JElement_t, JCollection_t, JDistance_t > collection_type
collection_type::ordinate_type ordinate_type
function_type::JExceptionHandler exceptionhandler_type
collection_type::reverse_iterator reverse_iterator
collection_type::value_type value_type
collection_type::distance_type distance_type
function_type::argument_type argument_type
function_type::result_type result_type
collection_type::const_iterator const_iterator
JHermiteSplineFunction1D()
Default contructor.
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
Template definition for functional collection with spline interpolation.
Functional map with spline interpolation.
collection_type::reverse_iterator reverse_iterator
JMap_t< JKey_t, JValue_t, JDistance_t > collection_type
JHermiteSplineMap()
Default constructor.
collection_type::abscissa_type abscissa_type
JHermiteSplineFunction1D_t buffer
collection_type::value_type value_type
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
collection_type::const_reverse_iterator const_reverse_iterator
JHermiteSplineFunction1D< JSplineElement2D< argument_type, data_type >, JMapCollection< JMap_t >::template collection_type, result_type > JHermiteSplineFunction1D_t
JResultType< ordinate_type >::result_type data_type
function_type::result_type result_type
virtual void do_compile() override
Function compilation.
function_type::JExceptionHandler exceptionhandler_type
function_type::argument_type argument_type
collection_type::iterator iterator
collection_type::ordinate_type ordinate_type
JFunction< JKey_t, JResult_t > function_type
collection_type::distance_type distance_type
collection_type::const_iterator const_iterator
static const JZero zero
Function object to assign zero value.
Definition: JZero.hh:105
double getDistance(const JFirst_t &first, const JSecond_t &second)
Get distance between objects.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
std::vector< event_type > data_type
Definition: JPerth.cc:82
Auxiliary classes and methods for multi-dimensional interpolations and histograms.
data_type w[N+1][M+1]
Definition: JPolint.hh:867
JElement_t::ordinate_type integrate(const JHermiteSplineFunction1D< JElement_t, JCollection_t, JResultPDF< typename JElement_t::ordinate_type >, JDistance_t > &input, typename JMappable< JElement_t >::map_type &output)
Conversion of data points to integral values.
double u[N+1]
Definition: JPolint.hh:865
data_type v[N+1][M+1]
Definition: JPolint.hh:866
virtual void do_compile() override
Function compilation.
Definition: JPolint.hh:860
int j
Definition: JPolint.hh:792
Definition: JSTDTypes.hh:14
Template definition of function object interface in one dimension.
Definition: JFunctional.hh:334
functional_type::result_type result_type
Definition: JFunctional.hh:339
functional_type::argument_type argument_type
Definition: JFunctional.hh:338
Template definition of function object interface in multidimensions.
Definition: JFunctional.hh:320
functional_type::argument_type argument_type
Definition: JFunctional.hh:322
functional_type::result_type result_type
Definition: JFunctional.hh:323
Exception handler for functional object.
Definition: JFunctional.hh:133
Template class to define the corresponding JCollection for a given template JMap.
Template interface definition for associative collection of elements.
void put(typename JClass< key_type > ::argument_type key, typename JClass< mapped_type >::argument_type value)
Put pair-wise element (key,value) into collection.
Data structure for result including value and first derivative of function.
Definition: JResult.hh:45
Data structure for result including value, first derivative and integrals of function.
Definition: JResult.hh:339
Auxiliary class to evaluate result type.
Definition: JFunctional.hh:392
Auxiliary data structure for handling std::ostream.