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 Member Functions | Protected Attributes | Private Member Functions | List of all members
Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type > Class Template Reference

Generates three-term recurrence using the Discretized Stieltjes procedure. More...

#include <Stokhos_DiscretizedStieltjesBasis.hpp>

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

Public Member Functions

 DiscretizedStieltjesBasis (const std::string &name, const ordinal_type &p, value_type(*weightFn)(const value_type &), const value_type &leftEndPt, const value_type &rightEndPt, bool normalize, GrowthPolicy growth=SLOW_GROWTH)
 Constructor. More...
 
 ~DiscretizedStieltjesBasis ()
 Destructor. More...
 
value_type eval_inner_product (const ordinal_type &order1, const ordinal_type &order2) const
 Evaluate inner product of two basis functions to test orthogonality. 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...
 
- 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 void getQuadPoints (ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
 Compute quadrature points, weights, and values of basis polynomials at given set of points points. 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 Member Functions

value_type evaluateRecurrence (const value_type &point, ordinal_type order, const Teuchos::Array< value_type > &alpha, const Teuchos::Array< value_type > &beta) const
 Evaluate 3-term recurrence formula using supplied coefficients. More...
 
value_type evaluateWeight (const value_type &x) const
 Evaluates the scaled weight function. More...
 
value_type expectedValue_tJ_nsquared (const ordinal_type &order, const Teuchos::Array< value_type > &alpha, const Teuchos::Array< value_type > &beta) const
 Approximates $\langle t\psi_k(t) \rangle$ where $k$ = order. More...
 
value_type expectedValue_J_nsquared (const ordinal_type &order, const Teuchos::Array< value_type > &alpha, const Teuchos::Array< value_type > &beta) const
 Approximates $\langle \psi_k(t) \rangle$ where $k$ = order. More...
 
 DiscretizedStieltjesBasis (const ordinal_type &p, const DiscretizedStieltjesBasis &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...
 
virtual void setup ()
 Setup basis after computing recurrence coefficients. 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

value_type scaleFactor
 Scale for the weight. More...
 
const value_type leftEndPt_
 Left end point of domain. More...
 
const value_type rightEndPt_
 Right end point of domain. More...
 
value_type(* weightFn_ )(const value_type &)
 Weight function. More...
 
Teuchos::Array< value_type > quad_points
 Quadrature points for discretized stieltjes procedure. More...
 
Teuchos::Array< value_type > quad_weights
 Quadrature points for discretized stieltjes procedure. More...
 
Teuchos::Array< Teuchos::Array
< value_type > > 
quad_values
 Quadrature values for discretized stieltjes procedure. 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

 DiscretizedStieltjesBasis (const DiscretizedStieltjesBasis &)
 
DiscretizedStieltjesBasisoperator= (const DiscretizedStieltjesBasis &b)
 

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...
 

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::DiscretizedStieltjesBasis< ordinal_type, value_type >

Generates three-term recurrence using the Discretized Stieltjes procedure.

Class generates an orthogonal polynomial basis for a given weight function and interval. The Discretized Stiltjes procedure described in "On the Calculation of Rys Polynomials and Quadratures", Robin P. Sagar, Vedene H. Smith is used to generate the recurrence coefficients.

Please be aware that this method is not fullproof and that it appears to encounter trouble with strongly singular weights since Gaussan quadrature is used to compute the relevant integrals. For 'nice' weight functions the method seems relatively robust.

Definition at line 31 of file Stokhos_DiscretizedStieltjesBasis.hpp.

Constructor & Destructor Documentation

template<typename ordinal_type , typename value_type >
Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::DiscretizedStieltjesBasis ( const std::string &  name,
const ordinal_type &  p,
value_type(*)(const value_type &)  weightFn,
const value_type &  leftEndPt,
const value_type &  rightEndPt,
bool  normalize,
Stokhos::GrowthPolicy  growth = SLOW_GROWTH 
)

Constructor.

Parameters
namename of the basis displayed when printing
porder of the basis
weightFnfunction pointer defining weight function
leftEndPtleft end point of the domain of the weight function
rightEndPtright end point of the domain of the weight function
normalizewhether polynomials should be given unit norm

Definition at line 12 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

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

Destructor.

Definition at line 59 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

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

Copy constructor with specified order.

Definition at line 36 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

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

Member Function Documentation

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::eval_inner_product ( const ordinal_type &  order1,
const ordinal_type &  order2 
) const

Evaluate inner product of two basis functions to test orthogonality.

Definition at line 164 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::RCP< Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type > > Stokhos::DiscretizedStieltjesBasis< 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 >.

Reimplemented in Stokhos::RysBasis< ordinal_type, value_type >.

Definition at line 208 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
bool Stokhos::DiscretizedStieltjesBasis< 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 66 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::evaluateRecurrence ( const value_type &  point,
ordinal_type  order,
const Teuchos::Array< value_type > &  alpha,
const Teuchos::Array< value_type > &  beta 
) const
protected

Evaluate 3-term recurrence formula using supplied coefficients.

Definition at line 183 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::evaluateWeight ( const value_type &  x) const
protected

Evaluates the scaled weight function.

Definition at line 114 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::expectedValue_tJ_nsquared ( const ordinal_type &  order,
const Teuchos::Array< value_type > &  alpha,
const Teuchos::Array< value_type > &  beta 
) const
protected

Approximates $\langle t\psi_k(t) \rangle$ where $k$ = order.

Definition at line 122 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::expectedValue_J_nsquared ( const ordinal_type &  order,
const Teuchos::Array< value_type > &  alpha,
const Teuchos::Array< value_type > &  beta 
) const
protected

Approximates $\langle \psi_k(t) \rangle$ where $k$ = order.

Definition at line 143 of file Stokhos_DiscretizedStieltjesBasisImp.hpp.

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

Member Data Documentation

template<typename ordinal_type , typename value_type >
value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::scaleFactor
mutableprotected

Scale for the weight.

Definition at line 124 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
const value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::leftEndPt_
protected

Left end point of domain.

Definition at line 127 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
const value_type Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::rightEndPt_
protected

Right end point of domain.

Definition at line 130 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
value_type(* Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::weightFn_)(const value_type &)
protected

Weight function.

Definition at line 133 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::Array<value_type> Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::quad_points
protected

Quadrature points for discretized stieltjes procedure.

Definition at line 136 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::Array<value_type> Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::quad_weights
protected

Quadrature points for discretized stieltjes procedure.

Definition at line 139 of file Stokhos_DiscretizedStieltjesBasis.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::Array<Teuchos::Array< value_type > > Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::quad_values
protected

Quadrature values for discretized stieltjes procedure.

Definition at line 142 of file Stokhos_DiscretizedStieltjesBasis.hpp.


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