quadrature_simpson_2D.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 void QSimpson::init_2D(const ElemType type_in,
29  unsigned int)
30 {
31 #if LIBMESH_DIM > 1
32 
33  //-----------------------------------------------------------------------
34  // 2D quadrature rules
35  switch (type_in)
36  {
37 
38 
39  //---------------------------------------------
40  // Quadrilateral quadrature rules
41  case QUAD4:
42  case QUADSHELL4:
43  case QUAD8:
44  case QUADSHELL8:
45  case QUAD9:
46  {
47  // We compute the 2D quadrature rule as a tensor
48  // product of the 1D quadrature rule.
49  QSimpson q1D(1);
50  q1D.init(EDGE2);
51  tensor_product_quad( q1D );
52  return;
53  }
54 
55 
56  //---------------------------------------------
57  // Triangle quadrature rules
58  case TRI3:
59  case TRISHELL3:
60  case TRI6:
61  {
62  // I'm not sure if you would call this Simpson's
63  // rule for triangles. What it *Really* is is
64  // four trapezoidal rules combined to give a six
65  // point rule. The points lie at the nodal locations
66  // of the TRI6, so you can get diagonal element
67  // stiffness matrix entries for quadratic elements.
68  // This rule should be able to integrate a little
69  // better than linears exactly.
70 
71  _points.resize(6);
72  _weights.resize(6);
73 
74  _points[0](0) = 0.;
75  _points[0](1) = 0.;
76 
77  _points[1](0) = 1.;
78  _points[1](1) = 0.;
79 
80  _points[2](0) = 0.;
81  _points[2](1) = 1.;
82 
83  _points[3](0) = 0.5;
84  _points[3](1) = 0.;
85 
86  _points[4](0) = 0.;
87  _points[4](1) = 0.5;
88 
89  _points[5](0) = 0.5;
90  _points[5](1) = 0.5;
91 
92  _weights[0] = Real(1)/24;
93  _weights[1] = Real(1)/24;
94  _weights[2] = Real(1)/24;
95  _weights[3] = 0.125; // 1./8.
96  _weights[4] = 0.125; // 1./8.
97  _weights[5] = 0.125; // 1./8.
98 
99  return;
100  }
101 
102 
103  //---------------------------------------------
104  // Unsupported type
105  default:
106  libmesh_error_msg("Element type not supported!:" << type_in);
107  }
108 #endif
109 }
110 
111 } // namespace libMesh
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:28
Implements Simpson&#39;s rule, i.e. nodal quadrature for quadratic elements.
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355
virtual void init_2D(const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void tensor_product_quad(const QBase &q1D)
Definition: quadrature.C:127