quadrature_gauss_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 
28 
29 
31  unsigned int p)
32 {
33  //----------------------------------------------------------------------
34  // 1D quadrature rules
35  switch(_order + 2*p)
36  {
37  case CONSTANT:
38  case FIRST:
39  {
40  _points.resize (1);
41  _weights.resize(1);
42 
43  _points[0](0) = 0.;
44 
45  _weights[0] = 2.;
46 
47  return;
48  }
49  case SECOND:
50  case THIRD:
51  {
52  _points.resize (2);
53  _weights.resize(2);
54 
55  _points[0](0) = Real(-5.7735026918962576450914878050196e-01L); // -sqrt(3)/3
56  _points[1] = -_points[0];
57 
58  _weights[0] = 1.;
59  _weights[1] = _weights[0];
60 
61  return;
62  }
63  case FOURTH:
64  case FIFTH:
65  {
66  _points.resize (3);
67  _weights.resize(3);
68 
69  _points[ 0](0) = Real(-7.7459666924148337703585307995648e-01L);
70  _points[ 1](0) = 0.;
71  _points[ 2] = -_points[0];
72 
73  _weights[ 0] = Real(5.5555555555555555555555555555556e-01L);
74  _weights[ 1] = Real(8.8888888888888888888888888888889e-01L);
75  _weights[ 2] = _weights[0];
76 
77  return;
78  }
79  case SIXTH:
80  case SEVENTH:
81  {
82  _points.resize (4);
83  _weights.resize(4);
84 
85  _points[ 0](0) = Real(-8.6113631159405257522394648889281e-01L);
86  _points[ 1](0) = Real(-3.3998104358485626480266575910324e-01L);
87  _points[ 2] = -_points[1];
88  _points[ 3] = -_points[0];
89 
90  _weights[ 0] = Real(3.4785484513745385737306394922200e-01L);
91  _weights[ 1] = Real(6.5214515486254614262693605077800e-01L);
92  _weights[ 2] = _weights[1];
93  _weights[ 3] = _weights[0];
94 
95  return;
96  }
97  case EIGHTH:
98  case NINTH:
99  {
100  _points.resize (5);
101  _weights.resize(5);
102 
103  _points[ 0](0) = Real(-9.0617984593866399279762687829939e-01L);
104  _points[ 1](0) = Real(-5.3846931010568309103631442070021e-01L);
105  _points[ 2](0) = 0.;
106  _points[ 3] = -_points[1];
107  _points[ 4] = -_points[0];
108 
109  _weights[ 0] = Real(2.3692688505618908751426404071992e-01L);
110  _weights[ 1] = Real(4.7862867049936646804129151483564e-01L);
111  _weights[ 2] = Real(5.6888888888888888888888888888889e-01L);
112  _weights[ 3] = _weights[1];
113  _weights[ 4] = _weights[0];
114 
115  return;
116  }
117  case TENTH:
118  case ELEVENTH:
119  {
120  _points.resize (6);
121  _weights.resize(6);
122 
123  _points[ 0](0) = Real(-9.3246951420315202781230155449399e-01L);
124  _points[ 1](0) = Real(-6.6120938646626451366139959501991e-01L);
125  _points[ 2](0) = Real(-2.3861918608319690863050172168071e-01L);
126  _points[ 3] = -_points[2];
127  _points[ 4] = -_points[1];
128  _points[ 5] = -_points[0];
129 
130  _weights[ 0] = Real(1.7132449237917034504029614217273e-01L);
131  _weights[ 1] = Real(3.6076157304813860756983351383772e-01L);
132  _weights[ 2] = Real(4.6791393457269104738987034398955e-01L);
133  _weights[ 3] = _weights[2];
134  _weights[ 4] = _weights[1];
135  _weights[ 5] = _weights[0];
136 
137  return;
138  }
139  case TWELFTH:
140  case THIRTEENTH:
141  {
142  _points.resize (7);
143  _weights.resize(7);
144 
145  _points[ 0](0) = Real(-9.4910791234275852452618968404785e-01L);
146  _points[ 1](0) = Real(-7.4153118559939443986386477328079e-01L);
147  _points[ 2](0) = Real(-4.0584515137739716690660641207696e-01L);
148  _points[ 3](0) = 0.;
149  _points[ 4] = -_points[2];
150  _points[ 5] = -_points[1];
151  _points[ 6] = -_points[0];
152 
153  _weights[ 0] = Real(1.2948496616886969327061143267908e-01L);
154  _weights[ 1] = Real(2.7970539148927666790146777142378e-01L);
155  _weights[ 2] = Real(3.8183005050511894495036977548898e-01L);
156  _weights[ 3] = Real(4.1795918367346938775510204081633e-01L);
157  _weights[ 4] = _weights[2];
158  _weights[ 5] = _weights[1];
159  _weights[ 6] = _weights[0];
160 
161  return;
162  }
163  case FOURTEENTH:
164  case FIFTEENTH:
165  {
166  _points.resize (8);
167  _weights.resize(8);
168 
169  _points[ 0](0) = Real(-9.6028985649753623168356086856947e-01L);
170  _points[ 1](0) = Real(-7.9666647741362673959155393647583e-01L);
171  _points[ 2](0) = Real(-5.2553240991632898581773904918925e-01L);
172  _points[ 3](0) = Real(-1.8343464249564980493947614236018e-01L);
173  _points[ 4] = -_points[3];
174  _points[ 5] = -_points[2];
175  _points[ 6] = -_points[1];
176  _points[ 7] = -_points[0];
177 
178  _weights[ 0] = Real(1.0122853629037625915253135430996e-01L);
179  _weights[ 1] = Real(2.2238103445337447054435599442624e-01L);
180  _weights[ 2] = Real(3.1370664587788728733796220198660e-01L);
181  _weights[ 3] = Real(3.6268378337836198296515044927720e-01L);
182  _weights[ 4] = _weights[3];
183  _weights[ 5] = _weights[2];
184  _weights[ 6] = _weights[1];
185  _weights[ 7] = _weights[0];
186 
187  return;
188  }
189  case SIXTEENTH:
190  case SEVENTEENTH:
191  {
192  _points.resize (9);
193  _weights.resize(9);
194 
195  _points[ 0](0) = Real(-9.6816023950762608983557620290367e-01L);
196  _points[ 1](0) = Real(-8.3603110732663579429942978806973e-01L);
197  _points[ 2](0) = Real(-6.1337143270059039730870203934147e-01L);
198  _points[ 3](0) = Real(-3.2425342340380892903853801464334e-01L);
199  _points[ 4](0) = 0.;
200  _points[ 5] = -_points[3];
201  _points[ 6] = -_points[2];
202  _points[ 7] = -_points[1];
203  _points[ 8] = -_points[0];
204 
205  _weights[ 0] = Real(8.1274388361574411971892158110524e-02L);
206  _weights[ 1] = Real(1.8064816069485740405847203124291e-01L);
207  _weights[ 2] = Real(2.6061069640293546231874286941863e-01L);
208  _weights[ 3] = Real(3.1234707704000284006863040658444e-01L);
209  _weights[ 4] = Real(3.3023935500125976316452506928697e-01L);
210  _weights[ 5] = _weights[3];
211  _weights[ 6] = _weights[2];
212  _weights[ 7] = _weights[1];
213  _weights[ 8] = _weights[0];
214 
215  return;
216  }
217  case EIGHTTEENTH:
218  case NINETEENTH:
219  {
220  _points.resize (10);
221  _weights.resize(10);
222 
223  _points[ 0](0) = Real(-9.7390652851717172007796401208445e-01L);
224  _points[ 1](0) = Real(-8.6506336668898451073209668842349e-01L);
225  _points[ 2](0) = Real(-6.7940956829902440623432736511487e-01L);
226  _points[ 3](0) = Real(-4.3339539412924719079926594316578e-01L);
227  _points[ 4](0) = Real(-1.4887433898163121088482600112972e-01L);
228  _points[ 5] = -_points[4];
229  _points[ 6] = -_points[3];
230  _points[ 7] = -_points[2];
231  _points[ 8] = -_points[1];
232  _points[ 9] = -_points[0];
233 
234  _weights[ 0] = Real(6.6671344308688137593568809893332e-02L);
235  _weights[ 1] = Real(1.4945134915058059314577633965770e-01L);
236  _weights[ 2] = Real(2.1908636251598204399553493422816e-01L);
237  _weights[ 3] = Real(2.6926671930999635509122692156947e-01L);
238  _weights[ 4] = Real(2.9552422471475287017389299465134e-01L);
239  _weights[ 5] = _weights[4];
240  _weights[ 6] = _weights[3];
241  _weights[ 7] = _weights[2];
242  _weights[ 8] = _weights[1];
243  _weights[ 9] = _weights[0];
244 
245  return;
246  }
247 
248  case TWENTIETH:
249  case TWENTYFIRST:
250  {
251  _points.resize (11);
252  _weights.resize(11);
253 
254  _points[ 0](0) = Real(-9.7822865814605699280393800112286e-01L);
255  _points[ 1](0) = Real(-8.8706259976809529907515776930393e-01L);
256  _points[ 2](0) = Real(-7.3015200557404932409341625203115e-01L);
257  _points[ 3](0) = Real(-5.1909612920681181592572566945861e-01L);
258  _points[ 4](0) = Real(-2.6954315595234497233153198540086e-01L);
259  _points[ 5](0) = 0.;
260  _points[ 6] = -_points[4];
261  _points[ 7] = -_points[3];
262  _points[ 8] = -_points[2];
263  _points[ 9] = -_points[1];
264  _points[10] = -_points[0];
265 
266  _weights[ 0] = Real(5.5668567116173666482753720442549e-02L);
267  _weights[ 1] = Real(1.2558036946490462463469429922394e-01L);
268  _weights[ 2] = Real(1.8629021092773425142609764143166e-01L);
269  _weights[ 3] = Real(2.3319376459199047991852370484318e-01L);
270  _weights[ 4] = Real(2.6280454451024666218068886989051e-01L);
271  _weights[ 5] = Real(2.7292508677790063071448352833634e-01L);
272  _weights[ 6] = _weights[4];
273  _weights[ 7] = _weights[3];
274  _weights[ 8] = _weights[2];
275  _weights[ 9] = _weights[1];
276  _weights[10] = _weights[0];
277 
278  return;
279  }
280 
281  case TWENTYSECOND:
282  case TWENTYTHIRD:
283  {
284  _points.resize (12);
285  _weights.resize(12);
286 
287  _points[ 0](0) = Real(-9.8156063424671925069054909014928e-01L);
288  _points[ 1](0) = Real(-9.0411725637047485667846586611910e-01L);
289  _points[ 2](0) = Real(-7.6990267419430468703689383321282e-01L);
290  _points[ 3](0) = Real(-5.8731795428661744729670241894053e-01L);
291  _points[ 4](0) = Real(-3.6783149899818019375269153664372e-01L);
292  _points[ 5](0) = Real(-1.2523340851146891547244136946385e-01L);
293  _points[ 6] = -_points[5];
294  _points[ 7] = -_points[4];
295  _points[ 8] = -_points[3];
296  _points[ 9] = -_points[2];
297  _points[10] = -_points[1];
298  _points[11] = -_points[0];
299 
300  _weights[ 0] = Real(4.7175336386511827194615961485017e-02L);
301  _weights[ 1] = Real(1.0693932599531843096025471819400e-01L);
302  _weights[ 2] = Real(1.6007832854334622633465252954336e-01L);
303  _weights[ 3] = Real(2.0316742672306592174906445580980e-01L);
304  _weights[ 4] = Real(2.3349253653835480876084989892483e-01L);
305  _weights[ 5] = Real(2.4914704581340278500056243604295e-01L);
306  _weights[ 6] = _weights[5];
307  _weights[ 7] = _weights[4];
308  _weights[ 8] = _weights[3];
309  _weights[ 9] = _weights[2];
310  _weights[10] = _weights[1];
311  _weights[11] = _weights[0];
312 
313  return;
314  }
315 
316  case TWENTYFOURTH:
317  case TWENTYFIFTH:
318  {
319  _points.resize (13);
320  _weights.resize(13);
321 
322  _points[ 0](0) = Real(-9.8418305471858814947282944880711e-01L);
323  _points[ 1](0) = Real(-9.1759839922297796520654783650072e-01L);
324  _points[ 2](0) = Real(-8.0157809073330991279420648958286e-01L);
325  _points[ 3](0) = Real(-6.4234933944034022064398460699552e-01L);
326  _points[ 4](0) = Real(-4.4849275103644685287791285212764e-01L);
327  _points[ 5](0) = Real(-2.3045831595513479406552812109799e-01L);
328  _points[ 6](0) = 0.;
329  _points[ 7] = -_points[5];
330  _points[ 8] = -_points[4];
331  _points[ 9] = -_points[3];
332  _points[10] = -_points[2];
333  _points[11] = -_points[1];
334  _points[12] = -_points[0];
335 
336  _weights[ 0] = Real(4.0484004765315879520021592200986e-02L);
337  _weights[ 1] = Real(9.2121499837728447914421775953797e-02L);
338  _weights[ 2] = Real(1.3887351021978723846360177686887e-01L);
339  _weights[ 3] = Real(1.7814598076194573828004669199610e-01L);
340  _weights[ 4] = Real(2.0781604753688850231252321930605e-01L);
341  _weights[ 5] = Real(2.2628318026289723841209018603978e-01L);
342  _weights[ 6] = Real(2.3255155323087391019458951526884e-01L);
343  _weights[ 7] = _weights[5];
344  _weights[ 8] = _weights[4];
345  _weights[ 9] = _weights[3];
346  _weights[10] = _weights[2];
347  _weights[11] = _weights[1];
348  _weights[12] = _weights[0];
349 
350  return;
351  }
352 
353  case TWENTYSIXTH:
354  case TWENTYSEVENTH:
355  {
356  _points.resize (14);
357  _weights.resize(14);
358 
359  _points[ 0](0) = Real(-9.8628380869681233884159726670405e-01L);
360  _points[ 1](0) = Real(-9.2843488366357351733639113937787e-01L);
361  _points[ 2](0) = Real(-8.2720131506976499318979474265039e-01L);
362  _points[ 3](0) = Real(-6.8729290481168547014801980301933e-01L);
363  _points[ 4](0) = Real(-5.1524863635815409196529071855119e-01L);
364  _points[ 5](0) = Real(-3.1911236892788976043567182416848e-01L);
365  _points[ 6](0) = Real(-1.0805494870734366206624465021983e-01L);
366  _points[ 7] = -_points[6];
367  _points[ 8] = -_points[5];
368  _points[ 9] = -_points[4];
369  _points[10] = -_points[3];
370  _points[11] = -_points[2];
371  _points[12] = -_points[1];
372  _points[13] = -_points[0];
373 
374  _weights[ 0] = Real(3.5119460331751863031832876138192e-02L);
375  _weights[ 1] = Real(8.0158087159760209805633277062854e-02L);
376  _weights[ 2] = Real(1.2151857068790318468941480907248e-01L);
377  _weights[ 3] = Real(1.5720316715819353456960193862384e-01L);
378  _weights[ 4] = Real(1.8553839747793781374171659012516e-01L);
379  _weights[ 5] = Real(2.0519846372129560396592406566122e-01L);
380  _weights[ 6] = Real(2.1526385346315779019587644331626e-01L);
381  _weights[ 7] = _weights[6];
382  _weights[ 8] = _weights[5];
383  _weights[ 9] = _weights[4];
384  _weights[10] = _weights[3];
385  _weights[11] = _weights[2];
386  _weights[12] = _weights[1];
387  _weights[13] = _weights[0];
388 
389  return;
390  }
391 
392  case TWENTYEIGHTH:
393  case TWENTYNINTH:
394  {
395  _points.resize (15);
396  _weights.resize(15);
397 
398  _points[ 0](0) = Real(-9.8799251802048542848956571858661e-01L);
399  _points[ 1](0) = Real(-9.3727339240070590430775894771021e-01L);
400  _points[ 2](0) = Real(-8.4820658341042721620064832077422e-01L);
401  _points[ 3](0) = Real(-7.2441773136017004741618605461394e-01L);
402  _points[ 4](0) = Real(-5.7097217260853884753722673725391e-01L);
403  _points[ 5](0) = Real(-3.9415134707756336989720737098105e-01L);
404  _points[ 6](0) = Real(-2.0119409399743452230062830339460e-01L);
405  _points[ 7](0) = 0.;
406  _points[ 8] = -_points[6];
407  _points[ 9] = -_points[5];
408  _points[10] = -_points[4];
409  _points[11] = -_points[3];
410  _points[12] = -_points[2];
411  _points[13] = -_points[1];
412  _points[14] = -_points[0];
413 
414  _weights[ 0] = Real(3.0753241996117268354628393577204e-02L);
415  _weights[ 1] = Real(7.0366047488108124709267416450667e-02L);
416  _weights[ 2] = Real(1.0715922046717193501186954668587e-01L);
417  _weights[ 3] = Real(1.3957067792615431444780479451103e-01L);
418  _weights[ 4] = Real(1.6626920581699393355320086048121e-01L);
419  _weights[ 5] = Real(1.8616100001556221102680056186642e-01L);
420  _weights[ 6] = Real(1.9843148532711157645611832644384e-01L);
421  _weights[ 7] = Real(2.0257824192556127288062019996752e-01L);
422  _weights[ 8] = _weights[6];
423  _weights[ 9] = _weights[5];
424  _weights[10] = _weights[4];
425  _weights[11] = _weights[3];
426  _weights[12] = _weights[2];
427  _weights[13] = _weights[1];
428  _weights[14] = _weights[0];
429 
430  return;
431  }
432 
433  case THIRTIETH:
434  case THIRTYFIRST:
435  {
436  _points.resize (16);
437  _weights.resize(16);
438 
439  _points[ 0](0) = Real(-9.8940093499164993259615417345033e-01L);
440  _points[ 1](0) = Real(-9.4457502307323257607798841553461e-01L);
441  _points[ 2](0) = Real(-8.6563120238783174388046789771239e-01L);
442  _points[ 3](0) = Real(-7.5540440835500303389510119484744e-01L);
443  _points[ 4](0) = Real(-6.1787624440264374844667176404879e-01L);
444  _points[ 5](0) = Real(-4.5801677765722738634241944298358e-01L);
445  _points[ 6](0) = Real(-2.8160355077925891323046050146050e-01L);
446  _points[ 7](0) = Real(-9.5012509837637440185319335424958e-02L);
447  _points[ 8] = -_points[7];
448  _points[ 9] = -_points[6];
449  _points[10] = -_points[5];
450  _points[11] = -_points[4];
451  _points[12] = -_points[3];
452  _points[13] = -_points[2];
453  _points[14] = -_points[1];
454  _points[15] = -_points[0];
455 
456  _weights[ 0] = Real(2.7152459411754094851780572456018e-02L);
457  _weights[ 1] = Real(6.2253523938647892862843836994378e-02L);
458  _weights[ 2] = Real(9.5158511682492784809925107602246e-02L);
459  _weights[ 3] = Real(1.2462897125553387205247628219202e-01L);
460  _weights[ 4] = Real(1.4959598881657673208150173054748e-01L);
461  _weights[ 5] = Real(1.6915651939500253818931207903033e-01L);
462  _weights[ 6] = Real(1.8260341504492358886676366796922e-01L);
463  _weights[ 7] = Real(1.8945061045506849628539672320828e-01L);
464  _weights[ 8] = _weights[7];
465  _weights[ 9] = _weights[6];
466  _weights[10] = _weights[5];
467  _weights[11] = _weights[4];
468  _weights[12] = _weights[3];
469  _weights[13] = _weights[2];
470  _weights[14] = _weights[1];
471  _weights[15] = _weights[0];
472 
473  return;
474  }
475 
476  case THIRTYSECOND:
477  case THIRTYTHIRD:
478  {
479  _points.resize (17);
480  _weights.resize(17);
481 
482  _points[ 0](0) = Real(-9.9057547531441733567543401994067e-01L);
483  _points[ 1](0) = Real(-9.5067552176876776122271695789580e-01L);
484  _points[ 2](0) = Real(-8.8023915372698590212295569448816e-01L);
485  _points[ 3](0) = Real(-7.8151400389680140692523005552048e-01L);
486  _points[ 4](0) = Real(-6.5767115921669076585030221664300e-01L);
487  _points[ 5](0) = Real(-5.1269053708647696788624656862955e-01L);
488  _points[ 6](0) = Real(-3.5123176345387631529718551709535e-01L);
489  _points[ 7](0) = Real(-1.7848418149584785585067749365407e-01L);
490  _points[ 8](0) = 0.;
491  _points[ 9] = -_points[7];
492  _points[10] = -_points[6];
493  _points[11] = -_points[5];
494  _points[12] = -_points[4];
495  _points[13] = -_points[3];
496  _points[14] = -_points[2];
497  _points[15] = -_points[1];
498  _points[16] = -_points[0];
499 
500  _weights[ 0] = Real(2.4148302868547931960110026287565e-02L);
501  _weights[ 1] = Real(5.5459529373987201129440165358245e-02L);
502  _weights[ 2] = Real(8.5036148317179180883535370191062e-02L);
503  _weights[ 3] = Real(1.1188384719340397109478838562636e-01L);
504  _weights[ 4] = Real(1.3513636846852547328631998170235e-01L);
505  _weights[ 5] = Real(1.5404576107681028808143159480196e-01L);
506  _weights[ 6] = Real(1.6800410215645004450997066378832e-01L);
507  _weights[ 7] = Real(1.7656270536699264632527099011320e-01L);
508  _weights[ 8] = Real(1.7944647035620652545826564426189e-01L);
509  _weights[ 9] = _weights[7];
510  _weights[10] = _weights[6];
511  _weights[11] = _weights[5];
512  _weights[12] = _weights[4];
513  _weights[13] = _weights[3];
514  _weights[14] = _weights[2];
515  _weights[15] = _weights[1];
516  _weights[16] = _weights[0];
517 
518  return;
519  }
520 
521  case THIRTYFOURTH:
522  case THIRTYFIFTH:
523  {
524  _points.resize (18);
525  _weights.resize(18);
526 
527  _points[ 0](0) = Real(-9.9156516842093094673001600470615e-01L);
528  _points[ 1](0) = Real(-9.5582394957139775518119589292978e-01L);
529  _points[ 2](0) = Real(-8.9260246649755573920606059112715e-01L);
530  _points[ 3](0) = Real(-8.0370495897252311568241745501459e-01L);
531  _points[ 4](0) = Real(-6.9168704306035320787489108128885e-01L);
532  _points[ 5](0) = Real(-5.5977083107394753460787154852533e-01L);
533  _points[ 6](0) = Real(-4.1175116146284264603593179383305e-01L);
534  _points[ 7](0) = Real(-2.5188622569150550958897285487791e-01L);
535  _points[ 8](0) = Real(-8.4775013041735301242261852935784e-02L);
536  _points[ 9] = -_points[8];
537  _points[10] = -_points[7];
538  _points[11] = -_points[6];
539  _points[12] = -_points[5];
540  _points[13] = -_points[4];
541  _points[14] = -_points[3];
542  _points[15] = -_points[2];
543  _points[16] = -_points[1];
544  _points[17] = -_points[0];
545 
546  _weights[ 0] = Real(2.1616013526483310313342710266452e-02L);
547  _weights[ 1] = Real(4.9714548894969796453334946202639e-02L);
548  _weights[ 2] = Real(7.6425730254889056529129677616637e-02L);
549  _weights[ 3] = Real(1.0094204410628716556281398492483e-01L);
550  _weights[ 4] = Real(1.2255520671147846018451912680020e-01L);
551  _weights[ 5] = Real(1.4064291467065065120473130375195e-01L);
552  _weights[ 6] = Real(1.5468467512626524492541800383637e-01L);
553  _weights[ 7] = Real(1.6427648374583272298605377646593e-01L);
554  _weights[ 8] = Real(1.6914238296314359184065647013499e-01L);
555  _weights[ 9] = _weights[8];
556  _weights[10] = _weights[7];
557  _weights[11] = _weights[6];
558  _weights[12] = _weights[5];
559  _weights[13] = _weights[4];
560  _weights[14] = _weights[3];
561  _weights[15] = _weights[2];
562  _weights[16] = _weights[1];
563  _weights[17] = _weights[0];
564 
565  return;
566  }
567 
568  case THIRTYSIXTH:
569  case THIRTYSEVENTH:
570  {
571  _points.resize (19);
572  _weights.resize(19);
573 
574  _points[ 0](0) = Real(-9.9240684384358440318901767025326e-01L);
575  _points[ 1](0) = Real(-9.6020815213483003085277884068765e-01L);
576  _points[ 2](0) = Real(-9.0315590361481790164266092853231e-01L);
577  _points[ 3](0) = Real(-8.2271465653714282497892248671271e-01L);
578  _points[ 4](0) = Real(-7.2096617733522937861709586082378e-01L);
579  _points[ 5](0) = Real(-6.0054530466168102346963816494624e-01L);
580  _points[ 6](0) = Real(-4.6457074137596094571726714810410e-01L);
581  _points[ 7](0) = Real(-3.1656409996362983199011732884984e-01L);
582  _points[ 8](0) = Real(-1.6035864564022537586809611574074e-01L);
583  _points[ 9](0) = 0.;
584  _points[10] = -_points[8];
585  _points[11] = -_points[7];
586  _points[12] = -_points[6];
587  _points[13] = -_points[5];
588  _points[14] = -_points[4];
589  _points[15] = -_points[3];
590  _points[16] = -_points[2];
591  _points[17] = -_points[1];
592  _points[18] = -_points[0];
593 
594  _weights[ 0] = Real(1.9461788229726477036312041464438e-02L);
595  _weights[ 1] = Real(4.4814226765699600332838157401994e-02L);
596  _weights[ 2] = Real(6.9044542737641226580708258006013e-02L);
597  _weights[ 3] = Real(9.1490021622449999464462094123840e-02L);
598  _weights[ 4] = Real(1.1156664554733399471602390168177e-01L);
599  _weights[ 5] = Real(1.2875396253933622767551578485688e-01L);
600  _weights[ 6] = Real(1.4260670217360661177574610944190e-01L);
601  _weights[ 7] = Real(1.5276604206585966677885540089766e-01L);
602  _weights[ 8] = Real(1.5896884339395434764995643946505e-01L);
603  _weights[ 9] = Real(1.6105444984878369597916362532092e-01L);
604  _weights[10] = _weights[8];
605  _weights[11] = _weights[7];
606  _weights[12] = _weights[6];
607  _weights[13] = _weights[5];
608  _weights[14] = _weights[4];
609  _weights[15] = _weights[3];
610  _weights[16] = _weights[2];
611  _weights[17] = _weights[1];
612  _weights[18] = _weights[0];
613 
614  return;
615  }
616 
617  case THIRTYEIGHTH:
618  case THIRTYNINTH:
619  {
620  _points.resize (20);
621  _weights.resize(20);
622 
623  _points[ 0](0) = Real(-9.9312859918509492478612238847132e-01L);
624  _points[ 1](0) = Real(-9.6397192727791379126766613119728e-01L);
625  _points[ 2](0) = Real(-9.1223442825132590586775244120330e-01L);
626  _points[ 3](0) = Real(-8.3911697182221882339452906170152e-01L);
627  _points[ 4](0) = Real(-7.4633190646015079261430507035564e-01L);
628  _points[ 5](0) = Real(-6.3605368072651502545283669622629e-01L);
629  _points[ 6](0) = Real(-5.1086700195082709800436405095525e-01L);
630  _points[ 7](0) = Real(-3.7370608871541956067254817702493e-01L);
631  _points[ 8](0) = Real(-2.2778585114164507808049619536857e-01L);
632  _points[ 9](0) = Real(-7.6526521133497333754640409398838e-02L);
633  _points[10] = -_points[9];
634  _points[11] = -_points[8];
635  _points[12] = -_points[7];
636  _points[13] = -_points[6];
637  _points[14] = -_points[5];
638  _points[15] = -_points[4];
639  _points[16] = -_points[3];
640  _points[17] = -_points[2];
641  _points[18] = -_points[1];
642  _points[19] = -_points[0];
643 
644  _weights[ 0] = Real(1.7614007139152118311861962351853e-02L);
645  _weights[ 1] = Real(4.0601429800386941331039952274932e-02L);
646  _weights[ 2] = Real(6.2672048334109063569506535187042e-02L);
647  _weights[ 3] = Real(8.3276741576704748724758143222046e-02L);
648  _weights[ 4] = Real(1.0193011981724043503675013548035e-01L);
649  _weights[ 5] = Real(1.1819453196151841731237737771138e-01L);
650  _weights[ 6] = Real(1.3168863844917662689849449974816e-01L);
651  _weights[ 7] = Real(1.4209610931838205132929832506716e-01L);
652  _weights[ 8] = Real(1.4917298647260374678782873700197e-01L);
653  _weights[ 9] = Real(1.5275338713072585069808433195510e-01L);
654  _weights[10] = _weights[9];
655  _weights[11] = _weights[8];
656  _weights[12] = _weights[7];
657  _weights[13] = _weights[6];
658  _weights[14] = _weights[5];
659  _weights[15] = _weights[4];
660  _weights[16] = _weights[3];
661  _weights[17] = _weights[2];
662  _weights[18] = _weights[1];
663  _weights[19] = _weights[0];
664 
665  return;
666  }
667 
668  case FORTIETH:
669  case FORTYFIRST:
670  {
671  _points.resize (21);
672  _weights.resize(21);
673 
674  _points[ 0](0) = Real(-9.9375217062038950026024203593794e-01L);
675  _points[ 1](0) = Real(-9.6722683856630629431662221490770e-01L);
676  _points[ 2](0) = Real(-9.2009933415040082879018713371497e-01L);
677  _points[ 3](0) = Real(-8.5336336458331728364725063858757e-01L);
678  _points[ 4](0) = Real(-7.6843996347567790861587785130623e-01L);
679  _points[ 5](0) = Real(-6.6713880419741231930596666999034e-01L);
680  _points[ 6](0) = Real(-5.5161883588721980705901879672431e-01L);
681  _points[ 7](0) = Real(-4.2434212020743878357366888854379e-01L);
682  _points[ 8](0) = Real(-2.8802131680240109660079251606460e-01L);
683  _points[ 9](0) = Real(-1.4556185416089509093703098233869e-01L);
684  _points[10](0) = 0.;
685  _points[11] = -_points[9];
686  _points[12] = -_points[8];
687  _points[13] = -_points[7];
688  _points[14] = -_points[6];
689  _points[15] = -_points[5];
690  _points[16] = -_points[4];
691  _points[17] = -_points[3];
692  _points[18] = -_points[2];
693  _points[19] = -_points[1];
694  _points[20] = -_points[0];
695 
696  _weights[ 0] = Real(1.6017228257774333324224616858471e-02L);
697  _weights[ 1] = Real(3.6953789770852493799950668299330e-02L);
698  _weights[ 2] = Real(5.7134425426857208283635826472448e-02L);
699  _weights[ 3] = Real(7.6100113628379302017051653300183e-02L);
700  _weights[ 4] = Real(9.3444423456033861553289741113932e-02L);
701  _weights[ 5] = Real(1.0879729916714837766347457807011e-01L);
702  _weights[ 6] = Real(1.2183141605372853419536717712572e-01L);
703  _weights[ 7] = Real(1.3226893863333746178105257449678e-01L);
704  _weights[ 8] = Real(1.3988739479107315472213342386758e-01L);
705  _weights[ 9] = Real(1.4452440398997005906382716655375e-01L);
706  _weights[10] = Real(1.4608113364969042719198514768337e-01L);
707  _weights[11] = _weights[9];
708  _weights[12] = _weights[8];
709  _weights[13] = _weights[7];
710  _weights[14] = _weights[6];
711  _weights[15] = _weights[5];
712  _weights[16] = _weights[4];
713  _weights[17] = _weights[3];
714  _weights[18] = _weights[2];
715  _weights[19] = _weights[1];
716  _weights[20] = _weights[0];
717 
718  return;
719  }
720 
721  case FORTYSECOND:
722  case FORTYTHIRD:
723  {
724  _points.resize (22);
725  _weights.resize(22);
726 
727  _points[ 0](0) = Real(-9.9429458548239929207303142116130e-01L);
728  _points[ 1](0) = Real(-9.7006049783542872712395098676527e-01L);
729  _points[ 2](0) = Real(-9.2695677218717400052069293925905e-01L);
730  _points[ 3](0) = Real(-8.6581257772030013653642563701938e-01L);
731  _points[ 4](0) = Real(-7.8781680597920816200427795540835e-01L);
732  _points[ 5](0) = Real(-6.9448726318668278005068983576226e-01L);
733  _points[ 6](0) = Real(-5.8764040350691159295887692763865e-01L);
734  _points[ 7](0) = Real(-4.6935583798675702640633071096641e-01L);
735  _points[ 8](0) = Real(-3.4193582089208422515814742042738e-01L);
736  _points[ 9](0) = Real(-2.0786042668822128547884653391955e-01L);
737  _points[10](0) = Real(-6.9739273319722221213841796118628e-02L);
738  _points[11] = -_points[10];
739  _points[12] = -_points[9];
740  _points[13] = -_points[8];
741  _points[14] = -_points[7];
742  _points[15] = -_points[6];
743  _points[16] = -_points[5];
744  _points[17] = -_points[4];
745  _points[18] = -_points[3];
746  _points[19] = -_points[2];
747  _points[20] = -_points[1];
748  _points[21] = -_points[0];
749 
750  _weights[ 0] = Real(1.4627995298272200684991098047185e-02L);
751  _weights[ 1] = Real(3.3774901584814154793302246865913e-02L);
752  _weights[ 2] = Real(5.2293335152683285940312051273211e-02L);
753  _weights[ 3] = Real(6.9796468424520488094961418930218e-02L);
754  _weights[ 4] = Real(8.5941606217067727414443681372703e-02L);
755  _weights[ 5] = Real(1.0041414444288096493207883783054e-01L);
756  _weights[ 6] = Real(1.1293229608053921839340060742175e-01L);
757  _weights[ 7] = Real(1.2325237681051242428556098615481e-01L);
758  _weights[ 8] = Real(1.3117350478706237073296499253031e-01L);
759  _weights[ 9] = Real(1.3654149834601517135257383123152e-01L);
760  _weights[10] = Real(1.3925187285563199337541024834181e-01L);
761  _weights[11] = _weights[10];
762  _weights[12] = _weights[9];
763  _weights[13] = _weights[8];
764  _weights[14] = _weights[7];
765  _weights[15] = _weights[6];
766  _weights[16] = _weights[5];
767  _weights[17] = _weights[4];
768  _weights[18] = _weights[3];
769  _weights[19] = _weights[2];
770  _weights[20] = _weights[1];
771  _weights[21] = _weights[0];
772 
773  return;
774  }
775 
776 
777  default:
778  libmesh_error_msg("Quadrature rule " << _order << " not supported!");
779  }
780 }
781 
782 } // 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