quadrature_trap_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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
22 
23 namespace libMesh
24 {
25 
26 
27 
28 
29 
30 void QTrap::init_3D(const ElemType type_in,
31  unsigned int)
32 {
33 #if LIBMESH_DIM == 3
34 
35  //-----------------------------------------------------------------------
36  // 3D quadrature rules
37  switch (type_in)
38  {
39  //---------------------------------------------
40  // Hex quadrature rules
41  case HEX8:
42  case HEX20:
43  case HEX27:
44  {
45  // We compute the 3D quadrature rule as a tensor
46  // product of the 1D quadrature rule.
47  QTrap q1D(1);
48  q1D.init(EDGE2);
49 
50  tensor_product_hex( q1D );
51 
52  return;
53  }
54 
55 
56 
57  //---------------------------------------------
58  // Tetrahedral quadrature rules
59  case TET4:
60  case TET10:
61  {
62  _points.resize(4);
63  _weights.resize(4);
64 
65  _points[0](0) = 0.;
66  _points[0](1) = 0.;
67  _points[0](2) = 0.;
68 
69  _points[1](0) = 1.;
70  _points[1](1) = 0.;
71  _points[1](2) = 0.;
72 
73  _points[2](0) = 0.;
74  _points[2](1) = 1.;
75  _points[2](2) = 0.;
76 
77  _points[3](0) = 0.;
78  _points[3](1) = 0.;
79  _points[3](2) = 1.;
80 
81 
82 
83  _weights[0] = .0416666666666666666666666666666666666666666667;
84  _weights[1] = _weights[0];
85  _weights[2] = _weights[0];
86  _weights[3] = _weights[0];
87 
88  return;
89  }
90 
91 
92 
93  //---------------------------------------------
94  // Prism quadrature rules
95  case PRISM6:
96  case PRISM15:
97  case PRISM18:
98  {
99  // We compute the 3D quadrature rule as a tensor
100  // product of the 1D quadrature rule and a 2D
101  // triangle quadrature rule
102 
103  QTrap q1D(1);
104  QTrap q2D(2);
105 
106  // Initialize
107  q1D.init(EDGE2);
108  q2D.init(TRI3);
109 
110  tensor_product_prism(q1D, q2D);
111 
112  return;
113  }
114 
115 
116  //---------------------------------------------
117  // Unsupported type
118  default:
119  libmesh_error_msg("ERROR: Unsupported type: " << type_in);
120  }
121 #endif
122 }
123 
124 } // namespace libMesh
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:28
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
Definition: quadrature.C:181
virtual void init_3D(const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
void tensor_product_hex(const QBase &q1D)
Definition: quadrature.C:154
Implements trapezoidal rule, i.e. nodal quadrature for linear elements.