Chi-Tech
quadrature_triangle.cc
Go to the documentation of this file.
2
4
5#include <cmath>
6#include <stdexcept>
7#include <cassert>
8
9//###################################################################
10/**Initializes quadratures for use on triangles.*/
13 Quadrature(order)
14{
15 double x=0.0,y=0.0,z=0.0;
16
17 auto& _points = qpoints_;
18 auto& _weights = weights_;
19
20 typedef chi_mesh::Vector3 Point;
21 typedef double Real;
22
23 switch (order)
24 {
27 {
28 x = 1.0/3.0;
29 y = 1.0/3.0;
30 qpoints_.emplace_back(x, y, z);
31 weights_.push_back(0.5);
32 break;
33 }
35 {
36 x = 1.0/6.0;
37 y = 1.0/6.0;
38 qpoints_.emplace_back(x, y, z);
39 weights_.push_back(1.0 / 6.0);
40
41 x = 4.0/6.0;
42 y = 1.0/6.0;
43 qpoints_.emplace_back(x, y, z);
44 weights_.push_back(1.0 / 6.0);
45
46 x = 1.0/6.0;
47 y = 4.0/6.0;
48 qpoints_.emplace_back(x, y, z);
49 weights_.push_back(1.0 / 6.0);
50 break;
51 }
53 {
54 // Exact for cubics
55 qpoints_.resize(4);
56 weights_.resize(4);
57
58 // This rule is formed from a tensor product of
59 // appropriately-scaled Gauss and Jacobi rules. (See
60 // also the QConical quadrature class, this is a
61 // hard-coded version of one of those rules.) For high
62 // orders these rules generally have too many points,
63 // but at extremely low order they are competitive and
64 // have the additional benefit of having all positive
65 // weights.
66 qpoints_[0](0) = (1.5505102572168219018027159252941e-01L);
67 qpoints_[0](1) = (1.7855872826361642311703513337422e-01L);
68 qpoints_[1](0) = (6.4494897427831780981972840747059e-01L);
69 qpoints_[1](1) = (7.5031110222608118177475598324603e-02L);
70 qpoints_[2](0) = (1.5505102572168219018027159252941e-01L);
71 qpoints_[2](1) = (6.6639024601470138670269327409637e-01L);
72 qpoints_[3](0) = (6.4494897427831780981972840747059e-01L);
73 qpoints_[3](1) = (2.8001991549907407200279599420481e-01L);
74
75 weights_[0] = (1.5902069087198858469718450103758e-01L);
76 weights_[1] = (9.0979309128011415302815498962418e-02L);
77 weights_[2] = (1.5902069087198858469718450103758e-01L);
78 weights_[3] = (9.0979309128011415302815498962418e-02L);
79 break;
80 }
81
82 /** A degree 4 rule with six points. This rule can be found in many places
83 including:
84 J.N. Lyness and D. Jespersen, Moderate degree symmetric
85 quadrature rules for the triangle, J. Inst. Math. Appl. 15 (1975),
86 19--32.
87 We used the code in:
88 L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
89 on triangles and tetrahedra" Journal of Computational Mathematics,
90 v. 27, no. 1, 2009, pp. 89-96.
91 to generate additional precision.*/
93 {
94 const unsigned int n_wts = 2;
95 const double wts[n_wts] =
96 {
97 double(1.1169079483900573284750350421656140e-01L),
98 double(5.4975871827660933819163162450105264e-02L)
99 };
100
101 const double a[n_wts] =
102 {
103 double(4.4594849091596488631832925388305199e-01L),
104 double(9.1576213509770743459571463402201508e-02L)
105 };
106
107 const double b[n_wts] = {0., 0.}; // not used
108 const unsigned int permutation_ids[n_wts] = {3, 3};
109
110 dunavant_rule2(wts, a, b, permutation_ids, n_wts); // 6 total points
111 break;
112 }
113
114 /**Exact for quintics
115 Can be found in "Quadrature on Simplices of Arbitrary
116 Dimension" by Walkington.*/
118 {
119 const unsigned int n_wts = 3;
120 const double wts[n_wts] =
121 {
122 double(9)/80,
123 double(31)/480 + std::sqrt(double(15.0L))/2400,
124 double(31)/480 - std::sqrt(double(15.0L))/2400
125 };
126
127 const double a[n_wts] =
128 {
129 0., // 'a' parameter not used for origin permutation
130 double(2)/7 + std::sqrt(double(15.0L))/21,
131 double(2)/7 - std::sqrt(double(15.0L))/21
132 };
133
134 const double b[n_wts] = {0., 0., 0.}; // not used
135 const unsigned int permutation_ids[n_wts] = {1, 3, 3};
136
137 dunavant_rule2(wts, a, b, permutation_ids, n_wts); // 7 total points
138
139 break;
140 }
141
142 // A degree 6 rule with 12 points. This rule can be found in many places
143 // including:
144 //
145 // J.N. Lyness and D. Jespersen, Moderate degree symmetric
146 // quadrature rules for the triangle, J. Inst. Math. Appl. 15 (1975),
147 // 19--32.
148 //
149 // We used the code in:
150 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
151 // on triangles and tetrahedra" Journal of Computational Mathematics,
152 // v. 27, no. 1, 2009, pp. 89-96.
153 // to generate additional precision.
154 //
155 // Note that the following 7th-order Ro3-invariant rule also has only 12 points,
156 // which technically makes it the superior rule. This one is here for completeness.
158 {
159 const unsigned int n_wts = 3;
160 const Real wts[n_wts] =
161 {
162 Real(5.8393137863189683012644805692789721e-02L),
163 Real(2.5422453185103408460468404553434492e-02L),
164 Real(4.1425537809186787596776728210221227e-02L)
165 };
166
167 const Real a[n_wts] =
168 {
169 Real(2.4928674517091042129163855310701908e-01L),
170 Real(6.3089014491502228340331602870819157e-02L),
171 Real(3.1035245103378440541660773395655215e-01L)
172 };
173
174 const Real b[n_wts] =
175 {
176 0.,
177 0.,
178 Real(6.3650249912139864723014259441204970e-01L)
179 };
180
181 const unsigned int permutation_ids[n_wts] = {3, 3, 6}; // 12 total points
182
183 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
184
185 return;
186 }
187
188
189 // A degree 7 rule with 12 points. This rule can be found in:
190 //
191 // K. Gatermann, The construction of symmetric cubature
192 // formulas for the square and the triangle, Computing 40
193 // (1988), 229--240.
194 //
195 // This rule, which is provably minimal in the number of
196 // integration points, is said to be 'Ro3 invariant' which
197 // means that a given set of barycentric coordinates
198 // (z1,z2,z3) implies the quadrature points (z1,z2),
199 // (z3,z1), (z2,z3) which are formed by taking the first
200 // two entries in cyclic permutations of the barycentric
201 // point. Barycentric coordinates are related in the
202 // sense that: z3 = 1 - z1 - z2.
203 //
204 // The 12-point sixth-order rule for triangles given in
205 // Flaherty's (http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps)
206 // lecture notes has been removed in favor of this rule
207 // which is higher-order (for the same number of
208 // quadrature points) and has a few more digits of
209 // precision in the points and weights. Some 10-point
210 // degree 6 rules exist for the triangle but they have
211 // quadrature points outside the region of integration.
213 {
214 _points.resize (12);
215 _weights.resize(12);
216
217 const unsigned int nrows=4;
218
219 // In each of the rows below, the first two entries are (z1, z2) which imply
220 // z3. The third entry is the weight for each of the points in the cyclic permutation.
221 // The original publication tabulated about 16 decimal digits for each point and weight
222 // parameter. The additional digits shown here were obtained using a code in the
223 // mp-quadrature library, https://github.com/jwpeterson/mp-quadrature
224 const Real rule_data[nrows][3] = {
225 {Real(6.2382265094402118173683000996350e-02L), Real(6.7517867073916085442557131050869e-02L), Real(2.6517028157436251428754180460739e-02L)}, // group A
226 {Real(5.5225456656926611737479190275645e-02L), Real(3.2150249385198182266630784919920e-01L), Real(4.3881408714446055036769903139288e-02L)}, // group B
227 {Real(3.4324302945097146469630642483938e-02L), Real(6.6094919618673565761198031019780e-01L), Real(2.8775042784981585738445496900219e-02L)}, // group C
228 {Real(5.1584233435359177925746338682643e-01L), Real(2.7771616697639178256958187139372e-01L), Real(6.7493187009802774462697086166421e-02L)} // group D
229 };
230
231 for (unsigned int i=0, offset=0; i<nrows; ++i)
232 {
233 _points[offset + 0] = Point(rule_data[i][0], rule_data[i][1]); // (z1,z2)
234 _points[offset + 1] = Point(1.-rule_data[i][0]-rule_data[i][1], rule_data[i][0]); // (z3,z1)
235 _points[offset + 2] = Point(rule_data[i][1], 1.-rule_data[i][0]-rule_data[i][1]); // (z2,z3)
236
237 // All these points get the same weight
238 _weights[offset + 0] = rule_data[i][2];
239 _weights[offset + 1] = rule_data[i][2];
240 _weights[offset + 2] = rule_data[i][2];
241
242 // Increment offset
243 offset += 3;
244 }
245
246 return;
247
248
249 // // The following is an inferior 7th-order Lyness-style rule with 15 points.
250 // // It's here only for completeness and the Ro3-invariant rule above should
251 // // be used instead!
252 // const unsigned int n_wts = 3;
253 // const Real wts[n_wts] =
254 // {
255 // Real(2.6538900895116205835977487499847719e-02L),
256 // Real(3.5426541846066783659206291623201826e-02L),
257 // Real(3.4637341039708446756138297960207647e-02L)
258 // };
259 //
260 // const Real a[n_wts] =
261 // {
262 // Real(6.4930513159164863078379776030396538e-02L),
263 // Real(2.8457558424917033519741605734978046e-01L),
264 // Real(3.1355918438493150795585190219862865e-01L)
265 // };
266 //
267 // const Real b[n_wts] =
268 // {
269 // 0.,
270 // Real(1.9838447668150671917987659863332941e-01L),
271 // Real(4.3863471792372471511798695971295936e-02L)
272 // };
273 //
274 // const unsigned int permutation_ids[n_wts] = {3, 6, 6}; // 15 total points
275 //
276 // dunavant_rule2(wts, a, b, permutation_ids, n_wts);
277 //
278 // return;
279 }
280
281
282
283
284 // Another Dunavant rule. This one has all positive weights. This rule has
285 // 16 points while a comparable conical product rule would have 5*5=25.
286 //
287 // It was copied 23rd June 2008 from:
288 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
289 //
290 // Additional precision obtained from the code in:
291 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
292 // on triangles and tetrahedra" Journal of Computational Mathematics,
293 // v. 27, no. 1, 2009, pp. 89-96.
295 {
296 const unsigned int n_wts = 5;
297 const Real wts[n_wts] =
298 {
299 Real(7.2157803838893584125545555244532310e-02L),
300 Real(4.7545817133642312396948052194292159e-02L),
301 Real(5.1608685267359125140895775146064515e-02L),
302 Real(1.6229248811599040155462964170890299e-02L),
303 Real(1.3615157087217497132422345036954462e-02L)
304 };
305
306 const Real a[n_wts] =
307 {
308 0.0, // 'a' parameter not used for origin permutation
309 Real(4.5929258829272315602881551449416932e-01L),
310 Real(1.7056930775176020662229350149146450e-01L),
311 Real(5.0547228317030975458423550596598947e-02L),
312 Real(2.6311282963463811342178578628464359e-01L),
313 };
314
315 const Real b[n_wts] =
316 {
317 0.,
318 0.,
319 0.,
320 0.,
321 Real(7.2849239295540428124100037917606196e-01L)
322 };
323
324 const unsigned int permutation_ids[n_wts] = {1, 3, 3, 3, 6}; // 16 total points
325
326 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
327
328 return;
329 }
330
331
332
333 // Another Dunavant rule. This one has all positive weights. This rule has 19
334 // points. The comparable conical product rule would have 25.
335 // It was copied 23rd June 2008 from:
336 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
337 //
338 // Additional precision obtained from the code in:
339 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
340 // on triangles and tetrahedra" Journal of Computational Mathematics,
341 // v. 27, no. 1, 2009, pp. 89-96.
343 {
344 const unsigned int n_wts = 6;
345 const Real wts[n_wts] =
346 {
347 Real(4.8567898141399416909620991253644315e-02L),
348 Real(1.5667350113569535268427415643604658e-02L),
349 Real(1.2788837829349015630839399279499912e-02L),
350 Real(3.8913770502387139658369678149701978e-02L),
351 Real(3.9823869463605126516445887132022637e-02L),
352 Real(2.1641769688644688644688644688644689e-02L)
353 };
354
355 const Real a[n_wts] =
356 {
357 0.0, // 'a' parameter not used for origin permutation
358 Real(4.8968251919873762778370692483619280e-01L),
359 Real(4.4729513394452709865106589966276365e-02L),
360 Real(4.3708959149293663726993036443535497e-01L),
361 Real(1.8820353561903273024096128046733557e-01L),
362 Real(2.2196298916076569567510252769319107e-01L)
363 };
364
365 const Real b[n_wts] =
366 {
367 0.,
368 0.,
369 0.,
370 0.,
371 0.,
372 Real(7.4119859878449802069007987352342383e-01L)
373 };
374
375 const unsigned int permutation_ids[n_wts] = {1, 3, 3, 3, 3, 6}; // 19 total points
376
377 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
378
379 return;
380 }
381
382
383 // Another Dunavant rule with all positive weights. This rule has 25
384 // points. The comparable conical product rule would have 36.
385 // It was copied 23rd June 2008 from:
386 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
387 //
388 // Additional precision obtained from the code in:
389 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
390 // on triangles and tetrahedra" Journal of Computational Mathematics,
391 // v. 27, no. 1, 2009, pp. 89-96.
393 {
394 const unsigned int n_wts = 6;
395 const Real wts[n_wts] =
396 {
397 Real(4.5408995191376790047643297550014267e-02L),
398 Real(1.8362978878233352358503035945683300e-02L),
399 Real(2.2660529717763967391302822369298659e-02L),
400 Real(3.6378958422710054302157588309680344e-02L),
401 Real(1.4163621265528742418368530791049552e-02L),
402 Real(4.7108334818664117299637354834434138e-03L)
403 };
404
405 const Real a[n_wts] =
406 {
407 0.0, // 'a' parameter not used for origin permutation
408 Real(4.8557763338365737736750753220812615e-01L),
409 Real(1.0948157548503705479545863134052284e-01L),
410 Real(3.0793983876412095016515502293063162e-01L),
411 Real(2.4667256063990269391727646541117681e-01L),
412 Real(6.6803251012200265773540212762024737e-02L)
413 };
414
415 const Real b[n_wts] =
416 {
417 0.,
418 0.,
419 0.,
420 Real(5.5035294182099909507816172659300821e-01L),
421 Real(7.2832390459741092000873505358107866e-01L),
422 Real(9.2365593358750027664630697761508843e-01L)
423 };
424
425 const unsigned int permutation_ids[n_wts] = {1, 3, 3, 6, 6, 6}; // 25 total points
426
427 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
428
429 return;
430 }
431
432
433 // Dunavant's 11th-order rule contains points outside the region of
434 // integration, and is thus unacceptable for our FEM calculations.
435 //
436 // This 30-point, 11th-order rule was obtained by me [JWP] using the code in
437 //
438 // Additional precision obtained from the code in:
439 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
440 // on triangles and tetrahedra" Journal of Computational Mathematics,
441 // v. 27, no. 1, 2009, pp. 89-96.
442 //
443 // Note: the 28-point 11th-order rule obtained by Zhang in the paper above
444 // does not appear to be unique. It is a solution in the sense that it
445 // minimizes the error in the least-squares minimization problem, but
446 // it involves too many unknowns and the Jacobian is therefore singular
447 // when attempting to improve the solution via Newton's method.
449 {
450 const unsigned int n_wts = 6;
451 const Real wts[n_wts] =
452 {
453 Real(3.6089021198604635216985338480426484e-02L),
454 Real(2.1607717807680420303346736867931050e-02L),
455 Real(3.1144524293927978774861144478241807e-03L),
456 Real(2.9086855161081509446654185084988077e-02L),
457 Real(8.4879241614917017182977532679947624e-03L),
458 Real(1.3795732078224796530729242858347546e-02L)
459 };
460
461 const Real a[n_wts] =
462 {
463 Real(3.9355079629947969884346551941969960e-01L),
464 Real(4.7979065808897448654107733982929214e-01L),
465 Real(5.1003445645828061436081405648347852e-03L),
466 Real(2.6597620190330158952732822450744488e-01L),
467 Real(2.8536418538696461608233522814483715e-01L),
468 Real(1.3723536747817085036455583801851025e-01L)
469 };
470
471 const Real b[n_wts] =
472 {
473 0.,
474 0.,
475 Real(5.6817155788572446538150614865768991e-02L),
476 Real(1.2539956353662088473247489775203396e-01L),
477 Real(1.2409970153698532116262152247041742e-02L),
478 Real(5.2792057988217708934207928630851643e-02L)
479 };
480
481 const unsigned int permutation_ids[n_wts] = {3, 3, 6, 6, 6, 6}; // 30 total points
482
483 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
484
485 return;
486 }
487
488
489
490
491 // Another Dunavant rule with all positive weights. This rule has 33
492 // points. The comparable conical product rule would have 36 (ELEVENTH) or 49 (TWELFTH).
493 //
494 // It was copied 23rd June 2008 from:
495 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
496 //
497 // Additional precision obtained from the code in:
498 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
499 // on triangles and tetrahedra" Journal of Computational Mathematics,
500 // v. 27, no. 1, 2009, pp. 89-96.
502 {
503 const unsigned int n_wts = 8;
504 const Real wts[n_wts] =
505 {
506 Real(3.0831305257795086169332418926151771e-03L),
507 Real(3.1429112108942550177135256546441273e-02L),
508 Real(1.7398056465354471494664198647499687e-02L),
509 Real(2.1846272269019201067728631278737487e-02L),
510 Real(1.2865533220227667708895461535782215e-02L),
511 Real(1.1178386601151722855919538351159995e-02L),
512 Real(8.6581155543294461858210504055170332e-03L),
513 Real(2.0185778883190464758914349626118386e-02L)
514 };
515
516 const Real a[n_wts] =
517 {
518 Real(2.1317350453210370246856975515728246e-02L),
519 Real(2.7121038501211592234595134039689474e-01L),
520 Real(1.2757614554158592467389632515428357e-01L),
521 Real(4.3972439229446027297973662348436108e-01L),
522 Real(4.8821738977380488256466206525881104e-01L),
523 Real(2.8132558098993954824813069297455275e-01L),
524 Real(1.1625191590759714124135414784260182e-01L),
525 Real(2.7571326968551419397479634607976398e-01L)
526 };
527
528 const Real b[n_wts] =
529 {
530 0.,
531 0.,
532 0.,
533 0.,
534 0.,
535 Real(6.9583608678780342214163552323607254e-01L),
536 Real(8.5801403354407263059053661662617818e-01L),
537 Real(6.0894323577978780685619243776371007e-01L)
538 };
539
540 const unsigned int permutation_ids[n_wts] = {3, 3, 3, 3, 3, 6, 6, 6}; // 33 total points
541
542 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
543
544 return;
545 }
546
547
548 // Another Dunavant rule with all positive weights. This rule has 37
549 // points. The comparable conical product rule would have 49 points.
550 //
551 // It was copied 23rd June 2008 from:
552 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
553 //
554 // A second rule with additional precision obtained from the code in:
555 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
556 // on triangles and tetrahedra" Journal of Computational Mathematics,
557 // v. 27, no. 1, 2009, pp. 89-96.
559 {
560 const unsigned int n_wts = 9;
561 const Real wts[n_wts] =
562 {
563 Real(3.3980018293415822140887212340442440e-02L),
564 Real(2.7800983765226664353628733005230734e-02L),
565 Real(2.9139242559599990702383541756669905e-02L),
566 Real(3.0261685517695859208964000161454122e-03L),
567 Real(1.1997200964447365386855399725479827e-02L),
568 Real(1.7320638070424185232993414255459110e-02L),
569 Real(7.4827005525828336316229285664517190e-03L),
570 Real(1.2089519905796909568722872786530380e-02L),
571 Real(4.7953405017716313612975450830554457e-03L)
572 };
573
574 const Real a[n_wts] =
575 {
576 0., // 'a' parameter not used for origin permutation
577 Real(4.2694141425980040602081253503137421e-01L),
578 Real(2.2137228629183290065481255470507908e-01L),
579 Real(2.1509681108843183869291313534052083e-02L),
580 Real(4.8907694645253934990068971909020439e-01L),
581 Real(3.0844176089211777465847185254124531e-01L),
582 Real(1.1092204280346339541286954522167452e-01L),
583 Real(1.6359740106785048023388790171095725e-01L),
584 Real(2.7251581777342966618005046435408685e-01L)
585 };
586
587 const Real b[n_wts] =
588 {
589 0.,
590 0.,
591 0.,
592 0.,
593 0.,
594 Real(6.2354599555367557081585435318623659e-01L),
595 Real(8.6470777029544277530254595089569318e-01L),
596 Real(7.4850711589995219517301859578870965e-01L),
597 Real(7.2235779312418796526062013230478405e-01L)
598 };
599
600 const unsigned int permutation_ids[n_wts] = {1, 3, 3, 3, 3, 6, 6, 6, 6}; // 37 total points
601
602 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
603
604 return;
605 }
606
607
608 // Another Dunavant rule. This rule has 42 points, while
609 // a comparable conical product rule would have 64.
610 //
611 // It was copied 23rd June 2008 from:
612 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
613 //
614 // Additional precision obtained from the code in:
615 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
616 // on triangles and tetrahedra" Journal of Computational Mathematics,
617 // v. 27, no. 1, 2009, pp. 89-96.
619 {
620 const unsigned int n_wts = 10;
621 const Real wts[n_wts] =
622 {
623 Real(1.0941790684714445320422472981662986e-02L),
624 Real(1.6394176772062675320655489369312672e-02L),
625 Real(2.5887052253645793157392455083198201e-02L),
626 Real(2.1081294368496508769115218662093065e-02L),
627 Real(7.2168498348883338008549607403266583e-03L),
628 Real(2.4617018012000408409130117545210774e-03L),
629 Real(1.2332876606281836981437622591818114e-02L),
630 Real(1.9285755393530341614244513905205430e-02L),
631 Real(7.2181540567669202480443459995079017e-03L),
632 Real(2.5051144192503358849300465412445582e-03L)
633 };
634
635 const Real a[n_wts] =
636 {
637 Real(4.8896391036217863867737602045239024e-01L),
638 Real(4.1764471934045392250944082218564344e-01L),
639 Real(2.7347752830883865975494428326269856e-01L),
640 Real(1.7720553241254343695661069046505908e-01L),
641 Real(6.1799883090872601267478828436935788e-02L),
642 Real(1.9390961248701048178250095054529511e-02L),
643 Real(1.7226668782135557837528960161365733e-01L),
644 Real(3.3686145979634500174405519708892539e-01L),
645 Real(2.9837288213625775297083151805961273e-01L),
646 Real(1.1897449769695684539818196192990548e-01L)
647 };
648
649 const Real b[n_wts] =
650 {
651 0.,
652 0.,
653 0.,
654 0.,
655 0.,
656 0.,
657 Real(7.7060855477499648258903327416742796e-01L),
658 Real(5.7022229084668317349769621336235426e-01L),
659 Real(6.8698016780808783735862715402031306e-01L),
660 Real(8.7975717137017112951457163697460183e-01L)
661 };
662
663 const unsigned int permutation_ids[n_wts]
664 = {3, 3, 3, 3, 3, 3, 6, 6, 6, 6}; // 42 total points
665
666 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
667
668 return;
669 }
670
671
672 // This 49-point rule was found by me [JWP] using the code in:
673 //
674 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
675 // on triangles and tetrahedra" Journal of Computational Mathematics,
676 // v. 27, no. 1, 2009, pp. 89-96.
677 //
678 // A 54-point, 15th-order rule is reported by
679 //
680 // Stephen Wandzura, Hong Xiao,
681 // Symmetric Quadrature Rules on a Triangle,
682 // Computers and Mathematics with Applications,
683 // Volume 45, Number 12, June 2003, pages 1829-1840.
684 //
685 // can be found here:
686 // http://people.scs.fsu.edu/~burkardt/f_src/wandzura/wandzura.f90
687 //
688 // but this 49-point rule is superior.
690 {
691 const unsigned int n_wts = 11;
692 const Real wts[n_wts] =
693 {
694 Real(2.4777380743035579804788826970198951e-02L),
695 Real(9.2433943023307730591540642828347660e-03L),
696 Real(2.2485768962175402793245929133296627e-03L),
697 Real(6.7052581900064143760518398833360903e-03L),
698 Real(1.9011381726930579256700190357527956e-02L),
699 Real(1.4605445387471889398286155981802858e-02L),
700 Real(1.5087322572773133722829435011138258e-02L),
701 Real(1.5630213780078803020711746273129099e-02L),
702 Real(6.1808086085778203192616856133701233e-03L),
703 Real(3.2209366452594664857296985751120513e-03L),
704 Real(5.8747373242569702667677969985668817e-03L)
705 };
706
707 const Real a[n_wts] =
708 {
709 0.0, // 'a' parameter not used for origin
710 Real(7.9031013655541635005816956762252155e-02L),
711 Real(1.8789501810770077611247984432284226e-02L),
712 Real(4.9250168823249670532514526605352905e-01L),
713 Real(4.0886316907744105975059040108092775e-01L),
714 Real(5.3877851064220142445952549348423733e-01L),
715 Real(2.0250549804829997692885033941362673e-01L),
716 Real(5.5349674918711643207148086558288110e-01L),
717 Real(7.8345022567320812359258882143250181e-01L),
718 Real(8.9514624528794883409864566727625002e-01L),
719 Real(3.2515745241110782862789881780746490e-01L)
720 };
721
722 const Real b[n_wts] =
723 {
724 0.,
725 0.,
726 0.,
727 0.,
728 0.,
729 Real(1.9412620368774630292701241080996842e-01L),
730 Real(9.8765911355712115933807754318089099e-02L),
731 Real(7.7663767064308164090246588765178087e-02L),
732 Real(2.1594628433980258573654682690950798e-02L),
733 Real(1.2563596287784997705599005477153617e-02L),
734 Real(1.5082654870922784345283124845552190e-02L)
735 };
736
737 const unsigned int permutation_ids[n_wts]
738 = {1, 3, 3, 3, 3, 6, 6, 6, 6, 6, 6}; // 49 total points
739
740 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
741
742 return;
743 }
744
745
746
747
748 // Dunavant's 16th-order rule contains points outside the region of
749 // integration, and is thus unacceptable for our FEM calculations.
750 //
751 // This 55-point, 16th-order rule was obtained by me [JWP] using the code in
752 //
753 // Additional precision obtained from the code in:
754 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
755 // on triangles and tetrahedra" Journal of Computational Mathematics,
756 // v. 27, no. 1, 2009, pp. 89-96.
757 //
758 // Note: the 55-point 16th-order rule obtained by Zhang in the paper above
759 // does not appear to be unique. It is a solution in the sense that it
760 // minimizes the error in the least-squares minimization problem, but
761 // it involves too many unknowns and the Jacobian is therefore singular
762 // when attempting to improve the solution via Newton's method.
764 {
765 const unsigned int n_wts = 12;
766 const Real wts[n_wts] =
767 {
768 Real(2.2668082505910087151996321171534230e-02L),
769 Real(8.4043060714818596159798961899306135e-03L),
770 Real(1.0850949634049747713966288634484161e-03L),
771 Real(7.2252773375423638869298219383808751e-03L),
772 Real(1.2997715227338366024036316182572871e-02L),
773 Real(2.0054466616677715883228810959112227e-02L),
774 Real(9.7299841600417010281624372720122710e-03L),
775 Real(1.1651974438298104227427176444311766e-02L),
776 Real(9.1291185550484450744725847363097389e-03L),
777 Real(3.5568614040947150231712567900113671e-03L),
778 Real(5.8355861686234326181790822005304303e-03L),
779 Real(4.7411314396804228041879331486234396e-03L)
780 };
781
782 const Real a[n_wts] =
783 {
784 0.0, // 'a' parameter not used for centroid weight
785 Real(8.5402539407933203673769900926355911e-02L),
786 Real(1.2425572001444092841183633409631260e-02L),
787 Real(4.9174838341891594024701017768490960e-01L),
788 Real(4.5669426695387464162068900231444462e-01L),
789 Real(4.8506759880447437974189793537259677e-01L),
790 Real(2.0622099278664205707909858461264083e-01L),
791 Real(3.2374950270039093446805340265853956e-01L),
792 Real(7.3834330556606586255186213302750029e-01L),
793 Real(9.1210673061680792565673823935174611e-01L),
794 Real(6.6129919222598721544966837350891531e-01L),
795 Real(1.7807138906021476039088828811346122e-01L)
796 };
797
798 const Real b[n_wts] =
799 {
800 0.0,
801 0.0,
802 0.0,
803 0.0,
804 0.0,
805 Real(3.2315912848634384647700266402091638e-01L),
806 Real(1.5341553679414688425981898952416987e-01L),
807 Real(7.4295478991330687632977899141707872e-02L),
808 Real(7.1278762832147862035977841733532020e-02L),
809 Real(1.6623223223705792825395256602140459e-02L),
810 Real(1.4160772533794791868984026749196156e-02L),
811 Real(1.4539694958941854654807449467759690e-02L)
812 };
813
814 const unsigned int permutation_ids[n_wts]
815 = {1, 3, 3, 3, 3, 6, 6, 6, 6, 6, 6, 6}; // 55 total points
816
817 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
818
819 return;
820 }
821
822
823 // Dunavant's 17th-order rule has 61 points, while a
824 // comparable conical product rule would have 81 (16th and 17th orders).
825 //
826 // It can be found here:
827 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
828 //
829 // Zhang reports an identical rule in:
830 // L. Zhang, T. Cui, and H. Liu. "A set of symmetric quadrature rules
831 // on triangles and tetrahedra" Journal of Computational Mathematics,
832 // v. 27, no. 1, 2009, pp. 89-96.
833 //
834 // Note: the 61-point 17th-order rule obtained by Dunavant and Zhang
835 // does not appear to be unique. It is a solution in the sense that it
836 // minimizes the error in the least-squares minimization problem, but
837 // it involves too many unknowns and the Jacobian is therefore singular
838 // when attempting to improve the solution via Newton's method.
839 //
840 // Therefore, we prefer the following 63-point rule which
841 // I [JWP] found. It appears to be more accurate than the
842 // rule reported by Dunavant and Zhang, even though it has
843 // a few more points.
845 {
846 const unsigned int n_wts = 12;
847 const Real wts[n_wts] =
848 {
849 Real(1.7464603792572004485690588092246146e-02L),
850 Real(5.9429003555801725246549713984660076e-03L),
851 Real(1.2490753345169579649319736639588729e-02L),
852 Real(1.5386987188875607593083456905596468e-02L),
853 Real(1.1185807311917706362674684312990270e-02L),
854 Real(1.0301845740670206831327304917180007e-02L),
855 Real(1.1767783072977049696840016810370464e-02L),
856 Real(3.8045312849431209558329128678945240e-03L),
857 Real(4.5139302178876351271037137230354382e-03L),
858 Real(2.2178812517580586419412547665472893e-03L),
859 Real(5.2216271537483672304731416553063103e-03L),
860 Real(9.8381136389470256422419930926212114e-04L)
861 };
862
863 const Real a[n_wts] =
864 {
865 Real(2.8796825754667362165337965123570514e-01L),
866 Real(4.9216175986208465345536805750663939e-01L),
867 Real(4.6252866763171173685916780827044612e-01L),
868 Real(1.6730292951631792248498303276090273e-01L),
869 Real(1.5816335500814652972296428532213019e-01L),
870 Real(1.6352252138387564873002458959679529e-01L),
871 Real(6.2447680488959768233910286168417367e-01L),
872 Real(8.7317249935244454285263604347964179e-01L),
873 Real(3.4428164322282694677972239461699271e-01L),
874 Real(9.1584484467813674010523309855340209e-02L),
875 Real(2.0172088013378989086826623852040632e-01L),
876 Real(9.6538762758254643474731509845084691e-01L)
877 };
878
879 const Real b[n_wts] =
880 {
881 0.0,
882 0.0,
883 0.0,
884 Real(3.4429160695501713926320695771253348e-01L),
885 Real(2.2541623431550639817203145525444726e-01L),
886 Real(8.0670083153531811694942222940484991e-02L),
887 Real(6.5967451375050925655738829747288190e-02L),
888 Real(4.5677879890996762665044366994439565e-02L),
889 Real(1.1528411723154215812386518751976084e-02L),
890 Real(9.3057714323900610398389176844165892e-03L),
891 Real(1.5916814107619812717966560404970160e-02L),
892 Real(1.0734733163764032541125434215228937e-02L)
893 };
894
895 const unsigned int permutation_ids[n_wts]
896 = {3, 3, 3, 6, 6, 6, 6, 6, 6, 6, 6, 6}; // 63 total points
897
898 dunavant_rule2(wts, a, b, permutation_ids, n_wts);
899
900 return;
901
902 // _points.resize (61);
903 // _weights.resize(61);
904
905 // // The raw data for the quadrature rule.
906 // const Real p[15][4] = {
907 // { 1./3., 0., 0., 0.033437199290803e+00 / 2.0}, // 1-perm
908 // {0.005658918886452e+00, 0.497170540556774e+00, 0., 0.005093415440507e+00 / 2.0}, // 3-perm
909 // {0.035647354750751e+00, 0.482176322624625e+00, 0., 0.014670864527638e+00 / 2.0}, // 3-perm
910 // {0.099520061958437e+00, 0.450239969020782e+00, 0., 0.024350878353672e+00 / 2.0}, // 3-perm
911 // {0.199467521245206e+00, 0.400266239377397e+00, 0., 0.031107550868969e+00 / 2.0}, // 3-perm
912 // {0.495717464058095e+00, 0.252141267970953e+00, 0., 0.031257111218620e+00 / 2.0}, // 3-perm
913 // {0.675905990683077e+00, 0.162047004658461e+00, 0., 0.024815654339665e+00 / 2.0}, // 3-perm
914 // {0.848248235478508e+00, 0.075875882260746e+00, 0., 0.014056073070557e+00 / 2.0}, // 3-perm
915 // {0.968690546064356e+00, 0.015654726967822e+00, 0., 0.003194676173779e+00 / 2.0}, // 3-perm
916 // {0.010186928826919e+00, 0.334319867363658e+00, 0.655493203809423e+00, 0.008119655318993e+00 / 2.0}, // 6-perm
917 // {0.135440871671036e+00, 0.292221537796944e+00, 0.572337590532020e+00, 0.026805742283163e+00 / 2.0}, // 6-perm
918 // {0.054423924290583e+00, 0.319574885423190e+00, 0.626001190286228e+00, 0.018459993210822e+00 / 2.0}, // 6-perm
919 // {0.012868560833637e+00, 0.190704224192292e+00, 0.796427214974071e+00, 0.008476868534328e+00 / 2.0}, // 6-perm
920 // {0.067165782413524e+00, 0.180483211648746e+00, 0.752351005937729e+00, 0.018292796770025e+00 / 2.0}, // 6-perm
921 // {0.014663182224828e+00, 0.080711313679564e+00, 0.904625504095608e+00, 0.006665632004165e+00 / 2.0} // 6-perm
922 // };
923
924
925 // // Now call the dunavant routine to generate _points and _weights
926 // dunavant_rule(p, 15);
927
928 // return;
929 }
930
931
932
933 // Dunavant's 18th-order rule contains points outside the region and is therefore unsuitable
934 // for our FEM calculations. His 19th-order rule has 73 points, compared with 100 points for
935 // a comparable-order conical product rule.
936 //
937 // It was copied 23rd June 2008 from:
938 // http://people.scs.fsu.edu/~burkardt/f_src/dunavant/dunavant.f90
941 {
942 _points.resize (73);
943 _weights.resize(73);
944
945 // The raw data for the quadrature rule.
946 const Real rule_data[17][4] = {
947 { 1./3., 0., 0., 0.032906331388919e+00 / 2.0}, // 1-perm
948 {0.020780025853987e+00, 0.489609987073006e+00, 0., 0.010330731891272e+00 / 2.0}, // 3-perm
949 {0.090926214604215e+00, 0.454536892697893e+00, 0., 0.022387247263016e+00 / 2.0}, // 3-perm
950 {0.197166638701138e+00, 0.401416680649431e+00, 0., 0.030266125869468e+00 / 2.0}, // 3-perm
951 {0.488896691193805e+00, 0.255551654403098e+00, 0., 0.030490967802198e+00 / 2.0}, // 3-perm
952 {0.645844115695741e+00, 0.177077942152130e+00, 0., 0.024159212741641e+00 / 2.0}, // 3-perm
953 {0.779877893544096e+00, 0.110061053227952e+00, 0., 0.016050803586801e+00 / 2.0}, // 3-perm
954 {0.888942751496321e+00, 0.055528624251840e+00, 0., 0.008084580261784e+00 / 2.0}, // 3-perm
955 {0.974756272445543e+00, 0.012621863777229e+00, 0., 0.002079362027485e+00 / 2.0}, // 3-perm
956 {0.003611417848412e+00, 0.395754787356943e+00, 0.600633794794645e+00, 0.003884876904981e+00 / 2.0}, // 6-perm
957 {0.134466754530780e+00, 0.307929983880436e+00, 0.557603261588784e+00, 0.025574160612022e+00 / 2.0}, // 6-perm
958 {0.014446025776115e+00, 0.264566948406520e+00, 0.720987025817365e+00, 0.008880903573338e+00 / 2.0}, // 6-perm
959 {0.046933578838178e+00, 0.358539352205951e+00, 0.594527068955871e+00, 0.016124546761731e+00 / 2.0}, // 6-perm
960 {0.002861120350567e+00, 0.157807405968595e+00, 0.839331473680839e+00, 0.002491941817491e+00 / 2.0}, // 6-perm
961 {0.223861424097916e+00, 0.075050596975911e+00, 0.701087978926173e+00, 0.018242840118951e+00 / 2.0}, // 6-perm
962 {0.034647074816760e+00, 0.142421601113383e+00, 0.822931324069857e+00, 0.010258563736199e+00 / 2.0}, // 6-perm
963 {0.010161119296278e+00, 0.065494628082938e+00, 0.924344252620784e+00, 0.003799928855302e+00 / 2.0} // 6-perm
964 };
965
966
967 // Now call the dunavant routine to generate _points and _weights
968 dunavant_rule(rule_data, 17);
969
970 return;
971 }
972
973
974 // 20th-order rule by Wandzura.
975 //
976 // Stephen Wandzura, Hong Xiao,
977 // Symmetric Quadrature Rules on a Triangle,
978 // Computers and Mathematics with Applications,
979 // Volume 45, Number 12, June 2003, pages 1829-1840.
980 //
981 // Wandzura's work extends the work of Dunavant by providing degree
982 // 5,10,15,20,25, and 30 rules with positive weights for the triangle.
983 //
984 // Copied on 3rd July 2008 from:
985 // http://people.scs.fsu.edu/~burkardt/f_src/wandzura/wandzura.f90
987 {
988 // The equivalent conical product rule would have 121 points
989 _points.resize (85);
990 _weights.resize(85);
991
992 // The raw data for the quadrature rule.
993 const Real rule_data[19][4] = {
994 {0.33333333333333e+00, 0.0, 0.0, 0.2761042699769952e-01 / 2.0}, // 1-perm
995 {0.00150064932443e+00, 0.49924967533779e+00, 0.0, 0.1779029547326740e-02 / 2.0}, // 3-perm
996 {0.09413975193895e+00, 0.45293012403052e+00, 0.0, 0.2011239811396117e-01 / 2.0}, // 3-perm
997 {0.20447212408953e+00, 0.39776393795524e+00, 0.0, 0.2681784725933157e-01 / 2.0}, // 3-perm
998 {0.47099959493443e+00, 0.26450020253279e+00, 0.0, 0.2452313380150201e-01 / 2.0}, // 3-perm
999 {0.57796207181585e+00, 0.21101896409208e+00, 0.0, 0.1639457841069539e-01 / 2.0}, // 3-perm
1000 {0.78452878565746e+00, 0.10773560717127e+00, 0.0, 0.1479590739864960e-01 / 2.0}, // 3-perm
1001 {0.92186182432439e+00, 0.03906908783780e+00, 0.0, 0.4579282277704251e-02 / 2.0}, // 3-perm
1002 {0.97765124054134e+00, 0.01117437972933e+00, 0.0, 0.1651826515576217e-02 / 2.0}, // 3-perm
1003 {0.00534961818734e+00, 0.06354966590835e+00, 0.93110071590431e+00, 0.2349170908575584e-02 / 2.0}, // 6-perm
1004 {0.00795481706620e+00, 0.15710691894071e+00, 0.83493826399309e+00, 0.4465925754181793e-02 / 2.0}, // 6-perm
1005 {0.01042239828126e+00, 0.39564211436437e+00, 0.59393548735436e+00, 0.6099566807907972e-02 / 2.0}, // 6-perm
1006 {0.01096441479612e+00, 0.27316757071291e+00, 0.71586801449097e+00, 0.6891081327188203e-02 / 2.0}, // 6-perm
1007 {0.03856671208546e+00, 0.10178538248502e+00, 0.85964790542952e+00, 0.7997475072478163e-02 / 2.0}, // 6-perm
1008 {0.03558050781722e+00, 0.44665854917641e+00, 0.51776094300637e+00, 0.7386134285336024e-02 / 2.0}, // 6-perm
1009 {0.04967081636276e+00, 0.19901079414950e+00, 0.75131838948773e+00, 0.1279933187864826e-01 / 2.0}, // 6-perm
1010 {0.05851972508433e+00, 0.32426118369228e+00, 0.61721909122339e+00, 0.1725807117569655e-01 / 2.0}, // 6-perm
1011 {0.12149778700439e+00, 0.20853136321013e+00, 0.66997084978547e+00, 0.1867294590293547e-01 / 2.0}, // 6-perm
1012 {0.14071084494394e+00, 0.32317056653626e+00, 0.53611858851980e+00, 0.2281822405839526e-01 / 2.0} // 6-perm
1013 };
1014
1015
1016 // Now call the dunavant routine to generate _points and _weights
1017 dunavant_rule(rule_data, 19);
1018
1019 return;
1020 }
1021
1022
1023
1024 // 25th-order rule by Wandzura.
1025 //
1026 // Stephen Wandzura, Hong Xiao,
1027 // Symmetric Quadrature Rules on a Triangle,
1028 // Computers and Mathematics with Applications,
1029 // Volume 45, Number 12, June 2003, pages 1829-1840.
1030 //
1031 // Wandzura's work extends the work of Dunavant by providing degree
1032 // 5,10,15,20,25, and 30 rules with positive weights for the triangle.
1033 //
1034 // Copied on 3rd July 2008 from:
1035 // http://people.scs.fsu.edu/~burkardt/f_src/wandzura/wandzura.f90
1036 // case TWENTYFIRST: // fall through to 121 point conical product rule below
1041 {
1042 // The equivalent conical product rule would have 169 points
1043 _points.resize (126);
1044 _weights.resize(126);
1045
1046 // The raw data for the quadrature rule.
1047 const Real rule_data[26][4] = {
1048 {0.02794648307317e+00, 0.48602675846341e+00, 0.0, 0.8005581880020417e-02 / 2.0}, // 3-perm
1049 {0.13117860132765e+00, 0.43441069933617e+00, 0.0, 0.1594707683239050e-01 / 2.0}, // 3-perm
1050 {0.22022172951207e+00, 0.38988913524396e+00, 0.0, 0.1310914123079553e-01 / 2.0}, // 3-perm
1051 {0.40311353196039e+00, 0.29844323401980e+00, 0.0, 0.1958300096563562e-01 / 2.0}, // 3-perm
1052 {0.53191165532526e+00, 0.23404417233737e+00, 0.0, 0.1647088544153727e-01 / 2.0}, // 3-perm
1053 {0.69706333078196e+00, 0.15146833460902e+00, 0.0, 0.8547279074092100e-02 / 2.0}, // 3-perm
1054 {0.77453221290801e+00, 0.11273389354599e+00, 0.0, 0.8161885857226492e-02 / 2.0}, // 3-perm
1055 {0.84456861581695e+00, 0.07771569209153e+00, 0.0, 0.6121146539983779e-02 / 2.0}, // 3-perm
1056 {0.93021381277141e+00, 0.03489309361430e+00, 0.0, 0.2908498264936665e-02 / 2.0}, // 3-perm
1057 {0.98548363075813e+00, 0.00725818462093e+00, 0.0, 0.6922752456619963e-03 / 2.0}, // 3-perm
1058 {0.00129235270444e+00, 0.22721445215336e+00, 0.77149319514219e+00, 0.1248289199277397e-02 / 2.0}, // 6-perm
1059 {0.00539970127212e+00, 0.43501055485357e+00, 0.55958974387431e+00, 0.3404752908803022e-02 / 2.0}, // 6-perm
1060 {0.00638400303398e+00, 0.32030959927220e+00, 0.67330639769382e+00, 0.3359654326064051e-02 / 2.0}, // 6-perm
1061 {0.00502821150199e+00, 0.09175032228001e+00, 0.90322146621800e+00, 0.1716156539496754e-02 / 2.0}, // 6-perm
1062 {0.00682675862178e+00, 0.03801083585872e+00, 0.95516240551949e+00, 0.1480856316715606e-02 / 2.0}, // 6-perm
1063 {0.01001619963993e+00, 0.15742521848531e+00, 0.83255858187476e+00, 0.3511312610728685e-02 / 2.0}, // 6-perm
1064 {0.02575781317339e+00, 0.23988965977853e+00, 0.73435252704808e+00, 0.7393550149706484e-02 / 2.0}, // 6-perm
1065 {0.03022789811992e+00, 0.36194311812606e+00, 0.60782898375402e+00, 0.7983087477376558e-02 / 2.0}, // 6-perm
1066 {0.03050499010716e+00, 0.08355196095483e+00, 0.88594304893801e+00, 0.4355962613158041e-02 / 2.0}, // 6-perm
1067 {0.04595654736257e+00, 0.14844322073242e+00, 0.80560023190501e+00, 0.7365056701417832e-02 / 2.0}, // 6-perm
1068 {0.06744280054028e+00, 0.28373970872753e+00, 0.64881749073219e+00, 0.1096357284641955e-01 / 2.0}, // 6-perm
1069 {0.07004509141591e+00, 0.40689937511879e+00, 0.52305553346530e+00, 0.1174996174354112e-01 / 2.0}, // 6-perm
1070 {0.08391152464012e+00, 0.19411398702489e+00, 0.72197448833499e+00, 0.1001560071379857e-01 / 2.0}, // 6-perm
1071 {0.12037553567715e+00, 0.32413434700070e+00, 0.55549011732214e+00, 0.1330964078762868e-01 / 2.0}, // 6-perm
1072 {0.14806689915737e+00, 0.22927748355598e+00, 0.62265561728665e+00, 0.1415444650522614e-01 / 2.0}, // 6-perm
1073 {0.19177186586733e+00, 0.32561812259598e+00, 0.48261001153669e+00, 0.1488137956116801e-01 / 2.0} // 6-perm
1074 };
1075
1076
1077 // Now call the dunavant routine to generate _points and _weights
1078 dunavant_rule(rule_data, 26);
1079
1080 return;
1081 }
1082
1083
1084
1085 // 30th-order rule by Wandzura.
1086 //
1087 // Stephen Wandzura, Hong Xiao,
1088 // Symmetric Quadrature Rules on a Triangle,
1089 // Computers and Mathematics with Applications,
1090 // Volume 45, Number 12, June 2003, pages 1829-1840.
1091 //
1092 // Wandzura's work extends the work of Dunavant by providing degree
1093 // 5,10,15,20,25, and 30 rules with positive weights for the triangle.
1094 //
1095 // Copied on 3rd July 2008 from:
1096 // http://people.scs.fsu.edu/~burkardt/f_src/wandzura/wandzura.f90
1102 {
1103 // The equivalent conical product rule would have 256 points
1104 _points.resize (175);
1105 _weights.resize(175);
1106
1107 // The raw data for the quadrature rule.
1108 const Real rule_data[36][4] = {
1109 {0.33333333333333e+00, 0.0, 0.0, 0.1557996020289920e-01 / 2.0}, // 1-perm
1110 {0.00733011643277e+00, 0.49633494178362e+00, 0.0, 0.3177233700534134e-02 / 2.0}, // 3-perm
1111 {0.08299567580296e+00, 0.45850216209852e+00, 0.0, 0.1048342663573077e-01 / 2.0}, // 3-perm
1112 {0.15098095612541e+00, 0.42450952193729e+00, 0.0, 0.1320945957774363e-01 / 2.0}, // 3-perm
1113 {0.23590585989217e+00, 0.38204707005392e+00, 0.0, 0.1497500696627150e-01 / 2.0}, // 3-perm
1114 {0.43802430840785e+00, 0.28098784579608e+00, 0.0, 0.1498790444338419e-01 / 2.0}, // 3-perm
1115 {0.54530204829193e+00, 0.22734897585403e+00, 0.0, 0.1333886474102166e-01 / 2.0}, // 3-perm
1116 {0.65088177698254e+00, 0.17455911150873e+00, 0.0, 0.1088917111390201e-01 / 2.0}, // 3-perm
1117 {0.75348314559713e+00, 0.12325842720144e+00, 0.0, 0.8189440660893461e-02 / 2.0}, // 3-perm
1118 {0.83983154221561e+00, 0.08008422889220e+00, 0.0, 0.5575387588607785e-02 / 2.0}, // 3-perm
1119 {0.90445106518420e+00, 0.04777446740790e+00, 0.0, 0.3191216473411976e-02 / 2.0}, // 3-perm
1120 {0.95655897063972e+00, 0.02172051468014e+00, 0.0, 0.1296715144327045e-02 / 2.0}, // 3-perm
1121 {0.99047064476913e+00, 0.00476467761544e+00, 0.0, 0.2982628261349172e-03 / 2.0}, // 3-perm
1122 {0.00092537119335e+00, 0.41529527091331e+00, 0.58377935789334e+00, 0.9989056850788964e-03 / 2.0}, // 6-perm
1123 {0.00138592585556e+00, 0.06118990978535e+00, 0.93742416435909e+00, 0.4628508491732533e-03 / 2.0}, // 6-perm
1124 {0.00368241545591e+00, 0.16490869013691e+00, 0.83140889440718e+00, 0.1234451336382413e-02 / 2.0}, // 6-perm
1125 {0.00390322342416e+00, 0.02503506223200e+00, 0.97106171434384e+00, 0.5707198522432062e-03 / 2.0}, // 6-perm
1126 {0.00323324815501e+00, 0.30606446515110e+00, 0.69070228669389e+00, 0.1126946125877624e-02 / 2.0}, // 6-perm
1127 {0.00646743211224e+00, 0.10707328373022e+00, 0.88645928415754e+00, 0.1747866949407337e-02 / 2.0}, // 6-perm
1128 {0.00324747549133e+00, 0.22995754934558e+00, 0.76679497516308e+00, 0.1182818815031657e-02 / 2.0}, // 6-perm
1129 {0.00867509080675e+00, 0.33703663330578e+00, 0.65428827588746e+00, 0.1990839294675034e-02 / 2.0}, // 6-perm
1130 {0.01559702646731e+00, 0.05625657618206e+00, 0.92814639735063e+00, 0.1900412795035980e-02 / 2.0}, // 6-perm
1131 {0.01797672125369e+00, 0.40245137521240e+00, 0.57957190353391e+00, 0.4498365808817451e-02 / 2.0}, // 6-perm
1132 {0.01712424535389e+00, 0.24365470201083e+00, 0.73922105263528e+00, 0.3478719460274719e-02 / 2.0}, // 6-perm
1133 {0.02288340534658e+00, 0.16538958561453e+00, 0.81172700903888e+00, 0.4102399036723953e-02 / 2.0}, // 6-perm
1134 {0.03273759728777e+00, 0.09930187449585e+00, 0.86796052821639e+00, 0.4021761549744162e-02 / 2.0}, // 6-perm
1135 {0.03382101234234e+00, 0.30847833306905e+00, 0.65770065458860e+00, 0.6033164660795066e-02 / 2.0}, // 6-perm
1136 {0.03554761446002e+00, 0.46066831859211e+00, 0.50378406694787e+00, 0.3946290302129598e-02 / 2.0}, // 6-perm
1137 {0.05053979030687e+00, 0.21881529945393e+00, 0.73064491023920e+00, 0.6644044537680268e-02 / 2.0}, // 6-perm
1138 {0.05701471491573e+00, 0.37920955156027e+00, 0.56377573352399e+00, 0.8254305856078458e-02 / 2.0}, // 6-perm
1139 {0.06415280642120e+00, 0.14296081941819e+00, 0.79288637416061e+00, 0.6496056633406411e-02 / 2.0}, // 6-perm
1140 {0.08050114828763e+00, 0.28373128210592e+00, 0.63576756960645e+00, 0.9252778144146602e-02 / 2.0}, // 6-perm
1141 {0.10436706813453e+00, 0.19673744100444e+00, 0.69889549086103e+00, 0.9164920726294280e-02 / 2.0}, // 6-perm
1142 {0.11384489442875e+00, 0.35588914121166e+00, 0.53026596435959e+00, 0.1156952462809767e-01 / 2.0}, // 6-perm
1143 {0.14536348771552e+00, 0.25981868535191e+00, 0.59481782693256e+00, 0.1176111646760917e-01 / 2.0}, // 6-perm
1144 {0.18994565282198e+00, 0.32192318123130e+00, 0.48813116594672e+00, 0.1382470218216540e-01 / 2.0} // 6-perm
1145 };
1146
1147
1148 // Now call the dunavant routine to generate _points and _weights
1149 dunavant_rule(rule_data, 36);
1150
1151 return;
1152 }
1153
1154 default:
1155 {
1156 chi_math::QuadratureConical conical(order);
1158 qpoints_.swap(conical.qpoints_);
1159 weights_.swap(conical.weights_);
1160 }
1161 }//switch order
1162}
1163
1164/*** The Dunavant rules are for triangles. This function takes
1165permutation points and weights in a specific format as input and
1166fills the _points and _weights vectors.*/
1167void chi_math::QuadratureTriangle::dunavant_rule(const double rule_data[][4],
1168 const unsigned int n_pts)
1169{
1170 // The input data array has 4 columns. The first 3 are the permutation points.
1171 // The last column is the weights for a given set of permutation points. A zero
1172 // in two of the first 3 columns implies the point is a 1-permutation (centroid).
1173 // A zero in one of the first 3 columns implies the point is a 3-permutation.
1174 // Otherwise each point is assumed to be a 6-permutation.
1175
1176 // Always insert into the points & weights vector relative to the offset
1177 unsigned int offset=0;
1178
1179
1180 for (unsigned int p=0; p<n_pts; ++p)
1181 {
1182
1183 // There must always be a non-zero entry to start the row
1184 assert( rule_data[p][0] != static_cast<double >(0.0) );
1185
1186 // A zero weight may imply you did not set up the raw data correctly
1187 assert ( rule_data[p][3] != static_cast<double>(0.0) );
1188
1189 // What kind of point is this?
1190 // One non-zero entry in first 3 cols ? 1-perm (centroid) point = 1
1191 // Two non-zero entries in first 3 cols ? 3-perm point = 3
1192 // Three non-zero entries ? 6-perm point = 6
1193 unsigned int pointtype=1;
1194
1195 if (rule_data[p][1] != static_cast<double>(0.0))
1196 {
1197 if (rule_data[p][2] != static_cast<double>(0.0))
1198 pointtype = 6;
1199 else
1200 pointtype = 3;
1201 }
1202
1203 switch (pointtype)
1204 {
1205 case 1:
1206 {
1207 // Be sure we have enough space to insert this point
1208 assert (offset + 0 < qpoints_.size());
1209
1210 // The point has only a single permutation (the centroid!)
1211 qpoints_[offset + 0] = chi_mesh::Vector3(rule_data[p][0], rule_data[p][0]);
1212
1213 // The weight is always the last entry in the row.
1214 weights_[offset + 0] = rule_data[p][3];
1215
1216 offset += 1;
1217 break;
1218 }
1219
1220 case 3:
1221 {
1222 // Be sure we have enough space to insert these points
1223 assert (offset + 2 < qpoints_.size());
1224
1225 // Here it's understood the second entry is to be used twice, and
1226 // thus there are three possible permutations.
1227 qpoints_[offset + 0] = chi_mesh::Vector3(rule_data[p][0], rule_data[p][1]);
1228 qpoints_[offset + 1] = chi_mesh::Vector3(rule_data[p][1], rule_data[p][0]);
1229 qpoints_[offset + 2] = chi_mesh::Vector3(rule_data[p][1], rule_data[p][1]);
1230
1231 for (unsigned int j=0; j<3; ++j)
1232 weights_[offset + j] = rule_data[p][3];
1233
1234 offset += 3;
1235 break;
1236 }
1237
1238 case 6:
1239 {
1240 // Be sure we have enough space to insert these points
1241 assert (offset + 5 < qpoints_.size());
1242
1243 // Three individual entries with six permutations.
1244 qpoints_[offset + 0] = chi_mesh::Vector3(rule_data[p][0], rule_data[p][1]);
1245 qpoints_[offset + 1] = chi_mesh::Vector3(rule_data[p][0], rule_data[p][2]);
1246 qpoints_[offset + 2] = chi_mesh::Vector3(rule_data[p][1], rule_data[p][0]);
1247 qpoints_[offset + 3] = chi_mesh::Vector3(rule_data[p][1], rule_data[p][2]);
1248 qpoints_[offset + 4] = chi_mesh::Vector3(rule_data[p][2], rule_data[p][0]);
1249 qpoints_[offset + 5] = chi_mesh::Vector3(rule_data[p][2], rule_data[p][1]);
1250
1251 for (unsigned int j=0; j<6; ++j)
1252 weights_[offset + j] = rule_data[p][3];
1253
1254 offset += 6;
1255 break;
1256 }
1257
1258 default:
1259 throw std::invalid_argument(std::string(__FUNCTION__) +
1260 ": Don't know what to do with this many permutation points!");
1261 }
1262 }
1263}
1264
1265// A number of different rules for triangles can be described by
1266// permutations of the following types of points:
1267// I: "1"-permutation, (1/3,1/3) (single point only)
1268// II: 3-permutation, (a,a,1-2a)
1269// III: 6-permutation, (a,b,1-a-b)
1270// The weights for a given set of permutations are all the same.
1272 const double * a,
1273 const double * b,
1274 const unsigned int * permutation_ids,
1275 unsigned int n_wts)
1276{
1277 auto& _points = qpoints_;
1278 auto& _weights = weights_;
1279
1280 typedef chi_mesh::Vector3 Point;
1281
1282 // Figure out how many total points by summing up the entries
1283 // in the permutation_ids array, and resize the _points and _weights
1284 // vectors appropriately.
1285 unsigned int total_pts = 0;
1286 for (unsigned int p=0; p<n_wts; ++p)
1287 total_pts += permutation_ids[p];
1288
1289 // Resize point and weight vectors appropriately.
1290 _points.resize(total_pts);
1291 _weights.resize(total_pts);
1292
1293 // Always insert into the points & weights vector relative to the offset
1294 unsigned int offset=0;
1295
1296 for (unsigned int p=0; p<n_wts; ++p)
1297 {
1298 switch (permutation_ids[p])
1299 {
1300 case 1:
1301 {
1302 // The point has only a single permutation (the centroid!)
1303 // So we don't even need to look in the a or b arrays.
1304 _points [offset + 0] = Point(double(1)/3, double(1)/3);
1305 _weights[offset + 0] = wts[p];
1306
1307 offset += 1;
1308 break;
1309 }
1310
1311
1312 case 3:
1313 {
1314 // For this type of rule, don't need to look in the b array.
1315 _points[offset + 0] = Point(a[p], a[p]); // (a,a)
1316 _points[offset + 1] = Point(a[p], 1-2*a[p]); // (a,1-2a)
1317 _points[offset + 2] = Point(1-2*a[p], a[p]); // (1-2a,a)
1318
1319 for (unsigned int j=0; j<3; ++j)
1320 _weights[offset + j] = wts[p];
1321
1322 offset += 3;
1323 break;
1324 }
1325
1326 case 6:
1327 {
1328 // This type of point uses all 3 arrays...
1329 _points[offset + 0] = Point(a[p], b[p]);
1330 _points[offset + 1] = Point(b[p], a[p]);
1331 _points[offset + 2] = Point(a[p], 1-a[p]-b[p]);
1332 _points[offset + 3] = Point(1-a[p]-b[p], a[p]);
1333 _points[offset + 4] = Point(b[p], 1-a[p]-b[p]);
1334 _points[offset + 5] = Point(1-a[p]-b[p], b[p]);
1335
1336 for (unsigned int j=0; j<6; ++j)
1337 _weights[offset + j] = wts[p];
1338
1339 offset += 6;
1340 break;
1341 }
1342
1343 default:
1344 throw std::invalid_argument(std::string(__FUNCTION__) +
1345 "Unknown permutation id: " + std::to_string(permutation_ids[p]) + "!");
1346 }
1347 }
1348
1349}
std::vector< chi_math::QuadraturePointXYZ > qpoints_
Definition: quadrature.h:37
std::vector< double > weights_
Definition: quadrature.h:38
void dunavant_rule2(const double *wts, const double *a, const double *b, const unsigned int *permutation_ids, const unsigned int n_wts)
void dunavant_rule(const double rule_data[][4], const unsigned int n_pts)
QuadratureTriangle(QuadratureOrder order)
QuadratureOrder
Definition: quadrature.h:12
VectorN< 3 > Vector3