libMesh
quadrature_monomial_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_monomial.h"
22 #include "libmesh/quadrature_gauss.h"
23 
24 namespace libMesh
25 {
26 
27 
28 void QMonomial::init_3D(const ElemType, unsigned int)
29 {
30 
31  switch (_type)
32  {
33  //---------------------------------------------
34  // Hex quadrature rules
35  case HEX8:
36  case HEX20:
37  case HEX27:
38  {
39  switch(get_order())
40  {
41 
42  // The CONSTANT/FIRST rule is the 1-point Gauss "product" rule...we fall
43  // through to the default case for this rule.
44 
45  case SECOND:
46  case THIRD:
47  {
48  // A degree 3, 6-point, "rotationally-symmetric" rule by
49  // Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
50  //
51  // Warning: this rule contains points on the boundary of the reference
52  // element, and therefore may be unsuitable for some problems. The alternative
53  // would be a 2x2x2 Gauss product rule.
54  const Real data[1][4] =
55  {
56  {1.0L, 0.0L, 0.0L, Real(4)/3}
57  };
58 
59  const unsigned int rule_id[1] = {
60  1 // (x,0,0) -> 6 permutations
61  };
62 
63  _points.resize(6);
64  _weights.resize(6);
65 
66  kim_rule(data, rule_id, 1);
67  return;
68  } // end case SECOND,THIRD
69 
70  case FOURTH:
71  case FIFTH:
72  {
73  // A degree 5, 13-point rule by Stroud,
74  // AH Stroud, "Some Fifth Degree Integration Formulas for Symmetric Regions II.",
75  // Numerische Mathematik 9, pp. 460-468 (1967).
76  //
77  // This rule is provably minimal in the number of points. The equations given for
78  // the n-cube on pg. 466 of the paper for mu/gamma and gamma are wrong, at least for
79  // the n=3 case. The analytical values given here were computed by me [JWP] in Maple.
80 
81  // Convenient intermediate values.
82  const Real sqrt19 = Real(std::sqrt(19.L));
83  const Real tp = Real(std::sqrt(71440 + 6802*sqrt19));
84 
85  // Point data for permutations.
86  const Real eta = 0;
87 
88  const Real lambda = Real(std::sqrt(1919.L/3285.L - 148.L*sqrt19/3285.L + 4.L*tp/3285.L));
89  // 8.8030440669930978047737818209860e-01L;
90 
91  const Real xi = Real(-std::sqrt(1121.L/3285.L + 74.L*sqrt19/3285.L - 2.L*tp/3285.L));
92  // -4.9584817142571115281421242364290e-01L;
93 
94  const Real mu = Real(std::sqrt(1121.L/3285.L + 74.L*sqrt19/3285.L + 2.L*tp/3285.L));
95  // 7.9562142216409541542982482567580e-01L;
96 
97  const Real gamma = Real(std::sqrt(1919.L/3285.L - 148.L*sqrt19/3285.L - 4.L*tp/3285.L));
98  // 2.5293711744842581347389255929324e-02L;
99 
100  // Weights: the centroid weight is given analytically. Weight B (resp C) goes
101  // with the {lambda,xi} (resp {gamma,mu}) permutation. The single-precision
102  // results reported by Stroud are given for reference.
103 
104  const Real A = Real(32)/19;
105  // Stroud: 0.21052632 * 8.0 = 1.684210560;
106 
107  const Real B = Real(1) / (Real(260072)/133225 - 1520*sqrt19/133225 + (133 - 37*sqrt19)*tp/133225);
108  // 5.4498735127757671684690782180890e-01L; // Stroud: 0.068123420 * 8.0 = 0.544987360;
109 
110  const Real C = Real(1) / (Real(260072)/133225 - 1520*sqrt19/133225 - (133 - 37*sqrt19)*tp/133225);
111  // 5.0764422766979170420572375713840e-01L; // Stroud: 0.063455527 * 8.0 = 0.507644216;
112 
113  _points.resize(13);
114  _weights.resize(13);
115 
116  unsigned int c=0;
117 
118  // Point with weight A (origin)
119  _points[c] = Point(eta, eta, eta);
120  _weights[c++] = A;
121 
122  // Points with weight B
123  _points[c] = Point(lambda, xi, xi);
124  _weights[c++] = B;
125  _points[c] = -_points[c-1];
126  _weights[c++] = B;
127 
128  _points[c] = Point(xi, lambda, xi);
129  _weights[c++] = B;
130  _points[c] = -_points[c-1];
131  _weights[c++] = B;
132 
133  _points[c] = Point(xi, xi, lambda);
134  _weights[c++] = B;
135  _points[c] = -_points[c-1];
136  _weights[c++] = B;
137 
138  // Points with weight C
139  _points[c] = Point(mu, mu, gamma);
140  _weights[c++] = C;
141  _points[c] = -_points[c-1];
142  _weights[c++] = C;
143 
144  _points[c] = Point(mu, gamma, mu);
145  _weights[c++] = C;
146  _points[c] = -_points[c-1];
147  _weights[c++] = C;
148 
149  _points[c] = Point(gamma, mu, mu);
150  _weights[c++] = C;
151  _points[c] = -_points[c-1];
152  _weights[c++] = C;
153 
154  return;
155 
156 
157  // // A degree 5, 14-point, "rotationally-symmetric" rule by
158  // // Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
159  // // Was also reported in Stroud's 1971 book.
160  // const Real data[2][4] =
161  // {
162  // {7.95822425754221463264548820476135e-01L, 0.00000000000000000000000000000000e+00L, 0.00000000000000000000000000000000e+00L, 8.86426592797783933518005540166204e-01L},
163  // {7.58786910639328146269034278112267e-01L, 7.58786910639328146269034278112267e-01L, 7.58786910639328146269034278112267e-01L, 3.35180055401662049861495844875346e-01L}
164  // };
165 
166  // const unsigned int rule_id[2] = {
167  // 1, // (x,0,0) -> 6 permutations
168  // 4 // (x,x,x) -> 8 permutations
169  // };
170 
171  // _points.resize(14);
172  // _weights.resize(14);
173 
174  // kim_rule(data, rule_id, 2);
175  // return;
176  } // end case FOURTH,FIFTH
177 
178 
179  case SIXTH:
180  case SEVENTH:
181  {
183  {
184  // A degree 7, 31-point, "rotationally-symmetric" rule by
185  // Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
186  // This rule contains a negative weight, so only use it if such type of
187  // rules are allowed.
188  const Real data[3][4] =
189  {
190  {Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(-1.27536231884057971014492753623188e+00L)},
191  {Real(5.85540043769119907612630781744060e-01L), Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(8.71111111111111111111111111111111e-01L)},
192  {Real(6.94470135991704766602025803883310e-01L), Real(9.37161638568208038511047377665396e-01L), Real(4.15659267604065126239606672567031e-01L), Real(1.68695652173913043478260869565217e-01L)}
193  };
194 
195  const unsigned int rule_id[3] = {
196  0, // (0,0,0) -> 1 permutation
197  1, // (x,0,0) -> 6 permutations
198  6 // (x,y,z) -> 24 permutations
199  };
200 
201  _points.resize(31);
202  _weights.resize(31);
203 
204  kim_rule(data, rule_id, 3);
205  return;
206  } // end if (allow_rules_with_negative_weights)
207 
208 
209  // A degree 7, 34-point, "fully-symmetric" rule, first published in
210  // P.C. Hammer and A.H. Stroud, "Numerical Evaluation of Multiple Integrals II",
211  // Mathematical Tables and Other Aids to Computation, vol 12., no 64, 1958, pp. 272-280
212  //
213  // This rule happens to fall under the same general
214  // construction as the Kim rules, so we've re-used
215  // that code here. Stroud gives 16 digits for his rule,
216  // and this is the most accurate version I've found.
217  //
218  // For comparison, a SEVENTH-order Gauss product rule
219  // (which integrates tri-7th order polynomials) would
220  // have 4^3=64 points.
221  const Real
222  r = Real(std::sqrt(6.L/7.L)),
223  s = Real(std::sqrt((960.L - 3.L*std::sqrt(28798.L)) / 2726.L)),
224  t = Real(std::sqrt((960.L + 3.L*std::sqrt(28798.L)) / 2726.L)),
225  B1 = Real(8624)/29160,
226  B2 = Real(2744)/29160,
227  B3 = 8*(774*t*t - 230)/(9720*(t*t-s*s)),
228  B4 = 8*(230 - 774*s*s)/(9720*(t*t-s*s));
229 
230  const Real data[4][4] =
231  {
232  {r, 0, 0, B1},
233  {r, r, 0, B2},
234  {s, s, s, B3},
235  {t, t, t, B4}
236  };
237 
238  const unsigned int rule_id[4] = {
239  1, // (x,0,0) -> 6 permutations
240  2, // (x,x,0) -> 12 permutations
241  4, // (x,x,x) -> 8 permutations
242  4 // (x,x,x) -> 8 permutations
243  };
244 
245  _points.resize(34);
246  _weights.resize(34);
247 
248  kim_rule(data, rule_id, 4);
249  return;
250 
251 
252  // // A degree 7, 38-point, "rotationally-symmetric" rule by
253  // // Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
254  // //
255  // // This rule is obviously inferior to the 34-point rule above...
256  // const Real data[3][4] =
257  //{
258  // {9.01687807821291289082811566285950e-01L, 0.00000000000000000000000000000000e+00L, 0.00000000000000000000000000000000e+00L, 2.95189738262622903181631100062774e-01L},
259  // {4.08372221499474674069588900002128e-01L, 4.08372221499474674069588900002128e-01L, 4.08372221499474674069588900002128e-01L, 4.04055417266200582425904380777126e-01L},
260  // {8.59523090201054193116477875786220e-01L, 8.59523090201054193116477875786220e-01L, 4.14735913727987720499709244748633e-01L, 1.24850759678944080062624098058597e-01L}
261  //};
262  //
263  // const unsigned int rule_id[3] = {
264  //1, // (x,0,0) -> 6 permutations
265  //4, // (x,x,x) -> 8 permutations
266  //5 // (x,x,z) -> 24 permutations
267  // };
268  //
269  // _points.resize(38);
270  // _weights.resize(38);
271  //
272  // kim_rule(data, rule_id, 3);
273  // return;
274  } // end case SIXTH,SEVENTH
275 
276  case EIGHTH:
277  {
278  // A degree 8, 47-point, "rotationally-symmetric" rule by
279  // Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
280  //
281  // A EIGHTH-order Gauss product rule (which integrates tri-8th order polynomials)
282  // would have 5^3=125 points.
283  const Real data[5][4] =
284  {
285  {Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(4.51903714875199690490763818699555e-01L)},
286  {Real(7.82460796435951590652813975429717e-01L), Real(0.00000000000000000000000000000000e+00L), Real(0.00000000000000000000000000000000e+00L), Real(2.99379177352338919703385618576171e-01L)},
287  {Real(4.88094669706366480526729301468686e-01L), Real(4.88094669706366480526729301468686e-01L), Real(4.88094669706366480526729301468686e-01L), Real(3.00876159371240019939698689791164e-01L)},
288  {Real(8.62218927661481188856422891110042e-01L), Real(8.62218927661481188856422891110042e-01L), Real(8.62218927661481188856422891110042e-01L), Real(4.94843255877038125738173175714853e-02L)},
289  {Real(2.81113909408341856058098281846420e-01L), Real(9.44196578292008195318687494773744e-01L), Real(6.97574833707236996779391729948984e-01L), Real(1.22872389222467338799199767122592e-01L)}
290  };
291 
292  const unsigned int rule_id[5] = {
293  0, // (0,0,0) -> 1 permutation
294  1, // (x,0,0) -> 6 permutations
295  4, // (x,x,x) -> 8 permutations
296  4, // (x,x,x) -> 8 permutations
297  6 // (x,y,z) -> 24 permutations
298  };
299 
300  _points.resize(47);
301  _weights.resize(47);
302 
303  kim_rule(data, rule_id, 5);
304  return;
305  } // end case EIGHTH
306 
307 
308  // By default: construct and use a Gauss quadrature rule
309  default:
310  {
311  // Break out and fall down into the default: case for the
312  // outer switch statement.
313  break;
314  }
315 
316  } // end switch(_order + 2*p)
317  } // end case HEX8/20/27
318 
319  libmesh_fallthrough();
320 
321  // By default: construct and use a Gauss quadrature rule
322  default:
323  {
324  QGauss gauss_rule(3, _order);
325  gauss_rule.init(_type, _p_level);
326 
327  // Swap points and weights with the about-to-be destroyed rule.
328  _points.swap (gauss_rule.get_points() );
329  _weights.swap(gauss_rule.get_weights());
330 
331  return;
332  }
333  } // end switch (_type)
334 }
335 
336 } // namespace libMesh
libMesh::HEX20
Definition: enum_elem_type.h:48
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::HEX8
Definition: enum_elem_type.h:47
libMesh::SIXTH
Definition: enum_order.h:47
libMesh::QMonomial::kim_rule
void kim_rule(const Real rule_data[][4], const unsigned int *rule_id, const unsigned int n_pts)
Rules from Kim and Song, Comm.
Definition: quadrature_monomial.C:205
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
B
Definition: assembly.h:38
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::QMonomial::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_monomial_3D.C:28
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::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::EIGHTH
Definition: enum_order.h:49
A
static PetscErrorCode Mat * A
Definition: petscdmlibmeshimpl.C:1026
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::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
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
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::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33