10 template <
typename ordinal_type,
typename value_type>
21 leftEndPt_(leftEndPt),
22 rightEndPt_(rightEndPt),
34 template <
typename ordinal_type,
typename value_type>
39 scaleFactor(basis.scaleFactor),
40 leftEndPt_(basis.leftEndPt_),
41 rightEndPt_(basis.rightEndPt_),
42 weightFn_(basis.weightFn_)
57 template <
typename ordinal_type,
typename value_type>
63 template <
typename ordinal_type,
typename value_type>
82 value_type oneNorm = expectedValue_J_nsquared(0, alpha, beta);
84 scaleFactor = 1/oneNorm;
92 value_type past_integral = expectedValue_J_nsquared(0, alpha, beta);
93 alpha[0] = expectedValue_tJ_nsquared(0, alpha, beta)/past_integral;
100 integral2 = expectedValue_J_nsquared(
n, alpha, beta);
101 alpha[
n] = expectedValue_tJ_nsquared(
n, alpha, beta)/integral2;
102 beta[
n] = integral2/past_integral;
103 past_integral = integral2;
111 template <
typename ordinal_type,
typename value_type>
116 return (x < leftEndPt_ || x > rightEndPt_) ? 0: scaleFactor*weightFn_(x);
119 template <
typename ordinal_type,
typename value_type>
130 quadIdx < static_cast<ordinal_type>(quad_points.size()); quadIdx++) {
131 value_type x = (rightEndPt_ - leftEndPt_)*.5*quad_points[quadIdx] +
132 (rightEndPt_ + leftEndPt_)*.5;
134 integral += x*val*val*evaluateWeight(x)*quad_weights[quadIdx];
137 return integral*(rightEndPt_ - leftEndPt_);
140 template <
typename ordinal_type,
typename value_type>
151 quadIdx < static_cast<ordinal_type>(quad_points.size()); quadIdx++){
152 value_type x = (rightEndPt_ - leftEndPt_)*.5*quad_points[quadIdx] +
153 (rightEndPt_ + leftEndPt_)*.5;
155 integral += val*val*evaluateWeight(x)*quad_weights[quadIdx];
158 return integral*(rightEndPt_ - leftEndPt_);
161 template <
typename ordinal_type,
typename value_type>
171 quadIdx < static_cast<ordinal_type>(quad_points.size()); quadIdx++){
172 value_type x = (rightEndPt_ - leftEndPt_)*.5*quad_points[quadIdx] +
173 (rightEndPt_ + leftEndPt_)*.5;
174 integral += this->evaluate(x,order1)*this->evaluate(x,order2)*evaluateWeight(x)*quad_weights[quadIdx];
177 return integral*(rightEndPt_ - leftEndPt_);
180 template <
typename ordinal_type,
typename value_type>
197 v2 = (x-alpha[i-1])*v1 - beta[i-1]*v0;
205 template <
typename ordinal_type,
typename value_type>
Teuchos::Array< value_type > delta
Recurrence coefficients.
value_type evaluateWeight(const value_type &x) const
Evaluates the scaled weight function.
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.
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship: for ...
Generates three-term recurrence using the Discretized Stieltjes procedure.
value_type eval_inner_product(const ordinal_type &order1, const ordinal_type &order2) const
Evaluate inner product of two basis functions to test orthogonality.
Teuchos::Array< value_type > beta
Recurrence coefficients.
GrowthPolicy
Enumerated type for determining Smolyak growth policies.
Teuchos::Array< value_type > quad_weights
Quadrature points for discretized stieltjes procedure.
Teuchos::Array< value_type > alpha
Recurrence coefficients.
Teuchos::Array< Teuchos::Array< value_type > > quad_values
Quadrature values for discretized stieltjes procedure.
value_type expectedValue_tJ_nsquared(const ordinal_type &order, const Teuchos::Array< value_type > &alpha, const Teuchos::Array< value_type > &beta) const
Approximates where = order.
value_type expectedValue_J_nsquared(const ordinal_type &order, const Teuchos::Array< value_type > &alpha, const Teuchos::Array< value_type > &beta) const
Approximates where = order.
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.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::Array< value_type > quad_points
Quadrature points for discretized stieltjes procedure.
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.
Legendre polynomial basis.
virtual void setup()
Setup basis after computing recurrence coefficients.
Teuchos::Array< value_type > gamma
Recurrence coefficients.
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.
~DiscretizedStieltjesBasis()
Destructor.