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