Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_DiscretizedStieltjesBasisImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 template <typename ordinal_type, typename value_type>
12 DiscretizedStieltjesBasis(const std::string& label,
13  const ordinal_type& p,
14  value_type (*weightFn)(const value_type&),
15  const value_type& leftEndPt,
16  const value_type& rightEndPt,
17  bool normalize,
18  Stokhos::GrowthPolicy growth) :
19  RecurrenceBasis<ordinal_type,value_type>(std::string("DiscretizedStieltjes -- ") + label, p, normalize, growth),
20  scaleFactor(1),
21  leftEndPt_(leftEndPt),
22  rightEndPt_(rightEndPt),
23  weightFn_(weightFn)
24 {
25  // Set up quadrature points for discretized stieltjes procedure
28  quadBasis->getQuadPoints(200*this->p, quad_points, quad_weights, quad_values);
29 
30  // Setup rest of recurrence basis
31  this->setup();
32 }
33 
34 template <typename ordinal_type, typename value_type>
37  const DiscretizedStieltjesBasis& basis) :
39  scaleFactor(basis.scaleFactor),
40  leftEndPt_(basis.leftEndPt_),
41  rightEndPt_(basis.rightEndPt_),
42  weightFn_(basis.weightFn_)
43 {
44  // Set up quadrature points for discretized stieltjes procedure
47  quadBasis->getQuadPoints(200*this->p, quad_points, quad_weights, quad_values);
48 
49  // Compute coefficients in 3-term recurrsion
50  computeRecurrenceCoefficients(p+1, this->alpha, this->beta, this->delta,
51  this->gamma);
52 
53  // Setup rest of recurrence basis
54  this->setup();
55 }
56 
57 template <typename ordinal_type, typename value_type>
60 {
61 }
62 
63 template <typename ordinal_type, typename value_type>
64 bool
70  Teuchos::Array<value_type>& gamma) const
71 {
72  //The Discretized Stieltjes polynomials are defined by a recurrance relation,
73  //P_n+1 = \gamma_n+1[(x-\alpha_n) P_n - \beta_n P_n-1].
74  //The alpha and beta coefficients are generated first using the
75  //discritized stilges procidure described in "On the Calculation of DiscretizedStieltjes Polynomials and Quadratures",
76  //Robin P. Sagar, Vedene H. Smith. The gamma coefficients are then optionally set so that each
77  //polynomial has norm 1. If normalization is not enabled then the gammas are set to 1.
78 
79  scaleFactor = 1;
80 
81  //First renormalize the weight function so that it has measure 1.
82  value_type oneNorm = expectedValue_J_nsquared(0, alpha, beta);
83  //future evaluations of the weight function will scale it by this factor.
84  scaleFactor = 1/oneNorm;
85 
86 
87  value_type integral2;
88  //NOTE!! This evaluation of 'expectedValue_J_nsquared(0)' is different
89  //from the one above since we rescaled the weight. Don't combine
90  //the two!!!
91 
92  value_type past_integral = expectedValue_J_nsquared(0, alpha, beta);
93  alpha[0] = expectedValue_tJ_nsquared(0, alpha, beta)/past_integral;
94  //beta[0] := \int_-c^c w(x) dx.
95  beta[0] = 1;
96  delta[0] = 1;
97  gamma[0] = 1;
98  //These formulas are from the above reference.
99  for (ordinal_type n = 1; n<k; n++){
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;
104  delta[n] = 1.0;
105  gamma[n] = 1.0;
106  }
107 
108  return false;
109 }
110 
111 template <typename ordinal_type, typename value_type>
112 value_type
114 evaluateWeight(const value_type& x) const
115 {
116  return (x < leftEndPt_ || x > rightEndPt_) ? 0: scaleFactor*weightFn_(x);
117 }
118 
119 template <typename ordinal_type, typename value_type>
123  const Teuchos::Array<value_type>& alpha,
124  const Teuchos::Array<value_type>& beta) const
125 {
126  //Impliments a gaussian quadrature routine to evaluate the integral,
127  // \int_-c^c J_n(x)^2w(x)dx. This is needed to compute the recurrance coefficients.
128  value_type integral = 0;
129  for(ordinal_type quadIdx = 0;
130  quadIdx < static_cast<ordinal_type>(quad_points.size()); quadIdx++) {
131  value_type x = (rightEndPt_ - leftEndPt_)*.5*quad_points[quadIdx] +
132  (rightEndPt_ + leftEndPt_)*.5;
133  value_type val = evaluateRecurrence(x,order,alpha,beta);
134  integral += x*val*val*evaluateWeight(x)*quad_weights[quadIdx];
135  }
136 
137  return integral*(rightEndPt_ - leftEndPt_);
138 }
139 
140 template <typename ordinal_type, typename value_type>
144  const Teuchos::Array<value_type>& alpha,
145  const Teuchos::Array<value_type>& beta) const
146 {
147  //Impliments a gaussian quadrature routineroutine to evaluate the integral,
148  // \int_-c^c J_n(x)^2w(x)dx. This is needed to compute the recurrance coefficients.
149  value_type integral = 0;
150  for(ordinal_type quadIdx = 0;
151  quadIdx < static_cast<ordinal_type>(quad_points.size()); quadIdx++){
152  value_type x = (rightEndPt_ - leftEndPt_)*.5*quad_points[quadIdx] +
153  (rightEndPt_ + leftEndPt_)*.5;
154  value_type val = evaluateRecurrence(x,order,alpha,beta);
155  integral += val*val*evaluateWeight(x)*quad_weights[quadIdx];
156  }
157 
158  return integral*(rightEndPt_ - leftEndPt_);
159 }
160 
161 template <typename ordinal_type, typename value_type>
162 value_type
164 eval_inner_product(const ordinal_type& order1, const ordinal_type& order2) const
165 {
166  //Impliments a gaussian quadrature routine to evaluate the integral,
167  // \int_-c^c J_n(x)J_m w(x)dx. This method is intended to allow the user to
168  // test for orthogonality and proper normalization.
169  value_type integral = 0;
170  for(ordinal_type quadIdx = 0;
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];
175  }
176 
177  return integral*(rightEndPt_ - leftEndPt_);
178 }
179 
180 template <typename ordinal_type, typename value_type>
181 value_type
184  ordinal_type k,
185  const Teuchos::Array<value_type>& alpha,
186  const Teuchos::Array<value_type>& beta) const
187 {
188  if (k == 0)
189  return value_type(1.0);
190  else if (k == 1)
191  return x-alpha[0];
192 
193  value_type v0 = value_type(1.0);
194  value_type v1 = x-alpha[0]*v0;
195  value_type v2 = value_type(0.0);
196  for (ordinal_type i=2; i<=k; i++) {
197  v2 = (x-alpha[i-1])*v1 - beta[i-1]*v0;
198  v0 = v1;
199  v1 = v2;
200  }
201 
202  return v2;
203 }
204 
205 template <typename ordinal_type, typename value_type>
209 {
211 }
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.
expr val()
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.
int n