libMesh
quadrature.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 // Local includes
20 #include "libmesh/elem.h"
21 #include "libmesh/quadrature.h"
22 #include "libmesh/int_range.h"
23 
24 namespace libMesh
25 {
26 
27 QBase::QBase(unsigned int d,
28  Order o) :
29  allow_rules_with_negative_weights(true),
30  _dim(d),
31  _order(o),
32  _type(INVALID_ELEM),
33  _p_level(0)
34 {}
35 
36 
37 void QBase::print_info(std::ostream & os) const
38 {
39  libmesh_assert(!_points.empty());
40  libmesh_assert(!_weights.empty());
41 
42  Real summed_weights=0;
43  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
44  for (auto qpoint: index_range(_points))
45  {
46  os << " Point " << qpoint << ":\n"
47  << " "
48  << _points[qpoint]
49  << "\n Weight:\n "
50  << " w=" << _weights[qpoint] << "\n" << std::endl;
51 
52  summed_weights += _weights[qpoint];
53  }
54  os << "Summed Weights: " << summed_weights << std::endl;
55 }
56 
57 
58 
59 void QBase::init(const ElemType t,
60  unsigned int p)
61 {
62  // check to see if we have already
63  // done the work for this quadrature rule
64  if (t == _type && p == _p_level)
65  return;
66  else
67  {
68  _type = t;
69  _p_level = p;
70  }
71 
72 
73 
74  switch(_dim)
75  {
76  case 0:
77  this->init_0D();
78 
79  return;
80 
81  case 1:
82  this->init_1D();
83 
84  return;
85 
86  case 2:
87  this->init_2D();
88 
89  return;
90 
91  case 3:
92  this->init_3D();
93 
94  return;
95 
96  default:
97  libmesh_error_msg("Invalid dimension _dim = " << _dim);
98  }
99 }
100 
101 
102 
103 void QBase::init (const Elem & elem,
104  const std::vector<Real> & /* vertex_distance_func */,
105  unsigned int p_level)
106 {
107  // dispatch generic implementation
108  this->init(elem.type(), p_level);
109 }
110 
111 
112 
113 void QBase::init_0D(const ElemType, unsigned int)
114 {
115  _points.resize(1);
116  _weights.resize(1);
117  _points[0] = Point(0.);
118  _weights[0] = 1.0;
119 }
120 
121 
122 
123 void QBase::init_2D (const ElemType, unsigned int)
124 {
125  libmesh_not_implemented();
126 }
127 
128 
129 
130 void QBase::init_3D (const ElemType, unsigned int)
131 {
132  libmesh_not_implemented();
133 }
134 
135 
136 
137 void QBase::scale(std::pair<Real, Real> old_range,
138  std::pair<Real, Real> new_range)
139 {
140  // Make sure we are in 1D
141  libmesh_assert_equal_to (_dim, 1);
142 
143  Real
144  h_new = new_range.second - new_range.first,
145  h_old = old_range.second - old_range.first;
146 
147  // Make sure that we have sane ranges
148  libmesh_assert_greater (h_new, 0.);
149  libmesh_assert_greater (h_old, 0.);
150 
151  // Make sure there are some points
152  libmesh_assert_greater (_points.size(), 0);
153 
154  // Compute the scale factor
155  Real scfact = h_new/h_old;
156 
157  // We're mapping from old_range -> new_range
158  for (auto i : index_range(_points))
159  {
160  _points[i](0) = new_range.first +
161  (_points[i](0) - old_range.first) * scfact;
162 
163  // Scale the weights
164  _weights[i] *= scfact;
165  }
166 }
167 
168 
169 
170 
172 {
173 
174  const unsigned int np = q1D.n_points();
175 
176  _points.resize(np * np);
177 
178  _weights.resize(np * np);
179 
180  unsigned int q=0;
181 
182  for (unsigned int j=0; j<np; j++)
183  for (unsigned int i=0; i<np; i++)
184  {
185  _points[q](0) = q1D.qp(i)(0);
186  _points[q](1) = q1D.qp(j)(0);
187 
188  _weights[q] = q1D.w(i)*q1D.w(j);
189 
190  q++;
191  }
192 }
193 
194 
195 
196 
197 
199 {
200  const unsigned int np = q1D.n_points();
201 
202  _points.resize(np * np * np);
203 
204  _weights.resize(np * np * np);
205 
206  unsigned int q=0;
207 
208  for (unsigned int k=0; k<np; k++)
209  for (unsigned int j=0; j<np; j++)
210  for (unsigned int i=0; i<np; i++)
211  {
212  _points[q](0) = q1D.qp(i)(0);
213  _points[q](1) = q1D.qp(j)(0);
214  _points[q](2) = q1D.qp(k)(0);
215 
216  _weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);
217 
218  q++;
219  }
220 }
221 
222 
223 
224 
225 void QBase::tensor_product_prism(const QBase & q1D, const QBase & q2D)
226 {
227  const unsigned int n_points1D = q1D.n_points();
228  const unsigned int n_points2D = q2D.n_points();
229 
230  _points.resize (n_points1D * n_points2D);
231  _weights.resize (n_points1D * n_points2D);
232 
233  unsigned int q=0;
234 
235  for (unsigned int j=0; j<n_points1D; j++)
236  for (unsigned int i=0; i<n_points2D; i++)
237  {
238  _points[q](0) = q2D.qp(i)(0);
239  _points[q](1) = q2D.qp(i)(1);
240  _points[q](2) = q1D.qp(j)(0);
241 
242  _weights[q] = q2D.w(i) * q1D.w(j);
243 
244  q++;
245  }
246 
247 }
248 
249 
250 
251 
252 std::ostream & operator << (std::ostream & os, const QBase & q)
253 {
254  q.print_info(os);
255  return os;
256 }
257 
258 } // namespace libMesh
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::QBase::_p_level
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:357
libMesh::QBase::tensor_product_hex
void tensor_product_hex(const QBase &q1D)
Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D.
Definition: quadrature.C:198
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::operator<<
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
Definition: fe_type.h:164
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::QBase::init_1D
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate val...
libMesh::QBase::w
Real w(const unsigned int i) const
Definition: quadrature.h:173
libMesh::QBase::init_3D
virtual void init_3D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:130
libMesh::QBase::init
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the data structures for a quadrature rule for an element of type type.
Definition: quadrature.C:59
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::QBase::QBase
QBase(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
Definition: quadrature.C:27
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
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::tensor_product_prism
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule.
Definition: quadrature.C:225
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::QBase::tensor_product_quad
void tensor_product_quad(const QBase &q1D)
Constructs a 2D rule from the tensor product of q1D with itself.
Definition: quadrature.C:171
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QBase::init_0D
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:113
libMesh::QBase::_dim
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:339
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::Elem::type
virtual ElemType type() const =0
libMesh::QBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information relevant to the quadrature rule, by default to libMesh::out.
Definition: quadrature.C:37
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::QBase::init_2D
virtual void init_2D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:123