libMesh
Public Member Functions | Static Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
libMesh::QNodal Class Reference

This class implements nodal quadrature rules for various element types. More...

#include <quadrature_nodal.h>

Inheritance diagram for libMesh::QNodal:
[legend]

Public Member Functions

 QNodal (unsigned int dim, Order order=FIRST)
 Constructor. More...
 
 QNodal (const QNodal &)=default
 Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class. More...
 
 QNodal (QNodal &&)=default
 
QNodaloperator= (const QNodal &)=default
 
QNodaloperator= (QNodal &&)=default
 
virtual ~QNodal ()=default
 
virtual QuadratureType type () const override
 
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)
 Initializes the data structures for a quadrature rule for an element of type type. 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
 
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 (const std::string &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 void print_info (std::ostream &out=libMesh::out)
 Prints the reference information, by default to libMesh::out. More...
 
static std::string get_info ()
 Gets a string containing the reference information. 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...
 

Protected Types

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

Protected Member Functions

virtual void init_0D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
 Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate 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)
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name)
 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...
 
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...
 

Private Member Functions

virtual void init_1D (const ElemType, unsigned int) override
 Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
virtual void init_2D (const ElemType, unsigned int) override
 Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 
virtual void init_3D (const ElemType, unsigned int) override
 Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate values. More...
 

Detailed Description

This class implements nodal quadrature rules for various element types.

For some element types, this corresponds to well-known quadratures such as the trapezoidal rule and Simpson's rule. For other element types (e.g. the serendipity QUAD8/HEX20/PRISM15) it implements a lower-order-accurate nodal quadrature which still generates a strictly positive definite (diagonal) mass matrix. Since the main purpose of this class is provide appropriate nodal quadrature rules, "order" arguments passed during initialization are ignored.

Author
John W. Peterson
Date
2019

Implements nodal quadrature for various element types.

Definition at line 44 of file quadrature_nodal.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 117 of file reference_counter.h.

Constructor & Destructor Documentation

◆ QNodal() [1/3]

libMesh::QNodal::QNodal ( unsigned int  dim,
Order  order = FIRST 
)
inlineexplicit

Constructor.

Declares the order of the quadrature rule. We explicitly call the init function in 1D since the other tensor-product rules require this one.

Note
The element type, EDGE2, will not be used internally, however if we called the function with INVALID_ELEM it would try to be smart and return, thinking it had already done the work.

Definition at line 58 of file quadrature_nodal.h.

59  :
60  QBase(dim,order)
61  {
62  if (_dim == 1)
63  init(EDGE2);
64  }

References libMesh::QBase::_dim, libMesh::EDGE2, and libMesh::QBase::init().

◆ QNodal() [2/3]

libMesh::QNodal::QNodal ( const QNodal )
default

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

◆ QNodal() [3/3]

libMesh::QNodal::QNodal ( QNodal &&  )
default

◆ ~QNodal()

virtual libMesh::QNodal::~QNodal ( )
virtualdefault

Member Function Documentation

◆ build() [1/2]

std::unique_ptr< QBase > libMesh::QBase::build ( const QuadratureType  qt,
const unsigned int  dim,
const Order  order = INVALID_ORDER 
)
staticinherited

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 52 of file quadrature_build.C.

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

References libMesh::QBase::_dim, libMesh::QBase::_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.

◆ build() [2/2]

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

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 41 of file quadrature_build.C.

44 {
45  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
46  _dim,
47  _order);
48 }

References libMesh::QBase::_dim, libMesh::QBase::_order, and libMesh::QBase::type().

Referenced by assemble_poisson(), libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), libMesh::InfFE< Dim, T_radial, T_map >::reinit(), QuadratureTest::test1DWeights(), QuadratureTest::test2DWeights(), QuadratureTest::test3DWeights(), QuadratureTest::testBuild(), QuadratureTest::testJacobi(), QuadratureTest::testMonomialQuadrature(), QuadratureTest::testTetQuadrature(), and QuadratureTest::testTriQuadrature().

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

107 {
108  _enable_print_counter = false;
109  return;
110 }

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ 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 100 of file reference_counter.C.

101 {
102  _enable_print_counter = true;
103  return;
104 }

References libMesh::ReferenceCounter::_enable_print_counter.

◆ get_dim()

unsigned int libMesh::QBase::get_dim ( ) const
inlineinherited

◆ get_elem_type()

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

Definition at line 116 of file quadrature.h.

116 { return _type; }

References libMesh::QBase::_type.

◆ 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.

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 & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }

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

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

◆ get_order()

Order libMesh::QBase::get_order ( ) const
inlineinherited
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.

Definition at line 211 of file quadrature.h.

211 { return static_cast<Order>(_order + 2 * _p_level); }

References libMesh::QBase::_order, and libMesh::QBase::_p_level.

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

◆ get_p_level()

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

Definition at line 121 of file quadrature.h.

121 { return _p_level; }

References libMesh::QBase::_p_level.

◆ get_points() [1/2]

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

Definition at line 147 of file quadrature.h.

147 { return _points; }

References libMesh::QBase::_points.

◆ get_points() [2/2]

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

◆ get_weights() [1/2]

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

Definition at line 159 of file quadrature.h.

159 { return _weights; }

References libMesh::QBase::_weights.

◆ get_weights() [2/2]

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

◆ increment_constructor_count()

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

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }

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

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

◆ increment_destructor_count()

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

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }

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

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

◆ init() [1/2]

void libMesh::QBase::init ( const Elem elem,
const std::vector< Real > &  vertex_distance_func,
unsigned int  p_level = 0 
)
virtualinherited

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

106 {
107  // dispatch generic implementation
108  this->init(elem.type(), p_level);
109 }

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

◆ init() [2/2]

void libMesh::QBase::init ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
virtualinherited

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

Definition at line 59 of file quadrature.C.

61 {
62  // check to see if we have already
63  // done the work for this quadrature rule
64  if (t == _type && p == _p_level)
65  return;
66  else
67  {
68  _type = t;
69  _p_level = p;
70  }
71 
72 
73 
74  switch(_dim)
75  {
76  case 0:
77  this->init_0D();
78 
79  return;
80 
81  case 1:
82  this->init_1D();
83 
84  return;
85 
86  case 2:
87  this->init_2D();
88 
89  return;
90 
91  case 3:
92  this->init_3D();
93 
94  return;
95 
96  default:
97  libmesh_error_msg("Invalid dimension _dim = " << _dim);
98  }
99 }

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

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

◆ init_0D()

void libMesh::QBase::init_0D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
protectedvirtualinherited

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.

Note
The arguments should no longer be used for anything in derived classes, they are only maintained for backwards compatibility and will eventually be removed.

Definition at line 113 of file quadrature.C.

114 {
115  _points.resize(1);
116  _weights.resize(1);
117  _points[0] = Point(0.);
118  _weights[0] = 1.0;
119 }

References libMesh::QBase::_points, and libMesh::QBase::_weights.

Referenced by libMesh::QBase::init().

◆ init_1D()

void libMesh::QNodal::init_1D ( const  type,
unsigned int  p_level 
)
overrideprivatevirtual

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.

Note
The arguments should no longer be used for anything in derived classes, they are only maintained for backwards compatibility and will eventually be removed.

Implements libMesh::QBase.

Definition at line 29 of file quadrature_nodal_1D.C.

30 {
31  switch (_type)
32  {
33  case EDGE2:
34  {
35  // Nodal quadrature on an Edge2 is QTrap
36  QTrap rule(/*dim=*/1, /*ignored*/_order);
37  rule.init(_type, /*ignored*/_p_level);
38  _points.swap (rule.get_points());
39  _weights.swap(rule.get_weights());
40  return;
41  }
42  case EDGE3:
43  {
44  // Nodal quadrature on an Edge3 is QSimpson
45  QSimpson rule(/*dim=*/1, /*ignored*/_order);
46  rule.init(_type, /*ignored*/_p_level);
47  _points.swap (rule.get_points());
48  _weights.swap(rule.get_weights());
49  return;
50  }
51  case EDGE4:
52  {
53  // The 4-point variant of Simpson's rule. The quadrature
54  // points are in the same order as the reference element
55  // nodes.
56  _points = {Point(-1,0.,0.), Point(+1,0.,0.),
57  Point(-Real(1)/3,0.,0.), Point(Real(1)/3,0.,0.)};
58  _weights = {0.25, 0.25, 0.75, 0.75};
59  return;
60  }
61  default:
62  libmesh_error_msg("Element type not supported:" << Utility::enum_to_string(_type));
63  }
64 }

References libMesh::QBase::_order, libMesh::QBase::_p_level, libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::Utility::enum_to_string(), libMesh::QBase::get_points(), libMesh::QBase::get_weights(), libMesh::QBase::init(), and libMesh::Real.

◆ init_2D()

void libMesh::QNodal::init_2D ( const  type,
unsigned int  p_level 
)
overrideprivatevirtual

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.

Note
The arguments should no longer be used for anything in derived classes, they are only maintained for backwards compatibility and will eventually be removed.

Reimplemented from libMesh::QBase.

Definition at line 29 of file quadrature_nodal_2D.C.

30 {
31 #if LIBMESH_DIM > 1
32 
33  switch (_type)
34  {
35 
36  case QUAD4:
37  case QUADSHELL4:
38  case TRI3:
39  case TRISHELL3:
40  {
41  QTrap rule(/*dim=*/2, /*ignored*/_order);
42  rule.init(_type, /*ignored*/_p_level);
43  _points.swap (rule.get_points());
44  _weights.swap(rule.get_weights());
45  return;
46  }
47 
48  case QUAD8:
49  case QUADSHELL8:
50  {
51  // A rule with 8 points which is exact for linears, and
52  // naturally produces a lumped approximation to the mass
53  // matrix. The quadrature points are numbered the same way as
54  // the reference element nodes.
55  _points =
56  {
57  Point(-1,-1), Point(+1,-1), Point(+1,+1), Point(-1,+1),
58  Point(0.,-1), Point(+1,0.), Point(0.,+1), Point(-1,0.)
59  };
60 
61  // The weights for the Quad8 nodal quadrature rule are
62  // obtained from the following specific steps. Other
63  // "serendipity" type rules are obtained similarly.
64  //
65  // 1.) Due to the symmetry of the bi-unit square domain, we
66  // first note that there are only two "classes" of node in the
67  // Quad8: vertices (with associated weight wv) and edges (with
68  // associated weight we).
69  //
70  // 2.) In order for such a nodal quadrature rule to be exact
71  // for constants, the weights must sum to the area of the
72  // reference element, and therefore we must have:
73  // 4*wv + 4*we = 4 --> wv + we = 1
74  //
75  // 3.) Due to symmetry (once again), such a rule is then
76  // automatically exact for the linear polynomials "x" and "y",
77  // regardless of the values of wv and we.
78  //
79  // 4.) We therefore still have two unknowns and one equation.
80  // Attempting to additionally make the nodal quadrature rule
81  // exact for the quadratic polynomials "x^2" and "y^2" leads
82  // to a rule with negative weights, namely:
83  // wv = -1/3
84  // we = 4/3
85  // Note: these are the same values one obtains by integrating
86  // the corresponding Quad8 Lagrange basis functions over the
87  // reference element.
88  //
89  // Since the weights appear on the diagonal of the nodal
90  // quadrature rule's approximation to the mass matrix, rules
91  // with negative weights yield an indefinite mass matrix,
92  // i.e. one with both positive and negative eigenvalues. Rules
93  // with negative weights can also produce a negative integral
94  // approximation to a strictly positive integrand, which may
95  // be unacceptable in some situations.
96  //
97  // 5.) Requiring all positive weights therefore restricts the
98  // nodal quadrature rule to only be exact for linears. But
99  // there is still one free parameter remaining, and thus we
100  // need some other criterion in order to complete the
101  // description of the rule.
102  //
103  // Here, we have decided to choose the quadrature weights in
104  // such a way that the resulting nodal quadrature mass matrix
105  // approximation is "proportional" to the true mass matrix's
106  // diagonal entries for the reference element. We therefore
107  // pose the following constrained optimization problem:
108  //
109  // { min_{wv, we} |diag(M) - C*diag(W)|^2
110  // (O) {
111  // { subject to wv + we = 1
112  //
113  // where:
114  // * M is the true mass matrix
115  // * W is the nodal quadrature approximation to the mass
116  // matrix. In this particular case:
117  // diag(W) = [wv,wv,wv,wv,we,we,we,we]
118  // * C = tr(M) / vol(E) is the ratio between the trace of the
119  // true mass matrix and the volume of the reference
120  // element. For all Lagrange finite elements, we have C<1.
121  //
122  // 6.) The optimization problem (O) is solved directly by
123  // substituting the algebraic constraint into the functional
124  // which is to be minimized, then setting the derivative with
125  // respect to the remaining parameter equal to zero and
126  // solving. In the Quad8 and Hex20 cases, there is only one
127  // free parameter, while in the Prism15 case there are two
128  // free parameters, so a 2x2 system of linear equations must
129  // be solved.
130  Real wv = Real(12) / 79;
131  Real we = Real(67) / 79;
132 
133  _weights = {wv, wv, wv, wv, we, we, we, we};
134 
135  return;
136  }
137 
138  case QUAD9:
139  case TRI6:
140  {
141  QSimpson rule(/*dim=*/2, /*ignored*/_order);
142  rule.init(_type, /*ignored*/_p_level);
143  _points.swap (rule.get_points());
144  _weights.swap(rule.get_weights());
145  return;
146  }
147 
148  default:
149  libmesh_error_msg("Element type not supported!:" << Utility::enum_to_string(_type));
150  }
151 #endif
152 }

References libMesh::QBase::_order, libMesh::QBase::_p_level, libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::Utility::enum_to_string(), libMesh::QBase::get_points(), libMesh::QBase::get_weights(), libMesh::QBase::init(), libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::Real, libMesh::TRI3, libMesh::TRI6, and libMesh::TRISHELL3.

◆ init_3D()

void libMesh::QNodal::init_3D ( const  type,
unsigned int  p_level 
)
overrideprivatevirtual

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.

Note
The arguments should no longer be used for anything in derived classes, they are only maintained for backwards compatibility and will eventually be removed.

Reimplemented from libMesh::QBase.

Definition at line 29 of file quadrature_nodal_3D.C.

30 {
31 #if LIBMESH_DIM == 3
32 
33  switch (_type)
34  {
35  case TET4:
36  case PRISM6:
37  case HEX8:
38  {
39  QTrap rule(/*dim=*/3, /*ignored*/_order);
40  rule.init(_type, /*ignored*/_p_level);
41  _points.swap (rule.get_points());
42  _weights.swap(rule.get_weights());
43  return;
44  }
45 
46  case PRISM15:
47  {
48  // A rule with 15 points which is exact for linears, and
49  // naturally produces a lumped approximation to the mass
50  // matrix. The quadrature points are numbered the same way as
51  // the reference element nodes.
52  _points =
53  {
54  Point(0.,0.,-1), Point(+1,0.,-1), Point(0.,+1,-1),
55  Point(0.,0.,+1), Point(+1,0.,+1), Point(0.,+1,+1),
56  Point(.5,0.,-1), Point(.5,.5,-1), Point(0.,.5,-1),
57  Point(0.,0.,0.), Point(+1,0.,0.), Point(0.,+1,0.),
58  Point(.5,0.,+1), Point(.5,.5,+1), Point(0.,.5,+1),
59  };
60 
61  // vertex (wv), tri edge (wt), and quad edge (wq) weights are
62  // obtained using the same approach that was used for the Quad8,
63  // see quadrature_nodal_2D.C for details.
64  Real wv = Real(1) / 34;
65  Real wt = Real(4) / 51;
66  Real wq = Real(2) / 17;
67 
68  _weights = {wv, wv, wv, wv, wv, wv,
69  wt, wt, wt,
70  wq, wq, wq,
71  wt, wt, wt};
72 
73  return;
74  }
75 
76  case HEX20:
77  {
78  // A rule with 20 points which is exact for linears, and
79  // naturally produces a lumped approximation to the mass
80  // matrix. The quadrature points are numbered the same way as
81  // the reference element nodes.
82  _points =
83  {
84  Point(-1,-1,-1), Point(+1,-1,-1), Point(+1,+1,-1), Point(-1,+1,-1),
85  Point(-1,-1,+1), Point(+1,-1,+1), Point(+1,+1,+1), Point(-1,+1,+1),
86  Point(0.,-1,-1), Point(+1,0.,-1), Point(0.,+1,-1), Point(-1,0.,-1),
87  Point(-1,-1,0.), Point(+1,-1,0.), Point(+1,+1,0.), Point(-1,+1,0.),
88  Point(0.,-1,+1), Point(+1,0.,+1), Point(0.,+1,+1), Point(-1,0.,+1)
89  };
90 
91  // vertex (wv), and edge (we) weights are obtained using the
92  // same approach that was used for the Quad8, see
93  // quadrature_nodal_2D.C for details.
94  Real wv = Real(7) / 31;
95  Real we = Real(16) / 31;
96 
97  _weights = {wv, wv, wv, wv, wv, wv, wv, wv,
98  we, we, we, we, we, we, we, we, we, we, we, we};
99 
100  return;
101  }
102 
103  case TET10:
104  case PRISM18:
105  case HEX27:
106  {
107  QSimpson rule(/*dim=*/3, /*ignored*/_order);
108  rule.init(_type, /*ignored*/_p_level);
109  _points.swap (rule.get_points());
110  _weights.swap(rule.get_weights());
111  return;
112  }
113 
114  default:
115  libmesh_error_msg("ERROR: Unsupported type: " << Utility::enum_to_string(_type));
116  }
117 #endif
118 }

References libMesh::QBase::_order, libMesh::QBase::_p_level, libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::Utility::enum_to_string(), libMesh::QBase::get_points(), libMesh::QBase::get_weights(), libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::QBase::init(), libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::Real, libMesh::TET10, and libMesh::TET4.

◆ n_objects()

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

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

Definition at line 83 of file reference_counter.h.

84  { return _n_objects; }

References libMesh::ReferenceCounter::_n_objects.

◆ n_points()

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

Definition at line 126 of file quadrature.h.

127  {
128  libmesh_assert (!_points.empty());
129  return cast_int<unsigned int>(_points.size());
130  }

References libMesh::QBase::_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(), AssemblyF0::boundary_assembly(), AssemblyA1::boundary_assembly(), AssemblyF1::boundary_assembly(), AssemblyF2::boundary_assembly(), A2::boundary_assembly(), AssemblyA2::boundary_assembly(), A3::boundary_assembly(), F0::boundary_assembly(), Output0::boundary_assembly(), libMesh::System::calculate_norm(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), SecondOrderScalarSystemSecondOrderTimeSolverBase::damping_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::damping_residual(), NavierSystem::element_constraint(), CoupledSystem::element_constraint(), PoissonSystem::element_postprocess(), LaplaceSystem::element_postprocess(), LaplaceQoI::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(), L2System::element_time_derivative(), CoupledSystem::element_time_derivative(), FirstOrderScalarSystemBase::element_time_derivative(), SecondOrderScalarSystemFirstOrderTimeSolverBase::element_time_derivative(), libMesh::RBEIMConstruction::enrich_RB_space(), B::interior_assembly(), A0::interior_assembly(), M0::interior_assembly(), A1::interior_assembly(), AssemblyA0::interior_assembly(), EIM_IP_assembly::interior_assembly(), AcousticsInnerProduct::interior_assembly(), AssemblyA1::interior_assembly(), A2::interior_assembly(), AssemblyA2::interior_assembly(), EIM_F::interior_assembly(), F0::interior_assembly(), OutputAssembly::interior_assembly(), InnerProductAssembly::interior_assembly(), AssemblyEIM::interior_assembly(), AssemblyF0::interior_assembly(), AssemblyF1::interior_assembly(), Ex6InnerProduct::interior_assembly(), Ex6EIMInnerProduct::interior_assembly(), LaplaceYoung::jacobian(), NavierSystem::mass_residual(), ElasticitySystem::mass_residual(), libMesh::FEMPhysics::mass_residual(), FirstOrderScalarSystemBase::mass_residual(), SecondOrderScalarSystemSecondOrderTimeSolverBase::mass_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::mass_residual(), libMesh::QBase::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(), libMesh::QBase::tensor_product_hex(), libMesh::QBase::tensor_product_prism(), libMesh::QBase::tensor_product_quad(), and libMesh::RBEIMConstruction::truth_solve().

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ print_info() [1/2]

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

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

Definition at line 37 of file quadrature.C.

38 {
39  libmesh_assert(!_points.empty());
40  libmesh_assert(!_weights.empty());
41 
42  Real summed_weights=0;
43  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
44  for (auto qpoint: index_range(_points))
45  {
46  os << " Point " << qpoint << ":\n"
47  << " "
48  << _points[qpoint]
49  << "\n Weight:\n "
50  << " w=" << _weights[qpoint] << "\n" << std::endl;
51 
52  summed_weights += _weights[qpoint];
53  }
54  os << "Summed Weights: " << summed_weights << std::endl;
55 }

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

Referenced by libMesh::operator<<().

◆ print_info() [2/2]

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

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

Definition at line 87 of file reference_counter.C.

88 {
90  out_stream << ReferenceCounter::get_info();
91 }

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

◆ qp()

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

Definition at line 164 of file quadrature.h.

165  {
166  libmesh_assert_less (i, _points.size());
167  return _points[i];
168  }

References libMesh::QBase::_points.

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

◆ scale()

void libMesh::QBase::scale ( std::pair< Real, Real old_range,
std::pair< Real, Real new_range 
)
inherited

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

139 {
140  // Make sure we are in 1D
141  libmesh_assert_equal_to (_dim, 1);
142 
143  Real
144  h_new = new_range.second - new_range.first,
145  h_old = old_range.second - old_range.first;
146 
147  // Make sure that we have sane ranges
148  libmesh_assert_greater (h_new, 0.);
149  libmesh_assert_greater (h_old, 0.);
150 
151  // Make sure there are some points
152  libmesh_assert_greater (_points.size(), 0);
153 
154  // Compute the scale factor
155  Real scfact = h_new/h_old;
156 
157  // We're mapping from old_range -> new_range
158  for (auto i : index_range(_points))
159  {
160  _points[i](0) = new_range.first +
161  (_points[i](0) - old_range.first) * scfact;
162 
163  // Scale the weights
164  _weights[i] *= scfact;
165  }
166 }

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

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

◆ shapes_need_reinit()

virtual bool libMesh::QBase::shapes_need_reinit ( )
inlinevirtualinherited
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 239 of file quadrature.h.

239 { return false; }

◆ tensor_product_hex()

void libMesh::QBase::tensor_product_hex ( const QBase q1D)
protectedinherited

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

199 {
200  const unsigned int np = q1D.n_points();
201 
202  _points.resize(np * np * np);
203 
204  _weights.resize(np * np * np);
205 
206  unsigned int q=0;
207 
208  for (unsigned int k=0; k<np; k++)
209  for (unsigned int j=0; j<np; j++)
210  for (unsigned int i=0; i<np; i++)
211  {
212  _points[q](0) = q1D.qp(i)(0);
213  _points[q](1) = q1D.qp(j)(0);
214  _points[q](2) = q1D.qp(k)(0);
215 
216  _weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);
217 
218  q++;
219  }
220 }

References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::n_points(), libMesh::QBase::qp(), and libMesh::QBase::w().

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

◆ tensor_product_prism()

void libMesh::QBase::tensor_product_prism ( const QBase q1D,
const QBase q2D 
)
protectedinherited

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

226 {
227  const unsigned int n_points1D = q1D.n_points();
228  const unsigned int n_points2D = q2D.n_points();
229 
230  _points.resize (n_points1D * n_points2D);
231  _weights.resize (n_points1D * n_points2D);
232 
233  unsigned int q=0;
234 
235  for (unsigned int j=0; j<n_points1D; j++)
236  for (unsigned int i=0; i<n_points2D; i++)
237  {
238  _points[q](0) = q2D.qp(i)(0);
239  _points[q](1) = q2D.qp(i)(1);
240  _points[q](2) = q1D.qp(j)(0);
241 
242  _weights[q] = q2D.w(i) * q1D.w(j);
243 
244  q++;
245  }
246 
247 }

References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::n_points(), libMesh::QBase::qp(), and libMesh::QBase::w().

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

◆ tensor_product_quad()

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

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

172 {
173 
174  const unsigned int np = q1D.n_points();
175 
176  _points.resize(np * np);
177 
178  _weights.resize(np * np);
179 
180  unsigned int q=0;
181 
182  for (unsigned int j=0; j<np; j++)
183  for (unsigned int i=0; i<np; i++)
184  {
185  _points[q](0) = q1D.qp(i)(0);
186  _points[q](1) = q1D.qp(j)(0);
187 
188  _weights[q] = q1D.w(i)*q1D.w(j);
189 
190  q++;
191  }
192 }

References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::n_points(), libMesh::QBase::qp(), and libMesh::QBase::w().

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

◆ type()

QuadratureType libMesh::QNodal::type ( ) const
overridevirtual
Returns
QNODAL.

Implements libMesh::QBase.

Definition at line 26 of file quadrature_nodal.C.

27 {
28  return QNODAL;
29 }

References libMesh::QNODAL.

◆ w()

Real libMesh::QBase::w ( const unsigned int  i) const
inlineinherited

Member Data Documentation

◆ _counts

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

◆ _dim

unsigned int libMesh::QBase::_dim
protectedinherited

◆ _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 141 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 135 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 130 of file reference_counter.h.

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

◆ _order

Order libMesh::QBase::_order
protectedinherited

◆ _p_level

unsigned int libMesh::QBase::_p_level
protectedinherited

◆ _points

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

◆ _type

ElemType libMesh::QBase::_type
protectedinherited

◆ _weights

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

◆ allow_rules_with_negative_weights

bool libMesh::QBase::allow_rules_with_negative_weights
inherited

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 254 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:
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::QBase::_p_level
unsigned int _p_level
The p-level of the element for which the current values have been computed.
Definition: quadrature.h:357
libMesh::QMONOMIAL
Definition: enum_quadrature_type.h:41
libMesh::QCLOUGH
Definition: enum_quadrature_type.h:44
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::EDGE4
Definition: enum_elem_type.h:37
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::ReferenceCounter::_counts
static Counts _counts
Actually holds the data.
Definition: reference_counter.h:122
libMesh::QBase::init_1D
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate val...
libMesh::ReferenceCounter::_n_objects
static Threads::atomic< unsigned int > _n_objects
The number of objects.
Definition: reference_counter.h:130
libMesh::QGAUSS
Definition: enum_quadrature_type.h:34
libMesh::ReferenceCounter::get_info
static std::string get_info()
Gets a string containing the reference information.
Definition: reference_counter.C:47
libMesh::QBase::init_3D
virtual void init_3D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:130
libMesh::QBase::init
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the data structures for a quadrature rule for an element of type type.
Definition: quadrature.C:59
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::QBase::build
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
Definition: quadrature_build.C:41
libMesh::FORTYTHIRD
Definition: enum_order.h:84
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::QSIMPSON
Definition: enum_quadrature_type.h:37
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::QBase::_order
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly.
Definition: quadrature.h:345
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::TWENTYTHIRD
Definition: enum_order.h:64
libMesh::QTRAP
Definition: enum_quadrature_type.h:38
libMesh::QNODAL
Definition: enum_quadrature_type.h:46
libMesh::Threads::spin_mtx
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:29
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::QBase::QBase
QBase(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
Definition: quadrature.C:27
libMesh::QBase::_type
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:351
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::QBase::_weights
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:369
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::QGAUSS_LOBATTO
Definition: enum_quadrature_type.h:43
libMesh::QBase::type
virtual QuadratureType type() const =0
libMesh::QJACOBI_1_0
Definition: enum_quadrature_type.h:35
libMesh::QUADSHELL8
Definition: enum_elem_type.h:73
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::QUADSHELL4
Definition: enum_elem_type.h:72
libMesh::QJACOBI_2_0
Definition: enum_quadrature_type.h:36
libMesh::THIRD
Definition: enum_order.h:44
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QCONICAL
Definition: enum_quadrature_type.h:42
libMesh::ReferenceCounter::_enable_print_counter
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called.
Definition: reference_counter.h:141
libMesh::PRISM18
Definition: enum_elem_type.h:52
libMesh::TRISHELL3
Definition: enum_elem_type.h:71
libMesh::QBase::init_0D
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:113
libMesh::QBase::_dim
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:339
libMesh::out
OStreamProxy out
libMesh::QBase::_points
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:363
libMesh::FIRST
Definition: enum_order.h:42
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::QGRUNDMANN_MOLLER
Definition: enum_quadrature_type.h:40
libMesh::QGRID
Definition: enum_quadrature_type.h:39
libMesh::QBase::init_2D
virtual void init_2D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...
Definition: quadrature.C:123