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