libMesh
quadrature_gauss_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 // Local includes
21 #include "libmesh/quadrature_gauss.h"
22 #include "libmesh/quadrature_conical.h"
23 #include "libmesh/quadrature_gm.h"
24 #include "libmesh/enum_to_string.h"
25 
26 namespace libMesh
27 {
28 
29 void QGauss::init_3D(const ElemType, unsigned int)
30 {
31 #if LIBMESH_DIM == 3
32 
33  //-----------------------------------------------------------------------
34  // 3D quadrature rules
35  switch (_type)
36  {
37  //---------------------------------------------
38  // Hex quadrature rules
39  case HEX8:
40  case HEX20:
41  case HEX27:
42  {
43  // We compute the 3D quadrature rule as a tensor
44  // product of the 1D quadrature rule.
45  QGauss q1D(1, _order);
46  q1D.init(EDGE2, _p_level);
47  tensor_product_hex( q1D );
48  return;
49  }
50 
51 
52 
53  //---------------------------------------------
54  // Tetrahedral quadrature rules
55  case TET4:
56  case TET10:
57  {
58  switch(get_order())
59  {
60  // Taken from pg. 222 of "The finite element method," vol. 1
61  // ed. 5 by Zienkiewicz & Taylor
62  case CONSTANT:
63  case FIRST:
64  {
65  // Exact for linears
66  _points.resize(1);
67  _weights.resize(1);
68 
69 
70  _points[0](0) = .25;
71  _points[0](1) = .25;
72  _points[0](2) = .25;
73 
74  _weights[0] = Real(1)/6;
75 
76  return;
77  }
78  case SECOND:
79  {
80  // Exact for quadratics
81  _points.resize(4);
82  _weights.resize(4);
83 
84 
85  const Real a = .585410196624969;
86  const Real b = .138196601125011;
87 
88  _points[0](0) = a;
89  _points[0](1) = b;
90  _points[0](2) = b;
91 
92  _points[1](0) = b;
93  _points[1](1) = a;
94  _points[1](2) = b;
95 
96  _points[2](0) = b;
97  _points[2](1) = b;
98  _points[2](2) = a;
99 
100  _points[3](0) = b;
101  _points[3](1) = b;
102  _points[3](2) = b;
103 
104 
105 
106  _weights[0] = Real(1)/24;
107  _weights[1] = _weights[0];
108  _weights[2] = _weights[0];
109  _weights[3] = _weights[0];
110 
111  return;
112  }
113 
114 
115 
116  // Can be found in the class notes
117  // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps
118  // by Flaherty.
119  //
120  // Caution: this rule has a negative weight and may be
121  // unsuitable for some problems.
122  // Exact for cubics.
123  //
124  // Note: Keast (see ref. elsewhere in this file) also gives
125  // a third-order rule with positive weights, but it contains points
126  // on the ref. elt. boundary, making it less suitable for FEM calculations.
127  case THIRD:
128  {
130  {
131  _points.resize(5);
132  _weights.resize(5);
133 
134 
135  _points[0](0) = .25;
136  _points[0](1) = .25;
137  _points[0](2) = .25;
138 
139  _points[1](0) = .5;
140  _points[1](1) = Real(1)/6;
141  _points[1](2) = Real(1)/6;
142 
143  _points[2](0) = Real(1)/6;
144  _points[2](1) = .5;
145  _points[2](2) = Real(1)/6;
146 
147  _points[3](0) = Real(1)/6;
148  _points[3](1) = Real(1)/6;
149  _points[3](2) = .5;
150 
151  _points[4](0) = Real(1)/6;
152  _points[4](1) = Real(1)/6;
153  _points[4](2) = Real(1)/6;
154 
155 
156  _weights[0] = Real(-2)/15;
157  _weights[1] = .075;
158  _weights[2] = _weights[1];
159  _weights[3] = _weights[1];
160  _weights[4] = _weights[1];
161 
162  return;
163  } // end if (allow_rules_with_negative_weights)
164  else
165  {
166  // If a rule with positive weights is required, the 2x2x2 conical
167  // product rule is third-order accurate and has less points than
168  // the next-available positive-weight rule at FIFTH order.
169  QConical conical_rule(3, _order);
170  conical_rule.init(_type, _p_level);
171 
172  // Swap points and weights with the about-to-be destroyed rule.
173  _points.swap (conical_rule.get_points() );
174  _weights.swap(conical_rule.get_weights());
175 
176  return;
177  }
178  // Note: if !allow_rules_with_negative_weights, fall through to next case.
179  }
180 
181 
182 
183  // Originally a Keast rule,
184  // Patrick Keast,
185  // Moderate Degree Tetrahedral Quadrature Formulas,
186  // Computer Methods in Applied Mechanics and Engineering,
187  // Volume 55, Number 3, May 1986, pages 339-348.
188  //
189  // Can also be found the class notes
190  // http://www.cs.rpi.edu/~flaherje/FEM/fem6.ps
191  // by Flaherty.
192  //
193  // Caution: this rule has a negative weight and may be
194  // unsuitable for some problems.
195  case FOURTH:
196  {
198  {
199  _points.resize(11);
200  _weights.resize(11);
201 
202  // The raw data for the quadrature rule.
203  const Real rule_data[3][4] = {
204  {0.250000000000000000e+00, 0., 0., -0.131555555555555556e-01}, // 1
205  {0.785714285714285714e+00, 0.714285714285714285e-01, 0., 0.762222222222222222e-02}, // 4
206  {0.399403576166799219e+00, 0., 0.100596423833200785e+00, 0.248888888888888889e-01} // 6
207  };
208 
209 
210  // Now call the keast routine to generate _points and _weights
211  keast_rule(rule_data, 3);
212 
213  return;
214  } // end if (allow_rules_with_negative_weights)
215  // Note: if !allow_rules_with_negative_weights, fall through to next case.
216  }
217 
218  libmesh_fallthrough();
219 
220 
221  // Walkington's fifth-order 14-point rule from
222  // "Quadrature on Simplices of Arbitrary Dimension"
223  //
224  // We originally had a Keast rule here, but this rule had
225  // more points than an equivalent rule by Walkington and
226  // also contained points on the boundary of the ref. elt,
227  // making it less suitable for FEM calculations.
228  case FIFTH:
229  {
230  _points.resize(14);
231  _weights.resize(14);
232 
233  // permutations of these points and suitably-modified versions of
234  // these points are the quadrature point locations
235  const Real a[3] = {0.31088591926330060980, // a1 from the paper
236  0.092735250310891226402, // a2 from the paper
237  0.045503704125649649492}; // a3 from the paper
238 
239  // weights. a[] and wt[] are the only floating-point inputs required
240  // for this rule.
241  const Real wt[3] = {0.018781320953002641800, // w1 from the paper
242  0.012248840519393658257, // w2 from the paper
243  0.0070910034628469110730}; // w3 from the paper
244 
245  // The first two sets of 4 points are formed in a similar manner
246  for (unsigned int i=0; i<2; ++i)
247  {
248  // Where we will insert values into _points and _weights
249  const unsigned int offset=4*i;
250 
251  // Stuff points and weights values into their arrays
252  const Real b = 1. - 3.*a[i];
253 
254  // Here are the permutations. Order of these is not important,
255  // all have the same weight
256  _points[offset + 0] = Point(a[i], a[i], a[i]);
257  _points[offset + 1] = Point(a[i], b, a[i]);
258  _points[offset + 2] = Point( b, a[i], a[i]);
259  _points[offset + 3] = Point(a[i], a[i], b);
260 
261  // These 4 points all have the same weights
262  for (unsigned int j=0; j<4; ++j)
263  _weights[offset + j] = wt[i];
264  } // end for
265 
266 
267  {
268  // The third set contains 6 points and is formed a little differently
269  const unsigned int offset = 8;
270  const Real b = 0.5*(1. - 2.*a[2]);
271 
272  // Here are the permutations. Order of these is not important,
273  // all have the same weight
274  _points[offset + 0] = Point(b , b, a[2]);
275  _points[offset + 1] = Point(b , a[2], a[2]);
276  _points[offset + 2] = Point(a[2], a[2], b);
277  _points[offset + 3] = Point(a[2], b, a[2]);
278  _points[offset + 4] = Point( b, a[2], b);
279  _points[offset + 5] = Point(a[2], b, b);
280 
281  // These 6 points all have the same weights
282  for (unsigned int j=0; j<6; ++j)
283  _weights[offset + j] = wt[2];
284  }
285 
286 
287  // Original rule by Keast, unsuitable because it has points on the
288  // reference element boundary!
289  // _points.resize(15);
290  // _weights.resize(15);
291 
292  // _points[0](0) = 0.25;
293  // _points[0](1) = 0.25;
294  // _points[0](2) = 0.25;
295 
296  // {
297  // const Real a = 0.;
298  // const Real b = Real(1)/3;
299 
300  // _points[1](0) = a;
301  // _points[1](1) = b;
302  // _points[1](2) = b;
303 
304  // _points[2](0) = b;
305  // _points[2](1) = a;
306  // _points[2](2) = b;
307 
308  // _points[3](0) = b;
309  // _points[3](1) = b;
310  // _points[3](2) = a;
311 
312  // _points[4](0) = b;
313  // _points[4](1) = b;
314  // _points[4](2) = b;
315  // }
316  // {
317  // const Real a = Real(8)/11;
318  // const Real b = Real(1)/11;
319 
320  // _points[5](0) = a;
321  // _points[5](1) = b;
322  // _points[5](2) = b;
323 
324  // _points[6](0) = b;
325  // _points[6](1) = a;
326  // _points[6](2) = b;
327 
328  // _points[7](0) = b;
329  // _points[7](1) = b;
330  // _points[7](2) = a;
331 
332  // _points[8](0) = b;
333  // _points[8](1) = b;
334  // _points[8](2) = b;
335  // }
336  // {
337  // const Real a = 0.066550153573664;
338  // const Real b = 0.433449846426336;
339 
340  // _points[9](0) = b;
341  // _points[9](1) = a;
342  // _points[9](2) = a;
343 
344  // _points[10](0) = a;
345  // _points[10](1) = a;
346  // _points[10](2) = b;
347 
348  // _points[11](0) = a;
349  // _points[11](1) = b;
350  // _points[11](2) = b;
351 
352  // _points[12](0) = b;
353  // _points[12](1) = a;
354  // _points[12](2) = b;
355 
356  // _points[13](0) = b;
357  // _points[13](1) = b;
358  // _points[13](2) = a;
359 
360  // _points[14](0) = a;
361  // _points[14](1) = b;
362  // _points[14](2) = a;
363  // }
364 
365  // _weights[0] = 0.030283678097089;
366  // _weights[1] = 0.006026785714286;
367  // _weights[2] = _weights[1];
368  // _weights[3] = _weights[1];
369  // _weights[4] = _weights[1];
370  // _weights[5] = 0.011645249086029;
371  // _weights[6] = _weights[5];
372  // _weights[7] = _weights[5];
373  // _weights[8] = _weights[5];
374  // _weights[9] = 0.010949141561386;
375  // _weights[10] = _weights[9];
376  // _weights[11] = _weights[9];
377  // _weights[12] = _weights[9];
378  // _weights[13] = _weights[9];
379  // _weights[14] = _weights[9];
380 
381  return;
382  }
383 
384 
385 
386 
387  // This rule is originally from Keast:
388  // Patrick Keast,
389  // Moderate Degree Tetrahedral Quadrature Formulas,
390  // Computer Methods in Applied Mechanics and Engineering,
391  // Volume 55, Number 3, May 1986, pages 339-348.
392  //
393  // It is accurate on 6th-degree polynomials and has 24 points
394  // vs. 64 for the comparable conical product rule.
395  //
396  // Values copied 24th June 2008 from:
397  // http://people.scs.fsu.edu/~burkardt/f_src/keast/keast.f90
398  case SIXTH:
399  {
400  _points.resize (24);
401  _weights.resize(24);
402 
403  // The raw data for the quadrature rule.
404  const Real rule_data[4][4] = {
405  {0.356191386222544953e+00 , 0.214602871259151684e+00 , 0., 0.00665379170969464506e+00}, // 4
406  {0.877978124396165982e+00 , 0.0406739585346113397e+00, 0., 0.00167953517588677620e+00}, // 4
407  {0.0329863295731730594e+00, 0.322337890142275646e+00 , 0., 0.00922619692394239843e+00}, // 4
408  {0.0636610018750175299e+00, 0.269672331458315867e+00 , 0.603005664791649076e+00, 0.00803571428571428248e+00} // 12
409  };
410 
411 
412  // Now call the keast routine to generate _points and _weights
413  keast_rule(rule_data, 4);
414 
415  return;
416  }
417 
418 
419  // Keast's 31 point, 7th-order rule contains points on the reference
420  // element boundary, so we've decided not to include it here.
421  //
422  // Keast's 8th-order rule has 45 points. and a negative
423  // weight, so if you've explicitly disallowed such rules
424  // you will fall through to the conical product rules
425  // below.
426  case SEVENTH:
427  case EIGHTH:
428  {
430  {
431  _points.resize (45);
432  _weights.resize(45);
433 
434  // The raw data for the quadrature rule.
435  const Real rule_data[7][4] = {
436  {0.250000000000000000e+00, 0., 0., -0.393270066412926145e-01}, // 1
437  {0.617587190300082967e+00, 0.127470936566639015e+00, 0., 0.408131605934270525e-02}, // 4
438  {0.903763508822103123e+00, 0.320788303926322960e-01, 0., 0.658086773304341943e-03}, // 4
439  {0.497770956432810185e-01, 0., 0.450222904356718978e+00, 0.438425882512284693e-02}, // 6
440  {0.183730447398549945e+00, 0., 0.316269552601450060e+00, 0.138300638425098166e-01}, // 6
441  {0.231901089397150906e+00, 0.229177878448171174e-01, 0.513280033360881072e+00, 0.424043742468372453e-02}, // 12
442  {0.379700484718286102e-01, 0.730313427807538396e+00, 0.193746475248804382e+00, 0.223873973961420164e-02} // 12
443  };
444 
445 
446  // Now call the keast routine to generate _points and _weights
447  keast_rule(rule_data, 7);
448 
449  return;
450  } // end if (allow_rules_with_negative_weights)
451  // Note: if !allow_rules_with_negative_weights, fall through to next case.
452  }
453 
454  libmesh_fallthrough();
455 
456 
457  // Fall back on Grundmann-Moller or Conical Product rules at high orders.
458  default:
459  {
461  {
462  // The Grundmann-Moller rules are defined to arbitrary order and
463  // can have significantly fewer evaluation points than conical product
464  // rules. If you allow rules with negative weights, the GM rules
465  // will be more efficient for degree up to 33 (for degree 34 and
466  // higher, CP is more efficient!) but may be more susceptible
467  // to round-off error. Safest is to disallow rules with negative
468  // weights, but this decision should be made on a case-by-case basis.
469  QGrundmann_Moller gm_rule(3, _order);
470  gm_rule.init(_type, _p_level);
471 
472  // Swap points and weights with the about-to-be destroyed rule.
473  _points.swap (gm_rule.get_points() );
474  _weights.swap(gm_rule.get_weights());
475 
476  return;
477  }
478 
479  else
480  {
481  // The following quadrature rules are generated as
482  // conical products. These tend to be non-optimal
483  // (use too many points, cluster points in certain
484  // regions of the domain) but they are quite easy to
485  // automatically generate using a 1D Gauss rule on
486  // [0,1] and two 1D Jacobi-Gauss rules on [0,1].
487  QConical conical_rule(3, _order);
488  conical_rule.init(_type, _p_level);
489 
490  // Swap points and weights with the about-to-be destroyed rule.
491  _points.swap (conical_rule.get_points() );
492  _weights.swap(conical_rule.get_weights());
493 
494  return;
495  }
496  }
497  }
498  } // end case TET4,TET10
499 
500 
501 
502  //---------------------------------------------
503  // Prism quadrature rules
504  case PRISM6:
505  case PRISM15:
506  case PRISM18:
507  {
508  // We compute the 3D quadrature rule as a tensor
509  // product of the 1D quadrature rule and a 2D
510  // triangle quadrature rule
511 
512  QGauss q1D(1,_order);
513  QGauss q2D(2,_order);
514 
515  // Initialize
516  q1D.init(EDGE2, _p_level);
517  q2D.init(TRI3, _p_level);
518 
519  tensor_product_prism(q1D, q2D);
520 
521  return;
522  }
523 
524 
525 
526  //---------------------------------------------
527  // Pyramid
528  case PYRAMID5:
529  case PYRAMID13:
530  case PYRAMID14:
531  {
532  // We compute the Pyramid rule as a conical product of a
533  // Jacobi rule with alpha==2 on the interval [0,1] two 1D
534  // Gauss rules with suitably modified points. The idea comes
535  // from: Stroud, A.H. "Approximate Calculation of Multiple
536  // Integrals."
537  QConical conical_rule(3, _order);
538  conical_rule.init(_type, _p_level);
539 
540  // Swap points and weights with the about-to-be destroyed rule.
541  _points.swap (conical_rule.get_points() );
542  _weights.swap(conical_rule.get_weights());
543 
544  return;
545 
546  }
547 
548 
549 
550  //---------------------------------------------
551  // Unsupported type
552  default:
553  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
554  }
555 #endif
556 }
557 
558 } // namespace libMesh
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::QBase::_p_level
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:357
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::QBase::tensor_product_hex
void tensor_product_hex(const QBase &q1D)
Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D.
Definition: quadrature.C:198
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::QConical
This class implements the so-called conical product quadrature rules for Tri and Tet elements.
Definition: quadrature_conical.h:43
libMesh::SIXTH
Definition: enum_order.h:47
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::QGauss::init_3D
virtual void init_3D(const ElemType, unsigned int) override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature_gauss_3D.C:29
libMesh::FIFTH
Definition: enum_order.h:46
libMesh::SECOND
Definition: enum_order.h:43
libMesh::QBase::init
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the data structures for a quadrature rule for an element of type type.
Definition: quadrature.C:59
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::QBase::_order
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly.
Definition: quadrature.h:345
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::QBase::get_points
const std::vector< Point > & get_points() const
Definition: quadrature.h:141
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::EIGHTH
Definition: enum_order.h:49
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::QBase::tensor_product_prism
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule.
Definition: quadrature.C:225
libMesh::FOURTH
Definition: enum_order.h:45
libMesh::QBase::get_weights
const std::vector< Real > & get_weights() const
Definition: quadrature.h:153
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::SEVENTH
Definition: enum_order.h:48
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::QGrundmann_Moller
This class implements the Grundmann-Moller quadrature rules for tetrahedra.
Definition: quadrature_gm.h:96
libMesh::QBase::allow_rules_with_negative_weights
bool allow_rules_with_negative_weights
Flag (default true) controlling the use of quadrature rules with negative weights.
Definition: quadrature.h:254
libMesh::THIRD
Definition: enum_order.h:44
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QBase::get_order
Order get_order() const
Definition: quadrature.h:211
libMesh::PRISM18
Definition: enum_elem_type.h:52
libMesh::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::FIRST
Definition: enum_order.h:42
libMesh::PYRAMID13
Definition: enum_elem_type.h:54
libMesh::PYRAMID14
Definition: enum_elem_type.h:55
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::QGauss::keast_rule
void keast_rule(const Real rule_data[][4], const unsigned int n_pts)
The Keast rules are for tets.
Definition: quadrature_gauss.C:38