libMesh
quadrature_conical.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 // libMesh includes
20 #include "libmesh/quadrature_conical.h"
21 #include "libmesh/quadrature_gauss.h"
22 #include "libmesh/quadrature_jacobi.h"
23 #include "libmesh/enum_quadrature_type.h"
24 
25 namespace libMesh
26 {
27 
28 // See also the files:
29 // quadrature_conical_2D.C
30 // quadrature_conical_3D.C
31 // for additional implementation.
32 
34 {
35  return QCONICAL;
36 }
37 
38 void QConical::init_1D(const ElemType, unsigned int)
39 {
40  QGauss gauss1D(1, get_order());
41 
42  // Swap points and weights with the about-to-be destroyed rule.
43  _points.swap(gauss1D.get_points());
44  _weights.swap(gauss1D.get_weights());
45 }
46 
47 
48 
49 // Builds and scales a Gauss rule and a Jacobi rule.
50 // Then combines them to compute points and weights
51 // of a 2D conical product rule.
53 {
54  // Be sure the underlying rule object was built with the same dimension as the
55  // rule we are about to construct.
56  libmesh_assert_equal_to (this->get_dim(), 2);
57 
58  QGauss gauss1D(1, get_order());
59  QJacobi jac1D(1, get_order(), 1, 0);
60 
61  // The Gauss rule needs to be scaled to [0,1]
62  std::pair<Real, Real> old_range(-1.0L, 1.0L);
63  std::pair<Real, Real> new_range( 0.0L, 1.0L);
64  gauss1D.scale(old_range,
65  new_range);
66 
67  // Now construct the points and weights for the conical product rule.
68 
69  // Both rules should have the same number of points.
70  libmesh_assert_equal_to (gauss1D.n_points(), jac1D.n_points());
71 
72  // Save the number of points as a convenient variable
73  const unsigned int np = gauss1D.n_points();
74 
75  // Both rules should be between x=0 and x=1
76  libmesh_assert_greater_equal (gauss1D.qp(0)(0), 0.0);
77  libmesh_assert_less_equal (gauss1D.qp(np-1)(0), 1.0);
78  libmesh_assert_greater_equal (jac1D.qp(0)(0), 0.0);
79  libmesh_assert_less_equal (jac1D.qp(np-1)(0), 1.0);
80 
81  // Resize the points and weights vectors
82  _points.resize(np * np);
83  _weights.resize(np * np);
84 
85  // Compute the conical product
86  unsigned int gp = 0;
87  for (unsigned int i=0; i<np; i++)
88  for (unsigned int j=0; j<np; j++)
89  {
90  _points[gp](0) = jac1D.qp(j)(0); //s[j];
91  _points[gp](1) = gauss1D.qp(i)(0) * (1.-jac1D.qp(j)(0)); //r[i]*(1.-s[j]);
92  _weights[gp] = gauss1D.w(i) * jac1D.w(j); //A[i]*B[j];
93  gp++;
94  }
95 }
96 
97 
98 
99 
100 // Builds and scales a Gauss rule and a Jacobi rule.
101 // Then combines them to compute points and weights
102 // of a 3D conical product rule for the Tet.
104 {
105  // Be sure the underlying rule object was built with the same dimension as the
106  // rule we are about to construct.
107  libmesh_assert_equal_to (this->get_dim(), 3);
108 
109  QGauss gauss1D(1, get_order());
110  QJacobi jacA1D(1, get_order(), /*alpha=*/1, /*beta=*/0);
111  QJacobi jacB1D(1, get_order(), /*alpha=*/2, /*beta=*/0);
112 
113  // The Gauss rule needs to be scaled to [0,1]
114  std::pair<Real, Real> old_range(-1.0L, 1.0L);
115  std::pair<Real, Real> new_range( 0.0L, 1.0L);
116  gauss1D.scale(old_range,
117  new_range);
118 
119  // Now construct the points and weights for the conical product rule.
120 
121  // All rules should have the same number of points
122  libmesh_assert_equal_to (gauss1D.n_points(), jacA1D.n_points());
123  libmesh_assert_equal_to (jacA1D.n_points(), jacB1D.n_points());
124 
125  // Save the number of points as a convenient variable
126  const unsigned int np = gauss1D.n_points();
127 
128  // All rules should be between x=0 and x=1
129  libmesh_assert_greater_equal (gauss1D.qp(0)(0), 0.0);
130  libmesh_assert_less_equal (gauss1D.qp(np-1)(0), 1.0);
131  libmesh_assert_greater_equal (jacA1D.qp(0)(0), 0.0);
132  libmesh_assert_less_equal (jacA1D.qp(np-1)(0), 1.0);
133  libmesh_assert_greater_equal (jacB1D.qp(0)(0), 0.0);
134  libmesh_assert_less_equal (jacB1D.qp(np-1)(0), 1.0);
135 
136  // Resize the points and weights vectors
137  _points.resize(np * np * np);
138  _weights.resize(np * np * np);
139 
140  // Compute the conical product
141  unsigned int gp = 0;
142  for (unsigned int i=0; i<np; i++)
143  for (unsigned int j=0; j<np; j++)
144  for (unsigned int k=0; k<np; k++)
145  {
146  _points[gp](0) = jacB1D.qp(k)(0); //t[k];
147  _points[gp](1) = jacA1D.qp(j)(0) * (1.-jacB1D.qp(k)(0)); //s[j]*(1.-t[k]);
148  _points[gp](2) = gauss1D.qp(i)(0) * (1.-jacA1D.qp(j)(0)) * (1.-jacB1D.qp(k)(0)); //r[i]*(1.-s[j])*(1.-t[k]);
149  _weights[gp] = gauss1D.w(i) * jacA1D.w(j) * jacB1D.w(k); //A[i]*B[j]*C[k];
150  gp++;
151  }
152 }
153 
154 
155 
156 
157 
158 // Builds and scales a Gauss rule and a Jacobi rule.
159 // Then combines them to compute points and weights
160 // of a 3D conical product rule for the Pyramid. The integral
161 // over the reference Pyramid can be written (in LaTeX notation) as:
162 //
163 // If := \int_0^1 dz \int_{-(1-z)}^{(1-z)} dy \int_{-(1-z)}^{(1-z)} f(x,y,z) dx (1)
164 //
165 // (Imagine a stack of infinitely thin squares which decrease in size as
166 // you approach the apex.) Under the transformation of variables:
167 //
168 // z=w
169 // y=(1-z)v
170 // x=(1-z)u,
171 //
172 // The Jacobian determinant of this transformation is |J|=(1-w)^2, and
173 // the integral itself is transformed to:
174 //
175 // If = \int_0^1 (1-w)^2 dw \int_{-1}^{1} dv \int_{-1}^{1} f((1-w)u, (1-w)v, w) du (2)
176 //
177 // The integral can now be approximated by the product of three 1D quadrature rules:
178 // A Jacobi rule with alpha==2, beta==0 in w, and Gauss rules in v and u. In this way
179 // we can obtain 3D rules to any order for which the 1D rules exist.
181 {
182  // Be sure the underlying rule object was built with the same dimension as the
183  // rule we are about to construct.
184  libmesh_assert_equal_to (this->get_dim(), 3);
185 
186  QGauss gauss1D(1, get_order());
187  QJacobi jac1D(1, get_order(), 2, 0);
188 
189  // These rules should have the same number of points
190  libmesh_assert_equal_to (gauss1D.n_points(), jac1D.n_points());
191 
192  // Save the number of points as a convenient variable
193  const unsigned int np = gauss1D.n_points();
194 
195  // Resize the points and weights vectors
196  _points.resize(np * np * np);
197  _weights.resize(np * np * np);
198 
199  // Compute the conical product
200  unsigned int q = 0;
201  for (unsigned int i=0; i<np; ++i)
202  for (unsigned int j=0; j<np; ++j)
203  for (unsigned int k=0; k<np; ++k, ++q)
204  {
205  const Real xi=gauss1D.qp(i)(0);
206  const Real yj=gauss1D.qp(j)(0);
207  const Real zk=jac1D.qp(k)(0);
208 
209  _points[q](0) = (1.-zk) * xi;
210  _points[q](1) = (1.-zk) * yj;
211  _points[q](2) = zk;
212  _weights[q] = gauss1D.w(i) * gauss1D.w(j) * jac1D.w(k);
213  }
214 }
215 
216 } // namespace libMesh
libMesh::QConical::conical_product_pyramid
void conical_product_pyramid()
Implementation of conical product rule for a Pyramid in 3D of order get_order().
Definition: quadrature_conical.C:180
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::QBase::get_dim
unsigned int get_dim() const
Definition: quadrature.h:135
libMesh::QConical::conical_product_tri
void conical_product_tri()
Implementation of conical product rule for a Tri in 2D of order get_order().
Definition: quadrature_conical.C:52
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::QuadratureType
QuadratureType
Defines an enum for currently available quadrature rules.
Definition: enum_quadrature_type.h:33
libMesh::QBase::w
Real w(const unsigned int i) const
Definition: quadrature.h:173
libMesh::QBase::get_points
const std::vector< Point > & get_points() const
Definition: quadrature.h:141
libMesh::QBase::scale
void scale(std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
Maps the points of a 1D quadrature rule defined by "old_range" to another 1D interval defined by "new...
Definition: quadrature.C:137
libMesh::QBase::get_weights
const std::vector< Real > & get_weights() const
Definition: quadrature.h:153
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::QConical::conical_product_tet
void conical_product_tet()
Implementation of conical product rule for a Tet in 3D of order get_order().
Definition: quadrature_conical.C:103
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::QConical::type
virtual QuadratureType type() const override
Definition: quadrature_conical.C:33
libMesh::QCONICAL
Definition: enum_quadrature_type.h:42
libMesh::QConical::init_1D
virtual void init_1D(const ElemType, unsigned int) override
In 1D, use a Gauss rule.
Definition: quadrature_conical.C:38
libMesh::QBase::qp
Point qp(const unsigned int i) const
Definition: quadrature.h:164
libMesh::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::QJacobi
This class implements two (for now) Jacobi-Gauss quadrature rules.
Definition: quadrature_jacobi.h:49
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33