libMesh::QBase Class Referenceabstract

Base class for all quadrature families and orders. More...

#include <quadrature.h>

Inheritance diagram for libMesh::QBase:

Public Member Functions

virtual ~QBase ()
 
virtual QuadratureType type () const =0
 
ElemType get_elem_type () const
 
unsigned int get_p_level () const
 
unsigned int n_points () const
 
unsigned int get_dim () const
 
const std::vector< Point > & get_points () const
 
std::vector< Point > & get_points ()
 
const std::vector< Real > & get_weights () const
 
std::vector< Real > & get_weights ()
 
Point qp (const unsigned int i) const
 
Real w (const unsigned int i) const
 
virtual void init (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
 
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0)
 
Order get_order () const
 
void print_info (std::ostream &os=libMesh::out) const
 
void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
 
virtual bool shapes_need_reinit ()
 

Static Public Member Functions

static UniquePtr< QBasebuild (const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
 
static UniquePtr< QBasebuild (const QuadratureType qt, const unsigned int dim, const Order order=INVALID_ORDER)
 
static std::string get_info ()
 
static void print_info (std::ostream &out=libMesh::out)
 
static unsigned int n_objects ()
 
static void enable_print_counter_info ()
 
static void disable_print_counter_info ()
 

Public Attributes

bool allow_rules_with_negative_weights
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 

Protected Member Functions

 QBase (const unsigned int _dim, const Order _order=INVALID_ORDER)
 
virtual void init_0D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
 
virtual void init_1D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
 
virtual void init_2D (const ElemType, unsigned int=0)
 
virtual void init_3D (const ElemType, unsigned int=0)
 
void tensor_product_quad (const QBase &q1D)
 
void tensor_product_hex (const QBase &q1D)
 
void tensor_product_prism (const QBase &q1D, const QBase &q2D)
 
void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

const unsigned int _dim
 
const Order _order
 
ElemType _type
 
unsigned int _p_level
 
std::vector< Point_points
 
std::vector< Real_weights
 

Static Protected Attributes

static Counts _counts
 
static Threads::atomic< unsigned int > _n_objects
 
static Threads::spin_mutex _mutex
 
static bool _enable_print_counter = true
 

Friends

std::ostream & operator<< (std::ostream &os, const QBase &q)
 

Detailed Description

Base class for all quadrature families and orders.

The QBase class provides the basic functionality from which various quadrature rules can be derived. It computes and stores the quadrature points (in reference element space) and associated weights.

Author
Benjamin S. Kirk
Date
2002

Definition at line 53 of file quadrature.h.

Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information. The log is identified by the class name.

Definition at line 110 of file reference_counter.h.

Constructor & Destructor Documentation

libMesh::QBase::QBase ( const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
)
inlineprotected

Constructor. Protected to prevent instantiation of this base class. Use the build() method instead.

Definition at line 348 of file quadrature.h.

349  :
351  _dim(d),
352  _order(o),
354  _p_level(0)
355 {
356 }
const unsigned int _dim
Definition: quadrature.h:309
bool allow_rules_with_negative_weights
Definition: quadrature.h:230
ElemType _type
Definition: quadrature.h:321
unsigned int _p_level
Definition: quadrature.h:327
const Order _order
Definition: quadrature.h:315
virtual libMesh::QBase::~QBase ( )
inlinevirtual

Destructor.

Definition at line 69 of file quadrature.h.

References build(), libMesh::INVALID_ORDER, libMesh::Quality::name(), and type().

69 {}

Member Function Documentation

UniquePtr< QBase > libMesh::QBase::build ( const std::string &  name,
const unsigned int  dim,
const Order  order = INVALID_ORDER 
)
static

Builds a specific quadrature rule based on the name string. This enables selection of the quadrature rule at run-time. The input parameter name must be mappable through the Utility::string_to_enum<>() function.

This function allocates memory, therefore a UniquePtr<QBase> is returned so that the user does not accidentally leak it.

Definition at line 42 of file quadrature_build.C.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), and ~QBase().

45 {
46  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
47  _dim,
48  _order);
49 }
const unsigned int _dim
Definition: quadrature.h:309
virtual QuadratureType type() const =0
static UniquePtr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
const Order _order
Definition: quadrature.h:315
UniquePtr< QBase > libMesh::QBase::build ( const QuadratureType  qt,
const unsigned int  dim,
const Order  order = INVALID_ORDER 
)
static

Builds a specific quadrature rule based on the QuadratureType. This enables selection of the quadrature rule at run-time.

This function allocates memory, therefore a UniquePtr<QBase> is returned so that the user does not accidentally leak it.

Definition at line 53 of file quadrature_build.C.

References libMesh::FIRST, libMesh::FORTYTHIRD, libMesh::out, libMesh::QCLOUGH, libMesh::QCONICAL, libMesh::QGAUSS, libMesh::QGAUSS_LOBATTO, libMesh::QGRID, libMesh::QGRUNDMANN_MOLLER, libMesh::QJACOBI_1_0, libMesh::QJACOBI_2_0, libMesh::QMONOMIAL, libMesh::QSIMPSON, libMesh::QTRAP, libMesh::THIRD, and libMesh::TWENTYTHIRD.

56 {
57  switch (_qt)
58  {
59 
60  case QCLOUGH:
61  {
62 #ifdef DEBUG
63  if (_order > TWENTYTHIRD)
64  {
65  libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
66  << " up to TWENTYTHIRD order." << std::endl;
67  }
68 #endif
69 
70  return UniquePtr<QBase>(new QClough(_dim, _order));
71  }
72 
73  case QGAUSS:
74  {
75 
76 #ifdef DEBUG
77  if (_order > FORTYTHIRD)
78  {
79  libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
80  << " up to FORTYTHIRD order." << std::endl;
81  }
82 #endif
83 
84  return UniquePtr<QBase>(new QGauss(_dim, _order));
85  }
86 
87  case QJACOBI_1_0:
88  {
89 
90 #ifdef DEBUG
91  if (_order > FORTYTHIRD)
92  {
93  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
94  << " up to FORTYTHIRD order." << std::endl;
95  }
96 
97  if (_dim > 1)
98  {
99  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
100  << " in 1D only." << std::endl;
101  }
102 #endif
103 
104  return UniquePtr<QBase>(new QJacobi(_dim, _order, 1, 0));
105  }
106 
107  case QJACOBI_2_0:
108  {
109 
110 #ifdef DEBUG
111  if (_order > FORTYTHIRD)
112  {
113  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
114  << " up to FORTYTHIRD order." << std::endl;
115  }
116 
117  if (_dim > 1)
118  {
119  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
120  << " in 1D only." << std::endl;
121  }
122 #endif
123 
124  return UniquePtr<QBase>(new QJacobi(_dim, _order, 2, 0));
125  }
126 
127  case QSIMPSON:
128  {
129 
130 #ifdef DEBUG
131  if (_order > THIRD)
132  {
133  libMesh::out << "WARNING: Simpson rule provides only" << std::endl
134  << " THIRD order!" << std::endl;
135  }
136 #endif
137 
138  return UniquePtr<QBase>(new QSimpson(_dim));
139  }
140 
141  case QTRAP:
142  {
143 
144 #ifdef DEBUG
145  if (_order > FIRST)
146  {
147  libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
148  << " FIRST order!" << std::endl;
149  }
150 #endif
151 
152  return UniquePtr<QBase>(new QTrap(_dim));
153  }
154 
155  case QGRID:
156  return UniquePtr<QBase>(new QGrid(_dim, _order));
157 
158  case QGRUNDMANN_MOLLER:
159  return UniquePtr<QBase>(new QGrundmann_Moller(_dim, _order));
160 
161  case QMONOMIAL:
162  return UniquePtr<QBase>(new QMonomial(_dim, _order));
163 
164  case QGAUSS_LOBATTO:
165  return UniquePtr<QBase>(new QGaussLobatto(_dim, _order));
166 
167  case QCONICAL:
168  return UniquePtr<QBase>(new QConical(_dim, _order));
169 
170  default:
171  libmesh_error_msg("ERROR: Bad qt=" << _qt);
172  }
173 
174 
175  libmesh_error_msg("We'll never get here!");
176  return UniquePtr<QBase>();
177 }
const unsigned int _dim
Definition: quadrature.h:309
OStreamProxy out(std::cout)
const Order _order
Definition: quadrature.h:315
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::ReferenceCounter::n_objects().

101 {
102  _enable_print_counter = true;
103  return;
104 }
unsigned int libMesh::QBase::get_dim ( ) const
inline
Returns
the spatial dimension of the quadrature rule.

Definition at line 122 of file quadrature.h.

References _dim.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), and libMesh::QConical::conical_product_tri().

122 { return _dim; }
const unsigned int _dim
Definition: quadrature.h:309
ElemType libMesh::QBase::get_elem_type ( ) const
inline
Returns
the element type we're currently using.

Definition at line 103 of file quadrature.h.

References _type.

103 { return _type; }
ElemType _type
Definition: quadrature.h:321
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
Order libMesh::QBase::get_order ( ) const
inline
Returns
the order of the quadrature rule.

Definition at line 189 of file quadrature.h.

References _order, _p_level, operator<<, libMesh::out, print_info(), and scale().

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().

189 { return static_cast<Order>(_order + _p_level); }
unsigned int _p_level
Definition: quadrature.h:327
const Order _order
Definition: quadrature.h:315
unsigned int libMesh::QBase::get_p_level ( ) const
inline
Returns
the p-refinement level we're currently using.

Definition at line 108 of file quadrature.h.

References _p_level.

108 { return _p_level; }
unsigned int _p_level
Definition: quadrature.h:327
const std::vector<Point>& libMesh::QBase::get_points ( ) const
inline
std::vector<Point>& libMesh::QBase::get_points ( )
inline
Returns
a std::vector containing the quadrature point locations in reference element space as a writeable reference.

Definition at line 134 of file quadrature.h.

References _points.

134 { return _points; }
std::vector< Point > _points
Definition: quadrature.h:333
const std::vector<Real>& libMesh::QBase::get_weights ( ) const
inline
std::vector<Real>& libMesh::QBase::get_weights ( )
inline
Returns
a writable references to a std::vector containing the quadrature weights.

Definition at line 146 of file quadrature.h.

References _weights.

146 { return _weights; }
std::vector< Real > _weights
Definition: quadrature.h:339
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 160 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

161 {
162  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
163  std::pair<unsigned int, unsigned int> & p = _counts[name];
164 
165  p.first++;
166 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 173 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

174 {
175  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
176  std::pair<unsigned int, unsigned int> & p = _counts[name];
177 
178  p.second++;
179 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::QBase::init ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
virtual

Initializes the data structures for a quadrature rule for an element of type type.

Definition at line 28 of file quadrature.C.

References _dim, _p_level, _type, init_0D(), init_1D(), init_2D(), and init_3D().

Referenced by libMesh::FESubdivision::attach_quadrature_rule(), init(), libMesh::QClough::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QJacobi::QJacobi(), libMesh::QSimpson::QSimpson(), libMesh::QTrap::QTrap(), and w().

30 {
31  // check to see if we have already
32  // done the work for this quadrature rule
33  if (t == _type && p == _p_level)
34  return;
35  else
36  {
37  _type = t;
38  _p_level = p;
39  }
40 
41 
42 
43  switch(_dim)
44  {
45  case 0:
46  this->init_0D(_type,_p_level);
47 
48  return;
49 
50  case 1:
51  this->init_1D(_type,_p_level);
52 
53  return;
54 
55  case 2:
56  this->init_2D(_type,_p_level);
57 
58  return;
59 
60  case 3:
61  this->init_3D(_type,_p_level);
62 
63  return;
64 
65  default:
66  libmesh_error_msg("Invalid dimension _dim = " << _dim);
67  }
68 }
const unsigned int _dim
Definition: quadrature.h:309
ElemType _type
Definition: quadrature.h:321
unsigned int _p_level
Definition: quadrature.h:327
virtual void init_2D(const ElemType, unsigned int=0)
Definition: quadrature.h:261
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:82
virtual void init_3D(const ElemType, unsigned int=0)
Definition: quadrature.h:277
void libMesh::QBase::init ( const Elem elem,
const std::vector< Real > &  vertex_distance_func,
unsigned int  p_level = 0 
)
virtual

Initializes the data structures for an element potentially "cut" by a signed distance function. The array vertex_distance_func contains vertex values of the signed distance function. If the signed distance function changes sign on the vertices, then the element is considered to be cut.) This interface can be extended by derived classes in order to subdivide the element and construct a composite quadrature rule.

Definition at line 72 of file quadrature.C.

References init(), and libMesh::Elem::type().

75 {
76  // dispatch generic implementation
77  this->init(elem.type(), p_level);
78 }
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Definition: quadrature.C:28
void libMesh::QBase::init_0D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
protectedvirtual

Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate values. Generally this is just one point with weight 1.

Definition at line 82 of file quadrature.C.

References _points, and _weights.

Referenced by init().

84 {
85  _points.resize(1);
86  _weights.resize(1);
87  _points[0] = Point(0.);
88  _weights[0] = 1.0;
89 }
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
virtual void libMesh::QBase::init_1D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
protectedpure virtual

Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. It is assumed that derived quadrature rules will at least define the init_1D function, therefore it is pure virtual.

Implemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QJacobi, libMesh::QGauss, libMesh::QGrid, libMesh::QSimpson, libMesh::QTrap, libMesh::QConical, libMesh::QGaussLobatto, and libMesh::QClough.

Referenced by init().

virtual void libMesh::QBase::init_2D ( const ElemType  ,
unsigned int  = 0 
)
inlineprotectedvirtual

Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. Should not be pure virtual since a derived quadrature rule may only be defined in 1D. If not redefined, gives an error (when DEBUG is defined) when called.

Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QConical, libMesh::QGauss, libMesh::QGrid, libMesh::QSimpson, libMesh::QTrap, libMesh::QGaussLobatto, and libMesh::QClough.

Definition at line 261 of file quadrature.h.

Referenced by init().

263  {
264 #ifdef DEBUG
265  libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 2D.");
266 #endif
267  }
virtual void libMesh::QBase::init_3D ( const ElemType  ,
unsigned int  = 0 
)
inlineprotectedvirtual

Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. Should not be pure virtual since a derived quadrature rule may only be defined in 1D. If not redefined, gives an error (when DEBUG is defined) when called.

Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QConical, libMesh::QGauss, libMesh::QGrid, libMesh::QSimpson, libMesh::QTrap, libMesh::QGaussLobatto, and libMesh::QClough.

Definition at line 277 of file quadrature.h.

References tensor_product_hex(), tensor_product_prism(), and tensor_product_quad().

Referenced by init().

279  {
280 #ifdef DEBUG
281  libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 3D.");
282 #endif
283  }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
unsigned int libMesh::QBase::n_points ( ) const
inline
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 88 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

Referenced by libMesh::LibMeshInit::LibMeshInit().

89 {
91 }
static std::string get_info()
void libMesh::QBase::print_info ( std::ostream &  os = libMesh::out) const
inline

Prints information relevant to the quadrature rule, by default to libMesh::out.

Definition at line 362 of file quadrature.h.

References _points, _weights, libMesh::libmesh_assert(), n_points(), and libMesh::Real.

Referenced by get_order(), and libMesh::operator<<().

363 {
364  libmesh_assert(!_points.empty());
365  libmesh_assert(!_weights.empty());
366 
367  Real summed_weights=0;
368  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
369  for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
370  {
371  os << " Point " << qpoint << ":\n"
372  << " "
373  << _points[qpoint]
374  << "\n Weight:\n "
375  << " w=" << _weights[qpoint] << "\n" << std::endl;
376 
377  summed_weights += _weights[qpoint];
378  }
379  os << "Summed Weights: " << summed_weights << std::endl;
380 }
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
libmesh_assert(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int n_points() const
Definition: quadrature.h:113
Point libMesh::QBase::qp ( const unsigned int  i) const
inline
Returns
the $ i^{th} $ quadrature point in reference element space.

Definition at line 151 of file quadrature.h.

References _points.

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), tensor_product_hex(), tensor_product_prism(), and tensor_product_quad().

152  {
153  libmesh_assert_less (i, _points.size());
154  return _points[i];
155  }
std::vector< Point > _points
Definition: quadrature.h:333
void libMesh::QBase::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_range" and scales the weights accordingly.

Definition at line 93 of file quadrature.C.

References _dim, _points, _weights, and libMesh::Real.

Referenced by libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), and get_order().

95 {
96  // Make sure we are in 1D
97  libmesh_assert_equal_to (_dim, 1);
98 
99  Real
100  h_new = new_range.second - new_range.first,
101  h_old = old_range.second - old_range.first;
102 
103  // Make sure that we have sane ranges
104  libmesh_assert_greater (h_new, 0.);
105  libmesh_assert_greater (h_old, 0.);
106 
107  // Make sure there are some points
108  libmesh_assert_greater (_points.size(), 0);
109 
110  // Compute the scale factor
111  Real scfact = h_new/h_old;
112 
113  // We're mapping from old_range -> new_range
114  for (std::size_t i=0; i<_points.size(); i++)
115  {
116  _points[i](0) = new_range.first +
117  (_points[i](0) - old_range.first) * scfact;
118 
119  // Scale the weights
120  _weights[i] *= scfact;
121  }
122 }
const unsigned int _dim
Definition: quadrature.h:309
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool libMesh::QBase::shapes_need_reinit ( )
inlinevirtual

Returns true if the shape functions need to be recalculated. This may be required if the number of quadrature points or their position changes. Returns false by default.

Definition at line 215 of file quadrature.h.

215 { return false; }
void libMesh::QBase::tensor_product_hex ( const QBase q1D)
protected

Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D. Used in the init_3D routines for hexahedral element types.

Definition at line 154 of file quadrature.C.

References _points, _weights, n_points(), qp(), and w().

Referenced by libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::init_3D(), libMesh::QGauss::init_3D(), and init_3D().

155 {
156  const unsigned int np = q1D.n_points();
157 
158  _points.resize(np * np * np);
159 
160  _weights.resize(np * np * np);
161 
162  unsigned int q=0;
163 
164  for (unsigned int k=0; k<np; k++)
165  for (unsigned int j=0; j<np; j++)
166  for (unsigned int i=0; i<np; i++)
167  {
168  _points[q](0) = q1D.qp(i)(0);
169  _points[q](1) = q1D.qp(j)(0);
170  _points[q](2) = q1D.qp(k)(0);
171 
172  _weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);
173 
174  q++;
175  }
176 }
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
void libMesh::QBase::tensor_product_prism ( const QBase q1D,
const QBase q2D 
)
protected

Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule. Used in the init_3D routines for prismatic element types.

Definition at line 181 of file quadrature.C.

References _points, _weights, n_points(), qp(), and w().

Referenced by libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGauss::init_3D(), libMesh::QGrid::init_3D(), and init_3D().

182 {
183  const unsigned int n_points1D = q1D.n_points();
184  const unsigned int n_points2D = q2D.n_points();
185 
186  _points.resize (n_points1D * n_points2D);
187  _weights.resize (n_points1D * n_points2D);
188 
189  unsigned int q=0;
190 
191  for (unsigned int j=0; j<n_points1D; j++)
192  for (unsigned int i=0; i<n_points2D; i++)
193  {
194  _points[q](0) = q2D.qp(i)(0);
195  _points[q](1) = q2D.qp(i)(1);
196  _points[q](2) = q1D.qp(j)(0);
197 
198  _weights[q] = q2D.w(i) * q1D.w(j);
199 
200  q++;
201  }
202 
203 }
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
void libMesh::QBase::tensor_product_quad ( const QBase q1D)
protected

Computes the tensor product of two 1D rules and returns a 2D rule. Used in the init_2D routines for quadrilateral element types.

Definition at line 127 of file quadrature.C.

References _points, _weights, n_points(), qp(), and w().

Referenced by libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QGauss::init_2D(), and init_3D().

128 {
129 
130  const unsigned int np = q1D.n_points();
131 
132  _points.resize(np * np);
133 
134  _weights.resize(np * np);
135 
136  unsigned int q=0;
137 
138  for (unsigned int j=0; j<np; j++)
139  for (unsigned int i=0; i<np; i++)
140  {
141  _points[q](0) = q1D.qp(i)(0);
142  _points[q](1) = q1D.qp(j)(0);
143 
144  _weights[q] = q1D.w(i)*q1D.w(j);
145 
146  q++;
147  }
148 }
std::vector< Point > _points
Definition: quadrature.h:333
std::vector< Real > _weights
Definition: quadrature.h:339
Real libMesh::QBase::w ( const unsigned int  i) const
inline
Returns
the $ i^{th} $ quadrature weight.

Definition at line 160 of file quadrature.h.

References _weights, init(), libMesh::INVALID_ELEM, and type().

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), tensor_product_hex(), tensor_product_prism(), and tensor_product_quad().

161  {
162  libmesh_assert_less (i, _weights.size());
163  return _weights[i];
164  }
std::vector< Real > _weights
Definition: quadrature.h:339

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const QBase q 
)
friend

Same as above, but allows you to use the stream syntax.

Definition at line 208 of file quadrature.C.

Referenced by get_order().

209 {
210  q.print_info(os);
211  return os;
212 }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
const unsigned int libMesh::QBase::_dim
protected

The spatial dimension of the quadrature rule.

Definition at line 309 of file quadrature.h.

Referenced by get_dim(), init(), libMesh::QGaussLobatto::QGaussLobatto(), and scale().

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 134 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 128 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 123 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

unsigned int libMesh::QBase::_p_level
protected

The p-level of the element for which the current values have been computed.

Definition at line 327 of file quadrature.h.

Referenced by get_order(), get_p_level(), and init().

bool libMesh::QBase::allow_rules_with_negative_weights

Flag (default true) controlling the use of quadrature rules with negative weights. Set this to false to require rules with all positive weights.

Rules with negative weights can be unsuitable for some problems. For example, it is possible for a rule with negative weights to obtain a negative result when integrating a positive function.

A particular example: if rules with negative weights are not allowed, a request for TET,THIRD (5 points) will return the TET,FIFTH (14 points) rule instead, nearly tripling the computational effort required!

Definition at line 230 of file quadrature.h.

Referenced by libMesh::QGrundmann_Moller::init_2D(), libMesh::QGauss::init_3D(), libMesh::QMonomial::init_3D(), and libMesh::QGrundmann_Moller::init_3D().


The documentation for this class was generated from the following files: