44 #include "Teuchos_Assert.hpp"
46 #include "Teuchos_TimeMonitor.hpp"
47 #include "Teuchos_Tuple.hpp"
49 template <
typename ordinal_type,
typename value_type,
typename node_type>
54 const Teuchos::RCP<
const Quadrature<ordinal_type, value_type> >& quad_,
55 const Teuchos::RCP<Teuchos::ParameterList>& params_) :
58 sz(this->basis->size()),
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()),
73 qv[qp*sz+i] = quad_values[qp][i];
74 sqv[qp*sz+i] = quad_values[qp][i]/norms[i];
77 Teuchos::RCP<Teuchos::ParameterList> params = params_;
78 if (params == Teuchos::null)
79 params = Teuchos::rcp(
new Teuchos::ParameterList);
80 use_quad_for_times = params->get(
"Use Quadrature for Times",
false);
81 use_quad_for_division = params->get(
"Use Quadrature for Division",
true);
84 template <
typename ordinal_type,
typename value_type,
typename node_type>
85 template <
typename FuncT>
107 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
108 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- Unary Polynomial Evaluation");
112 blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.
coeff(), 1, 0.0,
118 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
119 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- Unary Function Evaluation");
125 fvals[qp] = func(avals[qp])*quad_weights[qp];
133 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
134 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- Unary Polynomial Integration");
138 blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
144 template <
typename ordinal_type,
typename value_type,
typename node_type>
145 template <
typename FuncT>
156 if (pa == 1 && pb == 1)
164 c[0] = func(a[0], b[0]);
169 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
170 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PP Binary Polynomial Evaluation");
174 blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.
coeff(), 1, 0.0,
176 blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.
coeff(), 1, 0.0,
182 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
183 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PP Binary Function Evaluation");
189 fvals[qp] = func(avals[qp], bvals[qp])*quad_weights[qp];
196 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
197 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PP Binary Polynomial Integration");
201 blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
207 template <
typename ordinal_type,
typename value_type,
typename node_type>
208 template <
typename FuncT>
226 c[0] = func(a, b[0]);
231 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
232 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- CP Binary Polynomial Evaluation");
236 blas.GEMV(Teuchos::TRANS, pb, nqp, 1.0, &qv[0], sz, b.
coeff(), 1, 0.0,
242 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
243 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- CP Binary Function Evaluation");
249 fvals[qp] = func(a, bvals[qp])*quad_weights[qp];
256 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
257 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- CP Binary Polynomial Integration");
261 blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
267 template <
typename ordinal_type,
typename value_type,
typename node_type>
268 template <
typename FuncT>
286 c[0] = func(a[0], b);
291 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
292 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PC Binary Polynomial Evaluation");
296 blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, a.
coeff(), 1, 0.0,
302 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
303 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PC Binary Function Evaluation");
309 fvals[qp] = func(avals[qp], b)*quad_weights[qp];
316 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
317 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- PC Binary Polynomial Integration");
321 blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
327 template <
typename ordinal_type,
typename value_type,
typename node_type>
328 template <
typename FuncT>
335 const int N = FuncT::N;
337 for (
int i=0; i<N; i++) {
338 if (na[i]->size() > 1) {
353 for (
int i=0; i<N; i++)
354 val[i] = (*na[i])[0];
359 if (N >= navals.size())
361 if (navals[N].size() != N) {
363 for (
int i=0; i<N; i++)
364 navals[N][i].resize(nqp);
368 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
369 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- N(" << N <<
")-ary Polynomial Evaluation");
373 for (
int i=0; i<N; i++) {
376 blas.GEMV(Teuchos::TRANS, pa, nqp, 1.0, &qv[0], sz, na[i]->coeff(), 1, 0.0,
377 &navals[N][i][0], 1);
383 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
384 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- N(" << N <<
")-ary Function Evaluation");
391 for (
int i=0; i<N; i++)
392 val[i] = navals[N][i][qp];
393 fvals[qp] = func(
val)*quad_weights[qp];
402 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
403 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::QuadExp -- N(" << N <<
")-ary Polynomial Integration");
407 blas.GEMV(Teuchos::NO_TRANS, pc, nqp, 1.0, &sqv[0], sz, &fvals[0], 1, 0.0,
413 template <
typename ordinal_type,
typename value_type,
typename node_type>
423 template <
typename ordinal_type,
typename value_type,
typename node_type>
433 template <
typename ordinal_type,
typename value_type,
typename node_type>
440 if (use_quad_for_times)
446 template <
typename ordinal_type,
typename value_type,
typename node_type>
453 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
454 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::OrthogPolyExpansionBase::divideEqual(OPA)");
464 if (use_quad_for_division)
471 template <
typename ordinal_type,
typename value_type,
typename node_type>
478 if (use_quad_for_times)
484 template <
typename ordinal_type,
typename value_type,
typename node_type>
494 template <
typename ordinal_type,
typename value_type,
typename node_type>
504 template <
typename ordinal_type,
typename value_type,
typename node_type>
511 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
512 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::OrthogPolyExpansionBase::divide(OPA,OPA)");
527 if (use_quad_for_division)
534 template <
typename ordinal_type,
typename value_type,
typename node_type>
541 if (use_quad_for_division)
547 template <
typename ordinal_type,
typename value_type,
typename node_type>
557 template <
typename ordinal_type,
typename value_type,
typename node_type>
566 template <
typename ordinal_type,
typename value_type,
typename node_type>
575 template <
typename ordinal_type,
typename value_type,
typename node_type>
584 template <
typename ordinal_type,
typename value_type,
typename node_type>
593 template <
typename ordinal_type,
typename value_type,
typename node_type>
602 template <
typename ordinal_type,
typename value_type,
typename node_type>
612 template <
typename ordinal_type,
typename value_type,
typename node_type>
622 template <
typename ordinal_type,
typename value_type,
typename node_type>
632 template <
typename ordinal_type,
typename value_type,
typename node_type>
641 template <
typename ordinal_type,
typename value_type,
typename node_type>
650 template <
typename ordinal_type,
typename value_type,
typename node_type>
659 template <
typename ordinal_type,
typename value_type,
typename node_type>
668 template <
typename ordinal_type,
typename value_type,
typename node_type>
677 template <
typename ordinal_type,
typename value_type,
typename node_type>
686 template <
typename ordinal_type,
typename value_type,
typename node_type>
695 template <
typename ordinal_type,
typename value_type,
typename node_type>
704 template <
typename ordinal_type,
typename value_type,
typename node_type>
713 template <
typename ordinal_type,
typename value_type,
typename node_type>
723 template <
typename ordinal_type,
typename value_type,
typename node_type>
733 template <
typename ordinal_type,
typename value_type,
typename node_type>
743 template <
typename ordinal_type,
typename value_type,
typename node_type>
752 template <
typename ordinal_type,
typename value_type,
typename node_type>
761 template <
typename ordinal_type,
typename value_type,
typename node_type>
770 template <
typename ordinal_type,
typename value_type,
typename node_type>
771 template <
typename ExprT1,
typename ExprT2>
779 if (pa > 1 && pb > 1) {
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() <<
")!");
797 k_it != this->Cijk->k_end(i_it); ++k_it) {
804 aa = a.higher_order_coeff(k);
810 aa = b.higher_order_coeff(k);
813 j_it != this->Cijk->j_end(k_it); ++j_it) {
821 bb = b.higher_order_coeff(
j);
827 bb = a.higher_order_coeff(
j);
834 return cc / norms[i];
837 return a.val() * b.val();
839 return a.higher_order_coeff(i)*b.val();
842 return a.val()*b.higher_order_coeff(i);
846 template <
typename ordinal_type,
typename value_type,
typename node_type>
847 template <
typename ExprT1,
typename ExprT2>
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() <<
")!");
863 k_it != this->Cijk->k_end(i_it); ++k_it) {
868 aa = a.higher_order_coeff(k);
870 j_it != this->Cijk->j_end(k_it); ++j_it) {
876 bb = b.higher_order_coeff(
j);
881 return cc / norms[i];
Class to store coefficients of a projection onto an orthogonal polynomial basis.
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved)
pointer coeff()
Return coefficient array.
ordinal_type size() const
Return size.
Abstract base class for multivariate orthogonal polynomials.
Base class for consolidating common expansion implementations.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void log10(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)
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 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 tanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void sqrt(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)
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 > ¶ms=Teuchos::null)
Constructor.
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 nary_op(const FuncT &func, 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)
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 sinh(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 asinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cbrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
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.
void atan(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)
value_type fast_compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
void acos(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)
value_type compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
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 asin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void atanh(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.
ikj_sparse_array::const_iterator i_iterator
Iterator for looping over i entries.
SparseArrayIterator< index_iterator, value_iterator >::value_reference value(const SparseArrayIterator< index_iterator, value_iterator > &it)
SparseArrayIterator< index_iterator, value_iterator >::value_type index(const SparseArrayIterator< index_iterator, value_iterator > &it)
constexpr KOKKOS_INLINE_FUNCTION std::enable_if< is_view_uq_pce< view_type >::value, typename CijkType< view_type >::type >::type cijk(const view_type &view)
KOKKOS_INLINE_FUNCTION bool is_constant(const T &x)
const IndexType const IndexType const IndexType const IndexType const ValueType const ValueType * x
Bi-directional iterator for traversing a sparse array.