Jpp  18.3.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JPolint.hh
Go to the documentation of this file.
1 #ifndef __JTOOLS__JPOLINT__
2 #define __JTOOLS__JPOLINT__
3 
4 #include <cmath>
5 #include <iterator>
6 #include <algorithm>
7 #include <sstream>
8 
9 #include "JLang/JException.hh"
10 #include "JLang/JAssert.hh"
11 #include "JLang/JBool.hh"
13 #include "JTools/JFunctional.hh"
14 #include "JTools/JDistance.hh"
15 #include "JTools/JResult.hh"
16 #include "JTools/JElement.hh"
17 #include "JTools/JMapCollection.hh"
18 #include "JTools/JQuadrature.hh"
19 
20 
21 /**
22  * \author mdejong
23  */
24 
25 namespace JTOOLS {}
26 namespace JPP { using namespace JTOOLS; }
27 
28 namespace JTOOLS {
29 
30  using JLANG::JException;
34 
35 
36  /**
37  * Template definition for functional collection with polynomial interpolation.
38  */
39  template<unsigned int N,
40  class JElement_t,
41  template<class, class> class JCollection_t,
42  class JResult_t,
43  class JDistance_t>
45 
46 
47  /**
48  * Template specialisation for functional collection with polynomial interpolation.
49  *
50  * Polynomial interpolation code is taken from reference:
51  * Numerical Recipes in C++, W.H. Press, S.A. Teukolsky, W.T. Vetterling and B.P. Flannery,
52  * Cambridge University Press.
53  */
54  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
56  JElement_t,
57  JCollection_t,
58  typename JResultType<typename JElement_t::ordinate_type>::result_type,
59  JDistance_t> :
60  public JCollection_t<JElement_t, JDistance_t>,
61  public virtual JFunction<typename JElement_t::abscissa_type,
62  typename JResultType<typename JElement_t::ordinate_type>::result_type>
63  {
64  public:
65 
66  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
67 
68  typedef typename collection_type::abscissa_type abscissa_type;
69  typedef typename collection_type::ordinate_type ordinate_type;
70  typedef typename collection_type::value_type value_type;
71  typedef typename collection_type::distance_type distance_type;
72 
73  typedef typename collection_type::const_iterator const_iterator;
74  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
75  typedef typename collection_type::iterator iterator;
76  typedef typename collection_type::reverse_iterator reverse_iterator;
77 
80 
84 
85 
86  /**
87  * Default constructor.
88  */
90  {}
91 
92 
93  /**
94  * Recursive interpolation method implementation.
95  *
96  * \param pX pointer to abscissa values
97  * \return function value
98  */
99  virtual result_type evaluate(const argument_type* pX) const override
100  {
101  const argument_type x = *pX;
102 
103  if (this->size() > 1u) {
104 
105  const_iterator p = this->lower_bound(x);
106 
107  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
108  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
109 
110  try {
111  return this->getExceptionHandler().action();
112  }
113  catch (const JException& error) {
114 
115  std::ostringstream os;
116 
117  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
118  << STREAM("?") << x << " <> "
119  << STREAM("?") << this->begin() ->getX() << ' '
120  << STREAM("?") << this->rbegin()->getX();
121 
122  throw JValueOutOfRange(os.str());
123  }
124  }
125 
126  ++pX; // next argument value
127 
128 
129  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
130 
131  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
132  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
133 
134 
135  int j = 0;
136 
137  for (int i = 0; i != n; ++p, ++i) {
138 
139  u[i] = this->getDistance(x, p->getX());
140  v[i] = function_type::getValue(p->getY(), pX);
141  w[i] = v[i];
142 
143  if (fabs(u[i]) < fabs(u[j])) {
144  j = i;
145  }
146  }
147 
148 
149  result_type y = v[j];
150 
151  --j;
152 
153  for (int m = 1; m != n; ++m) {
154 
155  for (int i = 0; i != n-m; ++i) {
156 
157  const double ho = u[ i ];
158  const double hp = u[i+m];
159  const double dx = ho - hp;
160 
161  v[i] = v[i+1];
162  v[i] -= w[ i ];
163  w[i] = v[ i ];
164 
165  v[i] *= ho/dx;
166  w[i] *= hp/dx;
167  }
168 
169  if (2*(j+1) < n - m)
170  y += v[j+1];
171  else
172  y += w[j--];
173  }
174 
175  return y;
176 
177  } else if (this->size() == 1u && this->getDistance(x, this->begin()->getX()) <= distance_type::precision) {
178 
179  return function_type::getValue(this->begin()->getY(), ++pX);
180 
181  } else {
182 
183  try {
184  return this->getExceptionHandler().action();
185  }
186  catch (const JException& error) {
187 
188  std::ostringstream os;
189 
190  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
191 
192  throw JFunctionalException(os.str());
193  }
194  }
195  }
196 
197  protected:
198  /**
199  * Function compilation.
200  */
201  virtual void do_compile() override
202  {}
203 
204 
205  private:
206  mutable double u[N+1];
207  mutable result_type v[N+1];
208  mutable result_type w[N+1];
209  };
210 
211 
212  /**
213  * Template specialisation for zero-degree polynomial interpolation method.\n
214  * The interpolation is based on a simple lookup table.
215  */
216  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
217  class JPolintFunction<0,
218  JElement_t,
219  JCollection_t,
220  typename JResultType<typename JElement_t::ordinate_type>::result_type,
221  JDistance_t> :
222  public JCollection_t<JElement_t, JDistance_t>,
223  public virtual JFunction<typename JElement_t::abscissa_type,
224  typename JResultType<typename JElement_t::ordinate_type>::result_type>
225  {
226  public:
227 
228  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
229 
230  typedef typename collection_type::abscissa_type abscissa_type;
231  typedef typename collection_type::ordinate_type ordinate_type;
232  typedef typename collection_type::value_type value_type;
233  typedef typename collection_type::distance_type distance_type;
234 
235  typedef typename collection_type::const_iterator const_iterator;
236  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
237  typedef typename collection_type::iterator iterator;
238  typedef typename collection_type::reverse_iterator reverse_iterator;
239 
242 
246 
247 
248  /**
249  * Default constructor.
250  */
252  {}
253 
254 
255  /**
256  * Recursive interpolation method implementation.
257  *
258  * \param pX pointer to abscissa values
259  * \return function value
260  */
261  virtual result_type evaluate(const argument_type* pX) const override
262  {
263  const argument_type x = *pX;
264 
265  if (this->size() > 1u) {
266 
267  const_iterator p = this->lower_bound(x);
268 
269  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
270  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
271 
272  try {
273  return this->getExceptionHandler().action();
274  }
275  catch (const JException& error) {
276 
277  std::ostringstream os;
278 
279  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
280  << STREAM("?") << x << " <> "
281  << STREAM("?") << this->begin() ->getX() << ' '
282  << STREAM("?") << this->rbegin()->getX();
283 
284  throw JValueOutOfRange(os.str());
285  }
286  }
287 
288  ++pX; // next argument value
289 
290 
291  const_iterator q = p--;
292 
293  if (q == this->begin() || this->getDistance(x, q->getX()) < this->getDistance(p->getX(), x))
294  return function_type::getValue(q->getY(), pX);
295  else
296  return function_type::getValue(p->getY(), pX);
297 
298  } else if (this->size() == 1u && this->getDistance(x, this->begin()->getX()) <= distance_type::precision) {
299 
300  return function_type::getValue(this->begin()->getY(), ++pX);
301 
302  } else {
303 
304  try {
305  return this->getExceptionHandler().action();
306  }
307  catch (const JException& error) {
308 
309  std::ostringstream os;
310 
311  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
312 
313  throw JFunctionalException(os.str());
314  }
315  }
316  }
317 
318  protected:
319  /**
320  * Function compilation.
321  */
322  virtual void do_compile() override
323  {}
324  };
325 
326 
327  /**
328  * Template specialisation for first-degree polynomial interpolation method.\n
329  * The interpolation is based on a simple linear interpolation.
330  */
331  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
332  class JPolintFunction<1,
333  JElement_t,
334  JCollection_t,
335  typename JResultType<typename JElement_t::ordinate_type>::result_type,
336  JDistance_t> :
337  public JCollection_t<JElement_t, JDistance_t>,
338  public virtual JFunction<typename JElement_t::abscissa_type,
339  typename JResultType<typename JElement_t::ordinate_type>::result_type>
340  {
341  public:
342 
343  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
344 
345  typedef typename collection_type::abscissa_type abscissa_type;
346  typedef typename collection_type::ordinate_type ordinate_type;
347  typedef typename collection_type::value_type value_type;
348  typedef typename collection_type::distance_type distance_type;
349 
350  typedef typename collection_type::const_iterator const_iterator;
351  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
352  typedef typename collection_type::iterator iterator;
353  typedef typename collection_type::reverse_iterator reverse_iterator;
354 
357 
361 
362 
363  /**
364  * Default constructor.
365  */
367  {}
368 
369 
370  /**
371  * Recursive interpolation method implementation.
372  *
373  * \param pX pointer to abscissa values
374  * \return function value
375  */
376  virtual result_type evaluate(const argument_type* pX) const override
377  {
378  const argument_type x = *pX;
379 
380  if (this->size() > 1u) {
381 
382  const_iterator p = this->lower_bound(x);
383 
384  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
385  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
386 
387  try {
388  return this->getExceptionHandler().action();
389  }
390  catch (const JException& error) {
391 
392  std::ostringstream os;
393 
394  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
395  << STREAM("?") << x << " <> "
396  << STREAM("?") << this->begin() ->getX() << ' '
397  << STREAM("?") << this->rbegin()->getX();
398 
399  throw JValueOutOfRange(os.str());
400  }
401  }
402 
403  ++pX; // next argument value
404 
405 
406  const_iterator q = p--;
407 
408  const double dx = this->getDistance(p->getX(), q->getX());
409  const double a = this->getDistance(x, q->getX()) / dx;
410  const double b = 1.0 - a;
411 
412  ya = function_type::getValue(p->getY(), pX);
413  yb = function_type::getValue(q->getY(), pX);
414 
415  ya *= a;
416  yb *= b;
417 
418  ya += yb;
419 
420  return ya;
421 
422  } else if (this->size() == 1u && this->getDistance(x, this->begin()->getX()) <= distance_type::precision) {
423 
424  return function_type::getValue(this->begin()->getY(), ++pX);
425 
426  } else {
427 
428  try {
429  return this->getExceptionHandler().action();
430  }
431  catch (const JException& error) {
432 
433  std::ostringstream os;
434 
435  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
436 
437  throw JFunctionalException(os.str());
438  }
439  }
440  }
441 
442  protected:
443  /**
444  * Function compilation.
445  */
446  virtual void do_compile() override
447  {}
448 
449 
450  private:
451  mutable result_type ya;
452  mutable result_type yb;
453  };
454 
455 
456  /**
457  * Template specialisation for polynomial interpolation method with returning JResultPDF data structure.
458  *
459  * Note that the data structure of the elements in the collection should have the additional methods:
460  * <pre>
461  * ordinate_type getIntegral() const;
462  * void setIntegral(ordinate_type v);
463  * </pre>
464  * to get and set the integral values, respectively.
465  */
466  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
468  JElement_t,
469  JCollection_t,
470  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
471  JDistance_t> :
472  public JCollection_t<JElement_t, JDistance_t>,
473  public virtual JFunction<typename JElement_t::abscissa_type,
474  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
475  {
476  public:
477 
478  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
479 
480  typedef typename collection_type::abscissa_type abscissa_type;
481  typedef typename collection_type::ordinate_type ordinate_type;
482  typedef typename collection_type::value_type value_type;
483  typedef typename collection_type::distance_type distance_type;
484 
485  typedef typename collection_type::const_iterator const_iterator;
486  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
487  typedef typename collection_type::iterator iterator;
488  typedef typename collection_type::reverse_iterator reverse_iterator;
489 
492 
496 
498 
499 
500  /**
501  * Default contructor.
502  */
504  {}
505 
506 
507  /**
508  * Recursive interpolation method implementation.
509  *
510  * \param pX pointer to abscissa values
511  * \return function value
512  */
513  virtual result_type evaluate(const argument_type* pX) const override
514  {
515  const argument_type x = *pX;
516 
517  if (this->size() <= 1u) {
518 
519  try {
520  return this->getExceptionHandler().action();
521  }
522  catch (const JException& error) {
523 
524  std::ostringstream os;
525 
526  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
527 
528  throw JFunctionalException(os.str());
529  }
530  }
531 
532  const_iterator p = this->lower_bound(x);
533 
534  if (p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) {
535 
536  try {
537 
538  result = this->getExceptionHandler().action();
539 
540  // overwrite integral values
541 
542  result.v = 0;
543  result.V = this->rbegin()->getIntegral();
544 
545  } catch(const JException& error) {
546 
547  std::ostringstream os;
548 
549  os << __FILE__ << ':' << __LINE__ << " abscissa out of range " << STREAM("?") << x << " < " << STREAM("?") << this->begin() ->getX();
550 
551  throw JValueOutOfRange(os.str());
552  }
553 
554  return result;
555 
556  } else if (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision) {
557 
558  try {
559 
560  result = this->getExceptionHandler().action();
561 
562  // overwrite integral values
563 
564  result.v = this->rbegin()->getIntegral();
565  result.V = this->rbegin()->getIntegral();
566 
567  } catch(const JException& error) {
568 
569  std::ostringstream os;
570 
571  os << __FILE__ << ':' << __LINE__ << " abscissa out of range " << STREAM("?") << x << " > " << STREAM("?") << this->rbegin()->getX();
572 
573  throw JValueOutOfRange(os.str());
574  }
575 
576  return result;
577  }
578 
579  ++pX; // next argument value
580 
581  {
582  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
583 
584  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
585  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
586 
587 
588  int j = 0;
589 
590  for (int i = 0; i != n; ++p, ++i) {
591 
592  u[i] = this->getDistance(x, p->getX());
593 
594  w[i][0] = v[i][0] = JFunctional<argument_type, data_type>::getValue(p->getY(), pX);
595  w[i][1] = v[i][1] = JMATH::zero;
596  w[i][2] = v[i][2] = p->getIntegral();
597 
598  if (fabs(u[i]) < fabs(u[j])) {
599  j = i;
600  }
601  }
602 
603 
604  result.f = v[j][0];
605  result.fp = v[j][1];
606  result.v = v[j][2];
607  result.V = this->rbegin()->getIntegral();
608 
609  --j;
610 
611  for (int m = 1; m != n; ++m) {
612 
613  for (int i = 0; i != n-m; ++i) {
614 
615  const double ho = u[ i ];
616  const double hp = u[i+m];
617  const double dx = ho - hp;
618 
619  for (int k = 0; k != 3; ++k) {
620  r[k] = (v[i+1][k] - w[i][k]) / dx;
621  }
622 
623  v[i][0] = ho * r[0];
624  w[i][0] = hp * r[0];
625  v[i][1] = ho * r[1] - r[0];
626  w[i][1] = hp * r[1] - r[0];
627  v[i][2] = ho * r[2];
628  w[i][2] = hp * r[2];
629  }
630 
631  if (2*(j+1) < n - m) {
632 
633  result.f += v[j+1][0];
634  result.fp += v[j+1][1];
635  result.v += v[j+1][2];
636 
637  } else {
638 
639  result.f += w[j][0];
640  result.fp += w[j][1];
641  result.v += w[j][2];
642 
643  --j;
644  }
645  }
646  }
647 
648  return result;
649  }
650 
651  protected:
652  /**
653  * Function compilation.
654  */
655  virtual void do_compile() override
656  {
658 
659  if (this->getSize() > 1) {
660 
661  const JGaussLegendre engine(N);
662 
663  this->begin()->setIntegral(V);
664 
665  for (iterator j = this->begin(), i = j++; j != this->end(); ++i, ++j) {
666 
667  const abscissa_type xmin = i->getX();
668  const abscissa_type xmax = j->getX();
669 
670  for (JGaussLegendre::const_iterator m = engine.begin(); m != engine.end(); ++m) {
671 
672  const abscissa_type x = 0.5 * (xmax + xmin + m->getX() * (xmax - xmin));
673  const ordinate_type v = 0.5 * (xmax - xmin) * m->getY() * get_value(this->evaluate(&x));
674 
675  V += v;
676  }
677 
678  j->setIntegral(V);
679  }
680  }
681  }
682 
683 
684  mutable double u[N+1];
685  mutable data_type v[N+1][3];
686  mutable data_type w[N+1][3];
687  mutable data_type r[3];
688 
690  };
691 
692 
693  /**
694  * Template definition of base class for polynomial interpolations with polynomial result.
695  */
696  template<unsigned int N,
697  class JElement_t,
698  template<class, class> class JCollection_t,
699  class JResult_t,
700  class JDistance_t>
702 
703 
704  /**
705  * Template base class for polynomial interpolations with polynomial result.
706  *
707  * This class partially implements the JFunctional interface.
708  */
709  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t, unsigned int M>
710  class JPolintCollection<N,
711  JElement_t,
712  JCollection_t,
713  JResultPolynome<M, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
714  JDistance_t> :
715  public JCollection_t<JElement_t, JDistance_t>,
716  public virtual JFunctional<>,
717  private JLANG::JAssert<M <= N>
718  {
719  public:
720 
721  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
722 
723  typedef typename collection_type::abscissa_type abscissa_type;
724  typedef typename collection_type::ordinate_type ordinate_type;
725  typedef typename collection_type::value_type value_type;
726  typedef typename collection_type::distance_type distance_type;
727 
728  typedef typename collection_type::const_iterator const_iterator;
729  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
730  typedef typename collection_type::iterator iterator;
731  typedef typename collection_type::reverse_iterator reverse_iterator;
732 
733  typedef typename JResultType<ordinate_type>::result_type data_type;
734  typedef JFunction<abscissa_type, JResultPolynome<M, data_type> > function_type;
735 
736  typedef typename function_type::argument_type argument_type;
737  typedef typename function_type::result_type result_type;
738 
739  using JFunctional<>::compile;
740 
741 
742  /**
743  * Default contructor.
744  */
745  JPolintCollection()
746  {}
747 
748 
749  /**
750  * Recursive interpolation method implementation.
751  *
752  * \param pX pointer to abscissa values
753  * \return function value
754  */
755  result_type evaluate(const argument_type* pX) const
756  {
757  const argument_type x = *pX;
758 
759  if (this->size() <= N) {
760 
761  std::ostringstream os;
762 
763  os << __FILE__ << ':' << __LINE__ << " not enough data " << STREAM("?") << x;
764 
765  throw JFunctionalException(os.str());
766  }
767 
768  const_iterator p = this->lower_bound(x);
769 
770  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
771  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
772 
773  std::ostringstream os;
774 
775  os << __FILE__ << ':' << __LINE__ << " abscissa out of range "
776  << STREAM("?") << x << " <> "
777  << STREAM("?") << this->begin() ->getX() << ' '
778  << STREAM("?") << this->rbegin()->getX();
779 
780  throw JValueOutOfRange(os.str());
781  }
782 
783  ++pX; // next argument value
784 
785 
786  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
787 
788  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
789  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
790 
791 
792  int j = 0;
793 
794  for (int i = 0; i != n; ++p, ++i) {
795 
796  u[i] = this->getDistance(x, p->getX());
797 
798  w[i][0] = v[i][0] = JFunctional<argument_type, data_type>::getValue(p->getY(), pX);
799 
800  for (unsigned int k = 1; k != M+1; ++k) {
801  w[i][k] = v[i][k] = JMATH::zero;
802  }
803 
804  if (fabs(u[i]) < fabs(u[j])) {
805  j = i;
806  }
807  }
808 
809 
810  for (unsigned int k = 0; k != M+1; ++k) {
811  result.y[k] = v[j][k];
812  }
813 
814  --j;
815 
816  for (int m = 1; m != n; ++m) {
817 
818  for (int i = 0; i != n-m; ++i) {
819 
820  const double ho = u[ i ];
821  const double hp = u[i+m];
822  const double dx = ho - hp;
823 
824  for (int k = 0; k != M+1; ++k) {
825  r[k] = (v[i+1][k] - w[i][k]) / dx;
826  }
827 
828  v[i][0] = ho * r[0];
829  w[i][0] = hp * r[0];
830 
831  for (int k = 1; k != M+1; ++k) {
832  v[i][k] = ho * r[k] - k * r[k-1];
833  w[i][k] = hp * r[k] - k * r[k-1];
834  }
835  }
836 
837  if (2*(j+1) < n - m) {
838 
839  for (int k = 0; k != M+1; ++k) {
840  result.y[k] += v[j+1][k];
841  }
842 
843  } else {
844 
845  for (int k = 0; k != M+1; ++k) {
846  result.y[k] += w[j][k];
847  }
848 
849  --j;
850  }
851  }
852 
853  return result;
854  }
855 
856  protected:
857  /**
858  * Function compilation.
859  */
860  virtual void do_compile() override
861  {}
862 
863 
864  private:
865  mutable double u[N+1];
866  mutable data_type v[N+1][M+1];
867  mutable data_type w[N+1][M+1];
868  mutable data_type r[M+1];
869 
870  mutable result_type result;
871  };
872 
873 
874  /**
875  * Template specialisation for polynomial interpolation method with returning JResultPolynome data structure.
876  */
877  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t, unsigned int M>
878  class JPolintFunction<N,
879  JElement_t,
880  JCollection_t,
881  JResultPolynome<M, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
882  JDistance_t> :
883  public JPolintCollection<N,
884  JElement_t,
885  JCollection_t,
886  JResultPolynome<M, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
887  JDistance_t>,
888  public virtual JFunction<typename JElement_t::abscissa_type,
889  JResultPolynome<N, typename JResultType<typename JElement_t::ordinate_type>::result_type> >
890  {
891  public:
892 
893  typedef JPolintCollection<N,
894  JElement_t,
895  JCollection_t,
896  JResultPolynome<M, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
897  JDistance_t> collection_type;
898 
899  typedef typename collection_type::abscissa_type abscissa_type;
900  typedef typename collection_type::ordinate_type ordinate_type;
901  typedef typename collection_type::value_type value_type;
902  typedef typename collection_type::distance_type distance_type;
903 
904  typedef typename collection_type::const_iterator const_iterator;
905  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
906  typedef typename collection_type::iterator iterator;
907  typedef typename collection_type::reverse_iterator reverse_iterator;
908 
909  typedef typename JResultType<ordinate_type>::result_type data_type;
910  typedef JFunction<abscissa_type, JResultPolynome<M, data_type> > function_type;
911 
912  typedef typename function_type::argument_type argument_type;
913  typedef typename function_type::result_type result_type;
914  typedef typename function_type::JExceptionHandler exceptionhandler_type;
915 
916 
917  /**
918  * Default contructor.
919  */
920  JPolintFunction()
921  {}
922 
923 
924  /**
925  * Recursive interpolation method implementation.
926  *
927  * \param pX pointer to abscissa values
928  * \return function value
929  */
930  virtual result_type evaluate(const argument_type* pX) const override
931  {
932  try {
933  return collection_type::evaluate(pX);
934  }
935  catch(const JException& error) {
936  return this->getExceptionHandler().action(error);
937  }
938  }
939  };
940 
941 
942  /**
943  * Template specialisation for polynomial interpolation method with returning JResultDerivative data structure.
944  */
945  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
946  class JPolintFunction<N,
947  JElement_t,
948  JCollection_t,
949  JResultDerivative<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
950  JDistance_t> :
951  public JPolintCollection<N,
952  JElement_t,
953  JCollection_t,
954  JResultPolynome<1, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
955  JDistance_t>,
956  public virtual JFunction<typename JElement_t::abscissa_type,
957  JResultDerivative<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
958  {
959  public:
960 
961  typedef JPolintCollection<N,
962  JElement_t,
963  JCollection_t,
964  JResultPolynome<1, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
965  JDistance_t> collection_type;
966 
967  typedef typename collection_type::abscissa_type abscissa_type;
968  typedef typename collection_type::ordinate_type ordinate_type;
969  typedef typename collection_type::value_type value_type;
970  typedef typename collection_type::distance_type distance_type;
971 
972  typedef typename collection_type::const_iterator const_iterator;
973  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
974  typedef typename collection_type::iterator iterator;
975  typedef typename collection_type::reverse_iterator reverse_iterator;
976 
977  typedef typename JResultType<ordinate_type>::result_type data_type;
978  typedef JFunction<abscissa_type, JResultDerivative<data_type> > function_type;
979 
980  typedef typename function_type::argument_type argument_type;
981  typedef typename function_type::result_type result_type;
982  typedef typename function_type::JExceptionHandler exceptionhandler_type;
983 
984  using JFunctional<>::compile;
985 
986 
987  /**
988  * Default contructor.
989  */
990  JPolintFunction()
991  {}
992 
993 
994  /**
995  * Recursive interpolation method implementation.
996  *
997  * \param pX pointer to abscissa values
998  * \return function value
999  */
1000  virtual result_type evaluate(const argument_type* pX) const override
1001  {
1002  try {
1003  return collection_type::evaluate(pX);
1004  }
1005  catch(const JException& error) {
1006  return this->getExceptionHandler().action(error);
1007  }
1008  }
1009  };
1010 
1011 
1012  /**
1013  * Template specialisation for polynomial interpolation method with returning JResultHesse data structure.
1014  */
1015  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
1016  class JPolintFunction<N,
1017  JElement_t,
1018  JCollection_t,
1019  JResultHesse<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
1020  JDistance_t> :
1021  public JPolintCollection<N,
1022  JElement_t,
1023  JCollection_t,
1024  JResultPolynome<2, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
1025  JDistance_t>,
1026  public virtual JFunction<typename JElement_t::abscissa_type,
1027  JResultHesse<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
1028  {
1029  public:
1030 
1031  typedef JPolintCollection<N,
1032  JElement_t,
1033  JCollection_t,
1034  JResultPolynome<2, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
1035  JDistance_t> collection_type;
1036 
1037  typedef typename collection_type::abscissa_type abscissa_type;
1038  typedef typename collection_type::ordinate_type ordinate_type;
1039  typedef typename collection_type::value_type value_type;
1040  typedef typename collection_type::distance_type distance_type;
1041 
1042  typedef typename collection_type::const_iterator const_iterator;
1043  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
1044  typedef typename collection_type::iterator iterator;
1045  typedef typename collection_type::reverse_iterator reverse_iterator;
1046 
1047  typedef typename JResultType<ordinate_type>::result_type data_type;
1048  typedef JFunction<abscissa_type, JResultHesse<data_type> > function_type;
1049 
1050  typedef typename function_type::argument_type argument_type;
1051  typedef typename function_type::result_type result_type;
1052  typedef typename function_type::JExceptionHandler exceptionhandler_type;
1053 
1054  using JFunctional<>::compile;
1055 
1056 
1057  /**
1058  * Default contructor.
1059  */
1060  JPolintFunction()
1061  {}
1062 
1063 
1064  /**
1065  * Recursive interpolation method implementation.
1066  *
1067  * \param pX pointer to abscissa values
1068  * \return function value
1069  */
1070  virtual result_type evaluate(const argument_type* pX) const override
1071  {
1072  try {
1073  return collection_type::evaluate(pX);
1074  }
1075  catch(const JException& error) {
1076  return this->getExceptionHandler().action(error);
1077  }
1078  }
1079  };
1080 
1081 
1082  /**
1083  * Template class for polynomial interpolation in 1D
1084  *
1085  * This class implements the JFunction1D interface.
1086  */
1087  template<unsigned int N,
1088  class JElement_t,
1089  template<class, class> class JCollection_t,
1090  class JResult_t = typename JElement_t::ordinate_type,
1091  class JDistance_t = JDistance<typename JElement_t::abscissa_type> >
1092  class JPolintFunction1D :
1093  public JPolintFunction<N, JElement_t, JCollection_t, JResult_t, JDistance_t>,
1094  public virtual JFunction1D<typename JElement_t::abscissa_type, JResult_t>
1095  {
1096  public:
1097 
1098  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
1099 
1100  typedef typename collection_type::abscissa_type abscissa_type;
1101  typedef typename collection_type::ordinate_type ordinate_type;
1102  typedef typename collection_type::value_type value_type;
1103  typedef typename collection_type::distance_type distance_type;
1104 
1105  typedef typename collection_type::const_iterator const_iterator;
1106  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
1107  typedef typename collection_type::iterator iterator;
1108  typedef typename collection_type::reverse_iterator reverse_iterator;
1109 
1110  typedef JFunction1D<abscissa_type, JResult_t> function_type;
1111 
1112  typedef typename function_type::argument_type argument_type;
1113  typedef typename function_type::result_type result_type;
1114  typedef typename function_type::JExceptionHandler exceptionhandler_type;
1115 
1116 
1117  /**
1118  * Default contructor.
1119  */
1120  JPolintFunction1D()
1121  {}
1122  };
1123 
1124 
1125  /**
1126  * \cond NEVER
1127  * Forward declarations.
1128  * \endcond
1129  */
1130  template<class JAbscissa_t, class JOrdinate_t>
1131  struct JElement2D;
1132 
1133  template<template<class, class, class> class JMap_t>
1134  struct JMapCollection;
1135 
1136 
1137  /**
1138  * Functional map with polynomial interpolation.
1139  */
1140  template<unsigned int N,
1141  class JKey_t,
1142  class JValue_t,
1143  template<class, class, class> class JMap_t,
1144  class JResult_t,
1145  class JDistance_t = JDistance<JKey_t> >
1146  class JPolintMap :
1147  public JPolintFunction<N,
1148  JElement2D<JKey_t, JValue_t>,
1149  JMapCollection<JMap_t>::template collection_type,
1150  JResult_t,
1151  JDistance_t>
1152  {
1153  public:
1154 
1155  typedef JElement2D<JKey_t, JValue_t> element_type;
1156  typedef JPolintFunction<N,
1157  element_type,
1158  JMapCollection<JMap_t>::template collection_type,
1159  JResult_t,
1160  JDistance_t> JPolintFunction_t;
1161 
1162  typedef typename JPolintFunction_t::abscissa_type abscissa_type;
1163  typedef typename JPolintFunction_t::ordinate_type ordinate_type;
1164  typedef typename JPolintFunction_t::value_type value_type;
1165  typedef typename JPolintFunction_t::distance_type distance_type;
1166 
1167  typedef typename JPolintFunction_t::const_iterator const_iterator;
1168  typedef typename JPolintFunction_t::const_reverse_iterator const_reverse_iterator;
1169  typedef typename JPolintFunction_t::iterator iterator;
1170  typedef typename JPolintFunction_t::reverse_iterator reverse_iterator;
1171 
1172  typedef typename JPolintFunction_t::argument_type argument_type;
1173  typedef typename JPolintFunction_t::result_type result_type;
1174  typedef typename JPolintFunction_t::JExceptionHandler exceptionhandler_type;
1175 
1176 
1177  /**
1178  * Default constructor.
1179  */
1180  JPolintMap()
1181  {}
1182  };
1183 
1184 
1185  /**
1186  * Conversion of data points to integral values.
1187  *
1188  * This method transfers the integration to the corresponding specialised function.
1189  *
1190  * \param input collection
1191  * \param output mappable collection
1192  * \return integral
1193  */
1194  template<unsigned int N,
1195  class JElement_t,
1196  template<class, class> class JCollection_t,
1197  class JResult_t,
1198  class JDistance_t>
1199  inline typename JElement_t::ordinate_type
1200  integrate(const JPolintFunction1D<N, JElement_t, JCollection_t, JResult_t, JDistance_t>& input,
1201  typename JMappable<JElement_t>::map_type& output)
1202  {
1203  return integrate(input, output, JLANG::JBool<N == 0 || N == 1>());
1204  }
1205 
1206 
1207  /**
1208  * Conversion of data points to integral values.
1209  *
1210  * The integration uses the Gauss-Legendre quadratures with the number of points set
1211  * to the degree of the input polynomial interpolating function.
1212  *
1213  * \param input collection
1214  * \param output mappable collection
1215  * \param option false
1216  * \return integral
1217  */
1218  template<unsigned int N,
1219  class JElement_t,
1220  template<class, class> class JCollection_t,
1221  class JResult_t,
1222  class JDistance_t>
1223  inline typename JElement_t::ordinate_type
1224  integrate(const JPolintFunction1D<N, JElement_t, JCollection_t, JResult_t, JDistance_t>& input,
1225  typename JMappable<JElement_t>::map_type& output,
1226  const JLANG::JBool<false>& option)
1227  {
1228  typedef typename JElement_t::abscissa_type abscissa_type;
1229  typedef typename JElement_t::ordinate_type ordinate_type;
1230  typedef typename JPolintFunction1D<N, JElement_t, JCollection_t, JResult_t, JDistance_t>::const_iterator const_iterator;
1231 
1232  ordinate_type V(JMATH::zero);
1233 
1234  if (input.getSize() > 1) {
1235 
1236  output.put(input.begin()->getX(), V);
1237 
1238  const JGaussLegendre engine(N);
1239 
1240  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
1241 
1242  const abscissa_type xmin = i->getX();
1243  const abscissa_type xmax = j->getX();
1244 
1245  for (JGaussLegendre::const_iterator m = engine.begin(); m != engine.end(); ++m) {
1246 
1247  const abscissa_type x = 0.5 * (xmax + xmin + m->getX() * (xmax - xmin));
1248  const ordinate_type v = 0.5 * (xmax - xmin) * m->getY() * get_value(input(x));
1249 
1250  V += v;
1251  }
1252 
1253  output.put(xmax, V);
1254  }
1255  }
1256 
1257  return V;
1258  }
1259 
1260 
1261  /**
1262  * Conversion of data points to integral values.
1263  *
1264  * The integration is based on the sum of ordinates of the input data points.
1265  *
1266  * \param input collection
1267  * \param output mappable collection
1268  * \param option true
1269  * \return integral
1270  */
1271  template<class JElement_t,
1272  template<class, class> class JCollection_t,
1273  class JResult_t,
1274  class JDistance_t>
1275  inline typename JElement_t::ordinate_type
1276  integrate(const JPolintFunction1D<0, JElement_t, JCollection_t, JResult_t, JDistance_t>& input,
1277  typename JMappable<JElement_t>::map_type& output,
1278  const JLANG::JBool<true>& option)
1279  {
1280  typedef typename JElement_t::ordinate_type ordinate_type;
1281  typedef typename JPolintFunction1D<0, JElement_t, JCollection_t, JResult_t, JDistance_t>::const_iterator const_iterator;
1282 
1283  ordinate_type V(JMATH::zero);
1284 
1285  if (input.getSize() > 1) {
1286 
1287  output.put(input.begin()->getX(), V);
1288 
1289  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
1290 
1291  V += input.getDistance(i->getX(), j->getX()) * j->getY();
1292 
1293  output.put(j->getX(), V);
1294  }
1295  }
1296 
1297  return V;
1298  }
1299 
1300 
1301  /**
1302  * Conversion of data points to integral values.
1303  *
1304  * The integration is based on the trapezoidal rule applied to the input data points.
1305  *
1306  * \param input collection
1307  * \param output mappable collection
1308  * \param option true
1309  * \return integral
1310  */
1311  template<class JElement_t,
1312  template<class, class> class JCollection_t,
1313  class JResult_t,
1314  class JDistance_t>
1315  inline typename JElement_t::ordinate_type
1316  integrate(const JPolintFunction1D<1, JElement_t, JCollection_t, JResult_t, JDistance_t>& input,
1317  typename JMappable<JElement_t>::map_type& output,
1318  const JLANG::JBool<true>& option)
1319  {
1320  typedef typename JElement_t::ordinate_type ordinate_type;
1321  typedef typename JPolintFunction1D<1, JElement_t, JCollection_t, JResult_t, JDistance_t>::const_iterator const_iterator;
1322 
1323  ordinate_type V(JMATH::zero);
1324 
1325  if (input.getSize() > 1) {
1326 
1327  output.put(input.begin()->getX(), V);
1328 
1329  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
1330 
1331  V += 0.5 * input.getDistance(i->getX(), j->getX()) * (i->getY() + j->getY());
1332 
1333  output.put(j->getX(), V);
1334  }
1335  }
1336 
1337  return V;
1338  }
1339 }
1340 
1341 #endif
const double xmax
Definition: JQuadrature.cc:24
General exception.
Definition: JException.hh:24
data_type w[N+1][M+1]
Definition: JPolint.hh:867
Exceptions.
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
Template definition for functional collection with polynomial interpolation.
Definition: JPolint.hh:44
Numerical integrator for .
Definition: JQuadrature.hh:111
Exception for a functional operation.
Definition: JException.hh:142
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
Definition: JPolint.hh:376
The elements in a collection are sorted according to their abscissa values and a given distance opera...
Data structure for result including value and N derivatives of function.
Definition: JResult.hh:533
This include file containes various data structures that can be used as specific return types for the...
static const JZero zero
Function object to assign zero value.
Definition: JZero.hh:105
data_type r[M+1]
Definition: JPolint.hh:868
Generation of compiler error.
Definition: JAssert.hh:17
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
Definition: JPolint.hh:261
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
size_t getSize(T(&array)[N])
Get size of c-array.
Definition: JLangToolkit.hh:32
double getDistance(const JFirst_t &first, const JSecond_t &second)
Get distance between objects.
const int n
Definition: JPolint.hh:786
static result_type getValue(const JFunctional &function, const argument_type *pX)
Recursive function value evaluation.
Definition: JFunctional.hh:107
Template definition of base class for polynomial interpolations with polynomial result.
Definition: JPolint.hh:701
then JCalibrateToT a
Definition: JTuneHV.sh:113
Data structure for result including value, first derivative and integrals of function.
Definition: JResult.hh:337
Auxiliary data structure for handling std::ostream.
Template definition of function object interface.
Definition: JFunctional.hh:32
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
Definition: JPolint.hh:513
Auxiliary class to evaluate result type.
Definition: JFunctional.hh:392
Exception handler for functional object.
Definition: JFunctional.hh:131
Template definition of function object interface in multidimensions.
Definition: JFunctional.hh:318
const double xmin
Definition: JQuadrature.cc:23
Exception for numerical precision error.
Definition: JException.hh:268
functional_type::result_type result_type
Definition: JFunctional.hh:323
Auxiliary classes for numerical integration.
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
virtual result_type evaluate(const argument_type *pX) const override
Recursive interpolation method implementation.
Definition: JPolint.hh:99
container_type::const_iterator const_iterator
Definition: JCollection.hh:91
int j
Definition: JPolint.hh:792
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:178
functional_type::argument_type argument_type
Definition: JFunctional.hh:322
data_type v[N+1][M+1]
Definition: JPolint.hh:866
double u[N+1]
Definition: JPolint.hh:865
JResultEvaluator< JResult_t >::result_type get_value(const JResult_t &value)
Helper method to recursively evaluate a to function value.
Definition: JResult.hh:998