|
Stokhos Development
|
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship: More...
#include <Stokhos_RecurrenceBasis.hpp>


Public Member Functions | |
| virtual | ~RecurrenceBasis () |
| Destructor. | |
Public Member Functions inherited from Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type > | |
| OneDOrthogPolyBasis () | |
| Default constructor. | |
| virtual | ~OneDOrthogPolyBasis () |
| Destructor. | |
| 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. | |
| virtual void | setSparseGridGrowthRule (LevelToOrderFnPtr ptr)=0 |
| Set sparse grid rule. | |
Implementation of Stokhos::OneDOrthogPolyBasis methods | |
| typedef OneDOrthogPolyBasis< ordinal_type, value_type >::LevelToOrderFnPtr | LevelToOrderFnPtr |
| Function pointer needed for level_to_order mappings. | |
| std::string | name |
| Name of basis. | |
| ordinal_type | p |
| Order of basis. | |
| bool | normalize |
| Normalize basis. | |
| GrowthPolicy | growth |
| Smolyak growth policy. | |
| value_type | quad_zero_tol |
| Tolerance for quadrature points near zero. | |
| LevelToOrderFnPtr | sparse_grid_growth_rule |
| Sparse grid growth rule (as determined by Pecos) | |
| Teuchos::Array< value_type > | alpha |
| Recurrence | |
| Teuchos::Array< value_type > | beta |
| Recurrence | |
| Teuchos::Array< value_type > | delta |
| Recurrence | |
| Teuchos::Array< value_type > | gamma |
| Recurrence | |
| Teuchos::Array< value_type > | norms |
| Norms. | |
| virtual ordinal_type | order () const |
| Return order of basis (largest monomial degree | |
| virtual ordinal_type | size () const |
| Return total size of basis (given by order() + 1). | |
| virtual const Teuchos::Array< value_type > & | norm_squared () const |
| Return array storing norm-squared of each basis polynomial. | |
| virtual const value_type & | norm_squared (ordinal_type i) const |
Return norm squared of basis polynomial i. | |
| virtual Teuchos::RCP< Stokhos::Dense3Tensor< ordinal_type, value_type > > | computeTripleProductTensor () const |
| Compute triple product tensor. | |
| virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > | computeSparseTripleProductTensor (ordinal_type order) const |
| Compute triple product tensor. | |
| virtual Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > | computeDerivDoubleProductTensor () const |
| Compute derivative double product tensor. | |
| virtual void | evaluateBases (const value_type &point, Teuchos::Array< value_type > &basis_pts) const |
Evaluate each basis polynomial at given point point. | |
| virtual value_type | evaluate (const value_type &point, ordinal_type order) const |
Evaluate basis polynomial given by order order at given point point. | |
| virtual void | print (std::ostream &os) const |
Print basis to stream os. | |
| virtual const std::string & | getName () const |
| Return string name of basis. | |
| 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. | |
| virtual ordinal_type | quadDegreeOfExactness (ordinal_type n) const |
| virtual ordinal_type | coefficientGrowth (ordinal_type n) const |
| Evaluate coefficient growth rule for Smolyak-type bases. | |
| virtual ordinal_type | pointGrowth (ordinal_type n) const |
| Evaluate point growth rule for Smolyak-type bases. | |
| virtual LevelToOrderFnPtr | getSparseGridGrowthRule () const |
| Get sparse grid level_to_order mapping function. | |
| virtual void | setSparseGridGrowthRule (LevelToOrderFnPtr ptr) |
| Set sparse grid rule. | |
| 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. | |
| 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. | |
| virtual void | setQuadZeroTol (value_type tol) |
| Set tolerance for zero in quad point generation. | |
| RecurrenceBasis (const std::string &name, ordinal_type p, bool normalize, GrowthPolicy growth=SLOW_GROWTH) | |
| Constructor to be called by derived classes. | |
| RecurrenceBasis (ordinal_type p, const RecurrenceBasis &basis) | |
| Copy constructor with specified order. | |
| 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. | |
| virtual void | setup () |
| Setup basis after computing recurrence coefficients. | |
| 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. | |
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. | |
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)
\]](form_136.png)
for 



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)
\]](form_140.png)
for 



|
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).
|
virtual |
Evaluate coefficient growth rule for Smolyak-type bases.
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
Reimplemented in Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, and Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >.
|
virtual |
Compute derivative double product tensor.
The 






This method is implemented by computing 
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
protectedpure virtual |
Compute recurrence coefficients.
Derived classes should implement this method to compute their recurrence coefficients. n is the number of coefficients to compute. Return value indicates whether coefficients correspond to normalized (i.e., orthonormal) polynomials.
Note: Owing to the description above, gamma should be an array of length n+1.
Implemented in Stokhos::HouseTriDiagPCEBasis< ordinal_type, value_type >, Stokhos::MonoProjPCEBasis< ordinal_type, value_type >, Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >, Stokhos::HermiteBasis< ordinal_type, value_type >, Stokhos::JacobiBasis< ordinal_type, value_type >, Stokhos::LanczosPCEBasis< ordinal_type, value_type >, Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >, Stokhos::LegendreBasis< ordinal_type, value_type >, Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >, and Stokhos::StieltjesPCEBasis< ordinal_type, value_type >.
|
virtual |
Compute triple product tensor.
The 








order.
This method is implemented by computing 
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
virtual |
Compute triple product tensor.
The 








order.
This method is implemented by computing 
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
virtual |
Evaluate basis polynomial given by order order at given point point.
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
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 >.
|
virtual |
Return string name of basis.
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
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::HouseTriDiagPCEBasis< ordinal_type, value_type >, Stokhos::MonoProjPCEBasis< ordinal_type, value_type >, Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >, Stokhos::LanczosPCEBasis< ordinal_type, value_type >, Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >, and Stokhos::StieltjesPCEBasis< ordinal_type, value_type >.
Referenced by Stokhos::HouseTriDiagPCEBasis< ordinal_type, value_type >::getQuadPoints(), Stokhos::MonoProjPCEBasis< ordinal_type, value_type >::getQuadPoints(), Stokhos::LanczosPCEBasis< ordinal_type, value_type >::getQuadPoints(), Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >::getQuadPoints(), Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >::getQuadPoints(), Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::getQuadPoints(), and Stokhos::LanczosPCEBasis< ordinal_type, value_type >::LanczosPCEBasis().
|
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 >.
References Stokhos::RecurrenceBasis< ordinal_type, value_type >::sparse_grid_growth_rule.
|
virtual |
Return array storing norm-squared of each basis polynomial.
Entry 



Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
virtual |
Return norm squared of basis polynomial i.
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
virtual |
Return order of basis (largest monomial degree 
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
virtual |
Evaluate point growth rule for Smolyak-type bases.
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
Reimplemented in Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, and Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >.
|
virtual |
Print basis to stream os.
Reimplemented from Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.
|
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 >.
|
protectedvirtual |
Setup basis after computing recurrence coefficients.
Derived classes should call this method after computing their recurrence coefficients in their constructor to finish setting up the basis.
Reimplemented in Stokhos::LanczosPCEBasis< ordinal_type, value_type >, Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >, and Stokhos::StieltjesPCEBasis< ordinal_type, value_type >.
Referenced by Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::DiscretizedStieltjesBasis(), Stokhos::DiscretizedStieltjesBasis< ordinal_type, value_type >::DiscretizedStieltjesBasis(), Stokhos::HermiteBasis< ordinal_type, value_type >::HermiteBasis(), Stokhos::HermiteBasis< ordinal_type, value_type >::HermiteBasis(), Stokhos::HouseTriDiagPCEBasis< ordinal_type, value_type >::HouseTriDiagPCEBasis(), Stokhos::JacobiBasis< ordinal_type, value_type >::JacobiBasis(), Stokhos::JacobiBasis< ordinal_type, value_type >::JacobiBasis(), Stokhos::LegendreBasis< ordinal_type, value_type >::LegendreBasis(), Stokhos::LegendreBasis< ordinal_type, value_type >::LegendreBasis(), Stokhos::MonoProjPCEBasis< ordinal_type, value_type >::MonoProjPCEBasis(), Stokhos::LanczosPCEBasis< ordinal_type, value_type >::setup(), Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >::setup(), Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >::setup(), and Stokhos::StieltjesPCEBasis< ordinal_type, value_type >::setup().
|
virtual |
Return total size of basis (given by order() + 1).
Implements Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >.