Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Sacado_PCE_OrthogPolyImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Sacado_DynamicArrayTraits.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
13 
14 namespace Sacado {
15 namespace PCE {
16 
17 template <typename T, typename Storage>
19 OrthogPoly() :
20  expansion_(),
22 {
24  expansion_ = const_expansion_;
25 }
26 
27 template <typename T, typename Storage>
30  expansion_(),
32 {
34  expansion_ = const_expansion_;
35 }
36 
37 template <typename T, typename Storage>
39 OrthogPoly(const Teuchos::RCP<expansion_type>& expansion) :
40  expansion_(expansion),
41  th(new Stokhos::OrthogPolyApprox<int,value_type,Storage>(expansion_->getBasis()))
42 {
44 }
45 
46 template <typename T, typename Storage>
49  ordinal_type sz) :
50  expansion_(expansion),
51  th(new Stokhos::OrthogPolyApprox<int,value_type,Storage>(expansion_->getBasis(), sz))
52 {
54 }
55 
56 template <typename T, typename Storage>
59  expansion_(x.expansion_),
60  th(x.th)
61 {
63 }
64 
65 template <typename T, typename Storage>
68 {
69 }
70 
71 template <typename T, typename Storage>
72 void
74 reset(const Teuchos::RCP<expansion_type>& expansion)
75 {
76  expansion_ = expansion;
77  th->reset(expansion_->getBasis());
78 }
79 
80 template <typename T, typename Storage>
81 void
84 {
85  expansion_ = expansion;
86  th->reset(expansion_->getBasis(), sz);
87 }
88 
89 template <typename T, typename Storage>
92 evaluate(const Teuchos::Array<typename OrthogPoly<T,Storage>::value_type>& point) const
93 {
94  return th->evaluate(point);
95 }
96 
97 template <typename T, typename Storage>
100 evaluate(
101  const Teuchos::Array<typename OrthogPoly<T,Storage>::value_type>& point,
102  const Teuchos::Array<typename OrthogPoly<T,Storage>::value_type>& bvals) const
103 {
104  return th->evaluate(point, bvals);
105 }
106 
107 template <typename T, typename Storage>
108 bool
110 isEqualTo(const OrthogPoly& x) const {
111  typedef IsEqual<value_type> IE;
112  if (x.size() != this->size()) return false;
113  // Allow expansions to be different if their size is 1 and one
114  // of them is a constant expansion
115  if (expansion_ != x.expansion_) {
116  if (x.size() != 1)
117  return false;
118  if ((expansion_ != const_expansion_) &&
119  (x.expansion_ != x.const_expansion_))
120  return false;
121  }
122  bool eq = true;
123  for (int i=0; i<this->size(); i++)
124  eq = eq && IE::eval(x.coeff(i), this->coeff(i));
125  return eq;
126 }
127 
128 template <typename T, typename Storage>
132 {
133  th.makeOwnCopy();
134  *th = v;
135  return *this;
136 }
137 
138 template <typename T, typename Storage>
142 {
143  expansion_ = x.expansion_;
144  th = x.th;
145  return *this;
146 }
147 
148 template <typename T, typename Storage>
151 operator+() const
152 {
153  return *this;
154 }
155 
156 template <typename T, typename Storage>
159 operator-() const
160 {
161  OrthogPoly<T,Storage> x(expansion_);
162  expansion_->unaryMinus(*(x.th), *th);
163  return x;
164 }
165 
166 template <typename T, typename Storage>
170 {
171  th.makeOwnCopy();
172  expansion_->plusEqual(*th, v);
173  return *this;
174 }
175 
176 template <typename T, typename Storage>
180 {
181  th.makeOwnCopy();
182  expansion_->minusEqual(*th, v);
183  return *this;
184 }
185 
186 template <typename T, typename Storage>
190 {
191  th.makeOwnCopy();
192  expansion_->timesEqual(*th, v);
193  return *this;
194 }
195 
196 template <typename T, typename Storage>
200 {
201  th.makeOwnCopy();
202  expansion_->divideEqual(*th, v);
203  return *this;
204 }
205 
206 template <typename T, typename Storage>
210 {
211  th.makeOwnCopy();
213  if (x.size() > size()) {
214  e = x.expansion();
215  reset(e, size());
216  }
217  e->plusEqual(*th, *x.th);
218  return *this;
219 }
220 
221 template <typename T, typename Storage>
225 {
226  th.makeOwnCopy();
228  if (x.size() > size()) {
229  e = x.expansion();
230  reset(e, size());
231  }
232  e->minusEqual(*th, *x.th);
233  return *this;
234 }
235 
236 template <typename T, typename Storage>
240 {
241  th.makeOwnCopy();
243  if (x.size() > size()) {
244  e = x.expansion();
245  reset(e, size());
246  }
247  e->timesEqual(*th, *x.th);
248  return *this;
249 }
250 
251 template <typename T, typename Storage>
255 {
256  th.makeOwnCopy();
258  if (x.size() > size()) {
259  e = x.expansion();
260  reset(e, size());
261  }
262  e->divideEqual(*th, *x.th);
263  return *this;
264 }
265 
266 template <typename T, typename Storage>
269  const OrthogPoly<T,Storage>& b)
270 {
271  // Get expansion
273  ordinal_type da = a.size();
274  ordinal_type db = b.size();
276  if (da == db || da > 1)
277  e = a.expansion();
278  else
279  e = b.expansion();
280 
281  OrthogPoly<T,Storage> c(e, 0);
282  e->plus(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
283  b.getOrthogPolyApprox());
284 
285  return c;
286 }
287 
288 template <typename T, typename Storage>
289 OrthogPoly<T,Storage>
291  const OrthogPoly<T,Storage>& b)
292 {
293  OrthogPoly<T,Storage> c(b.expansion(), 0);
294  b.expansion()->plus(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
295  return c;
296 }
297 
298 template <typename T, typename Storage>
299 OrthogPoly<T,Storage>
301  const typename OrthogPoly<T,Storage>::value_type& b)
302 {
303  OrthogPoly<T,Storage> c(a.expansion(), 0);
304  a.expansion()->plus(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
305  return c;
306 }
307 
308 template <typename T, typename Storage>
309 OrthogPoly<T,Storage>
311  const OrthogPoly<T,Storage>& b)
312 {
313  // Get expansion
315  ordinal_type da = a.size();
316  ordinal_type db = b.size();
318  if (da == db || da > 1)
319  e = a.expansion();
320  else
321  e = b.expansion();
322 
323  OrthogPoly<T,Storage> c(e, 0);
324  e->minus(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
325  b.getOrthogPolyApprox());
326 
327  return c;
328 }
329 
330 template <typename T, typename Storage>
331 OrthogPoly<T,Storage>
333  const OrthogPoly<T,Storage>& b)
334 {
335  OrthogPoly<T,Storage> c(b.expansion(), 0);
336  b.expansion()->minus(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
337  return c;
338 }
339 
340 template <typename T, typename Storage>
341 OrthogPoly<T,Storage>
343  const typename OrthogPoly<T,Storage>::value_type& b)
344 {
345  OrthogPoly<T,Storage> c(a.expansion(), 0);
346  a.expansion()->minus(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
347  return c;
348 }
349 
350 template <typename T, typename Storage>
351 OrthogPoly<T,Storage>
353  const OrthogPoly<T,Storage>& b)
354 {
355  // Get expansion
357  ordinal_type da = a.size();
358  ordinal_type db = b.size();
360  if (da == db || da > 1)
361  e = a.expansion();
362  else
363  e = b.expansion();
364 
365  OrthogPoly<T,Storage> c(e, 0);
366  e->times(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
367  b.getOrthogPolyApprox());
368 
369  return c;
370 }
371 
372 template <typename T, typename Storage>
373 OrthogPoly<T,Storage>
375  const OrthogPoly<T,Storage>& b)
376 {
377  OrthogPoly<T,Storage> c(b.expansion(), 0);
378  b.expansion()->times(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
379  return c;
380 }
381 
382 template <typename T, typename Storage>
383 OrthogPoly<T,Storage>
385  const typename OrthogPoly<T,Storage>::value_type& b)
386 {
387  OrthogPoly<T,Storage> c(a.expansion(), 0);
388  a.expansion()->times(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
389  return c;
390 }
391 
392 template <typename T, typename Storage>
393 OrthogPoly<T,Storage>
395  const OrthogPoly<T,Storage>& b)
396 {
397  // Get expansion
399  ordinal_type da = a.size();
400  ordinal_type db = b.size();
402  if (da == db || da > 1)
403  e = a.expansion();
404  else
405  e = b.expansion();
406 
407  OrthogPoly<T,Storage> c(e, 0);
408  e->divide(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
409  b.getOrthogPolyApprox());
410 
411  return c;
412 }
413 
414 template <typename T, typename Storage>
415 OrthogPoly<T,Storage>
417  const OrthogPoly<T,Storage>& b)
418 {
419  OrthogPoly<T,Storage> c(b.expansion(), 0);
420  b.expansion()->divide(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
421  return c;
422 }
423 
424 template <typename T, typename Storage>
425 OrthogPoly<T,Storage>
427  const typename OrthogPoly<T,Storage>::value_type& b)
428 {
429  OrthogPoly<T,Storage> c(a.expansion(), 0);
430  a.expansion()->divide(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
431  return c;
432 }
433 
434 template <typename T, typename Storage>
435 OrthogPoly<T,Storage>
437 {
438  OrthogPoly<T,Storage> c(a.expansion(), 0);
439  a.expansion()->exp(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
440  return c;
441 }
442 
443 template <typename T, typename Storage>
444 OrthogPoly<T,Storage>
446 {
448  OrthogPoly<T,Storage> c(a.expansion(), 0);
449  {
450  TEUCHOS_FUNC_TIME_MONITOR("OPA LOG");
451  a.expansion()->log(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
452  }
453 
455  return d;
456 }
457 
458 template <typename T, typename Storage>
459 void
461 {
462  OrthogPoly<T,Storage> d(a.expansion(), 0);
463  a.expansion()->log(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
464 }
465 
466 template <typename T, typename Storage>
467 OrthogPoly<T,Storage>
469 {
470  OrthogPoly<T,Storage> c(a.expansion(), 0);
471  a.expansion()->log10(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
472  return c;
473 }
474 
475 template <typename T, typename Storage>
476 OrthogPoly<T,Storage>
478 {
479  OrthogPoly<T,Storage> c(a.expansion(), 0);
480  a.expansion()->sqrt(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
481  return c;
482 }
483 
484 template <typename T, typename Storage>
485 OrthogPoly<T,Storage>
487 {
488  OrthogPoly<T,Storage> c(a.expansion(), 0);
489  a.expansion()->cbrt(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
490  return c;
491 }
492 
493 template <typename T, typename Storage>
494 OrthogPoly<T,Storage>
496  const OrthogPoly<T,Storage>& b)
497 {
498  // Get expansion
500  ordinal_type da = a.size();
501  ordinal_type db = b.size();
503  if (da == db || da > 1)
504  e = a.expansion();
505  else
506  e = b.expansion();
507 
508  OrthogPoly<T,Storage> c(e, 0);
509  e->pow(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
510  b.getOrthogPolyApprox());
511 
512  return c;
513 }
514 
515 template <typename T, typename Storage>
516 OrthogPoly<T,Storage>
517 pow(const T& a,
518  const OrthogPoly<T,Storage>& b)
519 {
520  OrthogPoly<T,Storage> c(b.expansion(), 0);
521  b.expansion()->pow(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
522  return c;
523 }
524 
525 template <typename T, typename Storage>
526 OrthogPoly<T,Storage>
528  const T& b)
529 {
530  OrthogPoly<T,Storage> c(a.expansion(), 0);
531  a.expansion()->pow(c.getOrthogPolyApprox(),a.getOrthogPolyApprox(), b);
532  return c;
533 }
534 
535 template <typename T, typename Storage>
536 OrthogPoly<T,Storage>
538 {
539  OrthogPoly<T,Storage> c(a.expansion(), 0);
540  a.expansion()->sin(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
541  return c;
542 }
543 
544 template <typename T, typename Storage>
545 OrthogPoly<T,Storage>
547 {
548  OrthogPoly<T,Storage> c(a.expansion(), 0);
549  a.expansion()->cos(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
550  return c;
551 }
552 
553 template <typename T, typename Storage>
554 OrthogPoly<T,Storage>
556 {
557  OrthogPoly<T,Storage> c(a.expansion(), 0);
558  a.expansion()->tan(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
559  return c;
560 }
561 
562 template <typename T, typename Storage>
563 OrthogPoly<T,Storage>
565 {
566  OrthogPoly<T,Storage> c(a.expansion(), 0);
567  a.expansion()->sinh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
568  return c;
569 }
570 
571 template <typename T, typename Storage>
572 OrthogPoly<T,Storage>
574 {
575  OrthogPoly<T,Storage> c(a.expansion(), 0);
576  a.expansion()->cosh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
577  return c;
578 }
579 
580 template <typename T, typename Storage>
581 OrthogPoly<T,Storage>
583 {
584  OrthogPoly<T,Storage> c(a.expansion(), 0);
585  a.expansion()->tanh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
586  return c;
587 }
588 
589 template <typename T, typename Storage>
590 OrthogPoly<T,Storage>
592 {
593  OrthogPoly<T,Storage> c(a.expansion(), 0);
594  a.expansion()->acos(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
595  return c;
596 }
597 
598 template <typename T, typename Storage>
599 OrthogPoly<T,Storage>
601 {
602  OrthogPoly<T,Storage> c(a.expansion(), 0);
603  a.expansion()->asin(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
604  return c;
605 }
606 
607 template <typename T, typename Storage>
608 OrthogPoly<T,Storage>
610 {
611  OrthogPoly<T,Storage> c(a.expansion(), 0);
612  a.expansion()->atan(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
613  return c;
614 }
615 
616 template <typename T, typename Storage>
617 OrthogPoly<T,Storage>
619 {
620  OrthogPoly<T,Storage> c(a.expansion(), 0);
621  a.expansion()->acosh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
622  return c;
623 }
624 
625 template <typename T, typename Storage>
626 OrthogPoly<T,Storage>
628 {
629  OrthogPoly<T,Storage> c(a.expansion(), 0);
630  a.expansion()->asinh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
631  return c;
632 }
633 
634 template <typename T, typename Storage>
635 OrthogPoly<T,Storage>
637 {
638  OrthogPoly<T,Storage> c(a.expansion(), 0);
639  a.expansion()->atanh(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
640  return c;
641 }
642 
643 template <typename T, typename Storage>
644 OrthogPoly<T,Storage>
646 {
647  OrthogPoly<T,Storage> c(a.expansion(), 0);
648  a.expansion()->fabs(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
649  return c;
650 }
651 
652 template <typename T, typename Storage>
653 OrthogPoly<T,Storage>
655 {
656  OrthogPoly<T,Storage> c(a.expansion(), 0);
657  a.expansion()->abs(c.getOrthogPolyApprox(), a.getOrthogPolyApprox());
658  return c;
659 }
660 
661 template <typename T, typename Storage>
662 OrthogPoly<T,Storage>
664  const OrthogPoly<T,Storage>& b)
665 {
666  // Get expansion
668  ordinal_type da = a.size();
669  ordinal_type db = b.size();
671  if (da == db || da > 1)
672  e = a.expansion();
673  else
674  e = b.expansion();
675 
676  OrthogPoly<T,Storage> c(e, 0);
677  e->max(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
678  b.getOrthogPolyApprox());
679 
680  return c;
681 }
682 
683 template <typename T, typename Storage>
684 OrthogPoly<T,Storage>
686  const OrthogPoly<T,Storage>& b)
687 {
688  OrthogPoly<T,Storage> c(b.expansion(), 0);
689  b.expansion()->max(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
690  return c;
691 }
692 
693 template <typename T, typename Storage>
694 OrthogPoly<T,Storage>
696  const typename OrthogPoly<T,Storage>::value_type& b)
697 {
698  OrthogPoly<T,Storage> c(a.expansion(), 0);
699  a.expansion()->max(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
700  return c;
701 }
702 
703 template <typename T, typename Storage>
704 OrthogPoly<T,Storage>
706  const OrthogPoly<T,Storage>& b)
707 {
708  // Get expansion
710  ordinal_type da = a.size();
711  ordinal_type db = b.size();
713  if (da == db || da > 1)
714  e = a.expansion();
715  else
716  e = b.expansion();
717 
718  OrthogPoly<T,Storage> c(e, 0);
719  e->min(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(),
720  b.getOrthogPolyApprox());
721 
722  return c;
723 }
724 
725 template <typename T, typename Storage>
726 OrthogPoly<T,Storage>
728  const OrthogPoly<T,Storage>& b)
729 {
730  OrthogPoly<T,Storage> c(b.expansion(), 0);
731  b.expansion()->min(c.getOrthogPolyApprox(), a, b.getOrthogPolyApprox());
732  return c;
733 }
734 
735 template <typename T, typename Storage>
736 OrthogPoly<T,Storage>
738  const typename OrthogPoly<T,Storage>::value_type& b)
739 {
740  OrthogPoly<T,Storage> c(a.expansion(), 0);
741  a.expansion()->min(c.getOrthogPolyApprox(), a.getOrthogPolyApprox(), b);
742  return c;
743 }
744 
745 template <typename T, typename Storage>
746 bool
748  const OrthogPoly<T,Storage>& b)
749 {
750  int n = std::max(a.size(), b.size());
751  for (int i=0; i<n; i++)
752  if (a.coeff(i) != b.coeff(i))
753  return false;
754  return true;
755 }
756 
757 template <typename T, typename Storage>
758 bool
760  const OrthogPoly<T,Storage>& b)
761 {
762  if (a != b.coeff(0))
763  return false;
764  for (int i=1; i<b.size(); i++)
765  if (b.coeff(i) != T(0.0))
766  return false;
767  return true;
768 }
769 
770 template <typename T, typename Storage>
771 bool
773  const typename OrthogPoly<T,Storage>::value_type& b)
774 {
775  if (a.coeff(0) != b)
776  return false;
777  for (int i=1; i<a.size(); i++)
778  if (a.coeff(i) != T(0.0))
779  return false;
780  return true;
781 }
782 
783 template <typename T, typename Storage>
784 bool
786  const OrthogPoly<T,Storage>& b)
787 {
788  return !(a == b);
789 }
790 
791 template <typename T, typename Storage>
792 bool
794  const OrthogPoly<T,Storage>& b)
795 {
796  return !(a == b);
797 }
798 
799 template <typename T, typename Storage>
800 bool
802  const typename OrthogPoly<T,Storage>::value_type& b)
803 {
804  return !(a == b);
805 }
806 
807 template <typename T, typename Storage>
808 bool
809 operator<=(const OrthogPoly<T,Storage>& a,
810  const OrthogPoly<T,Storage>& b)
811 {
812  return a.two_norm() <= b.two_norm();
813 }
814 
815 template <typename T, typename Storage>
816 bool
818  const OrthogPoly<T,Storage>& b)
819 {
820  return a <= b.two_norm();
821 }
822 
823 template <typename T, typename Storage>
824 bool
825 operator<=(const OrthogPoly<T,Storage>& a,
826  const typename OrthogPoly<T,Storage>::value_type& b)
827 {
828  return a.two_norm() <= b;
829 }
830 
831 template <typename T, typename Storage>
832 bool
834  const OrthogPoly<T,Storage>& b)
835 {
836  return a.two_norm() >= b.two_norm();
837 }
838 
839 template <typename T, typename Storage>
840 bool
842  const OrthogPoly<T,Storage>& b)
843 {
844  return a >= b.two_norm();
845 }
846 
847 template <typename T, typename Storage>
848 bool
850  const typename OrthogPoly<T,Storage>::value_type& b)
851 {
852  return a.two_norm() >= b;
853 }
854 
855 template <typename T, typename Storage>
856 bool
857 operator<(const OrthogPoly<T,Storage>& a,
858  const OrthogPoly<T,Storage>& b)
859 {
860  return a.two_norm() < b.two_norm();
861 }
862 
863 template <typename T, typename Storage>
864 bool
866  const OrthogPoly<T,Storage>& b)
867 {
868  return a < b.two_norm();
869 }
870 
871 template <typename T, typename Storage>
872 bool
873 operator<(const OrthogPoly<T,Storage>& a,
874  const typename OrthogPoly<T,Storage>::value_type& b)
875 {
876  return a.two_norm() < b;
877 }
878 
879 template <typename T, typename Storage>
880 bool
882  const OrthogPoly<T,Storage>& b)
883 {
884  return a.two_norm() > b.two_norm();
885 }
886 
887 template <typename T, typename Storage>
888 bool
890  const OrthogPoly<T,Storage>& b)
891 {
892  return a > b.two_norm();
893 }
894 
895 template <typename T, typename Storage>
896 bool
898  const typename OrthogPoly<T,Storage>::value_type& b)
899 {
900  return a.two_norm() > b;
901 }
902 
903 template <typename T, typename Storage>
905  bool is_zero = true;
906  for (int i=0; i<x.size(); i++)
907  is_zero = is_zero && (x.coeff(i) == 0.0);
908  return !is_zero;
909 }
910 
911 template <typename T, typename Storage>
912 inline bool
914 {
915  return toBool(x1) && toBool(x2);
916 }
917 
918 template <typename T, typename Storage>
919 inline bool
921  const OrthogPoly<T,Storage>& x2)
922 {
923  return a && toBool(x2);
924 }
925 
926 template <typename T, typename Storage>
927 inline bool
929  const typename OrthogPoly<T,Storage>::value_type& b)
930 {
931  return toBool(x1) && b;
932 }
933 
934 template <typename T, typename Storage>
935 inline bool
937 {
938  return toBool(x1) || toBool(x2);
939 }
940 
941 template <typename T, typename Storage>
942 inline bool
944  const OrthogPoly<T,Storage>& x2)
945 {
946  return a || toBool(x2);
947 }
948 
949 template <typename T, typename Storage>
950 inline bool
952  const typename OrthogPoly<T,Storage>::value_type& b)
953 {
954  return toBool(x1) || b;
955 }
956 
957 template <typename T, typename Storage>
958 std::ostream&
959 operator << (std::ostream& os, const OrthogPoly<T,Storage>& a)
960 {
962 
963  os << "[ ";
964 
965  for (ordinal_type i=0; i<a.size(); i++) {
966  os << a.coeff(i) << " ";
967  }
968 
969  os << "]\n";
970  return os;
971 }
972 
973 template <typename T, typename Storage>
974 std::istream&
975 operator >> (std::istream& is, OrthogPoly<T,Storage>& a)
976 {
978 
979  // Read in the opening "["
980  char bracket;
981  is >> bracket;
982 
983  for (ordinal_type i=0; i<a.size(); i++) {
984  is >> a.fastAccessCoeff(i);
985  }
986 
987  // Read in the closing "]"
988 
989  is >> bracket;
990  return is;
991 }
992 
993 } // namespace PCE
994 } // namespace Sacado
OrthogPoly< T, Storage > exp(const OrthogPoly< T, Storage > &a)
bool operator||(const OrthogPoly< T, Storage > &x1, const OrthogPoly< T, Storage > &x2)
OrthogPoly< T, Storage > log(const OrthogPoly< T, Storage > &a)
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
OrthogPoly< T, Storage > sin(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > sqrt(const OrthogPoly< T, Storage > &a)
bool operator>=(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > operator-(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
bool operator&&(const OrthogPoly< T, Storage > &x1, const OrthogPoly< T, Storage > &x2)
OrthogPoly< T, Storage > pow(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > atan(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > cbrt(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > acos(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > atanh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > cosh(const OrthogPoly< T, Storage > &a)
std::istream & operator>>(std::istream &is, OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > sinh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > tan(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > asin(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > operator+(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > operator/(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > max(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
bool operator!=(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
OrthogPoly< T, Storage > cos(const OrthogPoly< T, Storage > &a)
bool operator==(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > acosh(const OrthogPoly< T, Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > max(const typename PCE< Storage >::value_type &a, const PCE< Storage > &b)
OrthogPoly< T, Storage > min(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
Orthogonal polynomial expansion class for constant (size 1) expansions.
OrthogPoly< T, Storage > log10(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > abs(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > fabs(const OrthogPoly< T, Storage > &a)
Class to store coefficients of a projection onto an orthogonal polynomial basis.
OrthogPoly< T, Storage > asinh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > tanh(const OrthogPoly< T, Storage > &a)
bool operator>(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > operator*(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
int n
bool toBool(const OrthogPoly< T, Storage > &x)