10 #ifndef STOKHOS_SMOLYAK_PSEUDO_SPECTRAL_OPERATOR_HPP
11 #define STOKHOS_SMOLYAK_PSEUDO_SPECTRAL_OPERATOR_HPP
25 template <
typename ordinal_t,
27 typename point_compare_type =
28 typename DefaultPointCompare<ordinal_t,value_t>::type>
49 template <
typename coeff_compare_type>
52 bool use_smolyak_apply =
true,
54 const point_compare_type& point_compare = point_compare_type());
93 it ==
points.end(), std::logic_error,
"Invalid term " <<
point);
94 return it->second.second;
113 bool trans =
false)
const;
128 bool trans =
false)
const;
Teuchos::SerialDenseMatrix< ordinal_type, value_type > pce2qp
Matrix mapping coefficients to points.
point_set_type::iterator set_iterator
virtual ~SmolyakPseudoSpectralOperator()
Destructor.
void apply_direct(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &A, const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Apply transformation operator using direct method.
std::map< point_type, std::pair< value_type, ordinal_type >, point_compare_type > point_set_type
point_map_type point_map
Map index to sparse grid term.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > qp2pce
Matrix mapping points to coefficients.
set_iterator set_end()
Iterator to end of point set.
Container storing a term in a generalized tensor product.
base_type::const_set_iterator const_set_iterator
ordinal_type index(const point_type &point) const
Get point index for given point.
base_type::point_map_type point_map_type
base_type::iterator iterator
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
An operator interface for building pseudo-spectral approximations.
iterator begin()
Iterator to begining of point set.
Teuchos::Array< Teuchos::RCP< operator_type > > operator_set_type
const_iterator begin() const
Iterator to begining of point set.
Teuchos::Array< Teuchos::Array< ordinal_type > > gather_maps
Gather maps for each operator for Smolyak apply.
Teuchos::Array< Teuchos::Array< ordinal_type > > scatter_maps
Scatter maps for each operator for Smolyak apply.
ordinal_type coeff_size() const
Number of coefficients.
ordinal_type point_size() const
Number of points.
ordinal_type coeff_sz
Number of coefficients.
base_type::point_set_type point_set_type
void transformQP2PCE(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans=false) const
Transform values at quadrature points to PCE coefficients.
iterator end()
Iterator to end of point set.
void transformPCE2QP_smolyak(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Transform PCE coefficients to values at quadrature points using Smolyak formula.
void scatter(const Teuchos::Array< ordinal_type > &map, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, bool trans, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result) const
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
Teuchos::BLAS< ordinal_type, value_type > blas
BLAS wrappers.
base_type::point_type point_type
base_type::const_iterator const_iterator
void gather(const Teuchos::Array< ordinal_type > &map, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, bool trans, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result) const
const_set_iterator set_end() const
Iterator to end of point set.
Teuchos::Array< value_type > smolyak_coeffs
Smolyak coefficients.
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
const_set_iterator set_begin() const
Iterator to begining of point set.
MultiIndex< ordinal_type > multiindex_type
point_map_type::iterator iterator
const_iterator end() const
Iterator to end of point set.
void transformQP2PCE_smolyak(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Transform values at quadrature points to PCE coefficients using Smolyak formula.
SmolyakPseudoSpectralOperator(const SmolyakBasis< ordinal_type, value_type, coeff_compare_type > &smolyak_basis, bool use_smolyak_apply=true, bool use_pst=true, const point_compare_type &point_compare=point_compare_type())
Constructor.
const point_type & point(ordinal_type n) const
Get point for given index.
An operator for building pseudo-spectral coefficients using tensor-product quadrature.
TensorProductPseudoSpectralOperator< ordinal_type, value_type > operator_type
set_iterator set_begin()
Iterator to begining of point set.
point_map_type::const_iterator const_iterator
base_type::set_iterator set_iterator
virtual void transformPCE2QP(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans=false) const
Transform PCE coefficients to quadrature values.
bool use_smolyak
Use Smolyak apply method.
PseudoSpectralOperator< ordinal_type, value_type, point_compare_type > base_type
point_set_type points
Smolyak sparse grid.
point_set_type::const_iterator const_set_iterator
operator_set_type operators
Tensor pseudospectral operators.