Stokhos Development
Loading...
Searching...
No Matches
Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type > Class Template Reference

Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product measure induced by these expansions. More...

#include <Stokhos_MonomialProjGramSchmidtPCEBasis2.hpp>

Inheritance diagram for Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >:
Collaboration diagram for Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >:

Public Member Functions

 MonomialProjGramSchmidtPCEBasis2 (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 &params=Teuchos::ParameterList())
 Constructor.
virtual ~MonomialProjGramSchmidtPCEBasis2 ()
 Destructor.
Implementation of Stokhos::OrthogPolyBasis methods
ordinal_type order () const
 Return order of basis.
ordinal_type dimension () const
 Return dimension of basis.
virtual ordinal_type size () const
 Return total size of basis.
virtual const Teuchos::Array< value_type > & norm_squared () const
 Return array storing norm-squared of each basis polynomial.
virtual const value_typenorm_squared (ordinal_type i) const
 Return norm squared of basis polynomial i.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor () const
 Compute triple product tensor.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor () const
 Compute linear triple product tensor where k = 0,1,..,d.
virtual value_type evaluateZero (ordinal_type i) const
 Evaluate basis polynomial i at zero.
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 std::string & getName () const
 Return string name of basis.
virtual void print (std::ostream &os) const
 Print basis to stream os.
Public Member Functions inherited from Stokhos::ReducedPCEBasis< ordinal_type, value_type >
 ReducedPCEBasis ()
 Default constructor.
virtual ~ReducedPCEBasis ()
 Destructor.
Public Member Functions inherited from Stokhos::OrthogPolyBasis< ordinal_type, value_type >
 OrthogPolyBasis ()
 Constructor.
virtual ~OrthogPolyBasis ()
 Destructor.

Implementation of Stokhos::ReducedPCEBasis methods

typedef Stokhos::CompletePolynomialBasisUtils< ordinal_type, value_typeCPBUtils
typedef Teuchos::SerialDenseVector< ordinal_type, value_typeSDV
typedef Teuchos::SerialDenseMatrix< ordinal_type, value_typeSDM
std::string name
 Name of basis.
Teuchos::ParameterList params
 Algorithm parameters.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > pce_basis
 Original pce basis.
ordinal_type pce_sz
 Size of original pce basis.
ordinal_type p
 Total order of basis.
ordinal_type d
 Total dimension of basis.
ordinal_type sz
 Total size of basis.
Teuchos::Array< Stokhos::MultiIndex< ordinal_type > > terms
 2-D array of basis terms
Teuchos::Array< ordinal_typenum_terms
 Number of terms up to each order.
Teuchos::Array< value_typenorms
 Norms.
SDM Q
 Values of transformed basis at quadrature points.
SDM Qp
 Coefficients of transformed basis in original basis.
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > reduced_quad
 Reduced quadrature object.
bool verbose
 Whether to print a bunch of stuff out.
value_type rank_threshold
 Rank threshold.
std::string orthogonalization_method
 Orthogonalization method.
Teuchos::BLAS< ordinal_type, value_typeblas
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.
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.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature () const
 Get reduced quadrature object.
ordinal_type buildQ (ordinal_type max_p, value_type threshold, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, 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 > &A_, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &F_)
 Build the reduced basis, parameterized by total order max_p.

Detailed Description

template<typename ordinal_type, typename value_type>
class Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >

Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product measure induced by these expansions.

Given the PCE expansions, first build a non-orthogonal monomial basis.
Orthogonalize this basis using Gram-Schmidt, then build a quadrature rule using the simplex method.

Constructor & Destructor Documentation

◆ MonomialProjGramSchmidtPCEBasis2()

template<typename ordinal_type, typename value_type>
Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::MonomialProjGramSchmidtPCEBasis2 ( 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 & params = Teuchos::ParameterList() )

Constructor.

Parameters
porder of the basis
pcepolynomial chaos expansions defining new measure
quadquadrature data for basis defining pce
Cijksparse triple product tensor for basis defining pce
sparse_toltolerance for dropping terms in sparse tensors

References buildQ(), Stokhos::ReducedQuadratureFactory< ordinal_type, value_type >::createReducedQuadrature(), d, name, norms, num_terms, orthogonalization_method, p, params, pce_basis, pce_sz, Q, Qp, rank_threshold, reduced_quad, size(), Stokhos::ProductContainer< coeff_type >::size(), sz, terms, and verbose.

Member Function Documentation

◆ buildQ()

template<typename ordinal_type, typename value_type>
ordinal_type Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::buildQ ( ordinal_type max_p,
value_type threshold,
const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > & pce,
const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > & quad,
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 > & A_,
Teuchos::SerialDenseMatrix< ordinal_type, value_type > & F_ )
protected

Build the reduced basis, parameterized by total order max_p.

Returns resulting size of reduced basis

References Stokhos::CompletePolynomialBasisUtils< ordinal_type, value_type >::compute_terms(), d, orthogonalization_method, Stokhos::ProductContainer< coeff_type >::size(), and verbose.

Referenced by MonomialProjGramSchmidtPCEBasis2().

◆ computeLinearTripleProductTensor()

template<typename ordinal_type, typename value_type>
Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::computeLinearTripleProductTensor ( ) const
virtual

Compute linear triple product tensor where k = 0,1,..,d.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

◆ computeTripleProductTensor()

template<typename ordinal_type, typename value_type>
Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::computeTripleProductTensor ( ) const
virtual

Compute triple product tensor.

The $(i,j,k)$ entry of the tensor $C_{ijk}$ is given by $C_{ijk} = \langle\Psi_i\Psi_j\Psi_k\rangle$ where $\Psi_l$ represents basis polynomial $l$ and $i,j,k=0,\dots,P$ where $P$ is size()-1.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

◆ dimension()

template<typename ordinal_type, typename value_type>
ordinal_type Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::dimension ( ) const
virtual

Return dimension of basis.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References d.

◆ evaluateBases()

template<typename ordinal_type, typename value_type>
void Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::evaluateBases ( const Teuchos::ArrayView< const value_type > & point,
Teuchos::Array< value_type > & basis_vals ) const
virtual

Evaluate basis polynomials at given point point.

Size of returned array is given by size(), and coefficients are ordered from order 0 up to size size()-1.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

◆ evaluateZero()

template<typename ordinal_type, typename value_type>
value_type Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::evaluateZero ( ordinal_type i) const
virtual

Evaluate basis polynomial i at zero.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

◆ getName()

template<typename ordinal_type, typename value_type>
const std::string & Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::getName ( ) const
virtual

Return string name of basis.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References name.

◆ getReducedQuadrature()

template<typename ordinal_type, typename value_type>
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::getReducedQuadrature ( ) const
virtual

Get reduced quadrature object.

Implements Stokhos::ReducedPCEBasis< ordinal_type, value_type >.

References reduced_quad.

◆ norm_squared() [1/2]

template<typename ordinal_type, typename value_type>
const Teuchos::Array< value_type > & Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::norm_squared ( ) const
virtual

Return array storing norm-squared of each basis polynomial.

Entry $l$ of returned array is given by $\langle\Psi_l^2\rangle$ for $l=0,\dots,P$ where $P$ is size()-1.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References norms.

◆ norm_squared() [2/2]

template<typename ordinal_type, typename value_type>
const value_type & Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::norm_squared ( ordinal_type i) const
virtual

Return norm squared of basis polynomial i.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References norms.

◆ order()

template<typename ordinal_type, typename value_type>
ordinal_type Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::order ( ) const
virtual

Return order of basis.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References p.

◆ print()

template<typename ordinal_type, typename value_type>
void Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::print ( std::ostream & os) const
virtual

Print basis to stream os.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References d, norms, p, Qp, and sz.

◆ size()

Return total size of basis.

Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.

References sz.

Referenced by MonomialProjGramSchmidtPCEBasis2().

◆ transformFromOriginalBasis()

template<typename ordinal_type, typename value_type>
void Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::transformFromOriginalBasis ( const value_type * in,
value_type * out,
ordinal_type ncol = 1,
bool transpose = false ) const
virtual

Transform coefficients from original basis to this basis.

Implements Stokhos::ReducedPCEBasis< ordinal_type, value_type >.

References pce_sz, Qp, and sz.

◆ transformToOriginalBasis()

template<typename ordinal_type, typename value_type>
void Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >::transformToOriginalBasis ( const value_type * in,
value_type * out,
ordinal_type ncol = 1,
bool transpose = false ) const
virtual

Transform coefficients to original basis from this basis.

Implements Stokhos::ReducedPCEBasis< ordinal_type, value_type >.

References pce_sz, Qp, and sz.


The documentation for this class was generated from the following files: