42 #ifndef STOKHOS_STIELTJESBASIS_HPP 43 #define STOKHOS_STIELTJESBASIS_HPP 45 #include "Teuchos_RCP.hpp" 46 #include "Teuchos_Array.hpp" 57 template <
typename ordinal_type,
typename value_type,
typename func_type>
74 const Teuchos::RCP<const func_type>&
func,
88 Teuchos::Array<value_type>& points,
89 Teuchos::Array<value_type>& weights,
90 Teuchos::Array< Teuchos::Array<value_type> >& values)
const;
114 Teuchos::Array<value_type>&
alpha,
115 Teuchos::Array<value_type>&
beta,
116 Teuchos::Array<value_type>&
delta,
117 Teuchos::Array<value_type>&
gamma)
const;
120 virtual void setup();
127 const Teuchos::Array<value_type>& weights,
128 const Teuchos::Array<value_type>& points,
129 Teuchos::Array<value_type>& a,
130 Teuchos::Array<value_type>& b,
131 Teuchos::Array<value_type>& nrm,
132 Teuchos::Array< Teuchos::Array<value_type> >&
phi_vals)
const;
140 const Teuchos::Array<value_type>& a,
141 const Teuchos::Array<value_type>& b,
142 const Teuchos::Array<value_type>& weights,
143 const Teuchos::Array<value_type>& points,
144 Teuchos::Array< Teuchos::Array<value_type> >&
phi_vals,
149 const Teuchos::Array<value_type>& a,
150 const Teuchos::Array<value_type>& b,
151 const Teuchos::Array<value_type>& points,
152 Teuchos::Array< Teuchos::Array<value_type> >& values)
const;
168 Teuchos::RCP<const func_type>
func;
171 Teuchos::RCP<const Stokhos::Quadrature<ordinal_type, value_type> >
quad;
183 mutable Teuchos::Array< Teuchos::Array<value_type> >
phi_vals;
Teuchos::Array< value_type > delta
Recurrence coefficients.
StieltjesBasis(ordinal_type p, const Teuchos::RCP< const func_type > &func, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, bool use_pce_quad_points, bool normalize=false)
Constructor.
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.
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship: for ...
StieltjesBasis & operator=(const StieltjesBasis &b)
const Teuchos::Array< Teuchos::Array< value_type > > & basis_values
Values of PCE basis functions at quadrature points.
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > quad
Quadrature object.
const Teuchos::Array< value_type > & pce_weights
PCE quadrature weights.
Teuchos::RCP< const func_type > func
PC expansion.
Teuchos::Array< Teuchos::Array< value_type > > phi_vals
Values of generated polynomials at PCE quadrature points.
Teuchos::Array< value_type > beta
Recurrence coefficients.
Teuchos::Array< value_type > alpha
Recurrence coefficients.
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a functional map...
virtual void setup()
Setup basis after computing recurrence coefficients.
~StieltjesBasis()
Destructor.
Abstract base class for quadrature methods.
Teuchos::Array< value_type > func_vals
Values of func at quadrature points.
Top-level namespace for Stokhos classes and functions.
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.
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 and .
ordinal_type p
Order of basis.
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.
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.
Teuchos::Array< value_type > gamma
Recurrence coefficients.
bool normalize
Normalize basis.
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.
bool use_pce_quad_points
Use underlying pce's quadrature data.