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