libMesh
quadrature_simpson_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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_simpson.h"
22 #include "libmesh/enum_to_string.h"
23 
24 namespace libMesh
25 {
26 
28 {
29 #if LIBMESH_DIM == 3
30 
31  //-----------------------------------------------------------------------
32  // 3D quadrature rules
33  switch (_type)
34  {
35  //---------------------------------------------
36  // Hex quadrature rules
37  case HEX8:
38  case HEX20:
39  case HEX27:
40  {
41  // We compute the 3D quadrature rule as a tensor
42  // product of the 1D quadrature rule.
43  QSimpson q1D(1);
44 
45  tensor_product_hex( q1D );
46 
47  return;
48  }
49 
50 
51 
52  //---------------------------------------------
53  // Tetrahedral quadrature rules
54  case TET4:
55  case TET10:
56  case TET14:
57  {
58  // This rule is created by combining 8 subtets
59  // which use the trapezoidal rule. The weights
60  // may seem a bit odd, but they are correct,
61  // and should add up to 1/6, the volume of the
62  // reference tet. The points of this rule are
63  // at the nodal points of the TET10, allowing
64  // you to generate diagonal element stiffness
65  // matrices when using quadratic elements.
66  // It should be able to integrate something
67  // better than linears, but I'm not sure how
68  // high.
69 
70  _points.resize(10);
71  _weights.resize(10);
72 
73  _points[0](0) = 0.; _points[5](0) = .5;
74  _points[0](1) = 0.; _points[5](1) = .5;
75  _points[0](2) = 0.; _points[5](2) = 0.;
76 
77  _points[1](0) = 1.; _points[6](0) = 0.;
78  _points[1](1) = 0.; _points[6](1) = .5;
79  _points[1](2) = 0.; _points[6](2) = 0.;
80 
81  _points[2](0) = 0.; _points[7](0) = 0.;
82  _points[2](1) = 1.; _points[7](1) = 0.;
83  _points[2](2) = 0.; _points[7](2) = .5;
84 
85  _points[3](0) = 0.; _points[8](0) = .5;
86  _points[3](1) = 0.; _points[8](1) = 0.;
87  _points[3](2) = 1.; _points[8](2) = .5;
88 
89  _points[4](0) = .5; _points[9](0) = 0.;
90  _points[4](1) = 0.; _points[9](1) = .5;
91  _points[4](2) = 0.; _points[9](2) = .5;
92 
93 
94  _weights[0] = Real(1)/192;
95  _weights[1] = _weights[0];
96  _weights[2] = _weights[0];
97  _weights[3] = _weights[0];
98 
99  _weights[4] = Real(14)/576;
100  _weights[5] = _weights[4];
101  _weights[6] = _weights[4];
102  _weights[7] = _weights[4];
103  _weights[8] = _weights[4];
104  _weights[9] = _weights[4];
105 
106  return;
107  }
108 
109 
110 
111  //---------------------------------------------
112  // Prism quadrature rules
113  case PRISM6:
114  case PRISM15:
115  case PRISM18:
116  case PRISM20:
117  case PRISM21:
118  {
119  // We compute the 3D quadrature rule as a tensor
120  // product of the 1D quadrature rule and a 2D
121  // triangle quadrature rule
122 
123  QSimpson q1D(1);
124  QSimpson q2D(2);
125 
126  // Initialize the 2D rule (1D is pre-initialized)
127  q2D.init(TRI3, _p_level, /*simple_type_only=*/true);
128 
129  tensor_product_prism(q1D, q2D);
130 
131  return;
132  }
133 
134 
135  //---------------------------------------------
136  // Pyramid quadrature rules
137  case PYRAMID5:
138  case PYRAMID13:
139  case PYRAMID14:
140  case PYRAMID18:
141  {
142  libmesh_error_msg_if(!allow_nodal_pyramid_quadrature,
143  "Nodal quadrature on Pyramid elements is not allowed by default since\n"
144  "the Jacobian of the inverse element map is not well-defined at the Pyramid apex.\n"
145  "Set the QBase::allow_nodal_pyramid_quadrature flag to true to ignore skip this check.");
146 
147  _points.resize(14);
148  _weights.resize(14);
149 
150  _points[0](0) = -1.;
151  _points[0](1) = -1.;
152  _points[0](2) = 0.;
153 
154  _points[1](0) = 1.;
155  _points[1](1) = -1.;
156  _points[1](2) = 0.;
157 
158  _points[2](0) = 1.;
159  _points[2](1) = 1.;
160  _points[2](2) = 0.;
161 
162  _points[3](0) = -1.;
163  _points[3](1) = 1.;
164  _points[3](2) = 0.;
165 
166  _points[4](0) = 0.;
167  _points[4](1) = 0.;
168  _points[4](2) = 1.;
169 
170  _points[5](0) = 0.;
171  _points[5](1) = -1.;
172  _points[5](2) = 0.;
173 
174  _points[6](0) = 1.;
175  _points[6](1) = 0.;
176  _points[6](2) = 0.;
177 
178  _points[7](0) = 0.;
179  _points[7](1) = 1.;
180  _points[7](2) = 0.;
181 
182  _points[8](0) = -1.;
183  _points[8](1) = 0.;
184  _points[8](2) = 0.;
185 
186  _points[9](0) = 0.;
187  _points[9](1) = -0.5;
188  _points[9](2) = 0.5;
189 
190  _points[10](0) = 0.5;
191  _points[10](1) = 0.;
192  _points[10](2) = 0.5;
193 
194  _points[11](0) = 0.;
195  _points[11](1) = 0.5;
196  _points[11](2) = 0.5;
197 
198  _points[12](0) = -0.5;
199  _points[12](1) = 0.;
200  _points[12](2) = 0.5;
201 
202  _points[13](0) = 0.;
203  _points[13](1) = 0.;
204  _points[13](2) = 0.;
205 
206  // These are of dubious value since we can't integrate on the
207  // vertex where the mapping Jacobian is ill-defined, and even
208  // if we could it looks like we'd need negative weight at that
209  // vertex to give us exact integrals of both z and z^2. So I
210  // punt and just use QTrap weights.
211  _weights[0] = 1/Real(4);
212  _weights[1] = _weights[0];
213  _weights[2] = _weights[0];
214  _weights[3] = _weights[0];
215  _weights[4] = 1/Real(3);
216  _weights[5] = 0;
217  _weights[6] = 0;
218  _weights[7] = 0;
219  _weights[8] = 0;
220  _weights[9] = 0;
221  _weights[10] = 0;
222  _weights[11] = 0;
223  _weights[12] = 0;
224  _weights[13] = 0;
225 
226  return;
227  }
228 
229 
230  //---------------------------------------------
231  // Unsupported type
232  default:
233  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
234  }
235 #endif
236 }
237 
238 } // namespace libMesh
This class implements Simpson quadrature.
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
The libMesh namespace provides an interface to certain functionality in the library.
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
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:310
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403
virtual void init_3D() override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
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:283
std::string enum_to_string(const T e)
bool allow_nodal_pyramid_quadrature
The flag&#39;s value defaults to false so that one does not accidentally use a nodal quadrature rule on P...
Definition: quadrature.h:314
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e.
Definition: quadrature.C:65