libMesh
quadrature_grid_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_grid.h"
22 #include "libmesh/enum_to_string.h"
23 
24 namespace libMesh
25 {
26 
27 
29 {
30 #if LIBMESH_DIM == 3
31 
32  switch (_type)
33  {
34  //---------------------------------------------
35  // Hex quadrature rules
36  case HEX8:
37  case HEX20:
38  case HEX27:
39  {
40  // We compute the 3D quadrature rule as a tensor
41  // product of the 1D quadrature rule.
42  //
43  // We ignore p_level in QGrid, since the "order" here is a
44  // user-requested number of grid points, nothing to do with
45  // the order of exactly-integrated polynomial spaces
46  QGrid q1D(1,_order);
47 
48  tensor_product_hex( q1D );
49 
50  return;
51  }
52 
53 
54 
55  //---------------------------------------------
56  // Tetrahedral quadrature rules
57  case TET4:
58  case TET10:
59  case TET14:
60  {
61  const unsigned int np = (_order+1)*(_order+2)*(_order+3)/6;
62  // Master tet has 1x1 triangle base, height 1, so volume = 1/6
63  const Real weight = Real(1)/Real(6)/np;
64  const Real dx = Real(1)/(_order+1);
65  _points.resize(np);
66  _weights.resize(np);
67 
68  unsigned int pt = 0;
69  for (int i = 0; i != _order + 1; ++i)
70  {
71  for (int j = 0; j != _order + 1 - i; ++j)
72  {
73  for (int k = 0; k != _order + 1 - i - j; ++k)
74  {
75  _points[pt](0) = (i+0.5)*dx;
76  _points[pt](1) = (j+0.5)*dx;
77  _points[pt](2) = (k+0.5)*dx;
78  _weights[pt] = weight;
79  pt++;
80  }
81  }
82  }
83  return;
84  }
85 
86 
87  // Prism quadrature rules
88  case PRISM6:
89  case PRISM15:
90  case PRISM18:
91  case PRISM20:
92  case PRISM21:
93  {
94  // We compute the 3D quadrature rule as a tensor
95  // product of the 1D quadrature rule and a 2D
96  // triangle quadrature rule
97  //
98  // We ignore p_level in QGrid, since the "order" here is a
99  // user-requested number of grid points, nothing to do with
100  // the order of exactly-integrated polynomial spaces
101  QGrid q1D(1,_order);
102  QGrid q2D(2,_order);
103 
104  // Initialize the 2D rule (1D is pre-initialized)
105  q2D.init(TRI3, 0, /*simple_type_only=*/true);
106 
107  tensor_product_prism(q1D, q2D);
108 
109  return;
110  }
111 
112 
113 
114  //---------------------------------------------
115  // Pyramid
116  case PYRAMID5:
117  case PYRAMID13:
118  case PYRAMID14:
119  case PYRAMID18:
120  {
121  const unsigned int np = (_order+1)*(_order+2)*(_order+3)/6;
122  _points.resize(np);
123  _weights.resize(np);
124  // Master pyramid has 2x2 base, height 1, so volume = 4/3
125  const Real weight = Real(4)/Real(3)/np;
126  const Real dx = Real(2)/(_order+1);
127  const Real dz = Real(1)/(_order+1);
128 
129  unsigned int pt = 0;
130  for (int k = 0; k != _order + 1; ++k)
131  {
132  for (int i = 0; i != _order + 1 - k; ++i)
133  {
134  for (int j = 0; j != _order + 1 - k - i; ++j)
135  {
136  _points[pt](0) = (i+0.5)*dx-1.0 +
137  (k+0.5)*dz;
138  _points[pt](1) = (j+0.5)*dx-1.0 +
139  (k+0.5)*dz;
140  _points[pt](2) = (k+0.5)*dz;
141  _weights[pt] = weight;
142  pt++;
143  }
144  }
145  }
146  return;
147  }
148 
149 
150 
151  //---------------------------------------------
152  // Unsupported type
153  default:
154  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
155  }
156 #endif
157 }
158 
159 } // namespace libMesh
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
virtual void init_3D() override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
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
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
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:437
This class creates quadrature points on a uniform grid, with order+1 points on an edge...
std::string enum_to_string(const T e)
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
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