queso-0.53.0
List of all members
QUESO::UniformLegendre1DQuadrature Class Reference

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

#include <1DQuadrature.h>

Inheritance diagram for QUESO::UniformLegendre1DQuadrature:
Inheritance graph
[legend]
Collaboration diagram for QUESO::UniformLegendre1DQuadrature:
Collaboration graph
[legend]

Public Member Functions

Constructor/Destructor methods
 UniformLegendre1DQuadrature (double minDomainValue, double maxDomainValue, unsigned int order, bool densityIsNormalized)
 Default constructor. More...
 
 ~UniformLegendre1DQuadrature ()
 Destructor. More...
 
Mathematical methods
void dumbRoutine () const
 A bogus method. More...
 
- Public Member Functions inherited from QUESO::Base1DQuadrature
 Base1DQuadrature (double minDomainValue, double maxDomainValue, unsigned int order)
 Default constructor. More...
 
virtual ~Base1DQuadrature ()
 Virtual destructor. 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...
 
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
 
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 162 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 118 of file 1DQuadrature.C.

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

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

Destructor.

Definition at line 415 of file 1DQuadrature.C.

416 {
417 }

Member Function Documentation

void QUESO::UniformLegendre1DQuadrature::dumbRoutine ( ) const
virtual

A bogus method.

Implements QUESO::Base1DQuadrature.

Definition at line 420 of file 1DQuadrature.C.

421 {
422  return;
423 }

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

Generated on Thu Jun 11 2015 13:52:35 for queso-0.53.0 by  doxygen 1.8.5