libMesh
quadrature_clough_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_clough.h"
22 #include "libmesh/quadrature_gauss.h"
23 #include "libmesh/enum_to_string.h"
24 
25 namespace libMesh
26 {
27 
28 
30 {
31 #if LIBMESH_DIM > 1
32  QGauss gauss_rule(2, _order);
33  gauss_rule.init(TRI6, _p_level, /*simple_type_only=*/true);
34 
35  //-----------------------------------------------------------------------
36  // 2D quadrature rules
37  switch (_type)
38  {
39 
40  //---------------------------------------------
41  // Triangle quadrature rules
42  case TRI3:
43  case TRI6:
44  case TRI7:
45  {
46  std::vector<Point> & gausspoints = gauss_rule.get_points();
47  std::vector<Real> & gaussweights = gauss_rule.get_weights();
48  std::size_t numgausspts = gausspoints.size();
49  _points.resize(numgausspts*3);
50  _weights.resize(numgausspts*3);
51  for (std::size_t i = 0; i != numgausspts; ++i)
52  {
53  _points[3*i](0) = gausspoints[i](0) +
54  gausspoints[i](1) / 3.;
55  _points[3*i](1) = gausspoints[i](1) / 3.;
56  _points[3*i+1](0) = gausspoints[i](1) / 3.;
57  _points[3*i+1](1) = gausspoints[i](0) +
58  gausspoints[i](1) / 3.;
59  _points[3*i+2](0) = 1./3. +
60  gausspoints[i](0) * 2./3. -
61  gausspoints[i](1) / 3.;
62  _points[3*i+2](1) = 1./3. -
63  gausspoints[i](0) / 3. +
64  gausspoints[i](1) * 2./3.;
65  _weights[3*i] = gaussweights[i] / 3.;
66  _weights[3*i+1] = _weights[3*i];
67  _weights[3*i+2] = _weights[3*i];
68  }
69  return;
70  }
71 
72 
73  //---------------------------------------------
74  // Unsupported type
75  default:
76  libmesh_error_msg("Element type not supported: " << Utility::enum_to_string(_type));
77  }
78 #endif
79 }
80 
81 } // namespace libMesh
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
const std::vector< Real > & get_weights() const
Definition: quadrature.h:168
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
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_2D() override
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
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
const std::vector< Point > & get_points() const
Definition: quadrature.h:156
This class implements specific orders of Gauss quadrature.
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