Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_PseudoSpectralOrthogPolyExpansionImp.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 "Teuchos_Assert.hpp"
12 #include "Teuchos_TimeMonitor.hpp"
13 #include "Teuchos_Tuple.hpp"
14 
15 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
21  const Teuchos::RCP<Teuchos::ParameterList>& params_) :
22  OrthogPolyExpansionBase<ordinal_type, value_type, node_type>(basis_, Cijk_, params_),
23  ps_op(ps_op_),
24  sz(ps_op->coeff_size()),
25  nqp(ps_op->point_size()),
26  avals(nqp),
27  bvals(nqp),
28  fvals(nqp)
29 {
31  if (params == Teuchos::null)
33  use_quad_for_times = params->get("Use Quadrature for Times", false);
34  use_quad_for_division = params->get("Use Quadrature for Division", true);
35 }
36 
37 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
38 template <typename FuncT>
39 void
41 unary_op(const FuncT& func,
44 {
45  ordinal_type pa = a.size();
46  ordinal_type pc;
47  if (a.size() == 1)
48  pc = 1;
49  else
50  pc = sz;
51  if (c.size() != pc)
52  c.resize(pc);
53 
54  if (pc == 1) {
55  c[0] = func(a[0]);
56  return;
57  }
58 
59  {
60 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
61  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Polynomial Evaluation");
62 #endif
63 
64  // Evaluate input
65  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
66  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
67 
68  }
69 
70  {
71 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
72  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Function Evaluation");
73 #endif
74 
75  // Evaluate function
76  for (ordinal_type qp=0; qp<nqp; qp++) {
77  fvals[qp] = func(avals[qp]);
78  }
79 
80  }
81 
82  {
83 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
84  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Polynomial Integration");
85 #endif
86 
87  // Integrate
88  SDV c_sdv(Teuchos::View, c.coeff(), pc);
89  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
90 
91  }
92 }
93 
94 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
95 template <typename FuncT>
96 void
98 binary_op(const FuncT& func,
102 {
103  ordinal_type pa = a.size();
104  ordinal_type pb = b.size();
105  ordinal_type pc;
106  if (pa == 1 && pb == 1)
107  pc = 1;
108  else
109  pc = sz;
110  if (c.size() != pc)
111  c.resize(pc);
112 
113  if (pc == 1) {
114  c[0] = func(a[0], b[0]);
115  return;
116  }
117 
118  {
119 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
120  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Polynomial Evaluation");
121 #endif
122 
123  // Evaluate input
124  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
125  SDV b_sdv(Teuchos::View, const_cast<value_type*>(b.coeff()), pb);
126  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
127  ps_op->transformPCE2QP(1.0, b_sdv, bvals, 0.0, false);
128 
129  }
130 
131  {
132 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
133  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Function Evaluation");
134 #endif
135 
136  // Evaluate function
137  for (ordinal_type qp=0; qp<nqp; qp++)
138  fvals[qp] = func(avals[qp], bvals[qp]);
139 
140  }
141 
142  {
143 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
144  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Polynomial Integration");
145 #endif
146 
147  // Integrate
148  SDV c_sdv(Teuchos::View, c.coeff(), pc);
149  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
150 
151  }
152 }
153 
154 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
155 template <typename FuncT>
156 void
158 binary_op(const FuncT& func,
160  const value_type& a,
162 {
163  ordinal_type pb = b.size();
164  ordinal_type pc;
165  if (pb == 1)
166  pc = 1;
167  else
168  pc = sz;
169  if (c.size() != pc)
170  c.resize(pc);
171 
172  if (pc == 1) {
173  c[0] = func(a, b[0]);
174  return;
175  }
176 
177  {
178 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
179  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Polynomial Evaluation");
180 #endif
181 
182  // Evaluate input
183  SDV b_sdv(Teuchos::View, const_cast<value_type*>(b.coeff()), pb);
184  ps_op->transformPCE2QP(1.0, b_sdv, bvals, 0.0, false);
185 
186  }
187 
188  {
189 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
190  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Function Evaluation");
191 #endif
192 
193  // Evaluate function
194  for (ordinal_type qp=0; qp<nqp; qp++)
195  fvals[qp] = func(a, bvals[qp]);
196 
197  }
198 
199  {
200 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
201  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Polynomial Integration");
202 #endif
203 
204  // Integrate
205  SDV c_sdv(Teuchos::View, c.coeff(), pc);
206  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
207 
208  }
209 }
210 
211 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
212 template <typename FuncT>
213 void
215 binary_op(const FuncT& func,
218  const value_type& b)
219 {
220  ordinal_type pa = a.size();
221  ordinal_type pc;
222  if (pa == 1)
223  pc = 1;
224  else
225  pc = sz;
226  if (c.size() != pc)
227  c.resize(pc);
228 
229  if (pc == 1) {
230  c[0] = func(a[0], b);
231  return;
232  }
233 
234  {
235 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
236  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Polynomial Evaluation");
237 #endif
238 
239  // Evaluate input
240  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
241  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
242 
243  }
244 
245  {
246 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
247  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Function Evaluation");
248 #endif
249 
250  // Evaluate function
251  for (ordinal_type qp=0; qp<nqp; qp++)
252  fvals[qp] = func(avals[qp], b);
253 
254  }
255 
256  {
257 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
258  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Polynomial Integration");
259 #endif
260 
261  // Integrate
262  SDV c_sdv(Teuchos::View, c.coeff(), pc);
263  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
264 
265  }
266 }
267 
268 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
269 template <typename FuncT>
270 void
272 nary_op(const FuncT& func,
275 {
276  const int N = FuncT::N;
277  bool is_constant = true;
278  for (int i=0; i<N; i++) {
279  if (na[i]->size() > 1) {
280  is_constant = false;
281  break;
282  }
283  }
284  ordinal_type pc;
285  if (is_constant)
286  pc = 1;
287  else
288  pc = sz;
289  if (c.size() != pc)
290  c.resize(pc);
291 
292  if (pc == 1) {
293  value_type val[N];
294  for (int i=0; i<N; i++)
295  val[i] = (*na[i])[0];
296  c[0] = func(val);
297  return;
298  }
299 
300  if (N >= navals.size())
301  navals.resize(N+1);
302  if (navals[N].size() != N) {
303  navals[N].resize(N);
304  for (int i=0; i<N; i++)
305  navals[N][i].resize(nqp);
306  }
307 
308  {
309 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
310  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Polynomial Evaluation");
311 #endif
312 
313  // Evaluate input
314  for (int i=0; i<N; i++) {
315  SDV sdv(Teuchos::View, const_cast<value_type*>(na[i]->coeff()),
316  na[i]->size());
317  ps_op->transformPCE2QP(1.0, sdv, navals[N][i], 0.0, false);
318  }
319 
320  }
321 
322  {
323 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
324  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Function Evaluation");
325 #endif
326 
327  // Evaluate function
328  value_type val[N];
329  for (ordinal_type qp=0; qp<nqp; qp++) {
330  for (int i=0; i<N; i++)
331  val[i] = navals[N][i][qp];
332  fvals[qp] = func(val);
333  }
334 
335  }
336 
337  {
338 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
339  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Polynomial Integration");
340 #endif
341 
342  // Integrate
343  SDV c_sdv(Teuchos::View, c.coeff(), pc);
344  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
345 
346  }
347 }
348 
349 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
350 void
354  const value_type& x)
355 {
357 }
358 
359 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
360 void
364  const value_type& x)
365 {
367 }
368 
369 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
370 void
375 {
376  if (use_quad_for_times)
377  binary_op(times_quad_func(), c, c, x);
378  else
380 }
381 
382 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
383 void
388 {
389 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
390  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divideEqual(OPA)");
391 #endif
392  if (x.size() == 1) {
393  ordinal_type p = c.size();
394  value_type* cc = c.coeff();
395  const value_type* xc = x.coeff();
396  for (ordinal_type i=0; i<p; i++)
397  cc[i] /= xc[0];
398  }
399  else {
400  if (use_quad_for_division)
401  binary_op(div_quad_func(), c, c, x);
402  else
404  }
405 }
406 
407 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
408 void
413 {
414  if (use_quad_for_times)
415  binary_op(times_quad_func(), c, a, b);
416  else
418 }
419 
420 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
421 void
424  const value_type& a,
426 {
428 }
429 
430 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
431 void
435  const value_type& b)
436 {
438 }
439 
440 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
441 void
446 {
447 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
448  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divide(OPA,OPA)");
449 #endif
450  if (b.size() == 1) {
451  ordinal_type pc = a.size();
452  if (c.size() != pc)
453  c.resize(pc);
454 
455  const value_type* ca = a.coeff();
456  const value_type* cb = b.coeff();
457  value_type* cc = c.coeff();
458 
459  for (ordinal_type i=0; i<pc; i++)
460  cc[i] = ca[i]/cb[0];
461  }
462  else {
463  if (use_quad_for_division)
464  binary_op(div_quad_func(), c, a, b);
465  else
467  }
468 }
469 
470 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
471 void
474  const value_type& a,
476 {
477  if (use_quad_for_division)
478  binary_op(div_quad_func(), c, a, b);
479  else
481 }
482 
483 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
484 void
488  const value_type& b)
489 {
491 }
492 
493 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
494 void
498 {
499  unary_op(exp_quad_func(), c, a);
500 }
501 
502 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
503 void
507 {
508  unary_op(log_quad_func(), c, a);
509 }
510 
511 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
512 void
516 {
517  unary_op(log10_quad_func(), c, a);
518 }
519 
520 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
521 void
525 {
526  unary_op(sqrt_quad_func(), c, a);
527 }
528 
529 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
530 void
534 {
535  unary_op(cbrt_quad_func(), c, a);
536 }
537 
538 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
539 void
544 {
545  binary_op(pow_quad_func(), c, a, b);
546 }
547 
548 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
549 void
552  const value_type& a,
554 {
555  binary_op(pow_quad_func(), c, a, b);
556 }
557 
558 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
559 void
563  const value_type& b)
564 {
565  binary_op(pow_quad_func(), c, a, b);
566 }
567 
568 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
569 void
573 {
574  unary_op(sin_quad_func(), s, a);
575 }
576 
577 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
578 void
582 {
583  unary_op(cos_quad_func(), c, a);
584 }
585 
586 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
587 void
591 {
592  unary_op(tan_quad_func(), t, a);
593 }
594 
595 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
596 void
600 {
601  unary_op(sinh_quad_func(), s, a);
602 }
603 
604 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
605 void
609 {
610  unary_op(cosh_quad_func(), c, a);
611 }
612 
613 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
614 void
618 {
619  unary_op(tanh_quad_func(), t, a);
620 }
621 
622 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
623 void
627 {
628  unary_op(acos_quad_func(), c, a);
629 }
630 
631 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
632 void
636 {
637  unary_op(asin_quad_func(), c, a);
638 }
639 
640 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
641 void
645 {
646  unary_op(atan_quad_func(), c, a);
647 }
648 
649 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
650 void
655 {
656  binary_op(atan2_quad_func(), c, a, b);
657 }
658 
659 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
660 void
663  const value_type& a,
665 {
666  binary_op(atan2_quad_func(), c, a, b);
667 }
668 
669 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
670 void
674  const value_type& b)
675 {
676  binary_op(atan2_quad_func(), c, a, b);
677 }
678 
679 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
680 void
684 {
685  unary_op(acosh_quad_func(), c, a);
686 }
687 
688 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
689 void
693 {
694  unary_op(asinh_quad_func(), c, a);
695 }
696 
697 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
698 void
702 {
703  unary_op(atanh_quad_func(), c, a);
704 }
705 
706 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
707 template <typename ExprT1, typename ExprT2>
710 compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
711 {
712  ordinal_type pa = a.size();
713  ordinal_type pb = b.size();
714 
715  if (pa > 1 && pb > 1) {
716  ordinal_type k_lim = pa;
717  ordinal_type j_lim = pb;
718  if (pb < pa) {
719  k_lim = pb;
720  j_lim = pa;
721  }
722  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
723 #ifdef STOKHOS_DEBUG
724  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
725  "Stokhos::PseudoSpectralOrthogPolyExpansion::compute_times_coeff()"
726  << ": Index " << i << " is out of range [0,"
727  << this->Cijk->num_i() << ")!");
728 #endif
729  value_type cc = value_type(0);
730  value_type aa, bb, cijk;
731  ordinal_type j, k;
732  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
733  k_it != this->Cijk->k_end(i_it); ++k_it) {
734  k = index(k_it);
735  if (k < k_lim) {
736  if (pa < pb) {
737  if (k == 0)
738  aa = a.val();
739  else
740  aa = a.higher_order_coeff(k);
741  }
742  else {
743  if (k == 0)
744  aa = b.val();
745  else
746  aa = b.higher_order_coeff(k);
747  }
748  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
749  j_it != this->Cijk->j_end(k_it); ++j_it) {
750  j = index(j_it);
751  cijk = value(j_it);
752  if (j < j_lim) {
753  if (pa < pb) {
754  if (j == 0)
755  bb = b.val();
756  else
757  bb = b.higher_order_coeff(j);
758  }
759  else {
760  if (j == 0)
761  bb = a.val();
762  else
763  bb = a.higher_order_coeff(j);
764  }
765  cc += cijk*aa*bb;
766  }
767  }
768  }
769  }
770  return cc / this->basis->norm_squared(i);
771  }
772  else if (i == 0)
773  return a.val() * b.val();
774  else if (pa > 1) {
775  return a.higher_order_coeff(i)*b.val();
776  }
777  else {
778  return a.val()*b.higher_order_coeff(i);
779  }
780 }
781 
782 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
783 template <typename ExprT1, typename ExprT2>
786 fast_compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
787 {
788  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
789 #ifdef STOKHOS_DEBUG
790  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
791  "Stokhos::PseudoSpectralOrthogPolyExpansion::fast_ompute_times_coeff()"
792  << ": Index " << i << " is out of range [0,"
793  << this->Cijk->num_i() << ")!");
794 #endif
795  value_type cc = value_type(0);
796  value_type aa, bb, cijk;
797  ordinal_type j, k;
798  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
799  k_it != this->Cijk->k_end(i_it); ++k_it) {
800  k = index(k_it);
801  if (k == 0)
802  aa = a.val();
803  else
804  aa = a.higher_order_coeff(k);
805  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
806  j_it != this->Cijk->j_end(k_it); ++j_it) {
807  j = index(j_it);
808  cijk = value(j_it);
809  if (j == 0)
810  bb = b.val();
811  else
812  bb = b.higher_order_coeff(j);
813  cc += cijk*aa*bb;
814  }
815  }
816 
817  return cc / this->basis->norm_squared(i);
818 }
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void log10(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
Teuchos::RCP< Teuchos::ParameterList > params
Parameter list.
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved)
void tanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
T & get(ParameterList &l, const std::string &name)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Base class for consolidating common expansion implementations.
void timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void nary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > **a)
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void acos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
pointer coeff()
Return coefficient array.
PseudoSpectralOrthogPolyExpansion(const Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > &basis, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::RCP< const PseudoSpectralOperator< ordinal_type, value_type, point_compare_type > > &ps_op, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor.
void sqrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void log(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Bi-directional iterator for traversing a sparse array.
Kokkos::Serial node_type
void atan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Abstract base class for multivariate orthogonal polynomials.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void binary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Nonlinear binary function.
value_type fast_compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void pow(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void atan2(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
bool use_quad_for_division
Use quadrature for division functions.
void tan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
KOKKOS_INLINE_FUNCTION bool is_constant(const T &x)
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void atanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< view_type >::value, typename CijkType< view_type >::type >::type cijk(const view_type &view)
void cos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Class to store coefficients of a projection onto an orthogonal polynomial basis.
expr val()
void cbrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
ordinal_type size() const
Return size.
void asinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void sinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void unary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Nonlinear unary function.
value_type compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
void asin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void acosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)