quadrature.h
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 #ifndef LIBMESH_QUADRATURE_H
21 #define LIBMESH_QUADRATURE_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
26 #include "libmesh/point.h"
27 #include "libmesh/enum_elem_type.h" // INVALID_ELEM
28 #include "libmesh/enum_order.h" // INVALID_ORDER
29 #include "libmesh/auto_ptr.h" // deprecated
30 
31 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
32 namespace libMesh
33 {
34 enum QuadratureType : int;
35 }
36 #else
38 #endif
39 
40 // C++ includes
41 #include <vector>
42 #include <string>
43 #include <utility>
44 #include <memory>
45 
46 namespace libMesh
47 {
48 
49 // forward declarations
50 class Elem;
51 
62 class QBase : public ReferenceCountedObject<QBase>
63 {
64 protected:
65 
70  QBase (unsigned int dim,
71  Order order=INVALID_ORDER);
72 
73 public:
74 
79  QBase (const QBase &) = default;
80  QBase (QBase &&) = default;
81  QBase & operator= (const QBase &) = default;
82  QBase & operator= (QBase &&) = default;
83  virtual ~QBase() = default;
84 
88  virtual QuadratureType type() const = 0;
89 
99  static std::unique_ptr<QBase> build (const std::string & name,
100  const unsigned int dim,
101  const Order order=INVALID_ORDER);
102 
110  static std::unique_ptr<QBase> build (const QuadratureType qt,
111  const unsigned int dim,
112  const Order order=INVALID_ORDER);
113 
117  ElemType get_elem_type() const { return _type; }
118 
122  unsigned int get_p_level() const { return _p_level; }
123 
127  unsigned int n_points() const
128  {
129  libmesh_assert (!_points.empty());
130  return cast_int<unsigned int>(_points.size());
131  }
132 
136  unsigned int get_dim() const { return _dim; }
137 
142  const std::vector<Point> & get_points() const { return _points; }
143 
148  std::vector<Point> & get_points() { return _points; }
149 
154  const std::vector<Real> & get_weights() const { return _weights; }
155 
160  std::vector<Real> & get_weights() { return _weights; }
161 
165  Point qp(const unsigned int i) const
166  {
167  libmesh_assert_less (i, _points.size());
168  return _points[i];
169  }
170 
174  Real w(const unsigned int i) const
175  {
176  libmesh_assert_less (i, _weights.size());
177  return _weights[i];
178  }
179 
184  virtual void init (const ElemType type=INVALID_ELEM,
185  unsigned int p_level=0);
186 
196  virtual void init (const Elem & elem,
197  const std::vector<Real> & vertex_distance_func,
198  unsigned int p_level=0);
199 
203  Order get_order() const { return static_cast<Order>(_order + _p_level); }
204 
209  void print_info(std::ostream & os=libMesh::out) const;
210 
216  void scale(std::pair<Real, Real> old_range,
217  std::pair<Real, Real> new_range);
218 
222  friend std::ostream & operator << (std::ostream & os, const QBase & q);
223 
231  virtual bool shapes_need_reinit() { return false; }
232 
247 
248 protected:
249 
250 
256  virtual void init_0D (const ElemType type=INVALID_ELEM,
257  unsigned int p_level=0);
258 
266  virtual void init_1D (const ElemType type=INVALID_ELEM,
267  unsigned int p_level=0) = 0;
268 
277  virtual void init_2D (const ElemType,
278  unsigned int = 0)
279  {
280 #ifdef DEBUG
281  libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 2D.");
282 #endif
283  }
284 
293  virtual void init_3D (const ElemType,
294  unsigned int = 0)
295  {
296 #ifdef DEBUG
297  libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 3D.");
298 #endif
299  }
300 
306  void tensor_product_quad (const QBase & q1D);
307 
313  void tensor_product_hex (const QBase & q1D);
314 
320  void tensor_product_prism (const QBase & q1D, const QBase & q2D);
321 
325  unsigned int _dim;
326 
332 
338 
343  unsigned int _p_level;
344 
349  std::vector<Point> _points;
350 
355  std::vector<Real> _weights;
356 };
357 
358 
359 
360 // ------------------------------------------------------------
361 // QBase class members
362 
363 inline
364 QBase::QBase(unsigned int d,
365  Order o) :
366  allow_rules_with_negative_weights(true),
367  _dim(d),
368  _order(o),
369  _type(INVALID_ELEM),
370  _p_level(0)
371 {
372 }
373 
374 
375 
376 
377 inline
378 void QBase::print_info(std::ostream & os) const
379 {
380  libmesh_assert(!_points.empty());
381  libmesh_assert(!_weights.empty());
382 
383  Real summed_weights=0;
384  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
385  for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
386  {
387  os << " Point " << qpoint << ":\n"
388  << " "
389  << _points[qpoint]
390  << "\n Weight:\n "
391  << " w=" << _weights[qpoint] << "\n" << std::endl;
392 
393  summed_weights += _weights[qpoint];
394  }
395  os << "Summed Weights: " << summed_weights << std::endl;
396 }
397 
398 } // namespace libMesh
399 
400 #endif // LIBMESH_QUADRATURE_H
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:28
virtual bool shapes_need_reinit()
Definition: quadrature.h:231
QBase & operator=(const QBase &)=default
std::vector< Point > & get_points()
Definition: quadrature.h:148
bool allow_rules_with_negative_weights
Definition: quadrature.h:246
Point qp(const unsigned int i) const
Definition: quadrature.h:165
ElemType get_elem_type() const
Definition: quadrature.h:117
The base class for all geometric element types.
Definition: elem.h:100
ElemType _type
Definition: quadrature.h:337
const std::vector< Real > & get_weights() const
Definition: quadrature.h:154
QBase(unsigned int dim, Order order=INVALID_ORDER)
Definition: quadrature.h:364
unsigned int _dim
Definition: quadrature.h:325
friend std::ostream & operator<<(std::ostream &os, const QBase &q)
Definition: quadrature.C:208
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355
virtual QuadratureType type() const =0
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
Definition: quadrature.C:181
unsigned int _p_level
Definition: quadrature.h:343
virtual void init_2D(const ElemType, unsigned int=0)
Definition: quadrature.h:277
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
void tensor_product_hex(const QBase &q1D)
Definition: quadrature.C:154
unsigned int get_p_level() const
Definition: quadrature.h:122
Order get_order() const
Definition: quadrature.h:203
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:82
unsigned int get_dim() const
Definition: quadrature.h:136
unsigned int n_points() const
Definition: quadrature.h:127
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
void print_info(std::ostream &os=libMesh::out) const
Definition: quadrature.h:378
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Point > & get_points() const
Definition: quadrature.h:142
virtual ~QBase()=default
void scale(std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
Definition: quadrature.C:93
void tensor_product_quad(const QBase &q1D)
Definition: quadrature.C:127
Real w(const unsigned int i) const
Definition: quadrature.h:174
OStreamProxy out(std::cout)
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual void init_3D(const ElemType, unsigned int=0)
Definition: quadrature.h:293
Base class for all quadrature families and orders.
Definition: quadrature.h:62
std::vector< Real > & get_weights()
Definition: quadrature.h:160