Jpp
 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 
8 #include "JLang/JException.hh"
9 #include "JLang/JBool.hh"
10 #include "JTools/JFunctional.hh"
11 #include "JTools/JDistance.hh"
12 #include "JTools/JResult.hh"
13 #include "JTools/JElement.hh"
14 #include "JTools/JMapCollection.hh"
15 #include "JTools/JQuadrature.hh"
16 
17 
18 /**
19  * \author mdejong
20  */
21 
22 namespace JTOOLS {}
23 namespace JPP { using namespace JTOOLS; }
24 
25 namespace JTOOLS {
26 
30 
31 
32  /**
33  * Template definition for functional collection with polynomial interpolation.
34  */
35  template<unsigned int N,
36  class JElement_t,
37  template<class, class> class JCollection_t,
38  class JResult_t,
39  class JDistance_t>
41 
42 
43  /**
44  * Template specialisation for functional collection with polynomial interpolation.
45  *
46  * Polynomial interpolation code is taken from reference:
47  * Numerical Recipes in C++, W.H. Press, S.A. Teukolsky, W.T. Vetterling and B.P. Flannery,
48  * Cambridge University Press.
49  */
50  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
51  class JPolintFunction<N,
52  JElement_t,
53  JCollection_t,
54  typename JResultType<typename JElement_t::ordinate_type>::result_type,
55  JDistance_t> :
56  public JCollection_t<JElement_t, JDistance_t>,
57  public JFunction<typename JElement_t::abscissa_type,
58  typename JResultType<typename JElement_t::ordinate_type>::result_type>
59  {
60  public:
61 
62  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
63 
64  typedef typename collection_type::abscissa_type abscissa_type;
65  typedef typename collection_type::ordinate_type ordinate_type;
66  typedef typename collection_type::value_type value_type;
67  typedef typename collection_type::distance_type distance_type;
68 
69  typedef typename collection_type::const_iterator const_iterator;
70  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
71  typedef typename collection_type::iterator iterator;
72  typedef typename collection_type::reverse_iterator reverse_iterator;
73 
76 
80 
81 
82  /**
83  * Default constructor.
84  */
86  {}
87 
88 
89  /**
90  * Recursive interpolation method implementation.
91  *
92  * \param pX pointer to abscissa values
93  * \return function value
94  */
95  virtual result_type evaluate(const argument_type* pX) const
96  {
97  if (this->empty()) {
98  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
99  }
100 
101  const argument_type x = *pX;
102 
103  const_iterator p = this->lower_bound(x);
104 
105  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
106  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
107  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
108  }
109 
110  ++pX; // next argument value
111 
112 
113  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
114 
115  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
116  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
117 
118 
119  int j = 0;
120 
121  for (int i = 0; i != n; ++p, ++i) {
122 
123  u[i] = this->getDistance(x, p->getX());
124  v[i] = function_type::getValue(p->getY(), pX);
125  w[i] = v[i];
126 
127  if (fabs(u[i]) < fabs(u[j])) {
128  j = i;
129  }
130  }
131 
132 
133  result_type y = v[j];
134 
135  --j;
136 
137  for (int m = 1; m != n; ++m) {
138 
139  for (int i = 0; i != n-m; ++i) {
140 
141  const double ho = u[ i ];
142  const double hp = u[i+m];
143  const double dx = ho - hp;
144 
145  v[i] = v[i+1];
146  v[i] -= w[ i ];
147  w[i] = v[ i ];
148 
149  v[i] *= ho/dx;
150  w[i] *= hp/dx;
151  }
152 
153  if (2*(j+1) < n - m)
154  y += v[j+1];
155  else
156  y += w[j--];
157  }
158 
159  return y;
160  }
161 
162  private:
163  /**
164  * Function compilation.
165  */
166  virtual void do_compile()
167  {}
168 
169 
170  mutable double u[N+1];
171  mutable result_type v[N+1];
172  mutable result_type w[N+1];
173  };
174 
175 
176  /**
177  * Template specialisation for zero-degree polynomial interpolation method.
178  * The interpolation is based on a simple lookup table.
179  */
180  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
181  class JPolintFunction<0,
182  JElement_t,
183  JCollection_t,
184  typename JResultType<typename JElement_t::ordinate_type>::result_type,
185  JDistance_t> :
186  public JCollection_t<JElement_t, JDistance_t>,
187  public JFunction<typename JElement_t::abscissa_type,
188  typename JResultType<typename JElement_t::ordinate_type>::result_type>
189  {
190  public:
191 
192  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
193 
194  typedef typename collection_type::abscissa_type abscissa_type;
195  typedef typename collection_type::ordinate_type ordinate_type;
196  typedef typename collection_type::value_type value_type;
197  typedef typename collection_type::distance_type distance_type;
198 
199  typedef typename collection_type::const_iterator const_iterator;
200  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
201  typedef typename collection_type::iterator iterator;
202  typedef typename collection_type::reverse_iterator reverse_iterator;
203 
206 
210 
211 
212  /**
213  * Default constructor.
214  */
216  {}
217 
218 
219  /**
220  * Recursive interpolation method implementation.
221  *
222  * \param pX pointer to abscissa values
223  * \return function value
224  */
225  virtual result_type evaluate(const argument_type* pX) const
226  {
227  if (this->empty()) {
228  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
229  }
230 
231  const argument_type x = *pX;
232 
233  const_iterator p = this->lower_bound(x);
234 
235  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
236  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
237 
238  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
239  }
240 
241  ++pX; // next argument value
242 
243 
244  const_iterator q = p--;
245 
246  if (q == this->begin() || this->getDistance(x, q->getX()) < this->getDistance(p->getX(), x))
247  return function_type::getValue(q->getY(), pX);
248  else
249  return function_type::getValue(p->getY(), pX);
250  }
251 
252  private:
253  /**
254  * Function compilation.
255  */
256  virtual void do_compile()
257  {}
258  };
259 
260 
261  /**
262  * Template specialisation for first-degree polynomial interpolation method.
263  * The interpolation is based on a simple linear interpolation.
264  */
265  template<class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
266  class JPolintFunction<1,
267  JElement_t,
268  JCollection_t,
269  typename JResultType<typename JElement_t::ordinate_type>::result_type,
270  JDistance_t> :
271  public JCollection_t<JElement_t, JDistance_t>,
272  public JFunction<typename JElement_t::abscissa_type,
273  typename JResultType<typename JElement_t::ordinate_type>::result_type>
274  {
275  public:
276 
277  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
278 
279  typedef typename collection_type::abscissa_type abscissa_type;
280  typedef typename collection_type::ordinate_type ordinate_type;
281  typedef typename collection_type::value_type value_type;
282  typedef typename collection_type::distance_type distance_type;
283 
284  typedef typename collection_type::const_iterator const_iterator;
285  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
286  typedef typename collection_type::iterator iterator;
287  typedef typename collection_type::reverse_iterator reverse_iterator;
288 
291 
295 
296 
297  /**
298  * Default constructor.
299  */
301  {}
302 
303 
304  /**
305  * Recursive interpolation method implementation.
306  *
307  * \param pX pointer to abscissa values
308  * \return function value
309  */
310  virtual result_type evaluate(const argument_type* pX) const
311  {
312  if (this->empty()) {
313  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
314  }
315 
316  const argument_type x = *pX;
317 
318  const_iterator p = this->lower_bound(x);
319 
320  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
321  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
322 
323  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
324  }
325 
326  ++pX; // next argument value
327 
328 
329  const_iterator q = p--;
330 
331  const double dx = this->getDistance(p->getX(), q->getX());
332  const double a = this->getDistance(x, q->getX()) / dx;
333  const double b = 1.0 - a;
334 
335  return a * function_type::getValue(p->getY(), pX) + b * function_type::getValue(q->getY(), pX);
336  }
337 
338  private:
339  /**
340  * Function compilation.
341  */
342  virtual void do_compile()
343  {}
344  };
345 
346 
347  /**
348  * Template specialisation for polynomial interpolation method with returning JResultHesse data structure.
349  */
350  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
352  JElement_t,
353  JCollection_t,
354  JResultHesse<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
355  JDistance_t> :
356  public JCollection_t<JElement_t, JDistance_t>,
357  public JFunction<typename JElement_t::abscissa_type,
358  JResultHesse<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
359  {
360  public:
361 
362  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
363 
364  typedef typename collection_type::abscissa_type abscissa_type;
365  typedef typename collection_type::ordinate_type ordinate_type;
366  typedef typename collection_type::value_type value_type;
367  typedef typename collection_type::distance_type distance_type;
368 
369  typedef typename collection_type::const_iterator const_iterator;
370  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
371  typedef typename collection_type::iterator iterator;
372  typedef typename collection_type::reverse_iterator reverse_iterator;
373 
376 
380 
382 
383 
384  /**
385  * Default contructor.
386  */
388  {}
389 
390 
391  /**
392  * Recursive interpolation method implementation.
393  *
394  * \param pX pointer to abscissa values
395  * \return function value
396  */
397  virtual result_type evaluate(const argument_type* pX) const
398  {
399  if (this->empty()) {
400  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
401  }
402 
403  const argument_type x = *pX;
404 
405  const_iterator p = this->lower_bound(x);
406 
407  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
408  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
409 
410  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
411  }
412 
413  ++pX; // next argument value
414 
415 
416  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
417 
418  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
419  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
420 
421 
422  int j = 0;
423 
424  for (int i = 0; i != n; ++p, ++i) {
425 
426  u [i] = this->getDistance(x, p->getX());
428  w [i] = v[i];
429  vp[i] = JMATH::zero;
430  wp[i] = JMATH::zero;
431 
432  if (fabs(u[i]) < fabs(u[j])) {
433  j = i;
434  }
435  }
436 
437 
438  result.f = v[j];
439  result.fp = JMATH::zero;
440 
441  --j;
442 
443  for (int m = 1; m != n; ++m) {
444 
445  for (int i = 0; i != n-m; ++i) {
446 
447  const double ho = u[ i ];
448  const double hp = u[i+m];
449  const double dx = ho - hp;
450 
451  const data_type r = (v [i+1] - w [i]) / dx;
452  const data_type rp = (vp[i+1] - wp[i]) / dx;
453 
454  v [i] = ho * r;
455  w [i] = hp * r;
456  vp[i] = ho * rp - r;
457  wp[i] = hp * rp - r;
458  }
459 
460  if (2*(j+1) < n - m) {
461 
462  result.f += v [j+1];
463  result.fp += vp[j+1];
464 
465  } else {
466 
467  result.f += w [j];
468  result.fp += wp[j];
469 
470  --j;
471  }
472  }
473 
474  return result;
475  }
476 
477  private:
478  /**
479  * Function compilation.
480  */
481  virtual void do_compile()
482  {}
483 
484 
485  mutable double u [N+1];
486  mutable data_type v [N+1];
487  mutable data_type w [N+1];
488  mutable data_type vp[N+1];
489  mutable data_type wp[N+1];
490 
492  };
493 
494 
495  /**
496  * Template specialisation for polynomial interpolation method with returning JResultPDF data structure.
497  *
498  * Note that the data structure of the elements in the collection should have the additional methods:
499  * <pre>
500  * ordinate_type getIntegral() const;
501  * void setIntegral(ordinate_type v);
502  * </pre>
503  * to get and set the integral values, respectively.
504  */
505  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
507  JElement_t,
508  JCollection_t,
509  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type>,
510  JDistance_t> :
511  public JCollection_t<JElement_t, JDistance_t>,
512  public JFunction<typename JElement_t::abscissa_type,
513  JResultPDF<typename JResultType<typename JElement_t::ordinate_type>::result_type> >
514  {
515  public:
516 
517  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
518 
519  typedef typename collection_type::abscissa_type abscissa_type;
520  typedef typename collection_type::ordinate_type ordinate_type;
521  typedef typename collection_type::value_type value_type;
522  typedef typename collection_type::distance_type distance_type;
523 
524  typedef typename collection_type::const_iterator const_iterator;
525  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
526  typedef typename collection_type::iterator iterator;
527  typedef typename collection_type::reverse_iterator reverse_iterator;
528 
531 
535 
537 
538 
539  /**
540  * Default contructor.
541  */
543  {}
544 
545 
546  /**
547  * Recursive interpolation method implementation.
548  *
549  * \param pX pointer to abscissa values
550  * \return function value
551  */
552  virtual result_type evaluate(const argument_type* pX) const
553  {
554  if (this->empty()) {
555  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
556  }
557 
558  const argument_type x = *pX;
559 
560  const_iterator p = this->lower_bound(x);
561 
562  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
563  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
564 
565  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
566  }
567 
568  ++pX; // next argument value
569 
570  {
571  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
572 
573  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
574  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
575 
576 
577  int j = 0;
578 
579  for (int i = 0; i != n; ++p, ++i) {
580 
581  u [i] = this->getDistance(x, p->getX());
583  w [i] = v[i];
584  vp[i] = JMATH::zero;
585  wp[i] = JMATH::zero;
586  vi[i] = p->getIntegral();
587  wi[i] = vi[i];
588 
589  if (fabs(u[i]) < fabs(u[j])) {
590  j = i;
591  }
592  }
593 
594 
595  result.f = v [j];
596  result.fp = JMATH::zero;
597  result.v = vi[j];
598  result.V = this->rbegin()->getIntegral();
599 
600  --j;
601 
602  for (int m = 1; m != n; ++m) {
603 
604  for (int i = 0; i != n-m; ++i) {
605 
606  const double ho = u[ i ];
607  const double hp = u[i+m];
608  const double dx = ho - hp;
609 
610  const data_type r = (v [i+1] - w [i]) / dx;
611  const data_type rp = (vp[i+1] - wp[i]) / dx;
612  const data_type ri = (vi[i+1] - wi[i]) / dx;
613 
614  v [i] = ho * r;
615  w [i] = hp * r;
616  vp[i] = ho * rp - r;
617  wp[i] = hp * rp - r;
618  vi[i] = ho * ri;
619  wi[i] = hp * ri;
620  }
621 
622  if (2*(j+1) < n - m) {
623 
624  result.f += v [j+1];
625  result.fp += vp[j+1];
626  result.v += vi[j+1];
627 
628  } else {
629 
630  result.f += w [j];
631  result.fp += wp[j];
632  result.v += wi[j];
633 
634  --j;
635  }
636  }
637  }
638 
639  return result;
640  }
641 
642  private:
643  /**
644  * Function compilation.
645  */
646  virtual void do_compile()
647  {
649 
650  if (this->getSize() > 1) {
651 
652  const JGaussLegendre engine(N);
653 
654  this->begin()->setIntegral(V);
655 
656  for (iterator j = this->begin(), i = j++; j != this->end(); ++i, ++j) {
657 
658  const abscissa_type xmin = i->getX();
659  const abscissa_type xmax = j->getX();
660 
661  for (JGaussLegendre::const_iterator m = engine.begin(); m != engine.end(); ++m) {
662 
663  const abscissa_type x = 0.5 * (xmax + xmin + m->getX() * (xmax - xmin));
664  const ordinate_type v = 0.5 * (xmax - xmin) * m->getY() * get_value(this->evaluate(&x));
665 
666  V += v;
667  }
668 
669  j->setIntegral(V);
670  }
671  }
672  }
673 
674 
675  mutable double u [N+1];
676  mutable data_type v [N+1];
677  mutable data_type w [N+1];
678  mutable data_type vp[N+1];
679  mutable data_type wp[N+1];
680  mutable data_type vi[N+1];
681  mutable data_type wi[N+1];
682 
684  };
685 
686 
687  /**
688  * Template specialisation for polynomial interpolation method with returning JResultPolynome data structure.
689  */
690  template<unsigned int N, class JElement_t, template<class, class> class JCollection_t, class JDistance_t>
692  JElement_t,
693  JCollection_t,
694  JResultPolynome<N, typename JResultType<typename JElement_t::ordinate_type>::result_type>,
695  JDistance_t> :
696  public JCollection_t<JElement_t, JDistance_t>,
697  public JFunction<typename JElement_t::abscissa_type,
698  JResultPolynome<N, typename JResultType<typename JElement_t::ordinate_type>::result_type> >
699  {
700  public:
701 
702  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
703 
704  typedef typename collection_type::abscissa_type abscissa_type;
705  typedef typename collection_type::ordinate_type ordinate_type;
706  typedef typename collection_type::value_type value_type;
707  typedef typename collection_type::distance_type distance_type;
708 
709  typedef typename collection_type::const_iterator const_iterator;
710  typedef typename collection_type::const_reverse_iterator const_reverse_iterator;
711  typedef typename collection_type::iterator iterator;
712  typedef typename collection_type::reverse_iterator reverse_iterator;
713 
716 
720 
722 
723 
724  /**
725  * Default contructor.
726  */
728  {}
729 
730 
731  /**
732  * Recursive interpolation method implementation.
733  *
734  * \param pX pointer to abscissa values
735  * \return function value
736  */
737  virtual result_type evaluate(const argument_type* pX) const
738  {
739  if (this->empty()) {
740  return this->getExceptionHandler().action(JEmptyCollection("JPolintFunction<>::evaluate() no data."));
741  }
742 
743  const argument_type x = *pX;
744 
745  const_iterator p = this->lower_bound(x);
746 
747  if ((p == this->begin() && this->getDistance(x, (p++)->getX()) > distance_type::precision) ||
748  (p == this->end() && this->getDistance((--p)->getX(), x) > distance_type::precision)) {
749 
750  return this->getExceptionHandler().action(JValueOutOfRange("JPolintFunction::evaluate() x out of range."));
751  }
752 
753  ++pX; // next argument value
754 
755 
756  const int n = std::min((int) (N + 1), (int) this->size()); // number of points to interpolate
757 
758  for (int i = n/2; i != 0 && p != this->end(); --i, ++p) {} // move p to begin of data
759  for (int i = n ; i != 0 && p != this->begin(); --i, --p) {}
760 
761 
762  for (unsigned int i = 0; i != N+1; ++i) {
763  for (unsigned int j = 0; j != N+1; ++j) {
764 
765  v[i][j] = JMATH::zero;
766  w[i][j] = JMATH::zero;
767  }
768 
769  result.y[i] = JMATH::zero;
770  }
771 
772 
773  int j = 0;
774 
775  for (int i = 0; i != n; ++p, ++i) {
776 
777  u[i] = this->getDistance(x, p->getX());
778  v[i][0] = JFunctional<argument_type, data_type>::getValue(p->getY(), pX);
779  w[i][0] = v[i][0];
780 
781  if (fabs(u[i]) < fabs(u[j])) {
782  j = i;
783  }
784  }
785 
786  result.y[0] = v[j][0];
787 
788  --j;
789 
790  for (int m = 1; m != n; ++m) {
791 
792  for (int i = 0; i != n-m; ++i) {
793 
794  const double ho = u[ i ];
795  const double hp = u[i+m];
796  const double dx = ho - hp;
797 
798  for (int k = 0; k != N+1; ++k) {
799  r[k] = (v[i+1][k] - w[i][k]) / dx;
800  }
801 
802  v[i][0] = ho * r[0];
803  w[i][0] = hp * r[0];
804 
805  for (int k = 1; k != N+1; ++k) {
806  v[i][k] = ho * r[k] - k * r[k-1];
807  w[i][k] = hp * r[k] - k * r[k-1];
808  }
809  }
810 
811  if (2*(j+1) < n - m) {
812 
813  for (int k = 0; k != N+1; ++k) {
814  result.y[k] += v[j+1][k];
815  }
816 
817  } else {
818 
819  for (int k = 0; k != N+1; ++k) {
820  result.y[k] += w[j][k];
821  }
822 
823  --j;
824  }
825  }
826 
827  return result;
828  }
829 
830  private:
831  /**
832  * Function compilation.
833  */
834  virtual void do_compile()
835  {}
836 
837 
838  mutable double u[N+1];
839  mutable data_type v[N+1][N+1];
840  mutable data_type w[N+1][N+1];
841  mutable data_type r[N+1];
842 
844  };
845 
846 
847  /**
848  * Template class for polynomial interpolation in 1D
849  *
850  * This class implements the JFunction1D interface.
851  */
852  template<unsigned int N,
853  class JElement_t,
854  template<class, class> class JCollection_t,
855  class JResult_t = typename JElement_t::ordinate_type,
858  public JPolintFunction<N, JElement_t, JCollection_t, JResult_t, JDistance_t>,
859  public JFunction1D<typename JElement_t::abscissa_type, JResult_t>
860  {
861  public:
862 
863  typedef JCollection_t<JElement_t, JDistance_t> collection_type;
864 
869 
874 
876 
880 
881 
882  /**
883  * Default contructor.
884  */
886  {}
887  };
888 
889 
890  /**
891  * Functional map with polynomial interpolation.
892  */
893  template<unsigned int N,
894  class JKey_t,
895  class JValue_t,
896  template<class, class, class> class JMap_t,
897  class JResult_t,
898  class JDistance_t = JDistance<JKey_t> >
899  class JPolintMap :
900  public JPolintFunction<N,
901  JElement2D<JKey_t, JValue_t>,
902  JMapCollection<JMap_t>::template collection_type,
903  JResult_t,
904  JDistance_t>
905  {
906  public:
907 
909  typedef JPolintFunction<N,
910  element_type,
911  JMapCollection<JMap_t>::template collection_type,
912  JResult_t,
913  JDistance_t> JPolintFunction_t;
914 
915  typedef typename JPolintFunction_t::abscissa_type abscissa_type;
916  typedef typename JPolintFunction_t::ordinate_type ordinate_type;
917  typedef typename JPolintFunction_t::value_type value_type;
918  typedef typename JPolintFunction_t::distance_type distance_type;
919 
920  typedef typename JPolintFunction_t::const_iterator const_iterator;
921  typedef typename JPolintFunction_t::const_reverse_iterator const_reverse_iterator;
922  typedef typename JPolintFunction_t::iterator iterator;
923  typedef typename JPolintFunction_t::reverse_iterator reverse_iterator;
924 
925  typedef typename JPolintFunction_t::argument_type argument_type;
926  typedef typename JPolintFunction_t::result_type result_type;
927  typedef typename JPolintFunction_t::JExceptionHandler exceptionhandler_type;
928 
929 
930  /**
931  * Default constructor.
932  */
934  {}
935  };
936 
937 
938  /**
939  * Conversion of data points to integral values.
940  *
941  * This method transfers the integration to the corresponding specialised function.
942  *
943  * \param input collection
944  * \param output mappable collection
945  * \return integral
946  */
947  template<unsigned int N,
948  class JElement_t,
949  template<class, class> class JCollection_t,
950  class JResult_t,
951  class JDistance_t>
952  inline typename JElement_t::ordinate_type
954  typename JMappable<JElement_t>::map_type& output)
955  {
956  return integrate(input, output, JLANG::JBool<N == 0 || N == 1>());
957  }
958 
959 
960  /**
961  * Conversion of data points to integral values.
962  *
963  * The integration uses the Gauss-Legendre quadratures with the number of points set
964  * to the degree of the input polynomial interpolating function.
965  *
966  * \param input collection
967  * \param output mappable collection
968  * \param option false
969  * \return integral
970  */
971  template<unsigned int N,
972  class JElement_t,
973  template<class, class> class JCollection_t,
974  class JResult_t,
975  class JDistance_t>
976  inline typename JElement_t::ordinate_type
978  typename JMappable<JElement_t>::map_type& output,
979  const JLANG::JBool<false>& option)
980  {
981  typedef typename JElement_t::abscissa_type abscissa_type;
982  typedef typename JElement_t::ordinate_type ordinate_type;
984 
985  ordinate_type V(JMATH::zero);
986 
987  if (input.getSize() > 1) {
988 
989  output.put(input.begin()->getX(), V);
990 
991  const JGaussLegendre engine(N);
992 
993  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
994 
995  const abscissa_type xmin = i->getX();
996  const abscissa_type xmax = j->getX();
997 
998  for (JGaussLegendre::const_iterator m = engine.begin(); m != engine.end(); ++m) {
999 
1000  const abscissa_type x = 0.5 * (xmax + xmin + m->getX() * (xmax - xmin));
1001  const ordinate_type v = 0.5 * (xmax - xmin) * m->getY() * get_value(input(x));
1002 
1003  V += v;
1004  }
1005 
1006  output.put(xmax, V);
1007  }
1008  }
1009 
1010  return V;
1011  }
1012 
1013 
1014  /**
1015  * Conversion of data points to integral values.
1016  *
1017  * The integration is based on the sum of ordinates of the input data points.
1018  *
1019  * \param input collection
1020  * \param output mappable collection
1021  * \param option true
1022  * \return integral
1023  */
1024  template<class JElement_t,
1025  template<class, class> class JCollection_t,
1026  class JResult_t,
1027  class JDistance_t>
1028  inline typename JElement_t::ordinate_type
1030  typename JMappable<JElement_t>::map_type& output,
1031  const JLANG::JBool<true>& option)
1032  {
1033  typedef typename JElement_t::ordinate_type ordinate_type;
1035 
1036  ordinate_type V(JMATH::zero);
1037 
1038  if (input.getSize() > 1) {
1039 
1040  output.put(input.begin()->getX(), V);
1041 
1042  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
1043 
1044  V += input.getDistance(i->getX(), j->getX()) * j->getY();
1045 
1046  output.put(j->getX(), V);
1047  }
1048  }
1049 
1050  return V;
1051  }
1052 
1053 
1054  /**
1055  * Conversion of data points to integral values.
1056  *
1057  * The integration is based on the trapezoidal rule applied to the input data points.
1058  *
1059  * \param input collection
1060  * \param output mappable collection
1061  * \param option true
1062  * \return integral
1063  */
1064  template<class JElement_t,
1065  template<class, class> class JCollection_t,
1066  class JResult_t,
1067  class JDistance_t>
1068  inline typename JElement_t::ordinate_type
1070  typename JMappable<JElement_t>::map_type& output,
1071  const JLANG::JBool<true>& option)
1072  {
1073  typedef typename JElement_t::ordinate_type ordinate_type;
1075 
1076  ordinate_type V(JMATH::zero);
1077 
1078  if (input.getSize() > 1) {
1079 
1080  output.put(input.begin()->getX(), V);
1081 
1082  for (const_iterator j = input.begin(), i = j++; j != input.end(); ++i, ++j) {
1083 
1084  V += 0.5 * input.getDistance(i->getX(), j->getX()) * (i->getY() + j->getY());
1085 
1086  output.put(j->getX(), V);
1087  }
1088  }
1089 
1090  return V;
1091  }
1092 }
1093 
1094 #endif
container_type::reverse_iterator reverse_iterator
Definition: JCollection.hh:92
JElement2D< JKey_t, JValue_t > element_type
Definition: JPolint.hh:908
Exceptions.
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
2D Element.
Definition: JElement.hh:44
Template definition for functional collection with polynomial interpolation.
Definition: JPolint.hh:40
function_type::JExceptionHandler exceptionhandler_type
Definition: JPolint.hh:879
Numerical integrator for W(x) = 1.
Definition: JQuadrature.hh:113
function_type::argument_type argument_type
Definition: JPolint.hh:877
JElement_t value_type
Definition: JCollection.hh:82
The elements in a collection are sorted according to their abscissa values and a given distance opera...
collection_type::value_type value_type
Definition: JPolint.hh:867
Template class for distance evaluation.
Definition: JDistance.hh:24
Data structure for result including value and N derivatives of function.
Definition: JResult.hh:407
This include file containes various data structures that can be used as specific return types for the...
JPolintFunction1D()
Default contructor.
Definition: JPolint.hh:885
functional_type::argument_type argument_type
Definition: JFunctional.hh:323
Template interface definition for associative collection of elements.
collection_type::const_iterator const_iterator
Definition: JPolint.hh:870
virtual result_type evaluate(const argument_type *pX) const
Recursive interpolation method implementation.
Definition: JPolint.hh:225
static const JZero zero
Function object to assign zero value.
Definition: JZero.hh:94
JPolintFunction_t::reverse_iterator reverse_iterator
Definition: JPolint.hh:923
size_t getSize(T(&array)[N])
Get size of c-array.
Definition: JLangToolkit.hh:32
JElement_t::abscissa_type abscissa_type
Definition: JCollection.hh:80
Template definition of function object interface in one dimension.
Definition: JFunctional.hh:317
double getDistance(const JFirst_t &first, const JSecond_t &second)
Get distance between objects.
Definition: JMathToolkit.hh:98
JPolintFunction_t::JExceptionHandler exceptionhandler_type
Definition: JPolint.hh:927
collection_type::ordinate_type ordinate_type
Definition: JPolint.hh:866
JPolintMap()
Default constructor.
Definition: JPolint.hh:933
JPolintFunction_t::iterator iterator
Definition: JPolint.hh:922
JPolintFunction_t::value_type value_type
Definition: JPolint.hh:917
JPolintFunction_t::ordinate_type ordinate_type
Definition: JPolint.hh:916
static result_type getValue(const JFunctional &function, const argument_type *pX)
Recursive function value evaluation.
Definition: JFunctional.hh:106
Template definition of function object interface in multidimensions.
Definition: JFunctional.hh:301
JPolintFunction_t::const_reverse_iterator const_reverse_iterator
Definition: JPolint.hh:921
Auxiliary template class for type bool.
Definition: JBool.hh:20
JPolintFunction< N, element_type, JMapCollection< JMap_t >::template collection_type, JResult_t, JDistance_t > JPolintFunction_t
Definition: JPolint.hh:913
container_type::iterator iterator
Definition: JCollection.hh:91
JPolintFunction_t::argument_type argument_type
Definition: JPolint.hh:925
Data structure for result including value, first derivative and integrals of function.
Definition: JResult.hh:211
Template definition of function object interface.
Definition: JFunctional.hh:31
functional_type::result_type result_type
Definition: JFunctional.hh:324
virtual result_type evaluate(const argument_type *pX) const
Recursive interpolation method implementation.
Definition: JPolint.hh:95
container_type::const_reverse_iterator const_reverse_iterator
Definition: JCollection.hh:90
Auxiliary class to evaluate result type.
Definition: JFunctional.hh:377
void put(typename JClass< key_type >::argument_type key, typename JClass< mapped_type >::argument_type value)
Put pair-wise element (key,value) into collection.
Template class to retreive underlying collection for the given template map.
collection_type::const_reverse_iterator const_reverse_iterator
Definition: JPolint.hh:871
JDistance_t distance_type
Definition: JCollection.hh:83
Exception handler for functional object.
Definition: JFunctional.hh:129
collection_type::abscissa_type abscissa_type
Definition: JPolint.hh:865
JCollection_t< JElement_t, JDistance_t > collection_type
Definition: JPolint.hh:863
collection_type::iterator iterator
Definition: JPolint.hh:872
virtual result_type evaluate(const argument_type *pX) const
Recursive interpolation method implementation.
Definition: JPolint.hh:310
Exception for numerical precision error.
Definition: JException.hh:234
Auxiliary classes for numerical integration.
JElement_t::ordinate_type ordinate_type
Definition: JCollection.hh:81
Functional map with polynomial interpolation.
Definition: JPolint.hh:899
Exception for an empty collection.
Definition: JException.hh:126
functional_type::argument_type argument_type
Definition: JFunctional.hh:307
container_type::const_iterator const_iterator
Definition: JCollection.hh:89
JPolintFunction_t::result_type result_type
Definition: JPolint.hh:926
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:144
JPolintFunction_t::abscissa_type abscissa_type
Definition: JPolint.hh:915
collection_type::reverse_iterator reverse_iterator
Definition: JPolint.hh:873
collection_type::distance_type distance_type
Definition: JPolint.hh:868
functional_type::result_type result_type
Definition: JFunctional.hh:308
Data structure for result including value and first derivative of function.
Definition: JResult.hh:40
JFunction1D< abscissa_type, JResult_t > function_type
Definition: JPolint.hh:875
function_type::result_type result_type
Definition: JPolint.hh:878
JResultEvaluator< JResult_t >::result_type get_value(const JResult_t &value)
Helper method to recursively evaluate a to function value.
Definition: JResult.hh:748
JPolintFunction_t::distance_type distance_type
Definition: JPolint.hh:918
JPolintFunction_t::const_iterator const_iterator
Definition: JPolint.hh:920
Template class for polynomial interpolation in 1D.
Definition: JPolint.hh:857
JElement_t::ordinate_type integrate(const JCollection< JElement_t, JDistance_t > &input, typename JMappable< JElement_t >::map_type &output)
Conversion of data points to integral values.
Definition: JCollection.hh:793