libMesh
quadrature_clough_2D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 
29 void QClough::init_2D(const ElemType, unsigned int)
30 {
31 #if LIBMESH_DIM > 1
32  QGauss gauss_rule(2, _order);
33  gauss_rule.init(TRI6, _p_level);
34 
35  //-----------------------------------------------------------------------
36  // 2D quadrature rules
37  switch (_type)
38  {
39 
40  //---------------------------------------------
41  // Triangle quadrature rules
42  case TRI3:
43  case TRI6:
44  {
45  std::vector<Point> & gausspoints = gauss_rule.get_points();
46  std::vector<Real> & gaussweights = gauss_rule.get_weights();
47  std::size_t numgausspts = gausspoints.size();
48  _points.resize(numgausspts*3);
49  _weights.resize(numgausspts*3);
50  for (std::size_t i = 0; i != numgausspts; ++i)
51  {
52  _points[3*i](0) = gausspoints[i](0) +
53  gausspoints[i](1) / 3.;
54  _points[3*i](1) = gausspoints[i](1) / 3.;
55  _points[3*i+1](0) = gausspoints[i](1) / 3.;
56  _points[3*i+1](1) = gausspoints[i](0) +
57  gausspoints[i](1) / 3.;
58  _points[3*i+2](0) = 1./3. +
59  gausspoints[i](0) * 2./3. -
60  gausspoints[i](1) / 3.;
61  _points[3*i+2](1) = 1./3. -
62  gausspoints[i](0) / 3. +
63  gausspoints[i](1) * 2./3.;
64  _weights[3*i] = gaussweights[i] / 3.;
65  _weights[3*i+1] = _weights[3*i];
66  _weights[3*i+2] = _weights[3*i];
67  }
68  return;
69  }
70 
71 
72  //---------------------------------------------
73  // Unsupported type
74  default:
75  libmesh_error_msg("Element type not supported: " << Utility::enum_to_string(_type));
76  }
77 #endif
78 }
79 
80 } // namespace libMesh
libMesh::QBase::_p_level
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:357
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::QBase::init
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the data structures for a quadrature rule for an element of type type.
Definition: quadrature.C:59
libMesh::QBase::_order
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly.
Definition: quadrature.h:345
libMesh::QBase::get_points
const std::vector< Point > & get_points() const
Definition: quadrature.h:141
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::QClough::init_2D
virtual void init_2D(const ElemType, unsigned int) override
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature_clough_2D.C:29
libMesh::QBase::get_weights
const std::vector< Real > & get_weights() const
Definition: quadrature.h:153
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33