quadrature_grid_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 void QGrid::init_2D(const ElemType type_in,
28  unsigned int)
29 {
30 #if LIBMESH_DIM > 1
31 
32  //-----------------------------------------------------------------------
33  // 2D quadrature rules
34 
35  // We ignore p - the grid rule is just for experimentation
36 
37  switch (type_in)
38  {
39 
40 
41  //---------------------------------------------
42  // Quadrilateral quadrature rules
43  case QUAD4:
44  case QUADSHELL4:
45  case QUAD8:
46  case QUADSHELL8:
47  case QUAD9:
48  {
49  // We compute the 2D quadrature rule as a tensor
50  // product of the 1D quadrature rule.
51  QGrid q1D(1,_order);
52  q1D.init(EDGE2);
53  tensor_product_quad( q1D );
54  return;
55  }
56 
57 
58  //---------------------------------------------
59  // Triangle quadrature rules
60  case TRI3:
61  case TRISHELL3:
62  case TRI6:
63  {
64  const unsigned int np = (_order + 1)*(_order + 2)/2;
65  const Real weight = Real(0.5)/np;
66  const Real dx = Real(1)/(_order+1);
67  _points.resize(np);
68  _weights.resize(np);
69 
70  unsigned int pt = 0;
71  for (int i = 0; i != _order + 1; ++i)
72  {
73  for (int j = 0; j != _order + 1 - i; ++j)
74  {
75  _points[pt](0) = (i+0.5)*dx;
76  _points[pt](1) = (j+0.5)*dx;
77  _weights[pt] = weight;
78  pt++;
79  }
80  }
81  return;
82  }
83 
84  //---------------------------------------------
85  // Unsupported type
86  default:
87  libmesh_error_msg("Element type not supported!:" << type_in);
88  }
89 #endif
90 }
91 
92 } // 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
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:233
virtual void init_2D(const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
Implements grid-based quadrature rules suitable for non-smooth functions.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void tensor_product_quad(const QBase &q1D)
Definition: quadrature.C:127