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 // $Id$
2 // $Source$
3 // @HEADER
4 // ***********************************************************************
5 //
6 // Stokhos Package
7 // Copyright (2009) Sandia Corporation
8 //
9 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
10 // license for use of this work by or on behalf of the U.S. Government.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
40 //
41 // ***********************************************************************
42 // @HEADER
43 
44 #include "Teuchos_Assert.hpp"
46 #include "Teuchos_TimeMonitor.hpp"
47 #include "Teuchos_Tuple.hpp"
48 
49 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
55  const Teuchos::RCP<Teuchos::ParameterList>& params_) :
56  OrthogPolyExpansionBase<ordinal_type, value_type, node_type>(basis_, Cijk_, params_),
57  ps_op(ps_op_),
58  sz(ps_op->coeff_size()),
59  nqp(ps_op->point_size()),
60  avals(nqp),
61  bvals(nqp),
62  fvals(nqp)
63 {
65  if (params == Teuchos::null)
67  use_quad_for_times = params->get("Use Quadrature for Times", false);
68  use_quad_for_division = params->get("Use Quadrature for Division", true);
69 }
70 
71 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
72 template <typename FuncT>
73 void
75 unary_op(const FuncT& func,
78 {
79  ordinal_type pa = a.size();
80  ordinal_type pc;
81  if (a.size() == 1)
82  pc = 1;
83  else
84  pc = sz;
85  if (c.size() != pc)
86  c.resize(pc);
87 
88  if (pc == 1) {
89  c[0] = func(a[0]);
90  return;
91  }
92 
93  {
94 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
95  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Polynomial Evaluation");
96 #endif
97 
98  // Evaluate input
99  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
100  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
101 
102  }
103 
104  {
105 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
106  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Function Evaluation");
107 #endif
108 
109  // Evaluate function
110  for (ordinal_type qp=0; qp<nqp; qp++) {
111  fvals[qp] = func(avals[qp]);
112  }
113 
114  }
115 
116  {
117 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
118  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- Unary Polynomial Integration");
119 #endif
120 
121  // Integrate
122  SDV c_sdv(Teuchos::View, c.coeff(), pc);
123  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
124 
125  }
126 }
127 
128 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
129 template <typename FuncT>
130 void
132 binary_op(const FuncT& func,
136 {
137  ordinal_type pa = a.size();
138  ordinal_type pb = b.size();
139  ordinal_type pc;
140  if (pa == 1 && pb == 1)
141  pc = 1;
142  else
143  pc = sz;
144  if (c.size() != pc)
145  c.resize(pc);
146 
147  if (pc == 1) {
148  c[0] = func(a[0], b[0]);
149  return;
150  }
151 
152  {
153 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
154  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Polynomial Evaluation");
155 #endif
156 
157  // Evaluate input
158  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
159  SDV b_sdv(Teuchos::View, const_cast<value_type*>(b.coeff()), pb);
160  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
161  ps_op->transformPCE2QP(1.0, b_sdv, bvals, 0.0, false);
162 
163  }
164 
165  {
166 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
167  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Function Evaluation");
168 #endif
169 
170  // Evaluate function
171  for (ordinal_type qp=0; qp<nqp; qp++)
172  fvals[qp] = func(avals[qp], bvals[qp]);
173 
174  }
175 
176  {
177 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
178  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PP Binary Polynomial Integration");
179 #endif
180 
181  // Integrate
182  SDV c_sdv(Teuchos::View, c.coeff(), pc);
183  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
184 
185  }
186 }
187 
188 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
189 template <typename FuncT>
190 void
192 binary_op(const FuncT& func,
194  const value_type& a,
196 {
197  ordinal_type pb = b.size();
198  ordinal_type pc;
199  if (pb == 1)
200  pc = 1;
201  else
202  pc = sz;
203  if (c.size() != pc)
204  c.resize(pc);
205 
206  if (pc == 1) {
207  c[0] = func(a, b[0]);
208  return;
209  }
210 
211  {
212 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
213  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Polynomial Evaluation");
214 #endif
215 
216  // Evaluate input
217  SDV b_sdv(Teuchos::View, const_cast<value_type*>(b.coeff()), pb);
218  ps_op->transformPCE2QP(1.0, b_sdv, bvals, 0.0, false);
219 
220  }
221 
222  {
223 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
224  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Function Evaluation");
225 #endif
226 
227  // Evaluate function
228  for (ordinal_type qp=0; qp<nqp; qp++)
229  fvals[qp] = func(a, bvals[qp]);
230 
231  }
232 
233  {
234 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
235  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- CP Binary Polynomial Integration");
236 #endif
237 
238  // Integrate
239  SDV c_sdv(Teuchos::View, c.coeff(), pc);
240  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
241 
242  }
243 }
244 
245 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
246 template <typename FuncT>
247 void
249 binary_op(const FuncT& func,
252  const value_type& b)
253 {
254  ordinal_type pa = a.size();
255  ordinal_type pc;
256  if (pa == 1)
257  pc = 1;
258  else
259  pc = sz;
260  if (c.size() != pc)
261  c.resize(pc);
262 
263  if (pc == 1) {
264  c[0] = func(a[0], b);
265  return;
266  }
267 
268  {
269 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
270  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Polynomial Evaluation");
271 #endif
272 
273  // Evaluate input
274  SDV a_sdv(Teuchos::View, const_cast<value_type*>(a.coeff()), pa);
275  ps_op->transformPCE2QP(1.0, a_sdv, avals, 0.0, false);
276 
277  }
278 
279  {
280 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
281  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Function Evaluation");
282 #endif
283 
284  // Evaluate function
285  for (ordinal_type qp=0; qp<nqp; qp++)
286  fvals[qp] = func(avals[qp], b);
287 
288  }
289 
290  {
291 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
292  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- PC Binary Polynomial Integration");
293 #endif
294 
295  // Integrate
296  SDV c_sdv(Teuchos::View, c.coeff(), pc);
297  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
298 
299  }
300 }
301 
302 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
303 template <typename FuncT>
304 void
306 nary_op(const FuncT& func,
309 {
310  const int N = FuncT::N;
311  bool is_constant = true;
312  for (int i=0; i<N; i++) {
313  if (na[i]->size() > 1) {
314  is_constant = false;
315  break;
316  }
317  }
318  ordinal_type pc;
319  if (is_constant)
320  pc = 1;
321  else
322  pc = sz;
323  if (c.size() != pc)
324  c.resize(pc);
325 
326  if (pc == 1) {
327  value_type val[N];
328  for (int i=0; i<N; i++)
329  val[i] = (*na[i])[0];
330  c[0] = func(val);
331  return;
332  }
333 
334  if (N >= navals.size())
335  navals.resize(N+1);
336  if (navals[N].size() != N) {
337  navals[N].resize(N);
338  for (int i=0; i<N; i++)
339  navals[N][i].resize(nqp);
340  }
341 
342  {
343 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
344  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Polynomial Evaluation");
345 #endif
346 
347  // Evaluate input
348  for (int i=0; i<N; i++) {
349  SDV sdv(Teuchos::View, const_cast<value_type*>(na[i]->coeff()),
350  na[i]->size());
351  ps_op->transformPCE2QP(1.0, sdv, navals[N][i], 0.0, false);
352  }
353 
354  }
355 
356  {
357 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
358  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Function Evaluation");
359 #endif
360 
361  // Evaluate function
362  value_type val[N];
363  for (ordinal_type qp=0; qp<nqp; qp++) {
364  for (int i=0; i<N; i++)
365  val[i] = navals[N][i][qp];
366  fvals[qp] = func(val);
367  }
368 
369  }
370 
371  {
372 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
373  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::PSExp -- N(" << N << ")-ary Polynomial Integration");
374 #endif
375 
376  // Integrate
377  SDV c_sdv(Teuchos::View, c.coeff(), pc);
378  ps_op->transformQP2PCE(1.0, fvals, c_sdv, 0.0, false);
379 
380  }
381 }
382 
383 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
384 void
388  const value_type& x)
389 {
391 }
392 
393 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
394 void
398  const value_type& x)
399 {
401 }
402 
403 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
404 void
409 {
410  if (use_quad_for_times)
411  binary_op(times_quad_func(), c, c, x);
412  else
414 }
415 
416 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
417 void
422 {
423 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
424  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divideEqual(OPA)");
425 #endif
426  if (x.size() == 1) {
427  ordinal_type p = c.size();
428  value_type* cc = c.coeff();
429  const value_type* xc = x.coeff();
430  for (ordinal_type i=0; i<p; i++)
431  cc[i] /= xc[0];
432  }
433  else {
434  if (use_quad_for_division)
435  binary_op(div_quad_func(), c, c, x);
436  else
438  }
439 }
440 
441 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
442 void
447 {
448  if (use_quad_for_times)
449  binary_op(times_quad_func(), c, a, b);
450  else
452 }
453 
454 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
455 void
458  const value_type& a,
460 {
462 }
463 
464 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
465 void
469  const value_type& b)
470 {
472 }
473 
474 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
475 void
480 {
481 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
482  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divide(OPA,OPA)");
483 #endif
484  if (b.size() == 1) {
485  ordinal_type pc = a.size();
486  if (c.size() != pc)
487  c.resize(pc);
488 
489  const value_type* ca = a.coeff();
490  const value_type* cb = b.coeff();
491  value_type* cc = c.coeff();
492 
493  for (ordinal_type i=0; i<pc; i++)
494  cc[i] = ca[i]/cb[0];
495  }
496  else {
497  if (use_quad_for_division)
498  binary_op(div_quad_func(), c, a, b);
499  else
501  }
502 }
503 
504 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
505 void
508  const value_type& a,
510 {
511  if (use_quad_for_division)
512  binary_op(div_quad_func(), c, a, b);
513  else
515 }
516 
517 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
518 void
522  const value_type& b)
523 {
525 }
526 
527 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
528 void
532 {
533  unary_op(exp_quad_func(), c, a);
534 }
535 
536 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
537 void
541 {
542  unary_op(log_quad_func(), c, a);
543 }
544 
545 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
546 void
550 {
551  unary_op(log10_quad_func(), c, a);
552 }
553 
554 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
555 void
559 {
560  unary_op(sqrt_quad_func(), c, a);
561 }
562 
563 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
564 void
568 {
569  unary_op(cbrt_quad_func(), c, a);
570 }
571 
572 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
573 void
578 {
579  binary_op(pow_quad_func(), c, a, b);
580 }
581 
582 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
583 void
586  const value_type& a,
588 {
589  binary_op(pow_quad_func(), c, a, b);
590 }
591 
592 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
593 void
597  const value_type& b)
598 {
599  binary_op(pow_quad_func(), c, a, b);
600 }
601 
602 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
603 void
607 {
608  unary_op(sin_quad_func(), s, a);
609 }
610 
611 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
612 void
616 {
617  unary_op(cos_quad_func(), c, a);
618 }
619 
620 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
621 void
625 {
626  unary_op(tan_quad_func(), t, a);
627 }
628 
629 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
630 void
634 {
635  unary_op(sinh_quad_func(), s, a);
636 }
637 
638 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
639 void
643 {
644  unary_op(cosh_quad_func(), c, a);
645 }
646 
647 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
648 void
652 {
653  unary_op(tanh_quad_func(), t, a);
654 }
655 
656 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
657 void
661 {
662  unary_op(acos_quad_func(), c, a);
663 }
664 
665 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
666 void
670 {
671  unary_op(asin_quad_func(), c, a);
672 }
673 
674 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
675 void
679 {
680  unary_op(atan_quad_func(), c, a);
681 }
682 
683 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
684 void
689 {
690  binary_op(atan2_quad_func(), c, a, b);
691 }
692 
693 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
694 void
697  const value_type& a,
699 {
700  binary_op(atan2_quad_func(), c, a, b);
701 }
702 
703 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
704 void
708  const value_type& b)
709 {
710  binary_op(atan2_quad_func(), c, a, b);
711 }
712 
713 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
714 void
718 {
719  unary_op(acosh_quad_func(), c, a);
720 }
721 
722 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
723 void
727 {
728  unary_op(asinh_quad_func(), c, a);
729 }
730 
731 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
732 void
736 {
737  unary_op(atanh_quad_func(), c, a);
738 }
739 
740 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
741 template <typename ExprT1, typename ExprT2>
744 compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
745 {
746  ordinal_type pa = a.size();
747  ordinal_type pb = b.size();
748 
749  if (pa > 1 && pb > 1) {
750  ordinal_type k_lim = pa;
751  ordinal_type j_lim = pb;
752  if (pb < pa) {
753  k_lim = pb;
754  j_lim = pa;
755  }
756  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
757 #ifdef STOKHOS_DEBUG
758  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
759  "Stokhos::PseudoSpectralOrthogPolyExpansion::compute_times_coeff()"
760  << ": Index " << i << " is out of range [0,"
761  << this->Cijk->num_i() << ")!");
762 #endif
763  value_type cc = value_type(0);
764  value_type aa, bb, cijk;
765  ordinal_type j, k;
766  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
767  k_it != this->Cijk->k_end(i_it); ++k_it) {
768  k = index(k_it);
769  if (k < k_lim) {
770  if (pa < pb) {
771  if (k == 0)
772  aa = a.val();
773  else
774  aa = a.higher_order_coeff(k);
775  }
776  else {
777  if (k == 0)
778  aa = b.val();
779  else
780  aa = b.higher_order_coeff(k);
781  }
782  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
783  j_it != this->Cijk->j_end(k_it); ++j_it) {
784  j = index(j_it);
785  cijk = value(j_it);
786  if (j < j_lim) {
787  if (pa < pb) {
788  if (j == 0)
789  bb = b.val();
790  else
791  bb = b.higher_order_coeff(j);
792  }
793  else {
794  if (j == 0)
795  bb = a.val();
796  else
797  bb = a.higher_order_coeff(j);
798  }
799  cc += cijk*aa*bb;
800  }
801  }
802  }
803  }
804  return cc / this->basis->norm_squared(i);
805  }
806  else if (i == 0)
807  return a.val() * b.val();
808  else if (pa > 1) {
809  return a.higher_order_coeff(i)*b.val();
810  }
811  else {
812  return a.val()*b.higher_order_coeff(i);
813  }
814 }
815 
816 template <typename ordinal_type, typename value_type, typename point_compare_type, typename node_type>
817 template <typename ExprT1, typename ExprT2>
820 fast_compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
821 {
822  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
823 #ifdef STOKHOS_DEBUG
824  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
825  "Stokhos::PseudoSpectralOrthogPolyExpansion::fast_ompute_times_coeff()"
826  << ": Index " << i << " is out of range [0,"
827  << this->Cijk->num_i() << ")!");
828 #endif
829  value_type cc = value_type(0);
830  value_type aa, bb, cijk;
831  ordinal_type j, k;
832  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
833  k_it != this->Cijk->k_end(i_it); ++k_it) {
834  k = index(k_it);
835  if (k == 0)
836  aa = a.val();
837  else
838  aa = a.higher_order_coeff(k);
839  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
840  j_it != this->Cijk->j_end(k_it); ++j_it) {
841  j = index(j_it);
842  cijk = value(j_it);
843  if (j == 0)
844  bb = b.val();
845  else
846  bb = b.higher_order_coeff(j);
847  cc += cijk*aa*bb;
848  }
849  }
850 
851  return cc / this->basis->norm_squared(i);
852 }
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)