Stokhos  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Public Types | Public Member Functions | List of all members
Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type > Class Template Referenceabstract

Abstract base class for 1-D orthogonal polynomials. More...

#include <Stokhos_OneDOrthogPolyBasis.hpp>

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

Public Types

typedef int(* LevelToOrderFnPtr )(int level, int growth)
 Function pointer needed for level_to_order mappings.
 

Public Member Functions

 OneDOrthogPolyBasis ()
 Default constructor.
 
virtual ~OneDOrthogPolyBasis ()
 Destructor.
 
virtual ordinal_type order () const =0
 Return order of basis (largest monomial degree $P$).
 
virtual ordinal_type size () const =0
 Return total size of basis (given by order() + 1).
 
virtual const Teuchos::Array
< value_type > & 
norm_squared () const =0
 Return array storing norm-squared of each basis polynomial. More...
 
virtual const value_type & norm_squared (ordinal_type i) const =0
 Return norm squared of basis polynomial i.
 
virtual Teuchos::RCP
< Stokhos::Dense3Tensor
< ordinal_type, value_type > > 
computeTripleProductTensor () const =0
 Compute triple product tensor. More...
 
virtual Teuchos::RCP
< Stokhos::Sparse3Tensor
< ordinal_type, value_type > > 
computeSparseTripleProductTensor (ordinal_type order) const =0
 Compute triple product tensor. More...
 
virtual Teuchos::RCP
< Teuchos::SerialDenseMatrix
< ordinal_type, value_type > > 
computeDerivDoubleProductTensor () const =0
 Compute derivative double product tensor. More...
 
virtual void evaluateBases (const value_type &point, Teuchos::Array< value_type > &basis_pts) const =0
 Evaluate each basis polynomial at given point point. More...
 
virtual value_type evaluate (const value_type &point, ordinal_type order) const =0
 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 =0
 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 =0
 Compute quadrature points, weights, and values of basis polynomials at given set of points points. More...
 
virtual ordinal_type quadDegreeOfExactness (ordinal_type n) const =0
 
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 ordinal_type coefficientGrowth (ordinal_type n) const =0
 Evaluate coefficient growth rule for Smolyak-type bases.
 
virtual ordinal_type pointGrowth (ordinal_type n) const =0
 Evaluate point growth rule for Smolyak-type bases.
 
virtual LevelToOrderFnPtr getSparseGridGrowthRule () const =0
 Get sparse grid level_to_order mapping function. More...
 
virtual void setSparseGridGrowthRule (LevelToOrderFnPtr ptr)=0
 Set sparse grid rule.
 

Detailed Description

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

Abstract base class for 1-D orthogonal polynomials.

This class provides an abstract interface for univariate orthogonal polynomials. Orthogonality is defined by the inner product

\[ (f,g) = \langle fg \rangle = \int_{-\infty}^{\infty} f(x)g(x) \rho(x) dx \]

where $\rho$ is the density function of the measure associated with the orthogonal polynomials. See Stokhos::RecurrenceBasis for a general implementation of this interface based on the three-term recurrence satisfied by these polynomials. Multivariate polynomials can be formed from a collection of univariate polynomials through tensor products (see Stokhos::CompletePolynomialBasis).

Like most classes in Stokhos, the class is templated on the ordinal and value types. Typically ordinal_type = int and value_type = double.

Member Function Documentation

template<typename ordinal_type, typename value_type>
virtual Teuchos::RCP<OneDOrthogPolyBasis<ordinal_type,value_type> > Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::cloneWithOrder ( ordinal_type  p) const
pure virtual
template<typename ordinal_type, typename value_type>
virtual Teuchos::RCP< Teuchos::SerialDenseMatrix<ordinal_type, value_type> > Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::computeDerivDoubleProductTensor ( ) const
pure 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.

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> > Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::computeSparseTripleProductTensor ( ordinal_type  order) const
pure 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.

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual Teuchos::RCP< Stokhos::Dense3Tensor<ordinal_type, value_type> > Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::computeTripleProductTensor ( ) const
pure 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.

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual void Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::evaluateBases ( const value_type &  point,
Teuchos::Array< value_type > &  basis_pts 
) const
pure 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().

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual void Stokhos::OneDOrthogPolyBasis< 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
pure 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.

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >, Stokhos::LanczosPCEBasis< ordinal_type, value_type >, Stokhos::LanczosProjPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesPCEBasis< ordinal_type, value_type >, Stokhos::HouseTriDiagPCEBasis< ordinal_type, value_type >, Stokhos::MonoProjPCEBasis< ordinal_type, value_type >, Stokhos::StieltjesBasis< ordinal_type, value_type, func_type >, Stokhos::ClenshawCurtisLegendreBasis< ordinal_type, value_type >, and Stokhos::GaussPattersonLegendreBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual LevelToOrderFnPtr Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::getSparseGridGrowthRule ( ) const
pure virtual

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

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

template<typename ordinal_type, typename value_type>
virtual const Teuchos::Array<value_type>& Stokhos::OneDOrthogPolyBasis< ordinal_type, value_type >::norm_squared ( ) const
pure 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().

Implemented in Stokhos::RecurrenceBasis< ordinal_type, value_type >, and Stokhos::PecosOneDOrthogPolyBasis< ordinal_type, value_type >.

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

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