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