42 #ifndef STOKHOS_PRODUCT_LANCZOS_PCE_BASIS_HPP 43 #define STOKHOS_PRODUCT_LANCZOS_PCE_BASIS_HPP 45 #include "Teuchos_RCP.hpp" 46 #include "Teuchos_Array.hpp" 47 #include "Teuchos_SerialDenseMatrix.hpp" 48 #include "Teuchos_SerialDenseVector.hpp" 49 #include "Teuchos_ParameterList.hpp" 70 template <
typename ordinal_type,
typename value_type>
87 const Teuchos::ParameterList&
params = Teuchos::ParameterList());
109 virtual const Teuchos::Array<value_type>&
norm_squared()
const;
122 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
127 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
139 const Teuchos::ArrayView<const value_type>& point,
140 Teuchos::Array<value_type>& basis_vals)
const;
143 virtual void print(std::ostream& os)
const;
146 virtual const std::string&
getName()
const;
189 bool transpose =
false)
const;
196 bool transpose =
false)
const;
199 virtual Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
223 typedef Teuchos::SerialDenseVector<ordinal_type,value_type>
SDV;
224 typedef Teuchos::SerialDenseMatrix<ordinal_type,value_type>
SDM;
245 Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
reduced_quad;
virtual void transformToOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients to original basis from this basis.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > SDM
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
virtual void print(std::ostream &os) const
Print basis to stream os.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
SDM Phi
Values of transformed basis at quadrature points.
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
virtual void evaluateBases(const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const
Evaluate basis polynomials at given point point.
virtual ordinal_type size() const
Return total size of basis.
SDM A
Transition matrix from reduced basis to original.
SDM Ainv
Projection matrix from original matrix to reduced.
virtual const std::string & getName() const
Return string name of basis.
Stokhos::OrthogPolyApprox< ordinal_type, value_type > tmp_pce
Temporary pce used in invariant subspace calculations.
std::string name
Name of basis.
ordinal_type order() const
Return order of basis.
virtual ~ProductLanczosPCEBasis()
Destructor.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
ProductLanczosPCEBasis(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::ParameterList ¶ms=Teuchos::ParameterList())
Constructor.
Abstract base class for quadrature methods.
Abstract base class for reduced basis strategies built from polynomial chaos expansions in some other...
virtual MultiIndex< ordinal_type > getMaxOrders() const
Return maximum order allowable for each coordinate basis.
Top-level namespace for Stokhos classes and functions.
ordinal_type dimension() const
Return dimension of basis.
virtual void transformFromOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients from original basis to this basis.
Teuchos::ParameterList params
Algorithm parameters.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > getCoordinateBases() const
Return coordinate bases.
Teuchos::RCP< Stokhos::CompletePolynomialBasis< ordinal_type, value_type > > tensor_lanczos_basis
Product Lanczos basis.
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > reduced_quad
Reduced quadrature object.
virtual ordinal_type index(const MultiIndex< ordinal_type > &term) const
Get index of the multivariate polynomial given orders of each coordinate.
virtual const MultiIndex< ordinal_type > & term(ordinal_type i) const
Get orders of each coordinate polynomial given an index i.
Abstract base class for 1-D orthogonal polynomials.
ordinal_type isInvariant(const Stokhos::OrthogPolyApprox< ordinal_type, value_type > &pce) const
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
Teuchos::SerialDenseVector< ordinal_type, value_type > SDV
ProductLanczosPCEBasis & operator=(const ProductLanczosPCEBasis &)
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.