quadrature_trap_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 QTrap::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 
48  // We compute the 2D quadrature rule as a tensor
49  // product of the 1D quadrature rule.
50  QTrap q1D(1);
51  q1D.init(EDGE2);
52 
53  tensor_product_quad( q1D );
54 
55  return;
56  }
57 
58 
59  //---------------------------------------------
60  // Triangle quadrature rules
61  case TRI3:
62  case TRISHELL3:
63  case TRI6:
64  {
65  _points.resize(3);
66  _weights.resize(3);
67 
68  _points[0](0) = 0.;
69  _points[0](1) = 0.;
70 
71  _points[1](0) = 1.;
72  _points[1](1) = 0.;
73 
74  _points[2](0) = 0.;
75  _points[2](1) = 1.;
76 
77 
78  _weights[0] = 1./6.;
79  _weights[1] = 1./6.;
80  _weights[2] = 1./6.;
81 
82  return;
83  }
84 
85 
86  //---------------------------------------------
87  // Unsupported type
88  default:
89  libmesh_error_msg("Element type not supported!:" << type_in);
90  }
91 #endif
92 }
93 
94 } // 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
virtual void init_2D(const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
void tensor_product_quad(const QBase &q1D)
Definition: quadrature.C:127
Implements trapezoidal rule, i.e. nodal quadrature for linear elements.