44#ifndef STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
45#define STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
50#include "Teuchos_RCP.hpp"
51#include "Teuchos_Array.hpp"
52#include "Teuchos_SerialDenseMatrix.hpp"
53#include "Teuchos_SerialDenseVector.hpp"
54#include "Teuchos_BLAS.hpp"
71 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
163 template <
typename FuncT>
164 void nary_op(
const FuncT& func,
168 template <
typename ExprT1,
typename ExprT2>
170 const ExprT2& b)
const;
172 template <
typename ExprT1,
typename ExprT2>
174 const ExprT2& b)
const;
190 Teuchos::RCP<const Quadrature<ordinal_type, value_type> >
quad;
202 Teuchos::BLAS<ordinal_type,value_type>
blas;
214 const Teuchos::Array<value_type>&
norms;
226 Teuchos::Array< Teuchos::Array< Teuchos::Array<value_type> > >
navals;
232 Teuchos::Array<value_type>
qv;
235 Teuchos::Array<value_type>
sqv;
240 template <
typename FuncT>
247 template <
typename FuncT>
255 template <
typename FuncT>
263 template <
typename FuncT>
298 return std::log10(a);
316 return std::pow(a,b);
376 return std::atan2(a,b);
382 return std::log(a+std::sqrt(a*a-
value_type(1.0)));
388 return std::log(a+std::sqrt(a*a+
value_type(1.0)));
Class to store coefficients of a projection onto an orthogonal polynomial basis.
Abstract base class for multivariate orthogonal polynomials.
Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > Cijk
Triple-product tensor.
OrthogPolyExpansionBase(const Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > &basis, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor.
Teuchos::RCP< Teuchos::ParameterList > params
Parameter list.
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > basis
Basis.
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
Short-hand for Cijk.
ordinal_type nqp
Number of Quad points.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::Array< value_type > sqv
Quad values scaled by norms.
Teuchos::Array< Teuchos::Array< Teuchos::Array< value_type > > > navals
Temporary array for values of n-ary arguments at quad points.
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)
Teuchos::Array< value_type > qv
Reshaped quad values into 1D array.
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)
Teuchos::RCP< const Quadrature< ordinal_type, value_type > > quad
Quadrature routine.
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)
Teuchos::BLAS< ordinal_type, value_type > blas
BLAS wrappers.
virtual ~QuadOrthogPolyExpansion()
Destructor.
OrthogPolyExpansionBase< ordinal_type, value_type, node_type >::Cijk_type Cijk_type
Short-hand for Cijk.
bool use_quad_for_times
Use quadrature for times functions.
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)
const Teuchos::Array< Teuchos::Array< value_type > > & quad_values
Values of basis at Quad points.
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)
Teuchos::Array< value_type > fvals
Temporary array for values of operation at quad points.
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.
Teuchos::Array< value_type > avals
Temporary array for values of first argument at quad points.
ordinal_type sz
Expansions size.
void atan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Teuchos::Array< value_type > & norms
Norms of basis vectors.
QuadOrthogPolyExpansion(const QuadOrthogPolyExpansion &)
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
Teuchos::Array< value_type > bvals
Temporary array for values of second argument at quad points.
bool use_quad_for_division
Use quadrature for division functions.
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)
const Teuchos::Array< Teuchos::Array< value_type > > & quad_points
Array of Quad points.
value_type compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
const Teuchos::Array< value_type > & quad_weights
Array of Quad weights.
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)
QuadOrthogPolyExpansion & operator=(const QuadOrthogPolyExpansion &b)
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)
Abstract base class for quadrature methods.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Top-level namespace for Stokhos classes and functions.
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a, const value_type &b) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a, const value_type &b) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a, const value_type &b) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a) const
value_type operator()(const value_type &a, const value_type &b) const