libMesh
quadrature_trap_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_trap.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  QTrap 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  _points.resize(4);
59  _weights.resize(4);
60 
61  _points[0](0) = 0.;
62  _points[0](1) = 0.;
63  _points[0](2) = 0.;
64 
65  _points[1](0) = 1.;
66  _points[1](1) = 0.;
67  _points[1](2) = 0.;
68 
69  _points[2](0) = 0.;
70  _points[2](1) = 1.;
71  _points[2](2) = 0.;
72 
73  _points[3](0) = 0.;
74  _points[3](1) = 0.;
75  _points[3](2) = 1.;
76 
77 
78 
79  _weights[0] = 1/Real(24);
80  _weights[1] = _weights[0];
81  _weights[2] = _weights[0];
82  _weights[3] = _weights[0];
83 
84  return;
85  }
86 
87 
88  //---------------------------------------------
89  // Pyramid quadrature rules
90  case PYRAMID5:
91  case PYRAMID13:
92  case PYRAMID14:
93  case PYRAMID18:
94  {
95  libmesh_error_msg_if(!allow_nodal_pyramid_quadrature,
96  "Nodal quadrature on Pyramid elements is not allowed by default since\n"
97  "the Jacobian of the inverse element map is not well-defined at the Pyramid apex.\n"
98  "Set the QBase::allow_nodal_pyramid_quadrature flag to true to ignore skip this check.");
99 
100  _points.resize(5);
101  _weights.resize(5);
102 
103  _points[0](0) = -1.;
104  _points[0](1) = -1.;
105  _points[0](2) = 0.;
106 
107  _points[1](0) = 1.;
108  _points[1](1) = -1.;
109  _points[1](2) = 0.;
110 
111  _points[2](0) = 1.;
112  _points[2](1) = 1.;
113  _points[2](2) = 0.;
114 
115  _points[3](0) = -1.;
116  _points[3](1) = 1.;
117  _points[3](2) = 0.;
118 
119  _points[4](0) = 0.;
120  _points[4](1) = 0.;
121  _points[4](2) = 1.;
122 
123 
124  // These are of dubious value since we can't integrate on the
125  // vertex where the mapping Jacobian is ill-defined, but if we
126  // could, this is what would give exact solutions for
127  // constants and linears on the master element.
128  _weights[0] = 1/Real(4);
129  _weights[1] = _weights[0];
130  _weights[2] = _weights[0];
131  _weights[3] = _weights[0];
132  _weights[4] = 1/Real(3);
133 
134  return;
135  }
136 
137  //---------------------------------------------
138  // Prism quadrature rules
139  case PRISM6:
140  case PRISM15:
141  case PRISM18:
142  case PRISM20:
143  case PRISM21:
144  {
145  // We compute the 3D quadrature rule as a tensor
146  // product of the 1D quadrature rule and a 2D
147  // triangle quadrature rule
148 
149  QTrap q1D(1);
150  QTrap q2D(2);
151 
152  // Initialize the 2D rule (1D is pre-initialized)
153  q2D.init(TRI3, _p_level, /*simple_type_only=*/true);
154 
155  tensor_product_prism(q1D, q2D);
156 
157  return;
158  }
159 
160 
161  //---------------------------------------------
162  // Unsupported type
163  default:
164  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
165  }
166 #endif
167 }
168 
169 } // namespace libMesh
virtual void init_3D() override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
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
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
This class implements trapezoidal quadrature.