49 template <
typename ordinal_type,
typename value_type,
typename node_type>
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()),
84 template <
typename ordinal_type,
typename value_type,
typename node_type>
85 template <
typename FuncT>
107 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
118 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
125 fvals[qp] = func(avals[qp])*quad_weights[qp];
133 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
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
182 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
189 fvals[qp] = func(avals[qp], bvals[qp])*quad_weights[qp];
196 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
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
242 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
249 fvals[qp] = func(a, bvals[qp])*quad_weights[qp];
256 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
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
302 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
309 fvals[qp] = func(avals[qp], b)*quad_weights[qp];
316 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
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
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
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
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
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
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) {
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>
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];
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 > ¶ms=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.
ordinal_type sz
Expansions size.
Bi-directional iterator for traversing a sparse array.
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)
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)