Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Stokhos::StieltjesPCEBasis< ordinal_type, value_type > Class Template Reference

Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial chaos expansion in another basis. More...

#include <Stokhos_StieltjesPCEBasis.hpp>

Inheritance diagram for Stokhos::StieltjesPCEBasis< ordinal_type, value_type >:
Inheritance graph
[legend]

Public Member Functions

 StieltjesPCEBasis (ordinal_type p, const Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, bool use_pce_quad_points, bool normalize=false, bool project_integrals=false, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk=Teuchos::null)
 Constructor. More...
 
 ~StieltjesPCEBasis ()
 Destructor. More...
 
void transformCoeffsFromStieltjes (const value_type *in, value_type *out) const
 Map expansion coefficients from this basis to original. More...
 
- Public Member Functions inherited from Stokhos::RecurrenceBasis< ordinal_type, value_type >
virtual ~RecurrenceBasis ()
 Destructor. More...
 
virtual void getRecurrenceCoefficients (Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &delta, Teuchos::Array< value_type > &gamma) const
 Return recurrence coefficients defined by above formula. More...
 
virtual void evaluateBasesAndDerivatives (const value_type &point, Teuchos::Array< value_type > &vals, Teuchos::Array< value_type > &derivs) const
 Evaluate basis polynomials and their derivatives at given point point. More...
 
virtual void setQuadZeroTol (value_type tol)
 Set tolerance for zero in quad point generation. More...
 
virtual ordinal_type order () const
 Return order of basis (largest monomial degree $P$). More...
 
virtual ordinal_type size () const
 Return total size of basis (given by order() + 1). More...
 
virtual const Teuchos::Array
< value_type > & 
norm_squared () const
 Return array storing norm-squared of each basis polynomial. More...
 
virtual const value_type & norm_squared (ordinal_type i) const
 Return norm squared of basis polynomial i. More...
 
virtual Teuchos::RCP
< Stokhos::Dense3Tensor
< ordinal_type, value_type > > 
computeTripleProductTensor () const
 Compute triple product tensor. More...
 
virtual Teuchos::RCP
< Stokhos::Sparse3Tensor
< ordinal_type, value_type > > 
computeSparseTripleProductTensor (ordinal_type order) const
 Compute triple product tensor. More...
 
virtual Teuchos::RCP
< Teuchos::SerialDenseMatrix
< ordinal_type, value_type > > 
computeDerivDoubleProductTensor () const
 Compute derivative double product tensor. More...
 
virtual void evaluateBases (const value_type &point, Teuchos::Array< value_type > &basis_pts) const
 Evaluate each basis polynomial at given point point. More...
 
virtual value_type evaluate (const value_type &point, ordinal_type order) const
 Evaluate basis polynomial given by order order at given point point. More...
 
virtual void print (std::ostream &os) const
 Print basis to stream os. More...
 
virtual const std::string & getName () const
 Return string name of basis. More...
 
virtual ordinal_type quadDegreeOfExactness (ordinal_type n) const
 
virtual ordinal_type coefficientGrowth (ordinal_type n) const
 Evaluate coefficient growth rule for Smolyak-type bases. More...
 
virtual ordinal_type pointGrowth (ordinal_type n) const
 Evaluate point growth rule for Smolyak-type bases. More...
 
virtual LevelToOrderFnPtr getSparseGridGrowthRule () const
 Get sparse grid level_to_order mapping function. More...
 
virtual void setSparseGridGrowthRule (LevelToOrderFnPtr ptr)
 Set sparse grid rule. More...
 
- Public Member Functions inherited from Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >
 OneDOrthogPolyBasis ()
 Default constructor. More...
 
virtual ~OneDOrthogPolyBasis ()
 Destructor. More...
 
virtual void setSparseGridGrowthRule (LevelToOrderFnPtr ptr)=0
 Set sparse grid rule. More...
 

Protected Types

typedef Stokhos::Sparse3Tensor
< ordinal_type, value_type > 
Cijk_type
 Short-hand for triple product. More...
 

Protected Member Functions

void stieltjes (ordinal_type nstart, ordinal_type nfinish, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &a, Teuchos::Array< value_type > &b, Teuchos::Array< value_type > &nrm, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals) const
 Compute 3-term recurrence using Stieljtes procedure. More...
 
void integrateBasisSquared (ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals, value_type &val1, value_type &val2) const
 Compute $\int\psi^2_k(t) d\lambda(t)$ and $\int t\psi^2_k(t) d\lambda(t)$. More...
 
void integrateBasisSquaredProj (ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals, value_type &val1, value_type &val2) const
 Compute $\int\psi^2_k(t) d\lambda(t)$ and $\int t\psi^2_k(t) d\lambda(t)$ by projecting onto original PCE basis. More...
 
void evaluateRecurrence (ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &values) const
 Evaluate polynomials via 3-term recurrence. More...
 
 StieltjesPCEBasis (ordinal_type p, const StieltjesPCEBasis &basis)
 Copy constructor with specified order. More...
 
- Protected Member Functions inherited from Stokhos::RecurrenceBasis< ordinal_type, value_type >
 RecurrenceBasis (const std::string &name, ordinal_type p, bool normalize, GrowthPolicy growth=SLOW_GROWTH)
 Constructor to be called by derived classes. More...
 
 RecurrenceBasis (ordinal_type p, const RecurrenceBasis &basis)
 Copy constructor with specified order. More...
 
void normalizeRecurrenceCoefficients (Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &delta, Teuchos::Array< value_type > &gamma) const
 Normalize coefficients. More...
 

Protected Attributes

Teuchos::RCP< const
Stokhos::OrthogPolyApprox
< ordinal_type, value_type > > 
pce
 PC expansion. More...
 
Teuchos::RCP< const
Stokhos::Quadrature
< ordinal_type, value_type > > 
quad
 Quadrature object. More...
 
const Teuchos::Array
< value_type > & 
pce_weights
 PCE quadrature weights. More...
 
const Teuchos::Array
< Teuchos::Array< value_type > > & 
basis_values
 Values of PCE basis functions at quadrature points. More...
 
Teuchos::Array< value_type > pce_vals
 Values of PCE at quadrature points. More...
 
Teuchos::Array< Teuchos::Array
< value_type > > 
phi_vals
 Values of generated polynomials at PCE quadrature points. More...
 
bool use_pce_quad_points
 Use underlying pce's quadrature data. More...
 
Teuchos::SerialDenseMatrix
< ordinal_type, value_type > 
fromStieltjesMat
 Matrix mapping coefficients in Stieltjes basis back to original basis. More...
 
bool project_integrals
 Project Stieltjes integrals to original PCE basis. More...
 
Teuchos::RCP< const
Stokhos::OrthogPolyBasis
< ordinal_type, value_type > > 
basis
 PCE basis (needed for integral projection method) More...
 
Teuchos::RCP< const Cijk_typeCijk
 Triple product tensor (needed for integral projection method) More...
 
Teuchos::Array< value_type > phi_pce_coeffs
 Array store PC expansion of generated orthogonal polynomials. More...
 
- Protected Attributes inherited from Stokhos::RecurrenceBasis< ordinal_type, value_type >
std::string name
 Name of basis. More...
 
ordinal_type p
 Order of basis. More...
 
bool normalize
 Normalize basis. More...
 
GrowthPolicy growth
 Smolyak growth policy. More...
 
value_type quad_zero_tol
 Tolerance for quadrature points near zero. More...
 
LevelToOrderFnPtr sparse_grid_growth_rule
 Sparse grid growth rule (as determined by Pecos) More...
 
Teuchos::Array< value_type > alpha
 Recurrence $\alpha$ coefficients. More...
 
Teuchos::Array< value_type > beta
 Recurrence $\beta$ coefficients. More...
 
Teuchos::Array< value_type > delta
 Recurrence $\delta$ coefficients. More...
 
Teuchos::Array< value_type > gamma
 Recurrence $\gamma$ coefficients. More...
 
Teuchos::Array< value_type > norms
 Norms. More...
 

Private Member Functions

 StieltjesPCEBasis (const StieltjesPCEBasis &)
 
StieltjesPCEBasisoperator= (const StieltjesPCEBasis &b)
 

Implementation of Stokhos::OneDOrthogPolyBasis methods

virtual void getQuadPoints (ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
 Get Gauss quadrature points, weights, and values of basis at points. More...
 
virtual Teuchos::RCP
< OneDOrthogPolyBasis
< ordinal_type, value_type > > 
cloneWithOrder (ordinal_type p) const
 Clone this object with the option of building a higher order basis. More...
 

Implementation of Stokhos::RecurrenceBasis methods

virtual bool computeRecurrenceCoefficients (ordinal_type n, Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &delta, Teuchos::Array< value_type > &gamma) const
 Compute recurrence coefficients. More...
 
virtual void setup ()
 Setup basis after computing recurrence coefficients. More...
 

Additional Inherited Members

- Public Types inherited from Stokhos::RecurrenceBasis< ordinal_type, value_type >
typedef OneDOrthogPolyBasis
< ordinal_type, value_type >
::LevelToOrderFnPtr 
LevelToOrderFnPtr
 Function pointer needed for level_to_order mappings. More...
 
- Public Types inherited from Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >
typedef int(* LevelToOrderFnPtr )(int level, int growth)
 Function pointer needed for level_to_order mappings. More...
 

Detailed Description

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

Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial chaos expansion in another basis.

Definition at line 30 of file Stokhos_StieltjesPCEBasis.hpp.

Member Typedef Documentation

template<typename ordinal_type, typename value_type>
typedef Stokhos::Sparse3Tensor<ordinal_type, value_type> Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::Cijk_type
protected

Short-hand for triple product.

Definition at line 190 of file Stokhos_StieltjesPCEBasis.hpp.

Constructor & Destructor Documentation

template<typename ordinal_type, typename value_type>
Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::StieltjesPCEBasis ( ordinal_type  p,
const Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &  pce,
const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &  quad,
bool  use_pce_quad_points,
bool  normalize = false,
bool  project_integrals = false,
const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &  Cijk = Teuchos::null 
)

Constructor.

Parameters
porder of the basis
pcepolynomial chaos expansion defining new density function
quadquadrature data for basis of PC expansion
use_pce_quad_pointswhether to use quad to define quadrature points for the new basis, or whether to use the Golub-Welsch system.
normalizewhether polynomials should be given unit norm

Definition at line 16 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type , typename value_type >
Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::~StieltjesPCEBasis ( )

Destructor.

Definition at line 44 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::StieltjesPCEBasis ( ordinal_type  p,
const StieltjesPCEBasis< ordinal_type, value_type > &  basis 
)
protected

Copy constructor with specified order.

Definition at line 318 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::StieltjesPCEBasis ( const StieltjesPCEBasis< ordinal_type, value_type > &  )
private

Member Function Documentation

template<typename ordinal_type, typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::getQuadPoints ( ordinal_type  quad_order,
Teuchos::Array< value_type > &  points,
Teuchos::Array< value_type > &  weights,
Teuchos::Array< Teuchos::Array< value_type > > &  values 
) const
virtual

Get Gauss quadrature points, weights, and values of basis at points.

Reimplemented from Stokhos::RecurrenceBasis< ordinal_type, value_type >.

Definition at line 51 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type >
Teuchos::RCP< Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type > > Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::cloneWithOrder ( ordinal_type  p) const
virtual

Clone this object with the option of building a higher order basis.

This method is following the Prototype pattern (see Design Pattern's textbook). The slight variation is that it allows the order of the polynomial to be modified, otherwise an exact copy is formed. The use case for this is creating basis functions for column indices in a spatially varying adaptive refinement context.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 311 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type , typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::transformCoeffsFromStieltjes ( const value_type *  in,
value_type *  out 
) const

Map expansion coefficients from this basis to original.

Definition at line 300 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
bool Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::computeRecurrenceCoefficients ( ordinal_type  n,
Teuchos::Array< value_type > &  alpha,
Teuchos::Array< value_type > &  beta,
Teuchos::Array< value_type > &  delta,
Teuchos::Array< value_type > &  gamma 
) const
protectedvirtual

Compute recurrence coefficients.

Implements Stokhos::RecurrenceBasis< ordinal_type, value_type >.

Definition at line 99 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type , typename value_type >
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::setup ( )
protectedvirtual

Setup basis after computing recurrence coefficients.

Reimplemented from Stokhos::RecurrenceBasis< ordinal_type, value_type >.

Definition at line 126 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::stieltjes ( ordinal_type  nstart,
ordinal_type  nfinish,
const Teuchos::Array< value_type > &  weights,
const Teuchos::Array< value_type > &  points,
Teuchos::Array< value_type > &  a,
Teuchos::Array< value_type > &  b,
Teuchos::Array< value_type > &  nrm,
Teuchos::Array< Teuchos::Array< value_type > > &  phi_vals 
) const
protected

Compute 3-term recurrence using Stieljtes procedure.

Definition at line 159 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::integrateBasisSquared ( ordinal_type  k,
const Teuchos::Array< value_type > &  a,
const Teuchos::Array< value_type > &  b,
const Teuchos::Array< value_type > &  weights,
const Teuchos::Array< value_type > &  points,
Teuchos::Array< Teuchos::Array< value_type > > &  phi_vals,
value_type &  val1,
value_type &  val2 
) const
protected

Compute $\int\psi^2_k(t) d\lambda(t)$ and $\int t\psi^2_k(t) d\lambda(t)$.

Definition at line 208 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::integrateBasisSquaredProj ( ordinal_type  k,
const Teuchos::Array< value_type > &  a,
const Teuchos::Array< value_type > &  b,
const Teuchos::Array< value_type > &  weights,
const Teuchos::Array< value_type > &  points,
Teuchos::Array< Teuchos::Array< value_type > > &  phi_vals,
value_type &  val1,
value_type &  val2 
) const
protected

Compute $\int\psi^2_k(t) d\lambda(t)$ and $\int t\psi^2_k(t) d\lambda(t)$ by projecting onto original PCE basis.

Definition at line 251 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
void Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::evaluateRecurrence ( ordinal_type  k,
const Teuchos::Array< value_type > &  a,
const Teuchos::Array< value_type > &  b,
const Teuchos::Array< value_type > &  points,
Teuchos::Array< Teuchos::Array< value_type > > &  values 
) const
protected

Evaluate polynomials via 3-term recurrence.

Definition at line 229 of file Stokhos_StieltjesPCEBasisImp.hpp.

template<typename ordinal_type, typename value_type>
StieltjesPCEBasis& Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::operator= ( const StieltjesPCEBasis< ordinal_type, value_type > &  b)
private

Member Data Documentation

template<typename ordinal_type, typename value_type>
Teuchos::RCP<const Stokhos::OrthogPolyApprox<ordinal_type, value_type> > Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::pce
protected

PC expansion.

Definition at line 160 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> > Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::quad
protected

Quadrature object.

Definition at line 163 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
const Teuchos::Array<value_type>& Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::pce_weights
protected

PCE quadrature weights.

Definition at line 166 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
const Teuchos::Array< Teuchos::Array<value_type> >& Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::basis_values
protected

Values of PCE basis functions at quadrature points.

Definition at line 169 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::pce_vals
mutableprotected

Values of PCE at quadrature points.

Definition at line 172 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array< Teuchos::Array<value_type> > Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::phi_vals
mutableprotected

Values of generated polynomials at PCE quadrature points.

Definition at line 175 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
bool Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::use_pce_quad_points
protected

Use underlying pce's quadrature data.

Definition at line 178 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::SerialDenseMatrix<ordinal_type, value_type> Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::fromStieltjesMat
protected

Matrix mapping coefficients in Stieltjes basis back to original basis.

Definition at line 181 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
bool Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::project_integrals
protected

Project Stieltjes integrals to original PCE basis.

Definition at line 184 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::RCP<const Stokhos::OrthogPolyBasis<ordinal_type, value_type> > Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::basis
protected

PCE basis (needed for integral projection method)

Definition at line 187 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::RCP<const Cijk_type> Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::Cijk
protected

Triple product tensor (needed for integral projection method)

Definition at line 193 of file Stokhos_StieltjesPCEBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::phi_pce_coeffs
mutableprotected

Array store PC expansion of generated orthogonal polynomials.

Definition at line 196 of file Stokhos_StieltjesPCEBasis.hpp.


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