queso-0.57.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
QUESO::UniformLegendre1DQuadrature Class Reference

Class for Legendre-Gauss quadrature rule for one-dimensional functions. More...

#include <1DQuadrature.h>

Inheritance diagram for QUESO::UniformLegendre1DQuadrature:
QUESO::Base1DQuadrature QUESO::BaseQuadrature

Public Member Functions

Constructor/Destructor methods
 UniformLegendre1DQuadrature (double minDomainValue, double maxDomainValue, unsigned int order, bool densityIsNormalized)
 Default constructor. More...
 
 ~UniformLegendre1DQuadrature ()
 Destructor. More...
 
- Public Member Functions inherited from QUESO::Base1DQuadrature
 Base1DQuadrature (double minDomainValue, double maxDomainValue, unsigned int order)
 Default constructor. More...
 
virtual ~Base1DQuadrature ()=0
 Pure virtual destructor, forcing this to be an abstract object. More...
 
double minDomainValue () const
 Returns the minimum value of the domain of the (one-dimensional) function. More...
 
double maxDomainValue () const
 Returns the maximum value of the domain of the (one-dimensional) function. More...
 
unsigned int order () const
 Returns the order of the quadrature rule. More...
 
const std::vector< double > & positions () const
 Array of the positions for the numerical integration. More...
 
- Public Member Functions inherited from QUESO::BaseQuadrature
 BaseQuadrature ()
 
virtual ~BaseQuadrature ()=0
 Pure virtual destructor, forcing this to be an abstract object. More...
 
const std::vector< double > & weights () const
 Array of the weights used in the numerical integration. More...
 

Additional Inherited Members

- Protected Attributes inherited from QUESO::Base1DQuadrature
double m_minDomainValue
 
double m_maxDomainValue
 
unsigned int m_order
 
std::vector< double > m_positions
 
- Protected Attributes inherited from QUESO::BaseQuadrature
std::vector< double > m_weights
 

Detailed Description

Class for Legendre-Gauss quadrature rule for one-dimensional functions.

In a general Gaussian quadrature rule, an definite integral of \( f(x)\) is first approximated over the interval [-1,1] by a polynomial approximable function \( g(x)\) and a known weighting function \( W(x)\):

\[\int_{-1}^1 f(x) \, dx = \int_{-1}^1 W(x) g(x) \, dx\]

Those are then approximated by a sum of function values at specified points \( x_i \) multiplied by some weights \( w_i \):

\[ \int_{-1}^1 W(x) g(x) \, dx \approx \sum_{i=1}^n w_i g(x_i) \]

In the case of Gauss-Legendre quadrature, the weighting function \( W(x) = 1 \), so we can approximate an integral of \( f(x) \) with:

\[ \int_{-1}^1 f(x)\,dx \approx \sum_{i=1}^n w_i f(x_i) \]

The abscissas for quadrature order \( n \) are given by the roots of the Legendre polynomials \( P_n(x)\), which occur symmetrically about 0. The weights are

\[ w_i = \frac{2}{(1-x_i^2)[P'_n(x_i)]^2}=\frac{2(1-x_i^2)}{(n+1)^2[P_{n+1}(x_i)]^2} \]

Several authors give a table of abscissas and weights:

\( n \)\( x_i \)\( w_i \)
2 \( \pm \frac{1}{3}\sqrt{3} \) \( 1 \)
3 \( 0 \) \( \frac{8}{9} \)
\( \pm \frac{1}{5} \sqrt{15} \) \( \frac{5}{9} \)
4 \( \pm \frac{1}{35}\sqrt{525-70\sqrt{30}} \) \( \frac{1}{36}(18+\sqrt{30})\)
\( \pm \frac{1}{35}\sqrt{525+70\sqrt{30}} \) \( \frac{1}{36}(18-\sqrt{30})\)
5 \( 0 \) \( \frac{128}{225}\)
\( \pm \frac{1}{21}\sqrt{245-14\sqrt{70}} \) \( \frac{1}{900}(322+13\sqrt{70})\)
\( \pm \frac{1}{21}\sqrt{245+14\sqrt{70}} \) \( \frac{1}{900}(322-13\sqrt{70})\)
See Also
Weisstein, Eric W. "Legendre-Gauss Quadrature." From MathWorld–A Wolfram Web Resource. http://mathworld.wolfram.com/Legendre-GaussQuadrature.html.

Definition at line 148 of file 1DQuadrature.h.

Constructor & Destructor Documentation

QUESO::UniformLegendre1DQuadrature::UniformLegendre1DQuadrature ( double  minDomainValue,
double  maxDomainValue,
unsigned int  order,
bool  densityIsNormalized 
)

Default constructor.

Constructs a Gaussian-Legendre quadrature of order order, in the interval [minDomainValue,maxDomainValue]. Valid values for the order of the quadrature rule are: 1-7, 10-12, 16. This method scales the abscissas (positions) of the quadrature from the interval [-1,1] to [minDomainValue,maxDomainValue], and the parameter densityIsNormalized determines whether the weights should be scaled as well.

Definition at line 98 of file 1DQuadrature.C.

References QUESO::Base1DQuadrature::m_maxDomainValue, QUESO::Base1DQuadrature::m_minDomainValue, QUESO::Base1DQuadrature::m_order, QUESO::Base1DQuadrature::m_positions, and QUESO::BaseQuadrature::m_weights.

103  :
105 {
106  m_positions.resize(m_order+1,0.); // Yes, '+1'
107  m_weights.resize (m_order+1,0.); // Yes, '+1'
108 
109  // http://www.holoborodko.com/pavel/?page_id=679
110  switch (m_order) { // eep2011
111  case 1:
112  m_weights [0] = 1.;
113  m_weights [1] = 1.;
114 
115  m_positions[0] = -1./sqrt(3.);
116  m_positions[1] = 1./sqrt(3.);
117  break;
118 
119  case 2:
120  m_weights [0] = 5./9.;
121  m_weights [1] = 8./9.;
122  m_weights [2] = 5./9.;
123 
124  m_positions[0] = -sqrt(.6);
125  m_positions[1] = 0.;
126  m_positions[2] = sqrt(.6);
127  break;
128 
129  case 3:
130  m_weights [0] = 0.5 - sqrt(30.)/36.;
131  m_weights [1] = 0.5 + sqrt(30.)/36.;
132  m_weights [2] = 0.5 + sqrt(30.)/36.;
133  m_weights [3] = 0.5 - sqrt(30.)/36.;
134 
135  m_positions[0] = -sqrt(3.+2.*sqrt(1.2))/sqrt(7.);
136  m_positions[1] = -sqrt(3.-2.*sqrt(1.2))/sqrt(7.);
137  m_positions[2] = sqrt(3.-2.*sqrt(1.2))/sqrt(7.);
138  m_positions[3] = sqrt(3.+2.*sqrt(1.2))/sqrt(7.);
139  break;
140 
141  case 4:
142  m_weights [0] = (322.-13.*sqrt(70.))/900.; // 0.236926885
143  m_weights [1] = (322.+13.*sqrt(70.))/900.; // 0.478628670
144  m_weights [2] = 128./225.; // 0.568888889
145  m_weights [3] = (322.+13.*sqrt(70.))/900.;
146  m_weights [4] = (322.-13.*sqrt(70.))/900.;
147 
148  m_positions[0] = -sqrt(5.+2.*sqrt(10./7.))/3.;
149  m_positions[1] = -sqrt(5.-2.*sqrt(10./7.))/3.;
150  m_positions[2] = 0.;
151  m_positions[3] = sqrt(5.-2.*sqrt(10./7.))/3.;
152  m_positions[4] = sqrt(5.+2.*sqrt(10./7.))/3.;
153  break;
154 
155  case 5:
156  m_weights [0] = 0.1713244923791703450402961;
157  m_weights [1] = 0.3607615730481386075698335;
158  m_weights [2] = 0.4679139345726910473898703;
159  m_weights [3] = 0.4679139345726910473898703;
160  m_weights [4] = 0.3607615730481386075698335;
161  m_weights [5] = 0.1713244923791703450402961;
162 
163  m_positions[0] = -0.9324695142031520278123016;
164  m_positions[1] = -0.6612093864662645136613996;
165  m_positions[2] = -0.2386191860831969086305017;
166  m_positions[3] = 0.2386191860831969086305017;
167  m_positions[4] = 0.6612093864662645136613996;
168  m_positions[5] = 0.9324695142031520278123016;
169  break;
170 
171  case 6:
172  m_weights [0] = 0.1294849661688696932706114;
173  m_weights [1] = 0.2797053914892766679014678;
174  m_weights [2] = 0.3818300505051189449503698;
175  m_weights [3] = 0.4179591836734693877551020;
176  m_weights [4] = 0.3818300505051189449503698;
177  m_weights [5] = 0.2797053914892766679014678;
178  m_weights [6] = 0.1294849661688696932706114;
179 
180  m_positions[0] = -0.9491079123427585245261897;
181  m_positions[1] = -0.7415311855993944398638648;
182  m_positions[2] = -0.4058451513773971669066064;
183  m_positions[3] = 0.;
184  m_positions[4] = 0.4058451513773971669066064;
185  m_positions[5] = 0.7415311855993944398638648;
186  m_positions[6] = 0.9491079123427585245261897;
187  break;
188 
189  case 7:
190  m_weights [0] = 0.10122854;
191  m_weights [1] = 0.22238103;
192  m_weights [2] = 0.31370665;
193  m_weights [3] = 0.36268378;
194  m_weights [4] = 0.36268378;
195  m_weights [5] = 0.31370665;
196  m_weights [6] = 0.22238103;
197  m_weights [7] = 0.10122854;
198 
199  m_positions[0] = -0.96028986;
200  m_positions[1] = -0.79666648;
201  m_positions[2] = -0.52553241;
202  m_positions[3] = -0.18343464;
203  m_positions[4] = 0.18343464;
204  m_positions[5] = 0.52553241;
205  m_positions[6] = 0.79666648;
206  m_positions[7] = 0.96028986;
207  break;
208 
209  case 10:
210  m_weights [ 0] = 0.0556685671161736664827537;
211  m_weights [ 1] = 0.1255803694649046246346943;
212  m_weights [ 2] = 0.1862902109277342514260976;
213  m_weights [ 3] = 0.2331937645919904799185237;
214  m_weights [ 4] = 0.2628045445102466621806889;
215  m_weights [ 5] = 0.2729250867779006307144835;
216  m_weights [ 6] = 0.2628045445102466621806889;
217  m_weights [ 7] = 0.2331937645919904799185237;
218  m_weights [ 8] = 0.1862902109277342514260976;
219  m_weights [ 9] = 0.1255803694649046246346943;
220  m_weights [10] = 0.0556685671161736664827537;
221 
222  m_positions[ 0] = -0.9782286581460569928039380;
223  m_positions[ 1] = -0.8870625997680952990751578;
224  m_positions[ 2] = -0.7301520055740493240934163;
225  m_positions[ 3] = -0.5190961292068118159257257;
226  m_positions[ 4] = -0.2695431559523449723315320;
227  m_positions[ 5] = 0.;
228  m_positions[ 6] = 0.2695431559523449723315320;
229  m_positions[ 7] = 0.5190961292068118159257257;
230  m_positions[ 8] = 0.7301520055740493240934163;
231  m_positions[ 9] = 0.8870625997680952990751578;
232  m_positions[10] = 0.9782286581460569928039380;
233  break;
234 
235  case 11:
236  m_weights [ 0] = 0.0471753363865118271946160;
237  m_weights [ 1] = 0.1069393259953184309602547;
238  m_weights [ 2] = 0.1600783285433462263346525;
239  m_weights [ 3] = 0.2031674267230659217490645;
240  m_weights [ 4] = 0.2334925365383548087608499;
241  m_weights [ 5] = 0.2491470458134027850005624;
242  m_weights [ 6] = 0.2491470458134027850005624;
243  m_weights [ 7] = 0.2334925365383548087608499;
244  m_weights [ 8] = 0.2031674267230659217490645;
245  m_weights [ 9] = 0.1600783285433462263346525;
246  m_weights [10] = 0.1069393259953184309602547;
247  m_weights [11] = 0.0471753363865118271946160;
248 
249  m_positions[ 0] = -0.9815606342467192506905491;
250  m_positions[ 1] = -0.9041172563704748566784659;
251  m_positions[ 2] = -0.7699026741943046870368938;
252  m_positions[ 3] = -0.5873179542866174472967024;
253  m_positions[ 4] = -0.3678314989981801937526915;
254  m_positions[ 5] = -0.1252334085114689154724414;
255  m_positions[ 6] = 0.1252334085114689154724414;
256  m_positions[ 7] = 0.3678314989981801937526915;
257  m_positions[ 8] = 0.5873179542866174472967024;
258  m_positions[ 9] = 0.7699026741943046870368938;
259  m_positions[10] = 0.9041172563704748566784659;
260  m_positions[11] = 0.9815606342467192506905491;
261  break;
262 
263  case 12:
264  m_weights [ 0] = 0.0404840047653158795200216;
265  m_weights [ 1] = 0.0921214998377284479144218;
266  m_weights [ 2] = 0.1388735102197872384636018;
267  m_weights [ 3] = 0.1781459807619457382800467;
268  m_weights [ 4] = 0.2078160475368885023125232;
269  m_weights [ 5] = 0.2262831802628972384120902;
270  m_weights [ 6] = 0.2325515532308739101945895;
271  m_weights [ 7] = 0.2262831802628972384120902;
272  m_weights [ 8] = 0.2078160475368885023125232;
273  m_weights [ 9] = 0.1781459807619457382800467;
274  m_weights [10] = 0.1388735102197872384636018;
275  m_weights [11] = 0.0921214998377284479144218;
276  m_weights [12] = 0.0404840047653158795200216;
277 
278  m_positions[ 0] = -0.9841830547185881494728294;
279  m_positions[ 1] = -0.9175983992229779652065478;
280  m_positions[ 2] = -0.8015780907333099127942065;
281  m_positions[ 3] = -0.6423493394403402206439846;
282  m_positions[ 4] = -0.4484927510364468528779129;
283  m_positions[ 5] = -0.2304583159551347940655281;
284  m_positions[ 6] = 0.;
285  m_positions[ 7] = 0.2304583159551347940655281;
286  m_positions[ 8] = 0.4484927510364468528779129;
287  m_positions[ 9] = 0.6423493394403402206439846;
288  m_positions[10] = 0.8015780907333099127942065;
289  m_positions[11] = 0.9175983992229779652065478;
290  m_positions[12] = 0.9841830547185881494728294;
291  break;
292 #if 0
293  case 13:
294  m_weights [ 0] = ;
295  m_weights [ 1] = ;
296  m_weights [ 2] = ;
297  m_weights [ 3] = ;
298  m_weights [ 4] = ;
299  m_weights [ 5] = ;
300  m_weights [ 6] = ;
301  m_weights [ 7] = ;
302  m_weights [ 8] = ;
303  m_weights [ 9] = ;
304  m_weights [10] = ;
305  m_weights [11] = ;
306  m_weights [12] = ;
307  m_weights [13] = ;
308 
309  m_positions[ 0] = -;
310  m_positions[ 1] = -;
311  m_positions[ 2] = -;
312  m_positions[ 3] = -;
313  m_positions[ 4] = -;
314  m_positions[ 5] = -;
315  m_positions[ 6] = -;
316  m_positions[ 7] = ;
317  m_positions[ 8] = ;
318  m_positions[ 9] = ;
319  m_positions[10] = ;
320  m_positions[11] = ;
321  m_positions[12] = ;
322  m_positions[13] = ;
323  break;
324 #endif
325 #if 0
326  14 ±0.1080549487073436620662447 0.2152638534631577901958764
327  ±0.3191123689278897604356718 0.2051984637212956039659241
328  ±0.5152486363581540919652907 0.1855383974779378137417166
329  ±0.6872929048116854701480198 0.1572031671581935345696019
330  ±0.8272013150697649931897947 0.1215185706879031846894148
331  ±0.9284348836635735173363911 0.0801580871597602098056333
332  ±0.9862838086968123388415973 0.0351194603317518630318329
333 #endif
334 
335  case 16:
336  m_weights [ 0] = 0.0241483028685479319601100;
337  m_weights [ 1] = 0.0554595293739872011294402;
338  m_weights [ 2] = 0.0850361483171791808835354;
339  m_weights [ 3] = 0.1118838471934039710947884;
340  m_weights [ 4] = 0.1351363684685254732863200;
341  m_weights [ 5] = 0.1540457610768102880814316;
342  m_weights [ 6] = 0.1680041021564500445099707;
343  m_weights [ 7] = 0.1765627053669926463252710;
344  m_weights [ 8] = 0.1794464703562065254582656;
345  m_weights [ 9] = 0.1765627053669926463252710;
346  m_weights [10] = 0.1680041021564500445099707;
347  m_weights [11] = 0.1540457610768102880814316;
348  m_weights [12] = 0.1351363684685254732863200;
349  m_weights [13] = 0.1118838471934039710947884;
350  m_weights [14] = 0.0850361483171791808835354;
351  m_weights [15] = 0.0554595293739872011294402;
352  m_weights [16] = 0.0241483028685479319601100;
353 
354  m_positions[ 0] = -0.9905754753144173356754340;
355  m_positions[ 1] = -0.9506755217687677612227170;
356  m_positions[ 2] = -0.8802391537269859021229557;
357  m_positions[ 3] = -0.7815140038968014069252301;
358  m_positions[ 4] = -0.6576711592166907658503022;
359  m_positions[ 5] = -0.5126905370864769678862466;
360  m_positions[ 6] = -0.3512317634538763152971855;
361  m_positions[ 7] = -0.1784841814958478558506775;
362  m_positions[ 8] = 0.;
363  m_positions[ 9] = 0.1784841814958478558506775;
364  m_positions[10] = 0.3512317634538763152971855;
365  m_positions[11] = 0.5126905370864769678862466;
366  m_positions[12] = 0.6576711592166907658503022;
367  m_positions[13] = 0.7815140038968014069252301;
368  m_positions[14] = 0.8802391537269859021229557;
369  m_positions[15] = 0.9506755217687677612227170;
370  m_positions[16] = 0.9905754753144173356754340;
371  break;
372 
373  default:
374  std::cerr << "In UniformLegendre1DQuadrature::constructor()"
375  << ": m_order = " << m_order
376  << std::endl;
377  queso_error_msg("order not supported");
378  break;
379  }
380 
381  // Scale positions from the interval [-1, 1] to the interval [min,max]
382  for (unsigned int j = 0; j < m_positions.size(); ++j) {
384  if (densityIsNormalized) {
385  // Since \rho is "1/\Delta", we just multiply by ".5"
386  m_weights[j] *= .5;
387  }
388  else {
389  // Since \rho is "1", we multiply by ".5 * \Delta"
391  }
392  }
393 }
double maxDomainValue() const
Returns the maximum value of the domain of the (one-dimensional) function.
Definition: 1DQuadrature.C:61
double minDomainValue() const
Returns the minimum value of the domain of the (one-dimensional) function.
Definition: 1DQuadrature.C:55
Base1DQuadrature(double minDomainValue, double maxDomainValue, unsigned int order)
Default constructor.
Definition: 1DQuadrature.C:34
std::vector< double > m_weights
std::vector< double > m_positions
Definition: 1DQuadrature.h:83
unsigned int order() const
Returns the order of the quadrature rule.
Definition: 1DQuadrature.C:67
QUESO::UniformLegendre1DQuadrature::~UniformLegendre1DQuadrature ( )

Destructor.

Definition at line 395 of file 1DQuadrature.C.

396 {
397 }

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

Generated on Tue Jun 5 2018 19:49:19 for queso-0.57.1 by  doxygen 1.8.5