libMesh
quadrature_nodal_3D.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_nodal.h"
22 #include "libmesh/quadrature_trap.h"
23 #include "libmesh/quadrature_simpson.h"
24 #include "libmesh/enum_to_string.h"
25 
26 namespace libMesh
27 {
28 
29 void QNodal::init_3D(const ElemType, unsigned int)
30 {
31 #if LIBMESH_DIM == 3
32 
33  switch (_type)
34  {
35  case TET4:
36  case PRISM6:
37  case HEX8:
38  {
39  QTrap rule(/*dim=*/3, /*ignored*/_order);
40  rule.init(_type, /*ignored*/_p_level);
41  _points.swap (rule.get_points());
42  _weights.swap(rule.get_weights());
43  return;
44  }
45 
46  case PRISM15:
47  {
48  // A rule with 15 points which is exact for linears, and
49  // naturally produces a lumped approximation to the mass
50  // matrix. The quadrature points are numbered the same way as
51  // the reference element nodes.
52  _points =
53  {
54  Point(0.,0.,-1), Point(+1,0.,-1), Point(0.,+1,-1),
55  Point(0.,0.,+1), Point(+1,0.,+1), Point(0.,+1,+1),
56  Point(.5,0.,-1), Point(.5,.5,-1), Point(0.,.5,-1),
57  Point(0.,0.,0.), Point(+1,0.,0.), Point(0.,+1,0.),
58  Point(.5,0.,+1), Point(.5,.5,+1), Point(0.,.5,+1),
59  };
60 
61  // vertex (wv), tri edge (wt), and quad edge (wq) weights are
62  // obtained using the same approach that was used for the Quad8,
63  // see quadrature_nodal_2D.C for details.
64  Real wv = Real(1) / 34;
65  Real wt = Real(4) / 51;
66  Real wq = Real(2) / 17;
67 
68  _weights = {wv, wv, wv, wv, wv, wv,
69  wt, wt, wt,
70  wq, wq, wq,
71  wt, wt, wt};
72 
73  return;
74  }
75 
76  case HEX20:
77  {
78  // A rule with 20 points which is exact for linears, and
79  // naturally produces a lumped approximation to the mass
80  // matrix. The quadrature points are numbered the same way as
81  // the reference element nodes.
82  _points =
83  {
84  Point(-1,-1,-1), Point(+1,-1,-1), Point(+1,+1,-1), Point(-1,+1,-1),
85  Point(-1,-1,+1), Point(+1,-1,+1), Point(+1,+1,+1), Point(-1,+1,+1),
86  Point(0.,-1,-1), Point(+1,0.,-1), Point(0.,+1,-1), Point(-1,0.,-1),
87  Point(-1,-1,0.), Point(+1,-1,0.), Point(+1,+1,0.), Point(-1,+1,0.),
88  Point(0.,-1,+1), Point(+1,0.,+1), Point(0.,+1,+1), Point(-1,0.,+1)
89  };
90 
91  // vertex (wv), and edge (we) weights are obtained using the
92  // same approach that was used for the Quad8, see
93  // quadrature_nodal_2D.C for details.
94  Real wv = Real(7) / 31;
95  Real we = Real(16) / 31;
96 
97  _weights = {wv, wv, wv, wv, wv, wv, wv, wv,
98  we, we, we, we, we, we, we, we, we, we, we, we};
99 
100  return;
101  }
102 
103  case TET10:
104  case PRISM18:
105  case HEX27:
106  {
107  QSimpson rule(/*dim=*/3, /*ignored*/_order);
108  rule.init(_type, /*ignored*/_p_level);
109  _points.swap (rule.get_points());
110  _weights.swap(rule.get_weights());
111  return;
112  }
113 
114  default:
115  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
116  }
117 #endif
118 }
119 
120 } // namespace libMesh
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::QSimpson
This class implements Simpson quadrature.
Definition: quadrature_simpson.h:38
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::HEX8
Definition: enum_elem_type.h:47
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::QTrap
This class implements trapezoidal quadrature.
Definition: quadrature_trap.h:38
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::TET4
Definition: enum_elem_type.h:45
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::QBase::_order
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly.
Definition: quadrature.h:345
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::QBase::get_points
const std::vector< Point > & get_points() const
Definition: quadrature.h:141
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::QNodal::init_3D
virtual void init_3D(const ElemType, unsigned int) override
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature_nodal_3D.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::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::PRISM18
Definition: enum_elem_type.h:52
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