Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_QuadOrthogPolyExpansionImp.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 node_type>
21  const Teuchos::RCP<Teuchos::ParameterList>& params_) :
22  OrthogPolyExpansionBase<ordinal_type, value_type, node_type>(basis_, Cijk_, params_),
23  quad(quad_),
24  sz(this->basis->size()),
25  blas(),
26  quad_points(quad->getQuadPoints()),
27  quad_weights(quad->getQuadWeights()),
28  quad_values(quad->getBasisAtQuadPoints()),
29  norms(this->basis->norm_squared()),
30  nqp(quad_points.size()),
31  avals(nqp),
32  bvals(nqp),
33  fvals(nqp),
34  qv(nqp*sz),
35  sqv(nqp*sz)
36 {
37  for (ordinal_type qp=0; qp<nqp; qp++)
38  for (ordinal_type i=0; i<sz; i++) {
39  qv[qp*sz+i] = quad_values[qp][i];
40  sqv[qp*sz+i] = quad_values[qp][i]/norms[i];
41  }
42 
44  if (params == Teuchos::null)
46  use_quad_for_times = params->get("Use Quadrature for Times", false);
47  use_quad_for_division = params->get("Use Quadrature for Division", true);
48 }
49 
50 template <typename ordinal_type, typename value_type, typename node_type>
51 template <typename FuncT>
52 void
54 unary_op(const FuncT& func,
57 {
58  ordinal_type pa = a.size();
59  ordinal_type pc;
60  if (a.size() == 1)
61  pc = 1;
62  else
63  pc = sz;
64  if (c.size() != pc)
65  c.resize(pc);
66 
67  if (pc == 1) {
68  c[0] = func(a[0]);
69  return;
70  }
71 
72  {
73 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
74  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Polynomial Evaluation");
75 #endif
76 
77  // Evaluate input
78  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
79  &avals[0], 1);
80 
81  }
82 
83  {
84 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
85  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Function Evaluation");
86 #endif
87 
88  // Evaluate function
89  for (ordinal_type qp=0; qp<nqp; qp++) {
90  if (quad_weights[qp] != value_type(0))
91  fvals[qp] = func(avals[qp])*quad_weights[qp];
92  else
93  fvals[qp] = value_type(0);
94  }
95 
96  }
97 
98  {
99 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
100  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Polynomial Integration");
101 #endif
102 
103  // Integrate
104  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
105  c.coeff(), 1);
106 
107  }
108 }
109 
110 template <typename ordinal_type, typename value_type, typename node_type>
111 template <typename FuncT>
112 void
114 binary_op(const FuncT& func,
118 {
119  ordinal_type pa = a.size();
120  ordinal_type pb = b.size();
121  ordinal_type pc;
122  if (pa == 1 && pb == 1)
123  pc = 1;
124  else
125  pc = sz;
126  if (c.size() != pc)
127  c.resize(pc);
128 
129  if (pc == 1) {
130  c[0] = func(a[0], b[0]);
131  return;
132  }
133 
134  {
135 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
136  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Polynomial Evaluation");
137 #endif
138 
139  // Evaluate input
140  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
141  &avals[0], 1);
142  blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.coeff(), 1, 0.0,
143  &bvals[0], 1);
144 
145  }
146 
147  {
148 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
149  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Function Evaluation");
150 #endif
151 
152  // Evaluate function
153  for (ordinal_type qp=0; qp<nqp; qp++)
154  if (quad_weights[qp] != value_type(0))
155  fvals[qp] = func(avals[qp], bvals[qp])*quad_weights[qp];
156  else
157  fvals[qp] = value_type(0);
158 
159  }
160 
161  {
162 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
163  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Polynomial Integration");
164 #endif
165 
166  // Integrate
167  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
168  c.coeff(), 1);
169 
170  }
171 }
172 
173 template <typename ordinal_type, typename value_type, typename node_type>
174 template <typename FuncT>
175 void
177 binary_op(const FuncT& func,
179  const value_type& a,
181 {
182  ordinal_type pb = b.size();
183  ordinal_type pc;
184  if (pb == 1)
185  pc = 1;
186  else
187  pc = sz;
188  if (c.size() != pc)
189  c.resize(pc);
190 
191  if (pc == 1) {
192  c[0] = func(a, b[0]);
193  return;
194  }
195 
196  {
197 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
198  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Polynomial Evaluation");
199 #endif
200 
201  // Evaluate input
202  blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.coeff(), 1, 0.0,
203  &bvals[0], 1);
204 
205  }
206 
207  {
208 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
209  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Function Evaluation");
210 #endif
211 
212  // Evaluate function
213  for (ordinal_type qp=0; qp<nqp; qp++)
214  if (quad_weights[qp] != value_type(0))
215  fvals[qp] = func(a, bvals[qp])*quad_weights[qp];
216  else
217  fvals[qp] = value_type(0);
218 
219  }
220 
221  {
222 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
223  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Polynomial Integration");
224 #endif
225 
226  // Integrate
227  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
228  c.coeff(), 1);
229 
230  }
231 }
232 
233 template <typename ordinal_type, typename value_type, typename node_type>
234 template <typename FuncT>
235 void
237 binary_op(const FuncT& func,
240  const value_type& b)
241 {
242  ordinal_type pa = a.size();
243  ordinal_type pc;
244  if (pa == 1)
245  pc = 1;
246  else
247  pc = sz;
248  if (c.size() != pc)
249  c.resize(pc);
250 
251  if (pc == 1) {
252  c[0] = func(a[0], b);
253  return;
254  }
255 
256  {
257 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
258  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Polynomial Evaluation");
259 #endif
260 
261  // Evaluate input
262  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
263  &avals[0], 1);
264 
265  }
266 
267  {
268 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
269  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Function Evaluation");
270 #endif
271 
272  // Evaluate function
273  for (ordinal_type qp=0; qp<nqp; qp++)
274  if (quad_weights[qp] != value_type(0))
275  fvals[qp] = func(avals[qp], b)*quad_weights[qp];
276  else
277  fvals[qp] = value_type(0);
278 
279  }
280 
281  {
282 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
283  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Polynomial Integration");
284 #endif
285 
286  // Integrate
287  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
288  c.coeff(), 1);
289 
290  }
291 }
292 
293 template <typename ordinal_type, typename value_type, typename node_type>
294 template <typename FuncT>
295 void
297 nary_op(const FuncT& func,
300 {
301  const int N = FuncT::N;
302  bool is_constant = true;
303  for (int i=0; i<N; i++) {
304  if (na[i]->size() > 1) {
305  is_constant = false;
306  break;
307  }
308  }
309  ordinal_type pc;
310  if (is_constant)
311  pc = 1;
312  else
313  pc = sz;
314  if (c.size() != pc)
315  c.resize(pc);
316 
317  if (pc == 1) {
318  value_type val[N];
319  for (int i=0; i<N; i++)
320  val[i] = (*na[i])[0];
321  c[0] = func(val);
322  return;
323  }
324 
325  if (N >= navals.size())
326  navals.resize(N+1);
327  if (navals[N].size() != N) {
328  navals[N].resize(N);
329  for (int i=0; i<N; i++)
330  navals[N][i].resize(nqp);
331  }
332 
333  {
334 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
335  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Polynomial Evaluation");
336 #endif
337 
338  // Evaluate input
339  for (int i=0; i<N; i++) {
340  //navals[i].resize(nqp);
341  ordinal_type pa = na[i]->size();
342  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, na[i]->coeff(), 1, 0.0,
343  &navals[N][i][0], 1);
344  }
345 
346  }
347 
348  {
349 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
350  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Function Evaluation");
351 #endif
352 
353  // Evaluate function
354  value_type val[N];
355  for (ordinal_type qp=0; qp<nqp; qp++) {
356  if (quad_weights[qp] != value_type(0)) {
357  for (int i=0; i<N; i++)
358  val[i] = navals[N][i][qp];
359  fvals[qp] = func(val)*quad_weights[qp];
360  }
361  else
362  fvals[qp] = value_type(0);
363  }
364 
365  }
366 
367  {
368 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
369  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Polynomial Integration");
370 #endif
371 
372  // Integrate
373  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
374  c.coeff(), 1);
375 
376  }
377 }
378 
379 template <typename ordinal_type, typename value_type, typename node_type>
380 void
384  const value_type& x)
385 {
387 }
388 
389 template <typename ordinal_type, typename value_type, typename node_type>
390 void
394  const value_type& x)
395 {
397 }
398 
399 template <typename ordinal_type, typename value_type, typename node_type>
400 void
405 {
406  if (use_quad_for_times)
407  binary_op(times_quad_func(), c, c, x);
408  else
410 }
411 
412 template <typename ordinal_type, typename value_type, typename node_type>
413 void
418 {
419 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
420  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divideEqual(OPA)");
421 #endif
422  if (x.size() == 1) {
423  ordinal_type p = c.size();
424  value_type* cc = c.coeff();
425  const value_type* xc = x.coeff();
426  for (ordinal_type i=0; i<p; i++)
427  cc[i] /= xc[0];
428  }
429  else {
430  if (use_quad_for_division)
431  binary_op(div_quad_func(), c, c, x);
432  else
434  }
435 }
436 
437 template <typename ordinal_type, typename value_type, typename node_type>
438 void
443 {
444  if (use_quad_for_times)
445  binary_op(times_quad_func(), c, a, b);
446  else
448 }
449 
450 template <typename ordinal_type, typename value_type, typename node_type>
451 void
454  const value_type& a,
456 {
458 }
459 
460 template <typename ordinal_type, typename value_type, typename node_type>
461 void
465  const value_type& b)
466 {
468 }
469 
470 template <typename ordinal_type, typename value_type, typename node_type>
471 void
476 {
477 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
478  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divide(OPA,OPA)");
479 #endif
480  if (b.size() == 1) {
481  ordinal_type pc = a.size();
482  if (c.size() != pc)
483  c.resize(pc);
484 
485  const value_type* ca = a.coeff();
486  const value_type* cb = b.coeff();
487  value_type* cc = c.coeff();
488 
489  for (ordinal_type i=0; i<pc; i++)
490  cc[i] = ca[i]/cb[0];
491  }
492  else {
493  if (use_quad_for_division)
494  binary_op(div_quad_func(), c, a, b);
495  else
497  }
498 }
499 
500 template <typename ordinal_type, typename value_type, typename node_type>
501 void
504  const value_type& a,
506 {
507  if (use_quad_for_division)
508  binary_op(div_quad_func(), c, a, b);
509  else
511 }
512 
513 template <typename ordinal_type, typename value_type, typename node_type>
514 void
518  const value_type& b)
519 {
521 }
522 
523 template <typename ordinal_type, typename value_type, typename node_type>
524 void
528 {
529  unary_op(exp_quad_func(), c, a);
530 }
531 
532 template <typename ordinal_type, typename value_type, typename node_type>
533 void
537 {
538  unary_op(log_quad_func(), c, a);
539 }
540 
541 template <typename ordinal_type, typename value_type, typename node_type>
542 void
546 {
547  unary_op(log10_quad_func(), c, a);
548 }
549 
550 template <typename ordinal_type, typename value_type, typename node_type>
551 void
555 {
556  unary_op(sqrt_quad_func(), c, a);
557 }
558 
559 template <typename ordinal_type, typename value_type, typename node_type>
560 void
564 {
565  unary_op(cbrt_quad_func(), c, a);
566 }
567 
568 template <typename ordinal_type, typename value_type, typename node_type>
569 void
574 {
575  binary_op(pow_quad_func(), c, a, b);
576 }
577 
578 template <typename ordinal_type, typename value_type, typename node_type>
579 void
582  const value_type& a,
584 {
585  binary_op(pow_quad_func(), c, a, b);
586 }
587 
588 template <typename ordinal_type, typename value_type, typename node_type>
589 void
593  const value_type& b)
594 {
595  binary_op(pow_quad_func(), c, a, b);
596 }
597 
598 template <typename ordinal_type, typename value_type, typename node_type>
599 void
603 {
604  unary_op(sin_quad_func(), s, a);
605 }
606 
607 template <typename ordinal_type, typename value_type, typename node_type>
608 void
612 {
613  unary_op(cos_quad_func(), c, a);
614 }
615 
616 template <typename ordinal_type, typename value_type, typename node_type>
617 void
621 {
622  unary_op(tan_quad_func(), t, a);
623 }
624 
625 template <typename ordinal_type, typename value_type, typename node_type>
626 void
630 {
631  unary_op(sinh_quad_func(), s, a);
632 }
633 
634 template <typename ordinal_type, typename value_type, typename node_type>
635 void
639 {
640  unary_op(cosh_quad_func(), c, a);
641 }
642 
643 template <typename ordinal_type, typename value_type, typename node_type>
644 void
648 {
649  unary_op(tanh_quad_func(), t, a);
650 }
651 
652 template <typename ordinal_type, typename value_type, typename node_type>
653 void
657 {
658  unary_op(acos_quad_func(), c, a);
659 }
660 
661 template <typename ordinal_type, typename value_type, typename node_type>
662 void
666 {
667  unary_op(asin_quad_func(), c, a);
668 }
669 
670 template <typename ordinal_type, typename value_type, typename node_type>
671 void
675 {
676  unary_op(atan_quad_func(), c, a);
677 }
678 
679 template <typename ordinal_type, typename value_type, typename node_type>
680 void
685 {
686  binary_op(atan2_quad_func(), c, a, b);
687 }
688 
689 template <typename ordinal_type, typename value_type, typename node_type>
690 void
693  const value_type& a,
695 {
696  binary_op(atan2_quad_func(), c, a, b);
697 }
698 
699 template <typename ordinal_type, typename value_type, typename node_type>
700 void
704  const value_type& b)
705 {
706  binary_op(atan2_quad_func(), c, a, b);
707 }
708 
709 template <typename ordinal_type, typename value_type, typename node_type>
710 void
714 {
715  unary_op(acosh_quad_func(), c, a);
716 }
717 
718 template <typename ordinal_type, typename value_type, typename node_type>
719 void
723 {
724  unary_op(asinh_quad_func(), c, a);
725 }
726 
727 template <typename ordinal_type, typename value_type, typename node_type>
728 void
732 {
733  unary_op(atanh_quad_func(), c, a);
734 }
735 
736 template <typename ordinal_type, typename value_type, typename node_type>
737 template <typename ExprT1, typename ExprT2>
740 compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
741 {
742  ordinal_type pa = a.size();
743  ordinal_type pb = b.size();
744 
745  if (pa > 1 && pb > 1) {
746  ordinal_type k_lim = pa;
747  ordinal_type j_lim = pb;
748  if (pb < pa) {
749  k_lim = pb;
750  j_lim = pa;
751  }
752  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
753 #ifdef STOKHOS_DEBUG
754  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
755  "Stokhos::QuadOrthogPolyExpansion::compute_times_coeff()"
756  << ": Index " << i << " is out of range [0,"
757  << this->Cijk->num_i() << ")!");
758 #endif
759  value_type cc = value_type(0);
760  value_type aa, bb, cijk;
761  ordinal_type j, k;
762  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
763  k_it != this->Cijk->k_end(i_it); ++k_it) {
764  k = index(k_it);
765  if (k < k_lim) {
766  if (pa < pb) {
767  if (k == 0)
768  aa = a.val();
769  else
770  aa = a.higher_order_coeff(k);
771  }
772  else {
773  if (k == 0)
774  aa = b.val();
775  else
776  aa = b.higher_order_coeff(k);
777  }
778  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
779  j_it != this->Cijk->j_end(k_it); ++j_it) {
780  j = index(j_it);
781  cijk = value(j_it);
782  if (j < j_lim) {
783  if (pa < pb) {
784  if (j == 0)
785  bb = b.val();
786  else
787  bb = b.higher_order_coeff(j);
788  }
789  else {
790  if (j == 0)
791  bb = a.val();
792  else
793  bb = a.higher_order_coeff(j);
794  }
795  cc += cijk*aa*bb;
796  }
797  }
798  }
799  }
800  return cc / norms[i];
801  }
802  else if (i == 0)
803  return a.val() * b.val();
804  else if (pa > 1) {
805  return a.higher_order_coeff(i)*b.val();
806  }
807  else {
808  return a.val()*b.higher_order_coeff(i);
809  }
810 }
811 
812 template <typename ordinal_type, typename value_type, typename node_type>
813 template <typename ExprT1, typename ExprT2>
816 fast_compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
817 {
818  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
819 #ifdef STOKHOS_DEBUG
820  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
821  "Stokhos::QuadOrthogPolyExpansion::fast_ompute_times_coeff()"
822  << ": Index " << i << " is out of range [0,"
823  << this->Cijk->num_i() << ")!");
824 #endif
825  value_type cc = value_type(0);
826  value_type aa, bb, cijk;
827  ordinal_type j, k;
828  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
829  k_it != this->Cijk->k_end(i_it); ++k_it) {
830  k = index(k_it);
831  if (k == 0)
832  aa = a.val();
833  else
834  aa = a.higher_order_coeff(k);
835  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
836  j_it != this->Cijk->j_end(k_it); ++j_it) {
837  j = index(j_it);
838  cijk = value(j_it);
839  if (j == 0)
840  bb = b.val();
841  else
842  bb = b.higher_order_coeff(j);
843  cc += cijk*aa*bb;
844  }
845  }
846 
847  return cc / norms[i];
848 }
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
const Teuchos::Array< value_type > & norms
Norms of basis vectors.
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.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
Teuchos::RCP< Teuchos::ParameterList > params
Parameter list.
void sin(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 sqrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved)
ordinal_type nqp
Number of Quad points.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
void acos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Teuchos::Array< Teuchos::Array< value_type > > & quad_values
Values of basis at Quad points.
QuadOrthogPolyExpansion(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 Quadrature< ordinal_type, value_type > > &quad, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
T & get(ParameterList &l, const std::string &name)
void atanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void log10(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cos(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 asin(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)
void atan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
pointer coeff()
Return coefficient array.
Bi-directional iterator for traversing a sparse array.
Kokkos::Serial node_type
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)
value_type compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
Abstract base class for multivariate orthogonal polynomials.
void sinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
value_type fast_compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
Abstract base class for quadrature methods.
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)
void nary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > **a)
bool use_quad_for_times
Use quadrature for times functions.
Teuchos::Array< value_type > sqv
Quad values scaled by norms.
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 log(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
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)
bool use_quad_for_division
Use quadrature for division functions.
Class to store coefficients of a projection onto an orthogonal polynomial basis.
void timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
expr val()
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void cbrt(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 tanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void tan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
ordinal_type size() const
Return size.
void cosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void asinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::Array< value_type > qv
Reshaped quad values into 1D array.
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.
void acosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)