42 #ifndef STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP 43 #define STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP 45 #include "Teuchos_RCP.hpp" 46 #include "Teuchos_Array.hpp" 47 #include "Teuchos_BLAS.hpp" 48 #include "Teuchos_SerialDenseMatrix.hpp" 49 #include "Teuchos_SerialDenseVector.hpp" 50 #include "Teuchos_ParameterList.hpp" 68 template <
typename ordinal_type,
typename value_type>
85 const Teuchos::ParameterList&
params = Teuchos::ParameterList());
107 virtual const Teuchos::Array<value_type>&
norm_squared()
const;
120 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
125 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
137 const Teuchos::ArrayView<const value_type>& point,
138 Teuchos::Array<value_type>& basis_vals)
const;
141 virtual void print(std::ostream& os)
const;
153 bool transpose =
false)
const;
160 bool transpose =
false)
const;
163 virtual Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
183 const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& A,
184 const Teuchos::SerialDenseMatrix<ordinal_type,value_type>& F,
185 const Teuchos::Array<value_type>& weights,
187 Teuchos::Array<ordinal_type>& num_terms_,
188 Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Qp_,
189 Teuchos::SerialDenseMatrix<ordinal_type,value_type>& Q_) = 0;
202 typedef Teuchos::SerialDenseVector<ordinal_type,value_type>
SDV;
203 typedef Teuchos::SerialDenseMatrix<ordinal_type,value_type>
SDM;
212 Teuchos::RCP<const Stokhos::OrthogPolyBasis<ordinal_type,value_type> >
pce_basis;
227 Teuchos::Array< Stokhos::MultiIndex<ordinal_type> >
terms;
242 Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
reduced_quad;
253 Teuchos::BLAS<ordinal_type,value_type>
blas;
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.
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.
ordinal_type order() const
Return order of basis.
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.
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)
Abstract base class for quadrature methods.
Abstract base class for reduced basis strategies built from polynomial chaos expansions in some other...
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
Top-level namespace for Stokhos classes and functions.
virtual void print(std::ostream &os) const
Print basis to stream os.
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 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 p
Total order of basis.
Teuchos::ParameterList params
Algorithm parameters.
ordinal_type dimension() const
Return dimension of basis.
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::SerialDenseVector< ordinal_type, value_type > SDV
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
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)
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
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 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.