libMesh
Public Member Functions | Static Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | Friends | List of all members
libMesh::QBase Class Referenceabstract

The QBase class provides the basic functionality from which various quadrature rules can be derived. More...

#include <quadrature.h>

Inheritance diagram for libMesh::QBase:
[legend]

Public Member Functions

 QBase (const QBase &)=default
 Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class. More...
 
 QBase (QBase &&)=default
 
QBaseoperator= (const QBase &)=default
 
QBaseoperator= (QBase &&)=default
 
virtual ~QBase ()=default
 
virtual std::unique_ptr< QBaseclone () const
 
virtual QuadratureType type () const =0
 
ElemType get_elem_type () const
 
unsigned int get_p_level () const
 
unsigned int n_points () const
 
unsigned int size () const
 Alias for n_points() to enable use in index_range. More...
 
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 Elem &e, unsigned int p_level=invalid_uint)
 Initializes the data structures for a quadrature rule for the element e. More...
 
virtual void init (const ElemType type=INVALID_ELEM, unsigned int p_level=0, bool simple_type_only=false)
 Initializes the data structures for a quadrature rule for an element of type type. More...
 
virtual void init (const QBase &other_rule)
 Initializes the data structures for a quadrature rule based on the element, element type, and p_level settings of other_rule. More...
 
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0)
 Initializes the data structures for an element potentially "cut" by a signed distance function. More...
 
Order get_order () const
 
Order get_base_order () const
 
void print_info (std::ostream &os=libMesh::out) const
 Prints information relevant to the quadrature rule, by default to libMesh::out. More...
 
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_range" and scales the weights accordingly. More...
 
virtual bool shapes_need_reinit ()
 

Static Public Member Functions

static std::unique_ptr< QBasebuild (std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
 Builds a specific quadrature rule based on the name string. More...
 
static std::unique_ptr< QBasebuild (const QuadratureType qt, const unsigned int dim, const Order order=INVALID_ORDER)
 Builds a specific quadrature rule based on the QuadratureType. More...
 
static std::string get_info ()
 Gets a string containing the reference information. More...
 
static void print_info (std::ostream &out_stream=libMesh::out)
 Prints the reference information, by default to libMesh::out. More...
 
static unsigned int n_objects ()
 Prints the number of outstanding (created, but not yet destroyed) objects. More...
 
static void enable_print_counter_info ()
 Methods to enable/disable the reference counter output from print_info() More...
 
static void disable_print_counter_info ()
 

Public Attributes

bool allow_rules_with_negative_weights
 Flag (default true) controlling the use of quadrature rules with negative weights. More...
 
bool allow_nodal_pyramid_quadrature
 The flag's value defaults to false so that one does not accidentally use a nodal quadrature rule on Pyramid elements, since evaluating the inverse element Jacobian (e.g. More...
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 Data structure to log the information. More...
 

Protected Member Functions

 QBase (unsigned int dim, Order order=INVALID_ORDER)
 Constructor. More...
 
virtual void init_0D ()
 Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
virtual void init_1D ()=0
 Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
virtual void init_2D ()
 Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
virtual void init_3D ()
 Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
void tensor_product_quad (const QBase &q1D)
 Constructs a 2D rule from the tensor product of q1D with itself. More...
 
void tensor_product_hex (const QBase &q1D)
 Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D. More...
 
void tensor_product_prism (const QBase &q1D, const QBase &q2D)
 Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule. More...
 
void increment_constructor_count (const std::string &name) noexcept
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name) noexcept
 Increments the destruction counter. More...
 

Protected Attributes

unsigned int _dim
 The spatial dimension of the quadrature rule. More...
 
Order _order
 The polynomial order which the quadrature rule is capable of integrating exactly. More...
 
ElemType _type
 The type of element for which the current values have been computed. More...
 
const Elem_elem
 The element for which the current values were computed, or nullptr if values were computed without a specific element. More...
 
unsigned int _p_level
 The p-level of the element for which the current values have been computed. More...
 
std::vector< Point_points
 The locations of the quadrature points in reference element space. More...
 
std::vector< Real_weights
 The quadrature weights. More...
 

Static Protected Attributes

static Counts _counts
 Actually holds the data. More...
 
static Threads::atomic< unsigned int_n_objects
 The number of objects. More...
 
static Threads::spin_mutex _mutex
 Mutual exclusion object to enable thread-safe reference counting. More...
 
static bool _enable_print_counter = true
 Flag to control whether reference count information is printed when print_info is called. More...
 

Friends

std::ostream & operator<< (std::ostream &os, const QBase &q)
 Same as above, but allows you to use the stream syntax. More...
 

Detailed Description

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 Base class for all quadrature families and orders.

Definition at line 61 of file quadrature.h.

Member Typedef Documentation

◆ Counts

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 119 of file reference_counter.h.

Constructor & Destructor Documentation

◆ QBase() [1/3]

libMesh::QBase::QBase ( unsigned int  dim,
Order  order = INVALID_ORDER 
)
protected

Constructor.

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

Definition at line 27 of file quadrature.C.

28  :
31  _dim(d),
32  _order(o),
34  _elem(nullptr),
35  _p_level(0)
36 {}
bool allow_rules_with_negative_weights
Flag (default true) controlling the use of quadrature rules with negative weights.
Definition: quadrature.h:301
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
const Elem * _elem
The element for which the current values were computed, or nullptr if values were computed without a ...
Definition: quadrature.h:397
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
bool allow_nodal_pyramid_quadrature
The flag&#39;s value defaults to false so that one does not accidentally use a nodal quadrature rule on P...
Definition: quadrature.h:314

◆ QBase() [2/3]

libMesh::QBase::QBase ( const QBase )
default

Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class.

◆ QBase() [3/3]

libMesh::QBase::QBase ( QBase &&  )
default

◆ ~QBase()

virtual libMesh::QBase::~QBase ( )
virtualdefault

Member Function Documentation

◆ build() [1/2]

std::unique_ptr< QBase > libMesh::QBase::build ( std::string_view  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 std::unique_ptr<QBase> is returned so that the user does not accidentally leak it.

Definition at line 43 of file quadrature_build.C.

References _dim, _order, and type().

Referenced by assemble_poisson(), libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), clone(), main(), libMesh::OverlapCoupling::OverlapCoupling(), libMesh::InfFE< Dim, T_radial, T_map >::reinit(), QuadratureTest::testBuild(), QuadratureTest::testJacobi(), and QuadratureTest::testPolynomials().

46 {
47  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
48  _dim,
49  _order);
50 }
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
virtual QuadratureType type() const =0
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.

◆ build() [2/2]

std::unique_ptr< 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 std::unique_ptr<QBase> is returned so that the user does not accidentally leak it.

Definition at line 54 of file quadrature_build.C.

References _dim, _order, 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::QNODAL, libMesh::QSIMPSON, libMesh::QTRAP, libMesh::THIRD, and libMesh::TWENTYTHIRD.

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

◆ clone()

std::unique_ptr< QBase > libMesh::QBase::clone ( ) const
virtual
Returns
A copy of this quadrature rule wrapped in a smart pointer.

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

Definition at line 38 of file quadrature.C.

References build(), get_dim(), get_order(), and type().

39 {
40  return QBase::build(this->type(), this->get_dim(), this->get_order());
41 }
virtual QuadratureType type() const =0
Order get_order() const
Definition: quadrature.h:249
unsigned int get_dim() const
Definition: quadrature.h:150
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = false;
103  return;
104 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 94 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

95 {
96  _enable_print_counter = true;
97  return;
98 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ get_base_order()

Order libMesh::QBase::get_base_order ( ) const
inline
Returns
The "base" order of the quadrature rule, independent of element.

This function should be used when comparing quadrature objects independently of their last initialization.

Definition at line 258 of file quadrature.h.

References _order.

258 { return static_cast<Order>(_order); }
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385

◆ get_dim()

unsigned int libMesh::QBase::get_dim ( ) const
inline
Returns
The spatial dimension of the quadrature rule.

Definition at line 150 of file quadrature.h.

References _dim.

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

150 { return _dim; }
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379

◆ get_elem_type()

ElemType libMesh::QBase::get_elem_type ( ) const
inline
Returns
The element type we're currently using.

Definition at line 121 of file quadrature.h.

References _type.

121 { return _type; }
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391

◆ get_info()

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 (const auto & [name, cd] : _counts)
59  oss << "| " << name << " reference count information:\n"
60  << "| Creations: " << cd.first << '\n'
61  << "| Destructions: " << cd.second << '\n';
62 
63  oss << " ---------------------------------------------------------------------------- \n";
64 
65  return oss.str();
66 
67 #else
68 
69  return "";
70 
71 #endif
72 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
static Counts _counts
Actually holds the data.

◆ get_order()

Order libMesh::QBase::get_order ( ) const
inline
Returns
The current "total" order of the quadrature rule which can vary element by element, depending on the Elem::p_level(), which gets passed to us during init().

Each additional power of p increases the quadrature order required to integrate the mass matrix by 2, hence the formula below.

Todo:
This function should also be used in all of the Order switch statements in the rules themselves.

Definition at line 249 of file quadrature.h.

References _order, and _p_level.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), clone(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_interiors(), libMesh::QGaussLobatto::init_1D(), libMesh::QGauss::init_1D(), libMesh::QConical::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QGrundmann_Moller::init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGrundmann_Moller::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QGauss::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGrundmann_Moller::init_3D(), and libMesh::RBParametrizedFunction::preevaluate_parametrized_function_on_mesh().

249 { return static_cast<Order>(_order + 2 * _p_level); }
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385

◆ get_p_level()

unsigned int libMesh::QBase::get_p_level ( ) const
inline
Returns
The p-refinement level we're currently using.

Definition at line 126 of file quadrature.h.

References _p_level.

126 { return _p_level; }
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403

◆ get_points() [1/2]

const std::vector<Point>& libMesh::QBase::get_points ( ) const
inline

◆ get_points() [2/2]

std::vector<Point>& libMesh::QBase::get_points ( )
inline
Returns
A std::vector containing the quadrature point locations in reference element space as a writable reference.

Definition at line 162 of file quadrature.h.

References _points.

162 { return _points; }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409

◆ get_weights() [1/2]

const std::vector<Real>& libMesh::QBase::get_weights ( ) const
inline

◆ get_weights() [2/2]

std::vector<Real>& libMesh::QBase::get_weights ( )
inline
Returns
A writable references to a std::vector containing the quadrature weights.

Definition at line 174 of file quadrature.h.

References _weights.

174 { return _weights; }
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectednoexceptinherited

Increments the construction counter.

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

Definition at line 183 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

184 {
185  libmesh_try
186  {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189  p.first++;
190  }
191  libmesh_catch (...)
192  {
193  auto stream = libMesh::err.get();
194  stream->exceptions(stream->goodbit); // stream must not throw
195  libMesh::err << "Encountered unrecoverable error while calling "
196  << "ReferenceCounter::increment_constructor_count() "
197  << "for a(n) " << name << " object." << std::endl;
198  std::terminate();
199  }
200 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectednoexceptinherited

Increments the destruction counter.

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

Definition at line 207 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

208 {
209  libmesh_try
210  {
211  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
212  std::pair<unsigned int, unsigned int> & p = _counts[name];
213  p.second++;
214  }
215  libmesh_catch (...)
216  {
217  auto stream = libMesh::err.get();
218  stream->exceptions(stream->goodbit); // stream must not throw
219  libMesh::err << "Encountered unrecoverable error while calling "
220  << "ReferenceCounter::increment_destructor_count() "
221  << "for a(n) " << name << " object." << std::endl;
222  std::terminate();
223  }
224 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ init() [1/4]

void libMesh::QBase::init ( const Elem e,
unsigned int  p_level = invalid_uint 
)
virtual

Initializes the data structures for a quadrature rule for the element e.

If p_level is specified it overrides the element p_level() elevation to use.

Definition at line 65 of file quadrature.C.

References _dim, _elem, _p_level, _type, libMesh::Elem::dim(), init_0D(), init_1D(), init_2D(), init_3D(), libMesh::invalid_uint, libMesh::Elem::p_level(), libMesh::Elem::runtime_topology(), and libMesh::Elem::type().

Referenced by init(), libMesh::QClough::init_1D(), libMesh::QConical::init_1D(), libMesh::QNodal::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QNodal::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGrid::init_3D(), libMesh::QNodal::init_3D(), libMesh::QMonomial::init_3D(), libMesh::OverlapCoupling::operator()(), libMesh::QClough::QClough(), libMesh::QConical::QConical(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QGrid::QGrid(), libMesh::QGrundmann_Moller::QGrundmann_Moller(), libMesh::QJacobi::QJacobi(), libMesh::QMonomial::QMonomial(), libMesh::QNodal::QNodal(), libMesh::QSimpson::QSimpson(), libMesh::QTrap::QTrap(), and libMesh::FE< Dim, LAGRANGE_VEC >::reinit_default_dual_shape_coeffs().

67 {
68  libmesh_assert_equal_to(elem.dim(), _dim);
69 
70  // Default to the element p_level() value
71  if (p == invalid_uint)
72  p = elem.p_level();
73 
74  ElemType t = elem.type();
75 
76  // check to see if we have already
77  // done the work for this quadrature rule
78  //
79  // If we have something like a Polygon subclass then we're going to
80  // need to recompute to be safe; even if we're using the same
81  // element, it might have been distorted enough that its subtriangle
82  // triangulation has been changed.
83  if (t == _type && p == _p_level && !elem.runtime_topology())
84  return;
85  else
86  {
87  _elem = &elem;
88  _type = t;
89  _p_level = p;
90  }
91 
92  switch(_elem->dim())
93  {
94  case 0:
95  this->init_0D();
96 
97  return;
98 
99  case 1:
100  this->init_1D();
101 
102  return;
103 
104  case 2:
105  this->init_2D();
106 
107  return;
108 
109  case 3:
110  this->init_3D();
111 
112  return;
113 
114  default:
115  libmesh_error_msg("Invalid dimension _dim = " << _dim);
116  }
117 }
ElemType
Defines an enum for geometric element types.
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
virtual void init_2D()
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:208
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
const Elem * _elem
The element for which the current values were computed, or nullptr if values were computed without a ...
Definition: quadrature.h:397
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403
virtual void init_0D()
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:198
virtual void init_1D()=0
Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate val...
virtual unsigned short dim() const =0
virtual void init_3D()
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:215

◆ init() [2/4]

void libMesh::QBase::init ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0,
bool  simple_type_only = false 
)
virtual

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

Some types, such as Polygon subclasses, might require more detailed element information and so might not be compatible with this API.

New code should use the Elem-based API for most use cases, but some code may initialize quadrature rules with simple ElemType values like triangles and edges for use in tensor product or conical product constructions; this code can set the simple_type_only flag to avoid being identified as deprecated.

Definition at line 121 of file quadrature.C.

References _dim, _elem, _p_level, _type, libMesh::C0POLYGON, libMesh::C0POLYHEDRON, libMesh::EDGE2, init_0D(), init_1D(), init_2D(), and init_3D().

124 {
125  // Some element types require data from a specific element, so can
126  // only be used with newer APIs.
127  if (t == C0POLYGON || t == C0POLYHEDRON)
128  libmesh_error_msg("Code (see stack trace) used an outdated quadrature function overload.\n"
129  "Quadrature rules on a C0Polygon are not defined by its ElemType alone.");
130 
131  // This API is dangerous to use on general meshes, which may include
132  // element types where the desired quadrature depends on the
133  // physical element, but we still want to be able to initialize
134  // based on only a type for certain simple cases
135  if (t != EDGE2 && !simple_type_only)
136  libmesh_deprecated();
137 
138  // check to see if we have already
139  // done the work for this quadrature rule
140  if (t == _type && p == _p_level)
141  return;
142  else
143  {
144  _elem = nullptr;
145  _type = t;
146  _p_level = p;
147  }
148 
149  switch(_dim)
150  {
151  case 0:
152  this->init_0D();
153 
154  return;
155 
156  case 1:
157  this->init_1D();
158 
159  return;
160 
161  case 2:
162  this->init_2D();
163 
164  return;
165 
166  case 3:
167  this->init_3D();
168 
169  return;
170 
171  default:
172  libmesh_error_msg("Invalid dimension _dim = " << _dim);
173  }
174 }
virtual void init_2D()
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:208
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
const Elem * _elem
The element for which the current values were computed, or nullptr if values were computed without a ...
Definition: quadrature.h:397
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:403
virtual void init_0D()
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:198
virtual void init_1D()=0
Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate val...
virtual void init_3D()
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:215

◆ init() [3/4]

void libMesh::QBase::init ( const QBase other_rule)
virtual

Initializes the data structures for a quadrature rule based on the element, element type, and p_level settings of other_rule.

Definition at line 178 of file quadrature.C.

References _elem, _p_level, _type, and init().

179 {
180  if (other_rule._elem)
181  this->init(*other_rule._elem, other_rule._p_level);
182  else
183  this->init(other_rule._type, other_rule._p_level, true);
184 }
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e.
Definition: quadrature.C:65

◆ init() [4/4]

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 188 of file quadrature.C.

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

191 {
192  // dispatch generic implementation
193  this->init(elem.type(), p_level);
194 }
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e.
Definition: quadrature.C:65

◆ init_0D()

void libMesh::QBase::init_0D ( )
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 198 of file quadrature.C.

References _points, and _weights.

Referenced by init().

199 {
200  _points.resize(1);
201  _weights.resize(1);
202  _points[0] = Point(0.);
203  _weights[0] = 1.0;
204 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

◆ init_1D()

virtual void libMesh::QBase::init_1D ( )
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::QNodal, libMesh::QConical, libMesh::QGrid, libMesh::QGauss, libMesh::QSimpson, libMesh::QTrap, libMesh::QClough, and libMesh::QGaussLobatto.

Referenced by init().

◆ init_2D()

void libMesh::QBase::init_2D ( )
protectedvirtual

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 overridden, throws an error.

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

Definition at line 208 of file quadrature.C.

Referenced by init().

209 {
210  libmesh_not_implemented();
211 }

◆ init_3D()

void libMesh::QBase::init_3D ( )
protectedvirtual

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 overridden, throws an error.

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

Definition at line 215 of file quadrature.C.

Referenced by init().

216 {
217  libmesh_not_implemented();
218 }

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 85 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

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

86  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects
The number of objects.

◆ n_points()

unsigned int libMesh::QBase::n_points ( ) const
inline
Returns
The number of points associated with the quadrature rule.

Definition at line 131 of file quadrature.h.

References _points, and libMesh::libmesh_assert().

Referenced by libMesh::ExactSolution::_compute_error(), assemble(), assemble_1D(), assemble_ellipticdg(), assemble_helmholtz(), assemble_poisson(), assemble_shell(), assemble_wave(), assembly_with_dg_fem_context(), AssemblyA0::boundary_assembly(), AssemblyA1::boundary_assembly(), AssemblyF0::boundary_assembly(), AssemblyF1::boundary_assembly(), AssemblyA2::boundary_assembly(), AssemblyF2::boundary_assembly(), A2::boundary_assembly(), A3::boundary_assembly(), F0::boundary_assembly(), Output0::boundary_assembly(), compute_enriched_soln(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::HDGProblem::create_identity_jacobian(), libMesh::HDGProblem::create_identity_residual(), SecondOrderScalarSystemSecondOrderTimeSolverBase::damping_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::damping_residual(), NavierSystem::element_constraint(), CoupledSystem::element_constraint(), PoissonSystem::element_postprocess(), LaplaceSystem::element_postprocess(), LaplaceQoI::element_qoi(), HeatSystem::element_qoi(), LaplaceQoI::element_qoi_derivative(), LaplaceSystem::element_qoi_derivative(), HeatSystem::element_qoi_derivative(), NavierSystem::element_time_derivative(), SolidSystem::element_time_derivative(), PoissonSystem::element_time_derivative(), LaplaceSystem::element_time_derivative(), CurlCurlSystem::element_time_derivative(), ElasticitySystem::element_time_derivative(), CoupledSystem::element_time_derivative(), HilbertSystem::element_time_derivative(), HeatSystem::element_time_derivative(), FirstOrderScalarSystemBase::element_time_derivative(), SecondOrderScalarSystemFirstOrderTimeSolverBase::element_time_derivative(), fe_assembly(), A0::interior_assembly(), B::interior_assembly(), M0::interior_assembly(), A1::interior_assembly(), AssemblyA0::interior_assembly(), AcousticsInnerProduct::interior_assembly(), AssemblyA1::interior_assembly(), A2::interior_assembly(), AssemblyA2::interior_assembly(), F0::interior_assembly(), OutputAssembly::interior_assembly(), EIM_IP_assembly::interior_assembly(), EIM_F::interior_assembly(), AssemblyEIM::interior_assembly(), InnerProductAssembly::interior_assembly(), AssemblyF0::interior_assembly(), AssemblyF1::interior_assembly(), Ex6InnerProduct::interior_assembly(), LaplaceYoung::jacobian(), NavierSystem::mass_residual(), ElasticitySystem::mass_residual(), libMesh::FEMPhysics::mass_residual(), FirstOrderScalarSystemBase::mass_residual(), SecondOrderScalarSystemSecondOrderTimeSolverBase::mass_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::mass_residual(), libMesh::FEAbstract::n_quadrature_points(), print_info(), LaplaceYoung::residual(), LaplaceSystem::side_constraint(), LaplaceSystem::side_postprocess(), CoupledSystemQoI::side_qoi(), CoupledSystemQoI::side_qoi_derivative(), LaplaceSystem::side_qoi_derivative(), SolidSystem::side_time_derivative(), CurlCurlSystem::side_time_derivative(), ElasticitySystem::side_time_derivative(), size(), tensor_product_hex(), tensor_product_prism(), tensor_product_quad(), and QuadratureTest::testPolynomial().

132  {
133  libmesh_assert (!_points.empty());
134  return cast_int<unsigned int>(_points.size());
135  }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
libmesh_assert(ctx)

◆ operator=() [1/2]

QBase& libMesh::QBase::operator= ( const QBase )
default

◆ operator=() [2/2]

QBase& libMesh::QBase::operator= ( QBase &&  )
default

◆ print_info() [1/2]

void libMesh::ReferenceCounter::print_info ( std::ostream &  out_stream = libMesh::out)
staticinherited

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

Definition at line 81 of file reference_counter.C.

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

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

82 {
84  out_stream << ReferenceCounter::get_info();
85 }
static std::string get_info()
Gets a string containing the reference information.
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ print_info() [2/2]

void libMesh::QBase::print_info ( std::ostream &  os = libMesh::out) const

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

Definition at line 43 of file quadrature.C.

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

Referenced by libMesh::operator<<().

44 {
45  libmesh_assert(!_points.empty());
46  libmesh_assert(!_weights.empty());
47 
48  Real summed_weights=0;
49  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
50  for (auto qpoint: index_range(_points))
51  {
52  os << " Point " << qpoint << ":\n"
53  << " "
54  << _points[qpoint]
55  << "\n Weight:\n "
56  << " w=" << _weights[qpoint] << "\n" << std::endl;
57 
58  summed_weights += _weights[qpoint];
59  }
60  os << "Summed Weights: " << summed_weights << std::endl;
61 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
libmesh_assert(ctx)
unsigned int n_points() const
Definition: quadrature.h:131
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ qp()

Point libMesh::QBase::qp ( const unsigned int  i) const
inline
Returns
The \( i^{th} \) quadrature point in reference element space.

Definition at line 179 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(), tensor_product_quad(), and QuadratureTest::testPolynomial().

180  {
181  libmesh_assert_less (i, _points.size());
182  return _points[i];
183  }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409

◆ scale()

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 222 of file quadrature.C.

References _dim, _points, _weights, libMesh::index_range(), and libMesh::Real.

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

224 {
225  // Make sure we are in 1D
226  libmesh_assert_equal_to (_dim, 1);
227 
228  Real
229  h_new = new_range.second - new_range.first,
230  h_old = old_range.second - old_range.first;
231 
232  // Make sure that we have sane ranges
233  libmesh_assert_greater (h_new, 0.);
234  libmesh_assert_greater (h_old, 0.);
235 
236  // Make sure there are some points
237  libmesh_assert_greater (_points.size(), 0);
238 
239  // Compute the scale factor
240  Real scfact = h_new/h_old;
241 
242  // We're mapping from old_range -> new_range
243  for (auto i : index_range(_points))
244  {
245  _points[i](0) = new_range.first +
246  (_points[i](0) - old_range.first) * scfact;
247 
248  // Scale the weights
249  _weights[i] *= scfact;
250  }
251 }
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117

◆ shapes_need_reinit()

virtual bool libMesh::QBase::shapes_need_reinit ( )
inlinevirtual
Returns
true if the shape functions need to be recalculated, false otherwise.

This may be required if the number of quadrature points or their position changes.

Definition at line 286 of file quadrature.h.

286 { return false; }

◆ size()

unsigned int libMesh::QBase::size ( ) const
inline

Alias for n_points() to enable use in index_range.

Returns
The number of points associated with the quadrature rule.

Definition at line 142 of file quadrature.h.

References n_points().

Referenced by libMesh::Elem::true_centroid().

143  {
144  return n_points();
145  }
unsigned int n_points() const
Definition: quadrature.h:131

◆ tensor_product_hex()

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 283 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::QGauss::init_3D(), and libMesh::QGrid::init_3D().

284 {
285  const unsigned int np = q1D.n_points();
286 
287  _points.resize(np * np * np);
288 
289  _weights.resize(np * np * np);
290 
291  unsigned int q=0;
292 
293  for (unsigned int k=0; k<np; k++)
294  for (unsigned int j=0; j<np; j++)
295  for (unsigned int i=0; i<np; i++)
296  {
297  _points[q](0) = q1D.qp(i)(0);
298  _points[q](1) = q1D.qp(j)(0);
299  _points[q](2) = q1D.qp(k)(0);
300 
301  _weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);
302 
303  q++;
304  }
305 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

◆ tensor_product_prism()

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 310 of file quadrature.C.

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

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

311 {
312  const unsigned int n_points1D = q1D.n_points();
313  const unsigned int n_points2D = q2D.n_points();
314 
315  _points.resize (n_points1D * n_points2D);
316  _weights.resize (n_points1D * n_points2D);
317 
318  unsigned int q=0;
319 
320  for (unsigned int j=0; j<n_points1D; j++)
321  for (unsigned int i=0; i<n_points2D; i++)
322  {
323  _points[q](0) = q2D.qp(i)(0);
324  _points[q](1) = q2D.qp(i)(1);
325  _points[q](2) = q1D.qp(j)(0);
326 
327  _weights[q] = q2D.w(i) * q1D.w(j);
328 
329  q++;
330  }
331 
332 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

◆ tensor_product_quad()

void libMesh::QBase::tensor_product_quad ( const QBase q1D)
protected

Constructs a 2D rule from the tensor product of q1D with itself.

Used in the init_2D() routines for quadrilateral element types.

Definition at line 256 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::QGauss::init_2D(), and libMesh::QGrid::init_2D().

257 {
258 
259  const unsigned int np = q1D.n_points();
260 
261  _points.resize(np * np);
262 
263  _weights.resize(np * np);
264 
265  unsigned int q=0;
266 
267  for (unsigned int j=0; j<np; j++)
268  for (unsigned int i=0; i<np; i++)
269  {
270  _points[q](0) = q1D.qp(i)(0);
271  _points[q](1) = q1D.qp(j)(0);
272 
273  _weights[q] = q1D.w(i)*q1D.w(j);
274 
275  q++;
276  }
277 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

◆ type()

virtual QuadratureType libMesh::QBase::type ( ) const
pure virtual

◆ w()

Real libMesh::QBase::w ( const unsigned int  i) const
inline
Returns
The \( i^{th} \) quadrature weight.

Definition at line 188 of file quadrature.h.

References _weights.

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

189  {
190  libmesh_assert_less (i, _weights.size());
191  return _weights[i];
192  }
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415

Friends And Related Function Documentation

◆ operator<<

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

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

Definition at line 337 of file quadrature.C.

338 {
339  q.print_info(os);
340  return os;
341 }

Member Data Documentation

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

Actually holds the data.

Definition at line 124 of file reference_counter.h.

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

◆ _dim

unsigned int libMesh::QBase::_dim
protected

◆ _elem

const Elem* libMesh::QBase::_elem
protected

The element for which the current values were computed, or nullptr if values were computed without a specific element.

Definition at line 397 of file quadrature.h.

Referenced by init(), libMesh::QGauss::init_2D(), libMesh::QNodal::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QNodal::init_3D().

◆ _enable_print_counter

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 143 of file reference_counter.h.

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

◆ _mutex

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

◆ _n_objects

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 132 of file reference_counter.h.

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

◆ _order

Order libMesh::QBase::_order
protected

◆ _p_level

unsigned int libMesh::QBase::_p_level
protected

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

Definition at line 403 of file quadrature.h.

Referenced by get_order(), get_p_level(), init(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGauss::init_3D().

◆ _points

std::vector<Point> libMesh::QBase::_points
protected

◆ _type

ElemType libMesh::QBase::_type
protected

◆ _weights

std::vector<Real> libMesh::QBase::_weights
protected

◆ allow_nodal_pyramid_quadrature

bool libMesh::QBase::allow_nodal_pyramid_quadrature

The flag's value defaults to false so that one does not accidentally use a nodal quadrature rule on Pyramid elements, since evaluating the inverse element Jacobian (e.g.

dphi) is not well-defined at the Pyramid apex because the element Jacobian is zero there.

We do not want to completely prevent someone from using a nodal quadrature rule on Pyramids, however, since there are legitimate use cases (lumped mass matrix) so the flag can be set to true to override this behavior.

Definition at line 314 of file quadrature.h.

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

◆ allow_rules_with_negative_weights

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 301 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: