quadrature_gauss_lobatto_1D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 
22 // Local includes
24 
25 namespace libMesh
26 {
27 
29  unsigned int p)
30 {
31  //----------------------------------------------------------------------
32  // 1D quadrature rules
33  switch(_order + 2*p)
34  {
35  // Since Gauss-Lobatto rules must include the endpoints of the
36  // domain, there is no 1-point rule. The two-point
37  // Gauss-Lobatto rule is equivalent to the trapezoidal rule.
38  case CONSTANT:
39  case FIRST:
40  {
41  _points.resize (2);
42  _weights.resize(2);
43 
44  _points[0](0) = -1.0L;
45  _points[1] = -_points[0];
46 
47  _weights[0] = 1.;
48  _weights[1] = _weights[0];
49 
50  return;
51  }
52 
53  // The three-point Gauss-Lobatto rule is equivalent to Simpsons' rule.
54  // It can integrate cubic polynomials exactly.
55  case SECOND:
56  case THIRD:
57  {
58  _points.resize (3);
59  _weights.resize(3);
60 
61  _points[0](0) = -1.0L;
62  _points[1] = 0.0L;
63  _points[2] = -_points[0];
64 
65  _weights[0] = Real(1)/3;
66  _weights[1] = Real(4)/3;
67  _weights[2] = _weights[0];
68  return;
69  }
70 
71  // The four-point Gauss-Lobatto rule can integrate 2*4-3 =
72  // 5th-order polynomials exactly.
73  case FOURTH:
74  case FIFTH:
75  {
76  _points.resize (4);
77  _weights.resize(4);
78 
79  _points[ 0](0) = -1.0L;
80  _points[ 1](0) = -std::sqrt(Real(1)/5);
81  _points[ 2] = -_points[1];
82  _points[ 3] = -_points[0];
83 
84  _weights[ 0] = Real(1)/6;
85  _weights[ 1] = Real(5)/6;
86  _weights[ 2] = _weights[1];
87  _weights[ 3] = _weights[0];
88 
89  return;
90  }
91 
92  // The five-point Gauss-Lobatto rule can integrate 2*5-3 =
93  // 7th-order polynomials exactly.
94  case SIXTH:
95  case SEVENTH:
96  {
97  _points.resize (5);
98  _weights.resize(5);
99 
100  _points[ 0](0) = -1.0L;
101  _points[ 1](0) = -std::sqrt(Real(3)/7);
102  _points[ 2](0) = 0.;
103  _points[ 3] = -_points[1];
104  _points[ 4] = -_points[0];
105 
106  _weights[ 0] = 0.1;
107  _weights[ 1] = Real(49)/90;
108  _weights[ 2] = Real(32)/45;
109  _weights[ 3] = _weights[1];
110  _weights[ 4] = _weights[0];
111 
112  return;
113  }
114 
115  // 2*6-3 = 9
116  case EIGHTH:
117  case NINTH:
118  {
119  _points.resize (6);
120  _weights.resize(6);
121 
122  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
123  _points[ 1](0) = Real(-7.6505532392946469285100297395934e-01L);
124  _points[ 2](0) = Real(-2.8523151648064509631415099404088e-01L);
125  _points[ 3] = -_points[2];
126  _points[ 4] = -_points[1];
127  _points[ 5] = -_points[0];
128 
129  _weights[ 0] = Real(6.6666666666666666666666666666667e-02L);
130  _weights[ 1] = Real(3.7847495629784698031661280821202e-01L);
131  _weights[ 2] = Real(5.5485837703548635301672052512131e-01L);
132  _weights[ 3] = _weights[2];
133  _weights[ 4] = _weights[1];
134  _weights[ 5] = _weights[0];
135 
136  return;
137  }
138 
139  // 2*7-3 = 11
140  case TENTH:
141  case ELEVENTH:
142  {
143  _points.resize (7);
144  _weights.resize(7);
145 
146  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
147  _points[ 1](0) = Real(-8.3022389627856692987203221396747e-01L);
148  _points[ 2](0) = Real(-4.6884879347071421380377188190877e-01L);
149  _points[ 3](0) = 0.;
150  _points[ 4] = -_points[2];
151  _points[ 5] = -_points[1];
152  _points[ 6] = -_points[0];
153 
154  _weights[ 0] = Real(4.7619047619047619047619047619048e-02L);
155  _weights[ 1] = Real(2.7682604736156594801070040629007e-01L);
156  _weights[ 2] = Real(4.3174538120986262341787102228136e-01L);
157  _weights[ 3] = Real(4.8761904761904761904761904761905e-01L);
158  _weights[ 4] = _weights[2];
159  _weights[ 5] = _weights[1];
160  _weights[ 6] = _weights[0];
161 
162  return;
163  }
164 
165  // 2*8-3 = 13
166  case TWELFTH:
167  case THIRTEENTH:
168  {
169  _points.resize (8);
170  _weights.resize(8);
171 
172  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
173  _points[ 1](0) = Real(-8.7174014850960661533744576122066e-01L);
174  _points[ 2](0) = Real(-5.9170018143314230214451073139795e-01L);
175  _points[ 3](0) = Real(-2.0929921790247886876865726034535e-01L);
176  _points[ 4] = -_points[3];
177  _points[ 5] = -_points[2];
178  _points[ 6] = -_points[1];
179  _points[ 7] = -_points[0];
180 
181  _weights[ 0] = Real(3.5714285714285714285714285714286e-02L);
182  _weights[ 1] = Real(2.1070422714350603938299206577576e-01L);
183  _weights[ 2] = Real(3.4112269248350436476424067710775e-01L);
184  _weights[ 3] = Real(4.1245879465870388156705297140221e-01L);
185  _weights[ 4] = _weights[3];
186  _weights[ 5] = _weights[2];
187  _weights[ 6] = _weights[1];
188  _weights[ 7] = _weights[0];
189 
190  return;
191  }
192 
193  // 2*9-3 = 15
194  case FOURTEENTH:
195  case FIFTEENTH:
196  {
197  _points.resize (9);
198  _weights.resize(9);
199 
200  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
201  _points[ 1](0) = Real(-8.9975799541146015731234524441834e-01L);
202  _points[ 2](0) = Real(-6.7718627951073775344588542709134e-01L);
203  _points[ 3](0) = Real(-3.6311746382617815871075206870866e-01L);
204  _points[ 4](0) = 0.;
205  _points[ 5] = -_points[3];
206  _points[ 6] = -_points[2];
207  _points[ 7] = -_points[1];
208  _points[ 8] = -_points[0];
209 
210  _weights[ 0] = Real(2.7777777777777777777777777777778e-02L);
211  _weights[ 1] = Real(1.6549536156080552504633972002921e-01L);
212  _weights[ 2] = Real(2.7453871250016173528070561857937e-01L);
213  _weights[ 3] = Real(3.4642851097304634511513153213972e-01L);
214  _weights[ 4] = Real(3.7151927437641723356009070294785e-01L);
215  _weights[ 5] = _weights[3];
216  _weights[ 6] = _weights[2];
217  _weights[ 7] = _weights[1];
218  _weights[ 8] = _weights[0];
219 
220  return;
221  }
222 
223  // 2*10-3 = 17
224  case SIXTEENTH:
225  case SEVENTEENTH:
226  {
227  _points.resize (10);
228  _weights.resize(10);
229 
230  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
231  _points[ 1](0) = Real(-9.1953390816645881382893266082234e-01L);
232  _points[ 2](0) = Real(-7.3877386510550507500310617485983e-01L);
233  _points[ 3](0) = Real(-4.7792494981044449566117509273126e-01L);
234  _points[ 4](0) = Real(-1.6527895766638702462621976595817e-01L);
235  _points[ 5] = -_points[4];
236  _points[ 6] = -_points[3];
237  _points[ 7] = -_points[2];
238  _points[ 8] = -_points[1];
239  _points[ 9] = -_points[0];
240 
241  _weights[ 0] = Real(2.2222222222222222222222222222222e-02L);
242  _weights[ 1] = Real(1.3330599085107011112622717075539e-01L);
243  _weights[ 2] = Real(2.2488934206312645211945782173105e-01L);
244  _weights[ 3] = Real(2.9204268367968375787558225737444e-01L);
245  _weights[ 4] = Real(3.2753976118389745665651052791689e-01L);
246  _weights[ 5] = _weights[4];
247  _weights[ 6] = _weights[3];
248  _weights[ 7] = _weights[2];
249  _weights[ 8] = _weights[1];
250  _weights[ 9] = _weights[0];
251 
252  return;
253  }
254 
255  // 2*11-3 = 19
256  case EIGHTTEENTH:
257  case NINETEENTH:
258  {
259  _points.resize (11);
260  _weights.resize(11);
261 
262  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
263  _points[ 1](0) = Real(-9.3400143040805913433227413609938e-01L);
264  _points[ 2](0) = Real(-7.8448347366314441862241781610846e-01L);
265  _points[ 3](0) = Real(-5.6523532699620500647096396947775e-01L);
266  _points[ 4](0) = Real(-2.9575813558693939143191151555906e-01L);
267  _points[ 5](0) = 0.;
268  _points[ 6] = -_points[4];
269  _points[ 7] = -_points[3];
270  _points[ 8] = -_points[2];
271  _points[ 9] = -_points[1];
272  _points[10] = -_points[0];
273 
274  _weights[ 0] = Real(1.8181818181818181818181818181818e-02L);
275  _weights[ 1] = Real(1.0961227326699486446140344958035e-01L);
276  _weights[ 2] = Real(1.8716988178030520410814152189943e-01L);
277  _weights[ 3] = Real(2.4804810426402831404008486642187e-01L);
278  _weights[ 4] = Real(2.8687912477900808867922240333154e-01L);
279  _weights[ 5] = Real(3.0021759545569069378593188116998e-01L);
280  _weights[ 6] = _weights[4];
281  _weights[ 7] = _weights[3];
282  _weights[ 8] = _weights[2];
283  _weights[ 9] = _weights[1];
284  _weights[10] = _weights[0];
285 
286  return;
287  }
288 
289  // 2*12-3 = 21
290  case TWENTIETH:
291  case TWENTYFIRST:
292  {
293  _points.resize (12);
294  _weights.resize(12);
295 
296  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
297  _points[ 1](0) = Real(-9.4489927222288222340758013830322e-01L);
298  _points[ 2](0) = Real(-8.1927932164400667834864158171690e-01L);
299  _points[ 3](0) = Real(-6.3287615303186067766240485444366e-01L);
300  _points[ 4](0) = Real(-3.9953094096534893226434979156697e-01L);
301  _points[ 5](0) = Real(-1.3655293285492755486406185573969e-01L);
302  _points[ 6] = -_points[5];
303  _points[ 7] = -_points[4];
304  _points[ 8] = -_points[3];
305  _points[ 9] = -_points[2];
306  _points[10] = -_points[1];
307  _points[11] = -_points[0];
308 
309  _weights[ 0] = Real(1.5151515151515151515151515151515e-02L);
310  _weights[ 1] = Real(9.1684517413196130668342594134079e-02L);
311  _weights[ 2] = Real(1.5797470556437011516467106270034e-01L);
312  _weights[ 3] = Real(2.1250841776102114535830207736687e-01L);
313  _weights[ 4] = Real(2.5127560319920128029324441214760e-01L);
314  _weights[ 5] = Real(2.7140524091069617700028833849960e-01L);
315  _weights[ 6] = _weights[5];
316  _weights[ 7] = _weights[4];
317  _weights[ 8] = _weights[3];
318  _weights[ 9] = _weights[2];
319  _weights[10] = _weights[1];
320  _weights[11] = _weights[0];
321 
322  return;
323  }
324 
325  // 2*13-3 = 23
326  case TWENTYSECOND:
327  case TWENTYTHIRD:
328  {
329  _points.resize (13);
330  _weights.resize(13);
331 
332  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
333  _points[ 1](0) = Real(-9.5330984664216391189690546475545e-01L);
334  _points[ 2](0) = Real(-8.4634756465187231686592560709875e-01L);
335  _points[ 3](0) = Real(-6.8618846908175742607275903956636e-01L);
336  _points[ 4](0) = Real(-4.8290982109133620174693723363693e-01L);
337  _points[ 5](0) = Real(-2.4928693010623999256867370037423e-01L);
338  _points[ 6](0) = 0.;
339  _points[ 7] = -_points[5];
340  _points[ 8] = -_points[4];
341  _points[ 9] = -_points[3];
342  _points[10] = -_points[2];
343  _points[11] = -_points[1];
344  _points[12] = -_points[0];
345 
346  _weights[ 0] = Real(1.2820512820512820512820512820513e-02L);
347  _weights[ 1] = Real(7.7801686746818927793588988333134e-02L);
348  _weights[ 2] = Real(1.3498192668960834911991476258937e-01L);
349  _weights[ 3] = Real(1.8364686520355009200749425874681e-01L);
350  _weights[ 4] = Real(2.2076779356611008608553400837940e-01L);
351  _weights[ 5] = Real(2.4401579030667635645857814836016e-01L);
352  _weights[ 6] = Real(2.5193084933344673604413864154124e-01L);
353  _weights[ 7] = _weights[5];
354  _weights[ 8] = _weights[4];
355  _weights[ 9] = _weights[3];
356  _weights[10] = _weights[2];
357  _weights[11] = _weights[1];
358  _weights[12] = _weights[0];
359 
360  return;
361  }
362 
363  // 2*14-3 = 25
364  case TWENTYFOURTH:
365  case TWENTYFIFTH:
366  {
367  _points.resize (14);
368  _weights.resize(14);
369 
370  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
371  _points[ 1](0) = Real(-9.5993504526726090135510016201542e-01L);
372  _points[ 2](0) = Real(-8.6780105383034725100022020290826e-01L);
373  _points[ 3](0) = Real(-7.2886859909132614058467240052088e-01L);
374  _points[ 4](0) = Real(-5.5063940292864705531662270585908e-01L);
375  _points[ 5](0) = Real(-3.4272401334271284504390340364167e-01L);
376  _points[ 6](0) = Real(-1.1633186888370386765877670973616e-01L);
377  _points[ 7] = -_points[6];
378  _points[ 8] = -_points[5];
379  _points[ 9] = -_points[4];
380  _points[10] = -_points[3];
381  _points[11] = -_points[2];
382  _points[12] = -_points[1];
383  _points[13] = -_points[0];
384 
385  _weights[ 0] = Real(1.0989010989010989010989010989011e-02L);
386  _weights[ 1] = Real(6.6837284497681284634070660746053e-02L);
387  _weights[ 2] = Real(1.1658665589871165154099667065465e-01L);
388  _weights[ 3] = Real(1.6002185176295214241282099798759e-01L);
389  _weights[ 4] = Real(1.9482614937341611864033177837588e-01L);
390  _weights[ 5] = Real(2.1912625300977075487116252395417e-01L);
391  _weights[ 6] = Real(2.3161279446845705888962835729264e-01L);
392  _weights[ 7] = _weights[6];
393  _weights[ 8] = _weights[5];
394  _weights[ 9] = _weights[4];
395  _weights[10] = _weights[3];
396  _weights[11] = _weights[2];
397  _weights[12] = _weights[1];
398  _weights[13] = _weights[0];
399 
400  return;
401  }
402 
403  // 2*15-3 = 27
404  case TWENTYSIXTH:
405  case TWENTYSEVENTH:
406  {
407  _points.resize (15);
408  _weights.resize(15);
409 
410  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
411  _points[ 1](0) = Real(-9.6524592650383857279585139206960e-01L);
412  _points[ 2](0) = Real(-8.8508204422297629882540163148223e-01L);
413  _points[ 3](0) = Real(-7.6351968995181520070411847597629e-01L);
414  _points[ 4](0) = Real(-6.0625320546984571112352993863673e-01L);
415  _points[ 5](0) = Real(-4.2063805471367248092189693873858e-01L);
416  _points[ 6](0) = Real(-2.1535395536379423822567944627292e-01L);
417  _points[ 7](0) = 0.;
418  _points[ 8] = -_points[6];
419  _points[ 9] = -_points[5];
420  _points[10] = -_points[4];
421  _points[11] = -_points[3];
422  _points[12] = -_points[2];
423  _points[13] = -_points[1];
424  _points[14] = -_points[0];
425 
426  _weights[ 0] = Real(9.5238095238095238095238095238095e-03L);
427  _weights[ 1] = Real(5.8029893028601249096880584025282e-02L);
428  _weights[ 2] = Real(1.0166007032571806760366617078880e-01L);
429  _weights[ 3] = Real(1.4051169980242810946044680564367e-01L);
430  _weights[ 4] = Real(1.7278964725360094905207709940835e-01L);
431  _weights[ 5] = Real(1.9698723596461335609250034650741e-01L);
432  _weights[ 6] = Real(2.1197358592682092012743007697722e-01L);
433  _weights[ 7] = Real(2.1704811634881564951495021425091e-01L);
434  _weights[ 8] = _weights[6];
435  _weights[ 9] = _weights[5];
436  _weights[10] = _weights[4];
437  _weights[11] = _weights[3];
438  _weights[12] = _weights[2];
439  _weights[13] = _weights[1];
440  _weights[14] = _weights[0];
441 
442  return;
443  }
444 
445  // 2*16-3 = 29
446  case TWENTYEIGHTH:
447  case TWENTYNINTH:
448  {
449  _points.resize (16);
450  _weights.resize(16);
451 
452  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
453  _points[ 1](0) = Real(-9.6956804627021793295224273836746e-01L);
454  _points[ 2](0) = Real(-8.9920053309347209299462826151985e-01L);
455  _points[ 3](0) = Real(-7.9200829186181506393108827096315e-01L);
456  _points[ 4](0) = Real(-6.5238870288249308946788321964058e-01L);
457  _points[ 5](0) = Real(-4.8605942188713761178189078584687e-01L);
458  _points[ 6](0) = Real(-2.9983046890076320809835345472230e-01L);
459  _points[ 7](0) = Real(-1.0132627352194944784303300504592e-01L);
460  _points[ 8] = -_points[7];
461  _points[ 9] = -_points[6];
462  _points[10] = -_points[5];
463  _points[11] = -_points[4];
464  _points[12] = -_points[3];
465  _points[13] = -_points[2];
466  _points[14] = -_points[1];
467  _points[15] = -_points[0];
468 
469  _weights[ 0] = Real(8.3333333333333333333333333333333e-03L);
470  _weights[ 1] = Real(5.0850361005919905403244919565455e-02L);
471  _weights[ 2] = Real(8.9393697325930800991052080166084e-02L);
472  _weights[ 3] = Real(1.2425538213251409834953633265731e-01L);
473  _weights[ 4] = Real(1.5402698080716428081564494048499e-01L);
474  _weights[ 5] = Real(1.7749191339170412530107566952836e-01L);
475  _weights[ 6] = Real(1.9369002382520358431691359885352e-01L);
476  _weights[ 7] = Real(2.0195830817822987148919912541094e-01L);
477  _weights[ 8] = _weights[7];
478  _weights[ 9] = _weights[6];
479  _weights[10] = _weights[5];
480  _weights[11] = _weights[4];
481  _weights[12] = _weights[3];
482  _weights[13] = _weights[2];
483  _weights[14] = _weights[1];
484  _weights[15] = _weights[0];
485 
486  return;
487  }
488 
489  // 2*17-3 = 31
490  case THIRTIETH:
491  case THIRTYFIRST:
492  {
493  _points.resize (17);
494  _weights.resize(17);
495 
496  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
497  _points[ 1](0) = Real(-9.7313217663141831415697950187372e-01L);
498  _points[ 2](0) = Real(-9.1087999591557359562380250639773e-01L);
499  _points[ 3](0) = Real(-8.1569625122177030710675055323753e-01L);
500  _points[ 4](0) = Real(-6.9102898062768470539491935737245e-01L);
501  _points[ 5](0) = Real(-5.4138539933010153912373340750406e-01L);
502  _points[ 6](0) = Real(-3.7217443356547704190723468073526e-01L);
503  _points[ 7](0) = Real(-1.8951197351831738830426301475311e-01L);
504  _points[ 8](0) = 0.;
505  _points[ 9] = -_points[7];
506  _points[10] = -_points[6];
507  _points[11] = -_points[5];
508  _points[12] = -_points[4];
509  _points[13] = -_points[3];
510  _points[14] = -_points[2];
511  _points[15] = -_points[1];
512  _points[16] = -_points[0];
513 
514  _weights[ 0] = Real(7.3529411764705882352941176470588e-03L);
515  _weights[ 1] = Real(4.4921940543254209647400954623212e-02L);
516  _weights[ 2] = Real(7.9198270503687119190264429952835e-02L);
517  _weights[ 3] = Real(1.1059290900702816137577270522008e-01L);
518  _weights[ 4] = Real(1.3798774620192655905620157495403e-01L);
519  _weights[ 5] = Real(1.6039466199762153951632836586475e-01L);
520  _weights[ 6] = Real(1.7700425351565787043694574536329e-01L);
521  _weights[ 7] = Real(1.8721633967761923589208848286062e-01L);
522  _weights[ 8] = Real(1.9066187475346943329940724702825e-01L);
523  _weights[ 9] = _weights[7];
524  _weights[10] = _weights[6];
525  _weights[11] = _weights[5];
526  _weights[12] = _weights[4];
527  _weights[13] = _weights[3];
528  _weights[14] = _weights[2];
529  _weights[15] = _weights[1];
530  _weights[16] = _weights[0];
531 
532  return;
533  }
534 
535  // 2*18-3 = 33
536  case THIRTYSECOND:
537  case THIRTYTHIRD:
538  {
539  _points.resize (18);
540  _weights.resize(18);
541 
542  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
543  _points[ 1](0) = Real(-9.7610555741219854286451892434170e-01L);
544  _points[ 2](0) = Real(-9.2064918534753387383785462543128e-01L);
545  _points[ 3](0) = Real(-8.3559353521809021371364636232794e-01L);
546  _points[ 4](0) = Real(-7.2367932928324268130621036530207e-01L);
547  _points[ 5](0) = Real(-5.8850483431866176117353589319356e-01L);
548  _points[ 6](0) = Real(-4.3441503691212397534228713674067e-01L);
549  _points[ 7](0) = Real(-2.6636265287828098416766533202560e-01L);
550  _points[ 8](0) = Real(-8.9749093484652111022645010088562e-02L);
551  _points[ 9] = -_points[8];
552  _points[10] = -_points[7];
553  _points[11] = -_points[6];
554  _points[12] = -_points[5];
555  _points[13] = -_points[4];
556  _points[14] = -_points[3];
557  _points[15] = -_points[2];
558  _points[16] = -_points[1];
559  _points[17] = -_points[0];
560 
561  _weights[ 0] = Real(6.5359477124183006535947712418301e-03L);
562  _weights[ 1] = Real(3.9970628810914066137599176410101e-02L);
563  _weights[ 2] = Real(7.0637166885633664999222960167786e-02L);
564  _weights[ 3] = Real(9.9016271717502802394423605318672e-02L);
565  _weights[ 4] = Real(1.2421053313296710026339635889675e-01L);
566  _weights[ 5] = Real(1.4541196157380226798300321049443e-01L);
567  _weights[ 6] = Real(1.6193951723760248926432670670023e-01L);
568  _weights[ 7] = Real(1.7326210948945622601061440382668e-01L);
569  _weights[ 8] = Real(1.7901586343970308229381880694353e-01L);
570  _weights[ 9] = _weights[8];
571  _weights[10] = _weights[7];
572  _weights[11] = _weights[6];
573  _weights[12] = _weights[5];
574  _weights[13] = _weights[4];
575  _weights[14] = _weights[3];
576  _weights[15] = _weights[2];
577  _weights[16] = _weights[1];
578  _weights[17] = _weights[0];
579 
580  return;
581  }
582 
583  // 2*19-3 = 35
584  case THIRTYFOURTH:
585  case THIRTYFIFTH:
586  {
587  _points.resize (19);
588  _weights.resize(19);
589 
590  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
591  _points[ 1](0) = Real(-9.7861176622208009515263406311022e-01L);
592  _points[ 2](0) = Real(-9.2890152815258624371794025879655e-01L);
593  _points[ 3](0) = Real(-8.5246057779664609308595597004106e-01L);
594  _points[ 4](0) = Real(-7.5149420255261301416363748963394e-01L);
595  _points[ 5](0) = Real(-6.2890813726522049776683230622873e-01L);
596  _points[ 6](0) = Real(-4.8822928568071350277790963762492e-01L);
597  _points[ 7](0) = Real(-3.3350484782449861029850010384493e-01L);
598  _points[ 8](0) = Real(-1.6918602340928157137515415344488e-01L);
599  _points[ 9](0) = 0.;
600  _points[10] = -_points[8];
601  _points[11] = -_points[7];
602  _points[12] = -_points[6];
603  _points[13] = -_points[5];
604  _points[14] = -_points[4];
605  _points[15] = -_points[3];
606  _points[16] = -_points[2];
607  _points[17] = -_points[1];
608  _points[18] = -_points[0];
609 
610  _weights[ 0] = Real(5.8479532163742690058479532163743e-03L);
611  _weights[ 1] = Real(3.5793365186176477115425569035122e-02L);
612  _weights[ 2] = Real(6.3381891762629736851695690418317e-02L);
613  _weights[ 3] = Real(8.9131757099207084448008790556153e-02L);
614  _weights[ 4] = Real(1.1231534147730504407091001546378e-01L);
615  _weights[ 5] = Real(1.3226728044875077692604673390973e-01L);
616  _weights[ 6] = Real(1.4841394259593888500968064366841e-01L);
617  _weights[ 7] = Real(1.6029092404406124197991096818359e-01L);
618  _weights[ 8] = Real(1.6755658452714286727013727774026e-01L);
619  _weights[ 9] = Real(1.7000191928482723464467271561652e-01L);
620  _weights[10] = _weights[8];
621  _weights[11] = _weights[7];
622  _weights[12] = _weights[6];
623  _weights[13] = _weights[5];
624  _weights[14] = _weights[4];
625  _weights[15] = _weights[3];
626  _weights[16] = _weights[2];
627  _weights[17] = _weights[1];
628  _weights[18] = _weights[0];
629 
630  return;
631  }
632 
633  // 2*20-3 = 37
634  case THIRTYSIXTH:
635  case THIRTYSEVENTH:
636  {
637  _points.resize (20);
638  _weights.resize(20);
639 
640  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
641  _points[ 1](0) = Real(-9.8074370489391417192544643858423e-01L);
642  _points[ 2](0) = Real(-9.3593449881266543571618158493063e-01L);
643  _points[ 3](0) = Real(-8.6687797808995014130984721461629e-01L);
644  _points[ 4](0) = Real(-7.7536826095205587041431752759469e-01L);
645  _points[ 5](0) = Real(-6.6377640229031128984640332297116e-01L);
646  _points[ 6](0) = Real(-5.3499286403188626164813596182898e-01L);
647  _points[ 7](0) = Real(-3.9235318371390929938647470381582e-01L);
648  _points[ 8](0) = Real(-2.3955170592298649518240135692709e-01L);
649  _points[ 9](0) = Real(-8.0545937238821837975944518159554e-02L);
650  _points[10] = -_points[9];
651  _points[11] = -_points[8];
652  _points[12] = -_points[7];
653  _points[13] = -_points[6];
654  _points[14] = -_points[5];
655  _points[15] = -_points[4];
656  _points[16] = -_points[3];
657  _points[17] = -_points[2];
658  _points[18] = -_points[1];
659  _points[19] = -_points[0];
660 
661  _weights[ 0] = Real(5.2631578947368421052631578947368e-03L);
662  _weights[ 1] = Real(3.2237123188488941491605028117294e-02L);
663  _weights[ 2] = Real(5.7181802127566826004753627173243e-02L);
664  _weights[ 3] = Real(8.0631763996119603144776846113721e-02L);
665  _weights[ 4] = Real(1.0199149969945081568378120573289e-01L);
666  _weights[ 5] = Real(1.2070922762867472509942970500239e-01L);
667  _weights[ 6] = Real(1.3630048235872418448978079298903e-01L);
668  _weights[ 7] = Real(1.4836155407091682581471301373397e-01L);
669  _weights[ 8] = Real(1.5658010264747548715816989679364e-01L);
670  _weights[ 9] = Real(1.6074328638784574900772672644908e-01L);
671  _weights[10] = _weights[9];
672  _weights[11] = _weights[8];
673  _weights[12] = _weights[7];
674  _weights[13] = _weights[6];
675  _weights[14] = _weights[5];
676  _weights[15] = _weights[4];
677  _weights[16] = _weights[3];
678  _weights[17] = _weights[2];
679  _weights[18] = _weights[1];
680  _weights[19] = _weights[0];
681 
682  return;
683  }
684 
685  // 2*21-3 = 39
686  case THIRTYEIGHTH:
687  case THIRTYNINTH:
688  {
689  _points.resize (21);
690  _weights.resize(21);
691 
692  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
693  _points[ 1](0) = Real(-9.8257229660454802823448127655541e-01L);
694  _points[ 2](0) = Real(-9.4197629695974553429610265066144e-01L);
695  _points[ 3](0) = Real(-8.7929475532359046445115359630494e-01L);
696  _points[ 4](0) = Real(-7.9600192607771240474431258966036e-01L);
697  _points[ 5](0) = Real(-6.9405102606222323262731639319467e-01L);
698  _points[ 6](0) = Real(-5.7583196026183068692702187033809e-01L);
699  _points[ 7](0) = Real(-4.4411578327900210119451634960735e-01L);
700  _points[ 8](0) = Real(-3.0198985650876488727535186785875e-01L);
701  _points[ 9](0) = Real(-1.5278551580218546600635832848567e-01L);
702  _points[10](0) = 0.;
703  _points[11] = -_points[9];
704  _points[12] = -_points[8];
705  _points[13] = -_points[7];
706  _points[14] = -_points[6];
707  _points[15] = -_points[5];
708  _points[16] = -_points[4];
709  _points[17] = -_points[3];
710  _points[18] = -_points[2];
711  _points[19] = -_points[1];
712  _points[20] = -_points[0];
713 
714  _weights[ 0] = Real(4.7619047619047619047619047619048e-03L);
715  _weights[ 1] = Real(2.9184840098505458609458543613171e-02L);
716  _weights[ 2] = Real(5.1843169000849625072722971852830e-02L);
717  _weights[ 3] = Real(7.3273918185074144252547861041894e-02L);
718  _weights[ 4] = Real(9.2985467957886065301137664149214e-02L);
719  _weights[ 5] = Real(1.1051708321912333526700048678439e-01L);
720  _weights[ 6] = Real(1.2545812119086894801515753570800e-01L);
721  _weights[ 7] = Real(1.3745846286004134358089961741515e-01L);
722  _weights[ 8] = Real(1.4623686244797745926727053063439e-01L);
723  _weights[ 9] = Real(1.5158757511168138445325068150529e-01L);
724  _weights[10] = Real(1.5338519033217494855158440506754e-01L);
725  _weights[11] = _weights[9];
726  _weights[12] = _weights[8];
727  _weights[13] = _weights[7];
728  _weights[14] = _weights[6];
729  _weights[15] = _weights[5];
730  _weights[16] = _weights[4];
731  _weights[17] = _weights[3];
732  _weights[18] = _weights[2];
733  _weights[19] = _weights[1];
734  _weights[20] = _weights[0];
735 
736  return;
737  }
738 
739  // 2*22-3 = 41
740  case FORTIETH:
741  case FORTYFIRST:
742  {
743  _points.resize (22);
744  _weights.resize(22);
745 
746  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
747  _points[ 1](0) = Real(-9.8415243845764617655228962221207e-01L);
748  _points[ 2](0) = Real(-9.4720428399922868052421376661573e-01L);
749  _points[ 3](0) = Real(-8.9006229019090447052965782577909e-01L);
750  _points[ 4](0) = Real(-8.1394892761192113604544184805614e-01L);
751  _points[ 5](0) = Real(-7.2048723996120215811988189639847e-01L);
752  _points[ 6](0) = Real(-6.1166943828425897122621160586993e-01L);
753  _points[ 7](0) = Real(-4.8981487518990234980875123568327e-01L);
754  _points[ 8](0) = Real(-3.5752071013891953806095728024018e-01L);
755  _points[ 9](0) = Real(-2.1760658515928504178795509346539e-01L);
756  _points[10](0) = Real(-7.3054540010898334761088790464107e-02L);
757  _points[11] = -_points[10];
758  _points[12] = -_points[9];
759  _points[13] = -_points[8];
760  _points[14] = -_points[7];
761  _points[15] = -_points[6];
762  _points[16] = -_points[5];
763  _points[17] = -_points[4];
764  _points[18] = -_points[3];
765  _points[19] = -_points[2];
766  _points[20] = -_points[1];
767  _points[21] = -_points[0];
768 
769  _weights[ 0] = Real(4.3290043290043290043290043290043e-03L);
770  _weights[ 1] = Real(2.6545747682501757911627904520543e-02L);
771  _weights[ 2] = Real(4.7214465293740752123775734864792e-02L);
772  _weights[ 3] = Real(6.6865605864553076012404194157097e-02L);
773  _weights[ 4] = Real(8.5090060391838447815711236095748e-02L);
774  _weights[ 5] = Real(1.0150057480164767437243730374960e-01L);
775  _weights[ 6] = Real(1.1574764465393906659003636772146e-01L);
776  _weights[ 7] = Real(1.2752769665343027553084445930883e-01L);
777  _weights[ 8] = Real(1.3658968861374142668617736220617e-01L);
778  _weights[ 9] = Real(1.4274049227136140033623599356679e-01L);
779  _weights[10] = Real(1.4584901944424179361642043947997e-01L);
780  _weights[11] = _weights[10];
781  _weights[12] = _weights[9];
782  _weights[13] = _weights[8];
783  _weights[14] = _weights[7];
784  _weights[15] = _weights[6];
785  _weights[16] = _weights[5];
786  _weights[17] = _weights[4];
787  _weights[18] = _weights[3];
788  _weights[19] = _weights[2];
789  _weights[20] = _weights[1];
790  _weights[21] = _weights[0];
791 
792  return;
793  }
794 
795  // 2*23-3 = 43
796  case FORTYSECOND:
797  case FORTYTHIRD:
798  {
799  _points.resize (23);
800  _weights.resize(23);
801 
802  _points[ 0](0) = Real(-1.0000000000000000000000000000000e+00L);
803  _points[ 1](0) = Real(-9.8552715587873257808146276673810e-01L);
804  _points[ 2](0) = Real(-9.5175795571071020413563967985143e-01L);
805  _points[ 3](0) = Real(-8.9945855804034501095016032034737e-01L);
806  _points[ 4](0) = Real(-8.2965109665128588622320061929000e-01L);
807  _points[ 5](0) = Real(-7.4369504117206068394516354306700e-01L);
808  _points[ 6](0) = Real(-6.4326364446013620847614553360277e-01L);
809  _points[ 7](0) = Real(-5.3031177113684416813011532015230e-01L);
810  _points[ 8](0) = Real(-4.0703793791447482919595048821510e-01L);
811  _points[ 9](0) = Real(-2.7584154894579306710687763267914e-01L);
812  _points[10](0) = Real(-1.3927620404066839859186261298277e-01L);
813  _points[11](0) = 0.;
814  _points[12] = -_points[10];
815  _points[13] = -_points[9];
816  _points[14] = -_points[8];
817  _points[15] = -_points[7];
818  _points[16] = -_points[6];
819  _points[17] = -_points[5];
820  _points[18] = -_points[4];
821  _points[19] = -_points[3];
822  _points[20] = -_points[2];
823  _points[21] = -_points[1];
824  _points[22] = -_points[0];
825 
826  _weights[ 0] = Real(3.9525691699604743083003952569170e-03L);
827  _weights[ 1] = Real(2.4248600771531736517399658937097e-02L);
828  _weights[ 2] = Real(4.3175871170241834748876465612042e-02L);
829  _weights[ 3] = Real(6.1252477129554206381382847440355e-02L);
830  _weights[ 4] = Real(7.8135449475569989741934255347965e-02L);
831  _weights[ 5] = Real(9.3497246163512341833500706906697e-02L);
832  _weights[ 6] = Real(1.0703910172433651153518362791547e-01L);
833  _weights[ 7] = Real(1.1849751066274913130212600472426e-01L);
834  _weights[ 8] = Real(1.2764947470175887663614855305567e-01L);
835  _weights[ 9] = Real(1.3431687263860381990156489770071e-01L);
836  _weights[10] = Real(1.3836993638580739452350273386294e-01L);
837  _weights[11] = Real(1.3972978001274736514015970647975e-01L);
838  _weights[12] = _weights[10];
839  _weights[13] = _weights[9];
840  _weights[14] = _weights[8];
841  _weights[15] = _weights[7];
842  _weights[16] = _weights[6];
843  _weights[17] = _weights[5];
844  _weights[18] = _weights[4];
845  _weights[19] = _weights[3];
846  _weights[20] = _weights[2];
847  _weights[21] = _weights[1];
848  _weights[22] = _weights[0];
849 
850  return;
851  }
852 
853  default:
854  libmesh_error_msg("Quadrature rule " << _order << " not supported!");
855  }
856 }
857 
858 } // namespace libMesh
virtual void init_1D(const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real