libMesh
quadrature_simpson_2D.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 > 1
30 
31  //-----------------------------------------------------------------------
32  // 2D quadrature rules
33  switch (_type)
34  {
35 
36 
37  //---------------------------------------------
38  // Quadrilateral quadrature rules
39  case QUAD4:
40  case QUADSHELL4:
41  case QUAD8:
42  case QUADSHELL8:
43  case QUAD9:
44  case QUADSHELL9:
45  {
46  // We compute the 2D quadrature rule as a tensor
47  // product of the 1D quadrature rule.
48  QSimpson q1D(1);
49  tensor_product_quad( q1D );
50  return;
51  }
52 
53 
54  //---------------------------------------------
55  // Triangle quadrature rules
56  case TRI3:
57  case TRISHELL3:
58  case TRI3SUBDIVISION:
59  case TRI6:
60  case TRI7:
61  {
62  // The integral of the vertex Lagrange basis functions of the
63  // TRI6 is equal to 0, while the integral of the edge Lagrange
64  // basis functions is 1/6, so attempting to derive the weights
65  // of a nodal quadrature rule using this approach shows that
66  // we only require three quadrature points to get a rule which
67  // is exact for quadratics.
68  //
69  // Unfortunately, it is not possible to derive a nodal
70  // quadrature rule on the TRI6 which is exact for cubics, so
71  // we instead choose the weights such that they are all
72  // strictly positive while still integrating linears
73  // exactly. This avoids issues with this nodal quadrature rule
74  // producing a singular elemental mass matrix. We use the
75  // "optimization" approach described in quadrature_nodal_2D.C
76  // to choose the weights.
77 
78  _points.resize(6);
79  _weights.resize(6);
80 
81  _points[0](0) = 0.;
82  _points[0](1) = 0.;
83 
84  _points[1](0) = 1.;
85  _points[1](1) = 0.;
86 
87  _points[2](0) = 0.;
88  _points[2](1) = 1.;
89 
90  _points[3](0) = 0.5;
91  _points[3](1) = 0.;
92 
93  _points[4](0) = 0.5;
94  _points[4](1) = 0.5;
95 
96  _points[5](0) = 0.;
97  _points[5](1) = 0.5;
98 
99  _weights[0] = Real(1)/38; // 0.0263157894736842
100  _weights[1] = Real(1)/38;
101  _weights[2] = Real(1)/38;
102  _weights[3] = Real(8)/57; // 0.140350877192982
103  _weights[4] = Real(8)/57;
104  _weights[5] = Real(8)/57;
105 
106  return;
107  }
108 
109 
110  //---------------------------------------------
111  // Unsupported type
112  default:
113  libmesh_error_msg("Element type not supported!:" << Utility::enum_to_string(_type));
114  }
115 #endif
116 }
117 
118 } // namespace libMesh
virtual void init_2D() override
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
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
std::string enum_to_string(const T e)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void tensor_product_quad(const QBase &q1D)
Constructs a 2D rule from the tensor product of q1D with itself.
Definition: quadrature.C:256