libMesh
quadrature_gm_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_gm.h"
22 #include "libmesh/enum_to_string.h"
23 
24 namespace libMesh
25 {
26 
27 void QGrundmann_Moller::init_3D(const ElemType, unsigned int)
28 {
29  // Nearly all GM rules contain negative weights, so if you are not
30  // allowing rules with negative weights, we cannot continue!
32  libmesh_error_msg("You requested a Grundmann-Moller rule but\n" \
33  << "are not allowing rules with negative weights!\n" \
34  << "Either select a different quadrature class or\n" \
35  << "set allow_rules_with_negative_weights==true.");
36 
37  switch (_type)
38  {
39  case TET4:
40  case TET10:
41  {
42  switch(get_order())
43  {
44  // We hard-code the first few orders based on output from
45  // the mp-quadrature library:
46  //
47  // https://code.google.com/p/mp-quadrature
48  //
49  // The points are ordered in such a way that the amount of
50  // round-off error in the quadrature calculations is
51  // (hopefully) minimized. These orderings were obtained
52  // via a simple permutation optimization strategy designed
53  // to produce the smallest overall average error while
54  // integrating all polynomials of degree <= d.
55  case CONSTANT:
56  case FIRST:
57  {
58  _points.resize(1);
59  _weights.resize(1);
60 
61  _points[0](0) = Real(1)/4;
62  _points[0](1) = Real(1)/4;
63  _points[0](2) = Real(1)/4;
64 
65  _weights[0] = Real(1)/6;
66  return;
67  }
68 
69  case SECOND:
70  case THIRD:
71  {
72  _points.resize(5);
73  _weights.resize(5);
74 
75  _points[0](0) = Real(1)/6; _points[0](1) = Real(1)/6; _points[0](2) = Real(1)/2; _weights[0] = Real(3)/40;
76  _points[1](0) = Real(1)/6; _points[1](1) = Real(1)/6; _points[1](2) = Real(1)/6; _weights[1] = Real(3)/40;
77  _points[2](0) = Real(1)/4; _points[2](1) = Real(1)/4; _points[2](2) = Real(1)/4; _weights[2] = -Real(2)/15;
78  _points[3](0) = Real(1)/2; _points[3](1) = Real(1)/6; _points[3](2) = Real(1)/6; _weights[3] = Real(3)/40;
79  _points[4](0) = Real(1)/6; _points[4](1) = Real(1)/2; _points[4](2) = Real(1)/6; _weights[4] = Real(3)/40;
80  return;
81  }
82 
83  case FOURTH:
84  case FIFTH:
85  {
86  _points.resize(15);
87  _weights.resize(15);
88 
89  _points[ 0](0) = Real(1)/8; _points[ 0](1) = Real(3)/8; _points[ 0](2) = Real(1)/8; _weights[ 0] = Real(16)/315;
90  _points[ 1](0) = Real(1)/8; _points[ 1](1) = Real(1)/8; _points[ 1](2) = Real(5)/8; _weights[ 1] = Real(16)/315;
91  _points[ 2](0) = Real(3)/8; _points[ 2](1) = Real(1)/8; _points[ 2](2) = Real(1)/8; _weights[ 2] = Real(16)/315;
92  _points[ 3](0) = Real(1)/6; _points[ 3](1) = Real(1)/2; _points[ 3](2) = Real(1)/6; _weights[ 3] = -Real(27)/280;
93  _points[ 4](0) = Real(3)/8; _points[ 4](1) = Real(1)/8; _points[ 4](2) = Real(3)/8; _weights[ 4] = Real(16)/315;
94  _points[ 5](0) = Real(1)/8; _points[ 5](1) = Real(3)/8; _points[ 5](2) = Real(3)/8; _weights[ 5] = Real(16)/315;
95  _points[ 6](0) = Real(1)/6; _points[ 6](1) = Real(1)/6; _points[ 6](2) = Real(1)/6; _weights[ 6] = -Real(27)/280;
96  _points[ 7](0) = Real(1)/6; _points[ 7](1) = Real(1)/6; _points[ 7](2) = Real(1)/2; _weights[ 7] = -Real(27)/280;
97  _points[ 8](0) = Real(1)/8; _points[ 8](1) = Real(1)/8; _points[ 8](2) = Real(1)/8; _weights[ 8] = Real(16)/315;
98  _points[ 9](0) = Real(1)/4; _points[ 9](1) = Real(1)/4; _points[ 9](2) = Real(1)/4; _weights[ 9] = Real(2)/45;
99  _points[10](0) = Real(1)/8; _points[10](1) = Real(5)/8; _points[10](2) = Real(1)/8; _weights[10] = Real(16)/315;
100  _points[11](0) = Real(1)/2; _points[11](1) = Real(1)/6; _points[11](2) = Real(1)/6; _weights[11] = -Real(27)/280;
101  _points[12](0) = Real(1)/8; _points[12](1) = Real(1)/8; _points[12](2) = Real(3)/8; _weights[12] = Real(16)/315;
102  _points[13](0) = Real(5)/8; _points[13](1) = Real(1)/8; _points[13](2) = Real(1)/8; _weights[13] = Real(16)/315;
103  _points[14](0) = Real(3)/8; _points[14](1) = Real(3)/8; _points[14](2) = Real(1)/8; _weights[14] = Real(16)/315;
104  return;
105  }
106 
107  case SIXTH:
108  case SEVENTH:
109  {
110  _points.resize(35);
111  _weights.resize(35);
112 
113  _points[ 0](0) = Real(3)/10; _points[ 0](1) = Real(1)/10; _points[ 0](2) = Real(3)/10; _weights[ 0] = Real(3125)/72576;
114  _points[ 1](0) = Real(1)/6; _points[ 1](1) = Real(1)/2; _points[ 1](2) = Real(1)/6; _weights[ 1] = Real(243)/4480;
115  _points[ 2](0) = Real(1)/6; _points[ 2](1) = Real(1)/6; _points[ 2](2) = Real(1)/2; _weights[ 2] = Real(243)/4480;
116  _points[ 3](0) = Real(1)/2; _points[ 3](1) = Real(1)/10; _points[ 3](2) = Real(1)/10; _weights[ 3] = Real(3125)/72576;
117  _points[ 4](0) = Real(3)/10; _points[ 4](1) = Real(1)/10; _points[ 4](2) = Real(1)/10; _weights[ 4] = Real(3125)/72576;
118  _points[ 5](0) = Real(3)/10; _points[ 5](1) = Real(3)/10; _points[ 5](2) = Real(1)/10; _weights[ 5] = Real(3125)/72576;
119  _points[ 6](0) = Real(1)/2; _points[ 6](1) = Real(1)/6; _points[ 6](2) = Real(1)/6; _weights[ 6] = Real(243)/4480;
120  _points[ 7](0) = Real(3)/10; _points[ 7](1) = Real(1)/10; _points[ 7](2) = Real(1)/2; _weights[ 7] = Real(3125)/72576;
121  _points[ 8](0) = Real(7)/10; _points[ 8](1) = Real(1)/10; _points[ 8](2) = Real(1)/10; _weights[ 8] = Real(3125)/72576;
122  _points[ 9](0) = Real(3)/8; _points[ 9](1) = Real(1)/8; _points[ 9](2) = Real(1)/8; _weights[ 9] = -Real(256)/2835;
123  _points[10](0) = Real(3)/10; _points[10](1) = Real(3)/10; _points[10](2) = Real(3)/10; _weights[10] = Real(3125)/72576;
124  _points[11](0) = Real(1)/10; _points[11](1) = Real(1)/2; _points[11](2) = Real(3)/10; _weights[11] = Real(3125)/72576;
125  _points[12](0) = Real(1)/10; _points[12](1) = Real(1)/10; _points[12](2) = Real(7)/10; _weights[12] = Real(3125)/72576;
126  _points[13](0) = Real(1)/2; _points[13](1) = Real(1)/10; _points[13](2) = Real(3)/10; _weights[13] = Real(3125)/72576;
127  _points[14](0) = Real(1)/10; _points[14](1) = Real(7)/10; _points[14](2) = Real(1)/10; _weights[14] = Real(3125)/72576;
128  _points[15](0) = Real(1)/10; _points[15](1) = Real(1)/2; _points[15](2) = Real(1)/10; _weights[15] = Real(3125)/72576;
129  _points[16](0) = Real(1)/6; _points[16](1) = Real(1)/6; _points[16](2) = Real(1)/6; _weights[16] = Real(243)/4480;
130  _points[17](0) = Real(3)/8; _points[17](1) = Real(1)/8; _points[17](2) = Real(3)/8; _weights[17] = -Real(256)/2835;
131  _points[18](0) = Real(1)/8; _points[18](1) = Real(1)/8; _points[18](2) = Real(5)/8; _weights[18] = -Real(256)/2835;
132  _points[19](0) = Real(1)/10; _points[19](1) = Real(1)/10; _points[19](2) = Real(3)/10; _weights[19] = Real(3125)/72576;
133  _points[20](0) = Real(1)/8; _points[20](1) = Real(3)/8; _points[20](2) = Real(3)/8; _weights[20] = -Real(256)/2835;
134  _points[21](0) = Real(5)/8; _points[21](1) = Real(1)/8; _points[21](2) = Real(1)/8; _weights[21] = -Real(256)/2835;
135  _points[22](0) = Real(1)/8; _points[22](1) = Real(5)/8; _points[22](2) = Real(1)/8; _weights[22] = -Real(256)/2835;
136  _points[23](0) = Real(1)/10; _points[23](1) = Real(3)/10; _points[23](2) = Real(1)/10; _weights[23] = Real(3125)/72576;
137  _points[24](0) = Real(1)/4; _points[24](1) = Real(1)/4; _points[24](2) = Real(1)/4; _weights[24] = -Real(8)/945;
138  _points[25](0) = Real(1)/8; _points[25](1) = Real(1)/8; _points[25](2) = Real(3)/8; _weights[25] = -Real(256)/2835;
139  _points[26](0) = Real(3)/8; _points[26](1) = Real(3)/8; _points[26](2) = Real(1)/8; _weights[26] = -Real(256)/2835;
140  _points[27](0) = Real(1)/8; _points[27](1) = Real(3)/8; _points[27](2) = Real(1)/8; _weights[27] = -Real(256)/2835;
141  _points[28](0) = Real(1)/10; _points[28](1) = Real(3)/10; _points[28](2) = Real(1)/2; _weights[28] = Real(3125)/72576;
142  _points[29](0) = Real(3)/10; _points[29](1) = Real(1)/2; _points[29](2) = Real(1)/10; _weights[29] = Real(3125)/72576;
143  _points[30](0) = Real(1)/10; _points[30](1) = Real(1)/10; _points[30](2) = Real(1)/2; _weights[30] = Real(3125)/72576;
144  _points[31](0) = Real(1)/2; _points[31](1) = Real(3)/10; _points[31](2) = Real(1)/10; _weights[31] = Real(3125)/72576;
145  _points[32](0) = Real(1)/8; _points[32](1) = Real(1)/8; _points[32](2) = Real(1)/8; _weights[32] = -Real(256)/2835;
146  _points[33](0) = Real(1)/10; _points[33](1) = Real(3)/10; _points[33](2) = Real(3)/10; _weights[33] = Real(3125)/72576;
147  _points[34](0) = Real(1)/10; _points[34](1) = Real(1)/10; _points[34](2) = Real(1)/10; _weights[34] = Real(3125)/72576;
148  return;
149  }
150 
151  default:
152  {
153  // Untested above _order=23 but should work...
154  gm_rule(get_order()/2, /*dim=*/3);
155  return;
156  }
157  } // end switch (order)
158  } // end case TET4, TET10
159 
160  default:
161  libmesh_error_msg("ERROR: Unsupported element type: " << Utility::enum_to_string(_type));
162  } // end switch (_type)
163 }
164 
165 } // namespace libMesh
libMesh::SIXTH
Definition: enum_order.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::FIFTH
Definition: enum_order.h:46
libMesh::SECOND
Definition: enum_order.h:43
libMesh::QGrundmann_Moller::gm_rule
void gm_rule(unsigned int s, unsigned int dim)
This routine is called from init_2D() and init_3D().
Definition: quadrature_gm.C:52
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::QGrundmann_Moller::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_gm_3D.C:27
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::FOURTH
Definition: enum_order.h:45
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::SEVENTH
Definition: enum_order.h:48
libMesh::QBase::allow_rules_with_negative_weights
bool allow_rules_with_negative_weights
Flag (default true) controlling the use of quadrature rules with negative weights.
Definition: quadrature.h:254
libMesh::THIRD
Definition: enum_order.h:44
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QBase::get_order
Order get_order() const
Definition: quadrature.h:211
libMesh::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::FIRST
Definition: enum_order.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33