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 // $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 node_type>
55  const Teuchos::RCP<Teuchos::ParameterList>& params_) :
56  OrthogPolyExpansionBase<ordinal_type, value_type, node_type>(basis_, Cijk_, params_),
57  quad(quad_),
58  sz(this->basis->size()),
59  blas(),
60  quad_points(quad->getQuadPoints()),
61  quad_weights(quad->getQuadWeights()),
62  quad_values(quad->getBasisAtQuadPoints()),
63  norms(this->basis->norm_squared()),
64  nqp(quad_points.size()),
65  avals(nqp),
66  bvals(nqp),
67  fvals(nqp),
68  qv(nqp*sz),
69  sqv(nqp*sz)
70 {
71  for (ordinal_type qp=0; qp<nqp; qp++)
72  for (ordinal_type i=0; i<sz; i++) {
73  qv[qp*sz+i] = quad_values[qp][i];
74  sqv[qp*sz+i] = quad_values[qp][i]/norms[i];
75  }
76 
78  if (params == Teuchos::null)
80  use_quad_for_times = params->get("Use Quadrature for Times", false);
81  use_quad_for_division = params->get("Use Quadrature for Division", true);
82 }
83 
84 template <typename ordinal_type, typename value_type, typename node_type>
85 template <typename FuncT>
86 void
88 unary_op(const FuncT& func,
91 {
92  ordinal_type pa = a.size();
93  ordinal_type pc;
94  if (a.size() == 1)
95  pc = 1;
96  else
97  pc = sz;
98  if (c.size() != pc)
99  c.resize(pc);
100 
101  if (pc == 1) {
102  c[0] = func(a[0]);
103  return;
104  }
105 
106  {
107 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
108  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Polynomial Evaluation");
109 #endif
110 
111  // Evaluate input
112  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
113  &avals[0], 1);
114 
115  }
116 
117  {
118 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
119  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Function Evaluation");
120 #endif
121 
122  // Evaluate function
123  for (ordinal_type qp=0; qp<nqp; qp++) {
124  if (quad_weights[qp] != value_type(0))
125  fvals[qp] = func(avals[qp])*quad_weights[qp];
126  else
127  fvals[qp] = value_type(0);
128  }
129 
130  }
131 
132  {
133 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
134  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- Unary Polynomial Integration");
135 #endif
136 
137  // Integrate
138  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
139  c.coeff(), 1);
140 
141  }
142 }
143 
144 template <typename ordinal_type, typename value_type, typename node_type>
145 template <typename FuncT>
146 void
148 binary_op(const FuncT& func,
152 {
153  ordinal_type pa = a.size();
154  ordinal_type pb = b.size();
155  ordinal_type pc;
156  if (pa == 1 && pb == 1)
157  pc = 1;
158  else
159  pc = sz;
160  if (c.size() != pc)
161  c.resize(pc);
162 
163  if (pc == 1) {
164  c[0] = func(a[0], b[0]);
165  return;
166  }
167 
168  {
169 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
170  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Polynomial Evaluation");
171 #endif
172 
173  // Evaluate input
174  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
175  &avals[0], 1);
176  blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.coeff(), 1, 0.0,
177  &bvals[0], 1);
178 
179  }
180 
181  {
182 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
183  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Function Evaluation");
184 #endif
185 
186  // Evaluate function
187  for (ordinal_type qp=0; qp<nqp; qp++)
188  if (quad_weights[qp] != value_type(0))
189  fvals[qp] = func(avals[qp], bvals[qp])*quad_weights[qp];
190  else
191  fvals[qp] = value_type(0);
192 
193  }
194 
195  {
196 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
197  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PP Binary Polynomial Integration");
198 #endif
199 
200  // Integrate
201  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
202  c.coeff(), 1);
203 
204  }
205 }
206 
207 template <typename ordinal_type, typename value_type, typename node_type>
208 template <typename FuncT>
209 void
211 binary_op(const FuncT& func,
213  const value_type& a,
215 {
216  ordinal_type pb = b.size();
217  ordinal_type pc;
218  if (pb == 1)
219  pc = 1;
220  else
221  pc = sz;
222  if (c.size() != pc)
223  c.resize(pc);
224 
225  if (pc == 1) {
226  c[0] = func(a, b[0]);
227  return;
228  }
229 
230  {
231 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
232  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Polynomial Evaluation");
233 #endif
234 
235  // Evaluate input
236  blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.coeff(), 1, 0.0,
237  &bvals[0], 1);
238 
239  }
240 
241  {
242 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
243  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Function Evaluation");
244 #endif
245 
246  // Evaluate function
247  for (ordinal_type qp=0; qp<nqp; qp++)
248  if (quad_weights[qp] != value_type(0))
249  fvals[qp] = func(a, bvals[qp])*quad_weights[qp];
250  else
251  fvals[qp] = value_type(0);
252 
253  }
254 
255  {
256 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
257  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- CP Binary Polynomial Integration");
258 #endif
259 
260  // Integrate
261  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
262  c.coeff(), 1);
263 
264  }
265 }
266 
267 template <typename ordinal_type, typename value_type, typename node_type>
268 template <typename FuncT>
269 void
271 binary_op(const FuncT& func,
274  const value_type& b)
275 {
276  ordinal_type pa = a.size();
277  ordinal_type pc;
278  if (pa == 1)
279  pc = 1;
280  else
281  pc = sz;
282  if (c.size() != pc)
283  c.resize(pc);
284 
285  if (pc == 1) {
286  c[0] = func(a[0], b);
287  return;
288  }
289 
290  {
291 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
292  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Polynomial Evaluation");
293 #endif
294 
295  // Evaluate input
296  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.coeff(), 1, 0.0,
297  &avals[0], 1);
298 
299  }
300 
301  {
302 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
303  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Function Evaluation");
304 #endif
305 
306  // Evaluate function
307  for (ordinal_type qp=0; qp<nqp; qp++)
308  if (quad_weights[qp] != value_type(0))
309  fvals[qp] = func(avals[qp], b)*quad_weights[qp];
310  else
311  fvals[qp] = value_type(0);
312 
313  }
314 
315  {
316 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
317  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- PC Binary Polynomial Integration");
318 #endif
319 
320  // Integrate
321  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
322  c.coeff(), 1);
323 
324  }
325 }
326 
327 template <typename ordinal_type, typename value_type, typename node_type>
328 template <typename FuncT>
329 void
331 nary_op(const FuncT& func,
334 {
335  const int N = FuncT::N;
336  bool is_constant = true;
337  for (int i=0; i<N; i++) {
338  if (na[i]->size() > 1) {
339  is_constant = false;
340  break;
341  }
342  }
343  ordinal_type pc;
344  if (is_constant)
345  pc = 1;
346  else
347  pc = sz;
348  if (c.size() != pc)
349  c.resize(pc);
350 
351  if (pc == 1) {
352  value_type val[N];
353  for (int i=0; i<N; i++)
354  val[i] = (*na[i])[0];
355  c[0] = func(val);
356  return;
357  }
358 
359  if (N >= navals.size())
360  navals.resize(N+1);
361  if (navals[N].size() != N) {
362  navals[N].resize(N);
363  for (int i=0; i<N; i++)
364  navals[N][i].resize(nqp);
365  }
366 
367  {
368 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
369  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Polynomial Evaluation");
370 #endif
371 
372  // Evaluate input
373  for (int i=0; i<N; i++) {
374  //navals[i].resize(nqp);
375  ordinal_type pa = na[i]->size();
376  blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, na[i]->coeff(), 1, 0.0,
377  &navals[N][i][0], 1);
378  }
379 
380  }
381 
382  {
383 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
384  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Function Evaluation");
385 #endif
386 
387  // Evaluate function
388  value_type val[N];
389  for (ordinal_type qp=0; qp<nqp; qp++) {
390  if (quad_weights[qp] != value_type(0)) {
391  for (int i=0; i<N; i++)
392  val[i] = navals[N][i][qp];
393  fvals[qp] = func(val)*quad_weights[qp];
394  }
395  else
396  fvals[qp] = value_type(0);
397  }
398 
399  }
400 
401  {
402 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
403  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::QuadExp -- N(" << N << ")-ary Polynomial Integration");
404 #endif
405 
406  // Integrate
407  blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
408  c.coeff(), 1);
409 
410  }
411 }
412 
413 template <typename ordinal_type, typename value_type, typename node_type>
414 void
418  const value_type& x)
419 {
421 }
422 
423 template <typename ordinal_type, typename value_type, typename node_type>
424 void
428  const value_type& x)
429 {
431 }
432 
433 template <typename ordinal_type, typename value_type, typename node_type>
434 void
439 {
440  if (use_quad_for_times)
441  binary_op(times_quad_func(), c, c, x);
442  else
444 }
445 
446 template <typename ordinal_type, typename value_type, typename node_type>
447 void
452 {
453 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
454  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divideEqual(OPA)");
455 #endif
456  if (x.size() == 1) {
457  ordinal_type p = c.size();
458  value_type* cc = c.coeff();
459  const value_type* xc = x.coeff();
460  for (ordinal_type i=0; i<p; i++)
461  cc[i] /= xc[0];
462  }
463  else {
464  if (use_quad_for_division)
465  binary_op(div_quad_func(), c, c, x);
466  else
468  }
469 }
470 
471 template <typename ordinal_type, typename value_type, typename node_type>
472 void
477 {
478  if (use_quad_for_times)
479  binary_op(times_quad_func(), c, a, b);
480  else
482 }
483 
484 template <typename ordinal_type, typename value_type, typename node_type>
485 void
488  const value_type& a,
490 {
492 }
493 
494 template <typename ordinal_type, typename value_type, typename node_type>
495 void
499  const value_type& b)
500 {
502 }
503 
504 template <typename ordinal_type, typename value_type, typename node_type>
505 void
510 {
511 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
512  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::OrthogPolyExpansionBase::divide(OPA,OPA)");
513 #endif
514  if (b.size() == 1) {
515  ordinal_type pc = a.size();
516  if (c.size() != pc)
517  c.resize(pc);
518 
519  const value_type* ca = a.coeff();
520  const value_type* cb = b.coeff();
521  value_type* cc = c.coeff();
522 
523  for (ordinal_type i=0; i<pc; i++)
524  cc[i] = ca[i]/cb[0];
525  }
526  else {
527  if (use_quad_for_division)
528  binary_op(div_quad_func(), c, a, b);
529  else
531  }
532 }
533 
534 template <typename ordinal_type, typename value_type, typename node_type>
535 void
538  const value_type& a,
540 {
541  if (use_quad_for_division)
542  binary_op(div_quad_func(), c, a, b);
543  else
545 }
546 
547 template <typename ordinal_type, typename value_type, typename node_type>
548 void
552  const value_type& b)
553 {
555 }
556 
557 template <typename ordinal_type, typename value_type, typename node_type>
558 void
562 {
563  unary_op(exp_quad_func(), c, a);
564 }
565 
566 template <typename ordinal_type, typename value_type, typename node_type>
567 void
571 {
572  unary_op(log_quad_func(), c, a);
573 }
574 
575 template <typename ordinal_type, typename value_type, typename node_type>
576 void
580 {
581  unary_op(log10_quad_func(), c, a);
582 }
583 
584 template <typename ordinal_type, typename value_type, typename node_type>
585 void
589 {
590  unary_op(sqrt_quad_func(), c, a);
591 }
592 
593 template <typename ordinal_type, typename value_type, typename node_type>
594 void
598 {
599  unary_op(cbrt_quad_func(), c, a);
600 }
601 
602 template <typename ordinal_type, typename value_type, typename node_type>
603 void
608 {
609  binary_op(pow_quad_func(), c, a, b);
610 }
611 
612 template <typename ordinal_type, typename value_type, typename node_type>
613 void
616  const value_type& a,
618 {
619  binary_op(pow_quad_func(), c, a, b);
620 }
621 
622 template <typename ordinal_type, typename value_type, typename node_type>
623 void
627  const value_type& b)
628 {
629  binary_op(pow_quad_func(), c, a, b);
630 }
631 
632 template <typename ordinal_type, typename value_type, typename node_type>
633 void
637 {
638  unary_op(sin_quad_func(), s, a);
639 }
640 
641 template <typename ordinal_type, typename value_type, typename node_type>
642 void
646 {
647  unary_op(cos_quad_func(), c, a);
648 }
649 
650 template <typename ordinal_type, typename value_type, typename node_type>
651 void
655 {
656  unary_op(tan_quad_func(), t, a);
657 }
658 
659 template <typename ordinal_type, typename value_type, typename node_type>
660 void
664 {
665  unary_op(sinh_quad_func(), s, a);
666 }
667 
668 template <typename ordinal_type, typename value_type, typename node_type>
669 void
673 {
674  unary_op(cosh_quad_func(), c, a);
675 }
676 
677 template <typename ordinal_type, typename value_type, typename node_type>
678 void
682 {
683  unary_op(tanh_quad_func(), t, a);
684 }
685 
686 template <typename ordinal_type, typename value_type, typename node_type>
687 void
691 {
692  unary_op(acos_quad_func(), c, a);
693 }
694 
695 template <typename ordinal_type, typename value_type, typename node_type>
696 void
700 {
701  unary_op(asin_quad_func(), c, a);
702 }
703 
704 template <typename ordinal_type, typename value_type, typename node_type>
705 void
709 {
710  unary_op(atan_quad_func(), c, a);
711 }
712 
713 template <typename ordinal_type, typename value_type, typename node_type>
714 void
719 {
720  binary_op(atan2_quad_func(), c, a, b);
721 }
722 
723 template <typename ordinal_type, typename value_type, typename node_type>
724 void
727  const value_type& a,
729 {
730  binary_op(atan2_quad_func(), c, a, b);
731 }
732 
733 template <typename ordinal_type, typename value_type, typename node_type>
734 void
738  const value_type& b)
739 {
740  binary_op(atan2_quad_func(), c, a, b);
741 }
742 
743 template <typename ordinal_type, typename value_type, typename node_type>
744 void
748 {
749  unary_op(acosh_quad_func(), c, a);
750 }
751 
752 template <typename ordinal_type, typename value_type, typename node_type>
753 void
757 {
758  unary_op(asinh_quad_func(), c, a);
759 }
760 
761 template <typename ordinal_type, typename value_type, typename node_type>
762 void
766 {
767  unary_op(atanh_quad_func(), c, a);
768 }
769 
770 template <typename ordinal_type, typename value_type, typename node_type>
771 template <typename ExprT1, typename ExprT2>
774 compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
775 {
776  ordinal_type pa = a.size();
777  ordinal_type pb = b.size();
778 
779  if (pa > 1 && pb > 1) {
780  ordinal_type k_lim = pa;
781  ordinal_type j_lim = pb;
782  if (pb < pa) {
783  k_lim = pb;
784  j_lim = pa;
785  }
786  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
787 #ifdef STOKHOS_DEBUG
788  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
789  "Stokhos::QuadOrthogPolyExpansion::compute_times_coeff()"
790  << ": Index " << i << " is out of range [0,"
791  << this->Cijk->num_i() << ")!");
792 #endif
793  value_type cc = value_type(0);
794  value_type aa, bb, cijk;
795  ordinal_type j, k;
796  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
797  k_it != this->Cijk->k_end(i_it); ++k_it) {
798  k = index(k_it);
799  if (k < k_lim) {
800  if (pa < pb) {
801  if (k == 0)
802  aa = a.val();
803  else
804  aa = a.higher_order_coeff(k);
805  }
806  else {
807  if (k == 0)
808  aa = b.val();
809  else
810  aa = b.higher_order_coeff(k);
811  }
812  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
813  j_it != this->Cijk->j_end(k_it); ++j_it) {
814  j = index(j_it);
815  cijk = value(j_it);
816  if (j < j_lim) {
817  if (pa < pb) {
818  if (j == 0)
819  bb = b.val();
820  else
821  bb = b.higher_order_coeff(j);
822  }
823  else {
824  if (j == 0)
825  bb = a.val();
826  else
827  bb = a.higher_order_coeff(j);
828  }
829  cc += cijk*aa*bb;
830  }
831  }
832  }
833  }
834  return cc / norms[i];
835  }
836  else if (i == 0)
837  return a.val() * b.val();
838  else if (pa > 1) {
839  return a.higher_order_coeff(i)*b.val();
840  }
841  else {
842  return a.val()*b.higher_order_coeff(i);
843  }
844 }
845 
846 template <typename ordinal_type, typename value_type, typename node_type>
847 template <typename ExprT1, typename ExprT2>
850 fast_compute_times_coeff(ordinal_type i, const ExprT1& a, const ExprT2& b) const
851 {
852  typename Cijk_type::i_iterator i_it = this->Cijk->find_i(i);
853 #ifdef STOKHOS_DEBUG
854  TEUCHOS_TEST_FOR_EXCEPTION(i_it == this->Cijk->i_end(), std::logic_error,
855  "Stokhos::QuadOrthogPolyExpansion::fast_ompute_times_coeff()"
856  << ": Index " << i << " is out of range [0,"
857  << this->Cijk->num_i() << ")!");
858 #endif
859  value_type cc = value_type(0);
860  value_type aa, bb, cijk;
861  ordinal_type j, k;
862  for (typename Cijk_type::ik_iterator k_it = this->Cijk->k_begin(i_it);
863  k_it != this->Cijk->k_end(i_it); ++k_it) {
864  k = index(k_it);
865  if (k == 0)
866  aa = a.val();
867  else
868  aa = a.higher_order_coeff(k);
869  for (typename Cijk_type::ikj_iterator j_it = this->Cijk->j_begin(k_it);
870  j_it != this->Cijk->j_end(k_it); ++j_it) {
871  j = index(j_it);
872  cijk = value(j_it);
873  if (j == 0)
874  bb = b.val();
875  else
876  bb = b.higher_order_coeff(j);
877  cc += cijk*aa*bb;
878  }
879  }
880 
881  return cc / norms[i];
882 }
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)