10 #ifndef STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP
11 #define STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP
36 template <
typename ordinal_type,
typename value_type>
109 virtual void print(std::ostream& os)
const;
121 bool transpose =
false)
const;
128 bool transpose =
false)
const;
std::string name
Name of basis.
value_type rank_threshold
Rank threshold.
ordinal_type sz
Total size of basis.
SDM Qp
Coefficients of transformed basis in original basis.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
Utilities for indexing a multi-variate complete polynomial basis.
bool verbose
Whether to print a bunch of stuff out.
virtual ~GSReducedPCEBasisBase()
Destructor.
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
SDM Q
Values of transformed basis at quadrature points.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > pce_basis
Original pce basis.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > SDM
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > reduced_quad
Reduced quadrature object.
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 const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
void setup(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad)
virtual void print(std::ostream &os) const
Print basis to stream os.
Abstract base class for quadrature methods.
Abstract base class for reduced basis strategies built from polynomial chaos expansions in some other...
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.
Stokhos::CompletePolynomialBasisUtils< ordinal_type, value_type > CPBUtils
Teuchos::Array< value_type > norms
Norms.
GSReducedPCEBasisBase(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::ParameterList ¶ms=Teuchos::ParameterList())
Constructor.
virtual ordinal_type size() const
Return total size of basis.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
std::string orthogonalization_method
Orthogonalization method.
Teuchos::Array< Stokhos::MultiIndex< ordinal_type > > terms
2-D array of basis terms
ordinal_type order() const
Return order of basis.
ordinal_type p
Total order of basis.
Teuchos::ParameterList params
Algorithm parameters.
Teuchos::SerialDenseVector< ordinal_type, value_type > SDV
ordinal_type dimension() const
Return dimension of basis.
ordinal_type pce_sz
Size of original pce basis.
Teuchos::Array< ordinal_type > num_terms
Number of terms up to each order.
GSReducedPCEBasisBase & operator=(const GSReducedPCEBasisBase &b)
ordinal_type d
Total dimension of basis.
Teuchos::BLAS< ordinal_type, value_type > blas
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
virtual ordinal_type buildReducedBasis(ordinal_type max_p, value_type threshold, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &A, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &F, const Teuchos::Array< value_type > &weights, Teuchos::Array< Stokhos::MultiIndex< ordinal_type > > &terms_, Teuchos::Array< ordinal_type > &num_terms_, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Qp_, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Q_)=0
Build the reduced basis, parameterized by total order max_p.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
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.