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