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::RecurrenceBasis< ordinal_type, value_type > Class Template Referenceabstract

Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship:

\[ \gamma_{k+1}\psi_{k+1}(x) = (\delta_k x - \alpha_k)\psi_k(x) - \beta_k\psi_{k-1}(x) \]

for $k=0,\dots,P$ where $\psi_{-1}(x) = 0$, $\psi_{0}(x) = 1/\gamma_0$, and $\beta_{0} = 1 = \int d\lambda$. More...

#include <Stokhos_RecurrenceBasis.hpp>

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

Public Member Functions

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...
 
- Public Member Functions inherited from Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >
 OneDOrthogPolyBasis ()
 Default constructor. More...
 
virtual ~OneDOrthogPolyBasis ()
 Destructor. More...
 
virtual Teuchos::RCP
< OneDOrthogPolyBasis
< ordinal_type, value_type > > 
cloneWithOrder (ordinal_type p) const =0
 Clone this object with the option of building a higher order basis. More...
 
virtual void setSparseGridGrowthRule (LevelToOrderFnPtr ptr)=0
 Set sparse grid rule. More...
 

Protected Member Functions

 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 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 =0
 Compute recurrence coefficients. 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

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

 RecurrenceBasis (const RecurrenceBasis &)
 
RecurrenceBasisoperator= (const RecurrenceBasis &b)
 

Implementation of Stokhos::OneDOrthogPolyBasis methods

typedef OneDOrthogPolyBasis
< ordinal_type, value_type >
::LevelToOrderFnPtr 
LevelToOrderFnPtr
 Function pointer needed for level_to_order mappings. 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...
 

Additional Inherited Members

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

Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship:

\[ \gamma_{k+1}\psi_{k+1}(x) = (\delta_k x - \alpha_k)\psi_k(x) - \beta_k\psi_{k-1}(x) \]

for $k=0,\dots,P$ where $\psi_{-1}(x) = 0$, $\psi_{0}(x) = 1/\gamma_0$, and $\beta_{0} = 1 = \int d\lambda$.

Derived classes implement the recurrence relationship by implementing computeRecurrenceCoefficients(). If normalize = true in the constructor, then the recurrence relationship becomes:

\[ \sqrt{\frac{\gamma_{k+1}\beta_{k+1}}{\delta_{k+1}\delta_k}} \psi_{k+1}(x) = (x - \alpha_k/\delta_k)\psi_k(x) - \sqrt{\frac{\gamma_k\beta_k}{\delta_k\delta_{k-1}}} \psi_{k-1}(x) \]

for $k=0,\dots,P$ where $\psi_{-1}(x) = 0$, $\psi_{0}(x) = 1/\sqrt{\beta_0}$, Note that a three term recurrence can always be defined with $\gamma_k = \delta_k = 1$ in which case the polynomials are monic. However typical normalizations of some polynomial families (see Stokhos::LegendreBasis) require the extra terms. Also, the quadrature rule (points and weights) is the same regardless if the polynomials are normalized. However the normalization can affect other algorithms.

Definition at line 53 of file Stokhos_RecurrenceBasis.hpp.

Member Typedef Documentation

template<typename ordinal_type, typename value_type>
typedef OneDOrthogPolyBasis<ordinal_type,value_type>::LevelToOrderFnPtr Stokhos::RecurrenceBasis< ordinal_type, value_type >::LevelToOrderFnPtr

Function pointer needed for level_to_order mappings.

Definition at line 177 of file Stokhos_RecurrenceBasis.hpp.

Constructor & Destructor Documentation

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

Destructor.

Definition at line 74 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type >
Stokhos::RecurrenceBasis< ordinal_type, value_type >::RecurrenceBasis ( const std::string &  name,
ordinal_type  p,
bool  normalize,
Stokhos::GrowthPolicy  growth_ = SLOW_GROWTH 
)
protected

Constructor to be called by derived classes.

name is the name for the basis that will be displayed when printing the basis, p is the order of the basis, normalize indicates whether the basis polynomials should have unit-norm, and quad_zero_tol is used to replace any quadrature point within this tolerance with zero (which can help with duplicate removal in sparse grid calculations).

Definition at line 15 of file Stokhos_RecurrenceBasisImp.hpp.

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

Copy constructor with specified order.

Definition at line 37 of file Stokhos_RecurrenceBasisImp.hpp.

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

Member Function Documentation

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

Return order of basis (largest monomial degree $P$).

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 81 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type >
ordinal_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::size ( ) const
virtual

Return total size of basis (given by order() + 1).

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 89 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type >
const Teuchos::Array< value_type > & Stokhos::RecurrenceBasis< 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 given by order().

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 97 of file Stokhos_RecurrenceBasisImp.hpp.

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

Return norm squared of basis polynomial i.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 105 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::RCP< Stokhos::Dense3Tensor< ordinal_type, value_type > > Stokhos::RecurrenceBasis< 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=0,\dots,P$ where $P$ is size()-1 and $k=0,\dots,p$ where $p$ is the supplied order.

This method is implemented by computing $C_{ijk}$ using Gaussian quadrature.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 113 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type >
Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > Stokhos::RecurrenceBasis< ordinal_type, value_type >::computeSparseTripleProductTensor ( ordinal_type  order) 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=0,\dots,P$ where $P$ is size()-1 and $k=0,\dots,p$ where $p$ is the supplied order.

This method is implemented by computing $C_{ijk}$ using Gaussian quadrature.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 143 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type >
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > Stokhos::RecurrenceBasis< ordinal_type, value_type >::computeDerivDoubleProductTensor ( ) const
virtual

Compute derivative double product tensor.

The $(i,j)$ entry of the tensor $B_{ij}$ is given by $B_{ij} = \langle\psi_i'\psi_j\rangle$ where $\psi_l$ represents basis polynomial $l$ and $i,j=0,\dots,P$ where $P$ is the order of the basis.

This method is implemented by computing $B_{ij}$ using Gaussian quadrature.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 176 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type>
void Stokhos::RecurrenceBasis< ordinal_type, value_type >::evaluateBases ( const value_type &  point,
Teuchos::Array< value_type > &  basis_pts 
) const
virtual

Evaluate each basis polynomial at given point point.

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

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 206 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type>
value_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::evaluate ( const value_type &  point,
ordinal_type  order 
) const
virtual

Evaluate basis polynomial given by order order at given point point.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 240 of file Stokhos_RecurrenceBasisImp.hpp.

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

Print basis to stream os.

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

Definition at line 269 of file Stokhos_RecurrenceBasisImp.hpp.

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

Return string name of basis.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 302 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type>
void Stokhos::RecurrenceBasis< 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

Compute quadrature points, weights, and values of basis polynomials at given set of points points.

quad_order specifies the order to which the quadrature should be accurate, not the number of quadrature points. The number of points is given by (quad_order + 1) / 2. Note however the passed arrays do NOT need to be sized correctly on input as they will be resized appropriately.

The quadrature points and weights are computed from the three-term recurrence by solving a tri-diagional symmetric eigenvalue problem (see Gene H. Golub and John H. Welsch, "Calculation of Gauss Quadrature Rules", Mathematics of Computation, Vol. 23, No. 106 (Apr., 1969), pp. 221-230).

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Reimplemented in Stokhos::LanczosPCEBasis< ordinal_type, value_type >, Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >, Stokhos::LanczosPCEBasis< OrdinalType, ValueType >, Stokhos::LanczosProjPCEBasis< OrdinalType, ValueType >, Stokhos::StieltjesPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesPCEBasis< OrdinalType, ValueType >, Stokhos::HouseTriDiagPCEBasis< ordinal_type, value_type >, Stokhos::MonoProjPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >, Stokhos::HouseTriDiagPCEBasis< OrdinalType, ValueType >, Stokhos::MonoProjPCEBasis< OrdinalType, ValueType >, Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, and Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >.

Definition at line 310 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type >
ordinal_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::quadDegreeOfExactness ( ordinal_type  n) const
virtual

Return polynomial degree of exactness for a given number of quadrature points.

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Reimplemented in Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, and Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >.

Definition at line 410 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type >
ordinal_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::coefficientGrowth ( ordinal_type  n) const
virtual
template<typename ordinal_type, typename value_type >
ordinal_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::pointGrowth ( ordinal_type  n) const
virtual
template<typename ordinal_type, typename value_type>
virtual LevelToOrderFnPtr Stokhos::RecurrenceBasis< ordinal_type, value_type >::getSparseGridGrowthRule ( ) const
inlinevirtual

Get sparse grid level_to_order mapping function.

Predefined functions are: webbur::level_to_order_linear_wn Symmetric Gaussian linear growth webbur::level_to_order_linear_nn Asymmetric Gaussian linear growth webbur::level_to_order_exp_cc Clenshaw-Curtis exponential growth webbur::level_to_order_exp_gp Gauss-Patterson exponential growth webbur::level_to_order_exp_hgk Genz-Keister exponential growth webbur::level_to_order_exp_f2 Fejer-2 exponential growth

Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.

Definition at line 189 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
virtual void Stokhos::RecurrenceBasis< ordinal_type, value_type >::setSparseGridGrowthRule ( LevelToOrderFnPtr  ptr)
inlinevirtual

Set sparse grid rule.

Definition at line 193 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type , typename value_type>
void Stokhos::RecurrenceBasis< ordinal_type, value_type >::getRecurrenceCoefficients ( Teuchos::Array< value_type > &  alpha,
Teuchos::Array< value_type > &  beta,
Teuchos::Array< value_type > &  delta,
Teuchos::Array< value_type > &  gamma 
) const
virtual

Return recurrence coefficients defined by above formula.

Definition at line 446 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type , typename value_type>
void Stokhos::RecurrenceBasis< ordinal_type, value_type >::evaluateBasesAndDerivatives ( const value_type &  point,
Teuchos::Array< value_type > &  vals,
Teuchos::Array< value_type > &  derivs 
) const
virtual

Evaluate basis polynomials and their derivatives at given point point.

Definition at line 224 of file Stokhos_RecurrenceBasisImp.hpp.

template<typename ordinal_type, typename value_type>
virtual void Stokhos::RecurrenceBasis< ordinal_type, value_type >::setQuadZeroTol ( value_type  tol)
inlinevirtual

Set tolerance for zero in quad point generation.

Definition at line 210 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
virtual bool Stokhos::RecurrenceBasis< 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
protectedpure virtual
template<typename ordinal_type , typename value_type >
void Stokhos::RecurrenceBasis< ordinal_type, value_type >::setup ( )
protectedvirtual
template<typename ordinal_type , typename value_type>
void Stokhos::RecurrenceBasis< ordinal_type, value_type >::normalizeRecurrenceCoefficients ( Teuchos::Array< value_type > &  alpha,
Teuchos::Array< value_type > &  beta,
Teuchos::Array< value_type > &  delta,
Teuchos::Array< value_type > &  gamma 
) const
protected

Normalize coefficients.

Definition at line 460 of file Stokhos_RecurrenceBasisImp.hpp.

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

Member Data Documentation

template<typename ordinal_type, typename value_type>
std::string Stokhos::RecurrenceBasis< ordinal_type, value_type >::name
protected

Name of basis.

Definition at line 272 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
ordinal_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::p
protected

Order of basis.

Definition at line 275 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
bool Stokhos::RecurrenceBasis< ordinal_type, value_type >::normalize
protected

Normalize basis.

Definition at line 278 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
GrowthPolicy Stokhos::RecurrenceBasis< ordinal_type, value_type >::growth
protected

Smolyak growth policy.

Definition at line 281 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
value_type Stokhos::RecurrenceBasis< ordinal_type, value_type >::quad_zero_tol
protected

Tolerance for quadrature points near zero.

Definition at line 284 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
LevelToOrderFnPtr Stokhos::RecurrenceBasis< ordinal_type, value_type >::sparse_grid_growth_rule
protected

Sparse grid growth rule (as determined by Pecos)

Definition at line 287 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::RecurrenceBasis< ordinal_type, value_type >::alpha
protected

Recurrence $\alpha$ coefficients.

Definition at line 290 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::RecurrenceBasis< ordinal_type, value_type >::beta
protected

Recurrence $\beta$ coefficients.

Definition at line 293 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::RecurrenceBasis< ordinal_type, value_type >::delta
protected

Recurrence $\delta$ coefficients.

Definition at line 296 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::RecurrenceBasis< ordinal_type, value_type >::gamma
protected

Recurrence $\gamma$ coefficients.

Definition at line 299 of file Stokhos_RecurrenceBasis.hpp.

template<typename ordinal_type, typename value_type>
Teuchos::Array<value_type> Stokhos::RecurrenceBasis< ordinal_type, value_type >::norms
protected

Norms.

Definition at line 302 of file Stokhos_RecurrenceBasis.hpp.


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