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::QGrundmann_Moller Class Referencefinal

This class implements the Grundmann-Moller quadrature rules for tetrahedra. More...

#include <quadrature_gm.h>

Inheritance diagram for libMesh::QGrundmann_Moller:
[legend]

Public Member Functions

 QGrundmann_Moller (unsigned int dim, Order order=INVALID_ORDER)
 Constructor. More...
 
 QGrundmann_Moller (const QGrundmann_Moller &)=default
 Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class. More...
 
 QGrundmann_Moller (QGrundmann_Moller &&)=default
 
QGrundmann_Molleroperator= (const QGrundmann_Moller &)=default
 
QGrundmann_Molleroperator= (QGrundmann_Moller &&)=default
 
virtual ~QGrundmann_Moller ()=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
 In 1D, use a Gauss rule. 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...
 
void gm_rule (unsigned int s, unsigned int dim)
 This routine is called from init_2D() and init_3D(). More...
 
void compose_all (unsigned int s, unsigned int p, std::vector< std::vector< unsigned int >> &result)
 Routine which generates p-compositions of a given order, s, as well as permutations thereof. More...
 

Detailed Description

This class implements the Grundmann-Moller quadrature rules for tetrahedra.

The GM rules are well-defined for simplices of arbitrary dimension and to any order, but the rules by Dunavant for two-dimensional simplices are in general superior. This is primarily due to the fact that the GM rules contain a significant proportion of negative weights, making them susceptible to round-off error at high-order.

The GM rules are interesting in 3D because they overlap with the conical product rules at higher order while having significantly fewer evaluation points, making them potentially much more efficient. The table below gives a comparison between the number of points in a conical product (CP) rule and the GM rule of equivalent order. The GM rules are defined to be exact for polynomials of degree d=2*s+1, s=0,1,2,3,... The table also gives the percentage of each GM rule's weights which are negative. The percentage of negative weights appears to approach 50, and the amplification factor (a measure of the effect of round-off) defined as

amp. factor = \( \frac{1}{V} \sum \|w_i\|, \)

where V is the volume of the reference element, grows like exp(C*s). (A rule with all positive weights has an amplification factor of 1.0 by definition.)

* s   degree  n_pts(conical)  n_pts(GM)   % neg wts  amp. factor
* ------------------------------------------------------------------------
* 0   1       1               1            0.00      1.00e+00
* 1   3       8               5           20.00      2.60e+00
* 2   5       27              15          26.67      5.63e+00
* 3   7       64              35          31.43      1.19e+01
* 4   9       125             70          34.29      2.54e+01
* 5   11      216             126         36.51      5.41e+01
* 6   13      343             210         38.10      1.16e+02
* 7   15      512             330         39.39      2.51e+02
* 8   17      729             495         40.40      5.45e+02
* 9   19      1000            715         41.26      1.19e+03
* 10  21      1331            1001        41.96      2.59e+03
* 11  23      1728            1365        42.56      5.68e+03
* 12  25      2197            1820        43.08      1.25e+04
* 13  27      2744            2380        43.53      2.75e+04
* 14  29      3375            3060        43.92      6.07e+04
* 15  31      4096            3876        44.27      1.34e+05
* 16  33      4913            4845        44.58      2.97e+05
* 17  35      5832            5985        44.86      6.59e+05 <= Conical rule has fewer points for degree >= 34
* 18  37      6859            7315        45.11      1.46e+06
* 19  39      8000            8855        45.34      3.25e+06
* 20  41      9261            10626       45.55      7.23e+06
* 21  43      10648           12650       45.74      1.61e+07
* 

Reference: Axel Grundmann and Michael M"{o}ller, "Invariant Integration Formulas for the N-Simplex by Combinatorial Methods," SIAM Journal on Numerical Analysis, Volume 15, Number 2, April 1978, pages 282-290.

Reference LGPL Fortran90 code by John Burkardt can be found here: http://people.scs.fsu.edu/~burkardt/f_src/gm_rules/gm_rules.html

Author
John W. Peterson
Date
2008

Implements the quadrature rules of Grundmann and Moller in 2D and 3D.

Definition at line 96 of file quadrature_gm.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

◆ QGrundmann_Moller() [1/3]

libMesh::QGrundmann_Moller::QGrundmann_Moller ( unsigned int  dim,
Order  order = INVALID_ORDER 
)
inline

Constructor.

Declares the order of the quadrature rule.

Definition at line 103 of file quadrature_gm.h.

104  :
105  QBase(dim,order)
106  {}

◆ QGrundmann_Moller() [2/3]

libMesh::QGrundmann_Moller::QGrundmann_Moller ( const QGrundmann_Moller )
default

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

◆ QGrundmann_Moller() [3/3]

libMesh::QGrundmann_Moller::QGrundmann_Moller ( QGrundmann_Moller &&  )
default

◆ ~QGrundmann_Moller()

virtual libMesh::QGrundmann_Moller::~QGrundmann_Moller ( )
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().

◆ compose_all()

void libMesh::QGrundmann_Moller::compose_all ( unsigned int  s,
unsigned int  p,
std::vector< std::vector< unsigned int >> &  result 
)
private

Routine which generates p-compositions of a given order, s, as well as permutations thereof.

This routine is called internally by the gm_rule() routine, you should not call this yourself!

Definition at line 150 of file quadrature_gm.C.

153 {
154  // Clear out results remaining from previous calls
155  result.clear();
156 
157  // Allocate storage for a workspace. The workspace will periodically
158  // be copied into the result container.
159  std::vector<unsigned int> workspace(p);
160 
161  // The first result is always (s,0,...,0)
162  workspace[0] = s;
163  result.push_back(workspace);
164 
165  // the value of the first non-zero entry
166  unsigned int head_value=s;
167 
168  // When head_index=-1, it refers to "off the front" of the array. Therefore,
169  // this needs to be a regular int rather than unsigned. I initially tried to
170  // do this with head_index unsigned and an else statement below, but then there
171  // is the special case: (1,0,...,0) which does not work correctly.
172  int head_index = -1;
173 
174  // At the end, all the entries will be in the final slot of workspace
175  while (workspace.back() != s)
176  {
177  // If the previous head value is still larger than 1, reset the index
178  // to "off the front" of the array
179  if (head_value > 1)
180  head_index = -1;
181 
182  // Either move the index onto the front of the array or on to
183  // the next value.
184  head_index++;
185 
186  // Get current value of the head entry
187  head_value = workspace[head_index];
188 
189  // Put a zero into the head_index of the array. If head_index==0,
190  // this will be overwritten in the next line with head_value-1.
191  workspace[head_index] = 0;
192 
193  // The initial entry gets the current head value, minus 1.
194  // If head_value > 1, the next loop iteration will start back
195  // at workspace[0] again.
196  libmesh_assert_greater (head_value, 0);
197  workspace[0] = head_value - 1;
198 
199  // Increment the head+1 value
200  workspace[head_index+1] += 1;
201 
202  // Save this composition in the results
203  result.push_back(workspace);
204  }
205 }

Referenced by gm_rule().

◆ 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(), init_1D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), init_2D(), libMesh::QGauss::init_3D(), libMesh::QMonomial::init_3D(), and 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

◆ gm_rule()

void libMesh::QGrundmann_Moller::gm_rule ( unsigned int  s,
unsigned int  dim 
)
private

This routine is called from init_2D() and init_3D().

It actually fills the _points and _weights vectors for a given rule index, s and dimension, dim.

Definition at line 52 of file quadrature_gm.C.

53 {
54  // Only dim==2 or dim==3 is allowed
55  libmesh_assert(dim==2 || dim==3);
56 
57  // A GM rule of index s can integrate polynomials of degree 2*s+1 exactly
58  const unsigned int degree = 2*s+1;
59 
60  // The number of points for rule of index s is
61  // (dim+1+s)! / (dim+1)! / s!
62  // In 3D, this is = 1/24 * Pi_{i=1}^4 (s+i)
63  // In 2D, this is = 1/6 * Pi_{i=1}^3 (s+i)
64  const unsigned int n_pts = dim==2 ? (s+3)*(s+2)*(s+1) / 6 : (s+4)*(s+3)*(s+2)*(s+1) / 24;
65  //libMesh::out << "n_pts=" << n_pts << std::endl;
66 
67  // Allocate space for points and weights
68  _points.resize(n_pts);
69  _weights.resize(n_pts);
70 
71  // (-1)^i -> This one flips sign at each iteration of the i-loop below.
72  int one_pm=1;
73 
74  // Where we store all the integer point compositions/permutations
75  std::vector<std::vector<unsigned int>> permutations;
76 
77  // Index into the vector where we should start adding the next round of points/weights
78  std::size_t offset=0;
79 
80  // Implement the GM formula 4.1 on page 286 of the paper
81  for (unsigned int i=0; i<=s; ++i)
82  {
83  // Get all the ordered compositions (and their permutations)
84  // of |beta| = s-i into dim+1 parts
85  compose_all(s-i, dim+1, permutations);
86  //libMesh::out << "n. permutations=" << permutations.size() << std::endl;
87 
88  for (auto p : index_range(permutations))
89  {
90  // We use the first dim entries of each permutation to
91  // construct an integration point.
92  for (unsigned int j=0; j<dim; ++j)
93  _points[offset+p](j) =
94  static_cast<Real>(2.*permutations[p][j] + 1.) /
95  static_cast<Real>( degree + dim - 2.*i );
96  }
97 
98  // Compute the weight for this i, being careful to avoid overflow.
99  // This technique is borrowed from Burkardt's code as well.
100  // Use once for each of the points obtained from the permutations array.
101  Real weight = one_pm;
102 
103  // This for loop needs to run for dim, degree, or dim+degree-i iterations,
104  // whichever is largest.
105  const unsigned int weight_loop_index =
106  std::max(dim, std::max(degree, degree+dim-i))+1;
107 
108  for (unsigned int j=1; j<weight_loop_index; ++j)
109  {
110  if (j <= degree) // Accumulate (d+n-2i)^d term
111  weight *= static_cast<Real>(degree+dim-2*i);
112 
113  if (j <= 2*s) // Accumulate 2^{-2s}
114  weight *= 0.5;
115 
116  if (j <= i) // Accumulate (i!)^{-1}
117  weight /= static_cast<Real>(j);
118 
119  if (j <= degree+dim-i) // Accumulate ( (d+n-i)! )^{-1}
120  weight /= static_cast<Real>(j);
121  }
122 
123  // This is the weight for each of the points computed previously
124  for (auto j : index_range(permutations))
125  _weights[offset+j] = weight;
126 
127  // Change sign for next iteration
128  one_pm = -one_pm;
129 
130  // Update offset for the next set of points
131  offset += permutations.size();
132  }
133 }

References libMesh::QBase::_points, libMesh::QBase::_weights, compose_all(), dim, libMesh::index_range(), libMesh::libmesh_assert(), libMesh::Real, and libMesh::MeshTools::weight().

Referenced by init_2D(), and init_3D().

◆ 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(), libMesh::QNodal::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QGrid::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QNodal::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGrid::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGauss::init_3D(), libMesh::QNodal::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QJacobi::QJacobi(), libMesh::QNodal::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::QGrundmann_Moller::init_1D ( const  ElemType,
unsigned int   
)
overrideprivatevirtual

In 1D, use a Gauss rule.

In 2D, the GM product rule is only defined for Tris. In 3D, the GM product rule is only defined for Tets.

Implements libMesh::QBase.

Definition at line 39 of file quadrature_gm.C.

40 {
41  QGauss gauss1D(1, get_order());
42 
43  // Swap points and weights with the about-to-be destroyed rule.
44  _points.swap(gauss1D.get_points());
45  _weights.swap(gauss1D.get_weights());
46 }

References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::get_order(), libMesh::QBase::get_points(), and libMesh::QBase::get_weights().

◆ init_2D()

void libMesh::QGrundmann_Moller::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 27 of file quadrature_gm_2D.C.

28 {
29  // Nearly all GM rules contain negative weights, so if you are not
30  // allowing rules with negative weights, we cannot continue!
32  libmesh_error_msg("You requested a Grundmann-Moller rule but\n" \
33  << "are not allowing rules with negative weights!\n" \
34  << "Either select a different quadrature class or\n" \
35  << "set allow_rules_with_negative_weights==true.");
36 
37  switch (_type)
38  {
39  case TRI3:
40  case TRISHELL3:
41  case TRI3SUBDIVISION:
42  case TRI6:
43  {
44  switch(get_order())
45  {
46 
47  default:
48  {
49  // Untested above _order=23 but should work...
50  gm_rule(get_order()/2, /*dim=*/2);
51  return;
52  }
53  } // end switch (order)
54  } // end case TRI3, TRI6
55 
56  default:
57  libmesh_error_msg("ERROR: Unsupported element type: " << Utility::enum_to_string(_type));
58  } // end switch (_type)
59 }

References libMesh::QBase::_type, libMesh::QBase::allow_rules_with_negative_weights, libMesh::Utility::enum_to_string(), libMesh::QBase::get_order(), gm_rule(), libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, and libMesh::TRISHELL3.

◆ init_3D()

void libMesh::QGrundmann_Moller::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 27 of file quadrature_gm_3D.C.

28 {
29  // Nearly all GM rules contain negative weights, so if you are not
30  // allowing rules with negative weights, we cannot continue!
32  libmesh_error_msg("You requested a Grundmann-Moller rule but\n" \
33  << "are not allowing rules with negative weights!\n" \
34  << "Either select a different quadrature class or\n" \
35  << "set allow_rules_with_negative_weights==true.");
36 
37  switch (_type)
38  {
39  case TET4:
40  case TET10:
41  {
42  switch(get_order())
43  {
44  // We hard-code the first few orders based on output from
45  // the mp-quadrature library:
46  //
47  // https://code.google.com/p/mp-quadrature
48  //
49  // The points are ordered in such a way that the amount of
50  // round-off error in the quadrature calculations is
51  // (hopefully) minimized. These orderings were obtained
52  // via a simple permutation optimization strategy designed
53  // to produce the smallest overall average error while
54  // integrating all polynomials of degree <= d.
55  case CONSTANT:
56  case FIRST:
57  {
58  _points.resize(1);
59  _weights.resize(1);
60 
61  _points[0](0) = Real(1)/4;
62  _points[0](1) = Real(1)/4;
63  _points[0](2) = Real(1)/4;
64 
65  _weights[0] = Real(1)/6;
66  return;
67  }
68 
69  case SECOND:
70  case THIRD:
71  {
72  _points.resize(5);
73  _weights.resize(5);
74 
75  _points[0](0) = Real(1)/6; _points[0](1) = Real(1)/6; _points[0](2) = Real(1)/2; _weights[0] = Real(3)/40;
76  _points[1](0) = Real(1)/6; _points[1](1) = Real(1)/6; _points[1](2) = Real(1)/6; _weights[1] = Real(3)/40;
77  _points[2](0) = Real(1)/4; _points[2](1) = Real(1)/4; _points[2](2) = Real(1)/4; _weights[2] = -Real(2)/15;
78  _points[3](0) = Real(1)/2; _points[3](1) = Real(1)/6; _points[3](2) = Real(1)/6; _weights[3] = Real(3)/40;
79  _points[4](0) = Real(1)/6; _points[4](1) = Real(1)/2; _points[4](2) = Real(1)/6; _weights[4] = Real(3)/40;
80  return;
81  }
82 
83  case FOURTH:
84  case FIFTH:
85  {
86  _points.resize(15);
87  _weights.resize(15);
88 
89  _points[ 0](0) = Real(1)/8; _points[ 0](1) = Real(3)/8; _points[ 0](2) = Real(1)/8; _weights[ 0] = Real(16)/315;
90  _points[ 1](0) = Real(1)/8; _points[ 1](1) = Real(1)/8; _points[ 1](2) = Real(5)/8; _weights[ 1] = Real(16)/315;
91  _points[ 2](0) = Real(3)/8; _points[ 2](1) = Real(1)/8; _points[ 2](2) = Real(1)/8; _weights[ 2] = Real(16)/315;
92  _points[ 3](0) = Real(1)/6; _points[ 3](1) = Real(1)/2; _points[ 3](2) = Real(1)/6; _weights[ 3] = -Real(27)/280;
93  _points[ 4](0) = Real(3)/8; _points[ 4](1) = Real(1)/8; _points[ 4](2) = Real(3)/8; _weights[ 4] = Real(16)/315;
94  _points[ 5](0) = Real(1)/8; _points[ 5](1) = Real(3)/8; _points[ 5](2) = Real(3)/8; _weights[ 5] = Real(16)/315;
95  _points[ 6](0) = Real(1)/6; _points[ 6](1) = Real(1)/6; _points[ 6](2) = Real(1)/6; _weights[ 6] = -Real(27)/280;
96  _points[ 7](0) = Real(1)/6; _points[ 7](1) = Real(1)/6; _points[ 7](2) = Real(1)/2; _weights[ 7] = -Real(27)/280;
97  _points[ 8](0) = Real(1)/8; _points[ 8](1) = Real(1)/8; _points[ 8](2) = Real(1)/8; _weights[ 8] = Real(16)/315;
98  _points[ 9](0) = Real(1)/4; _points[ 9](1) = Real(1)/4; _points[ 9](2) = Real(1)/4; _weights[ 9] = Real(2)/45;
99  _points[10](0) = Real(1)/8; _points[10](1) = Real(5)/8; _points[10](2) = Real(1)/8; _weights[10] = Real(16)/315;
100  _points[11](0) = Real(1)/2; _points[11](1) = Real(1)/6; _points[11](2) = Real(1)/6; _weights[11] = -Real(27)/280;
101  _points[12](0) = Real(1)/8; _points[12](1) = Real(1)/8; _points[12](2) = Real(3)/8; _weights[12] = Real(16)/315;
102  _points[13](0) = Real(5)/8; _points[13](1) = Real(1)/8; _points[13](2) = Real(1)/8; _weights[13] = Real(16)/315;
103  _points[14](0) = Real(3)/8; _points[14](1) = Real(3)/8; _points[14](2) = Real(1)/8; _weights[14] = Real(16)/315;
104  return;
105  }
106 
107  case SIXTH:
108  case SEVENTH:
109  {
110  _points.resize(35);
111  _weights.resize(35);
112 
113  _points[ 0](0) = Real(3)/10; _points[ 0](1) = Real(1)/10; _points[ 0](2) = Real(3)/10; _weights[ 0] = Real(3125)/72576;
114  _points[ 1](0) = Real(1)/6; _points[ 1](1) = Real(1)/2; _points[ 1](2) = Real(1)/6; _weights[ 1] = Real(243)/4480;
115  _points[ 2](0) = Real(1)/6; _points[ 2](1) = Real(1)/6; _points[ 2](2) = Real(1)/2; _weights[ 2] = Real(243)/4480;
116  _points[ 3](0) = Real(1)/2; _points[ 3](1) = Real(1)/10; _points[ 3](2) = Real(1)/10; _weights[ 3] = Real(3125)/72576;
117  _points[ 4](0) = Real(3)/10; _points[ 4](1) = Real(1)/10; _points[ 4](2) = Real(1)/10; _weights[ 4] = Real(3125)/72576;
118  _points[ 5](0) = Real(3)/10; _points[ 5](1) = Real(3)/10; _points[ 5](2) = Real(1)/10; _weights[ 5] = Real(3125)/72576;
119  _points[ 6](0) = Real(1)/2; _points[ 6](1) = Real(1)/6; _points[ 6](2) = Real(1)/6; _weights[ 6] = Real(243)/4480;
120  _points[ 7](0) = Real(3)/10; _points[ 7](1) = Real(1)/10; _points[ 7](2) = Real(1)/2; _weights[ 7] = Real(3125)/72576;
121  _points[ 8](0) = Real(7)/10; _points[ 8](1) = Real(1)/10; _points[ 8](2) = Real(1)/10; _weights[ 8] = Real(3125)/72576;
122  _points[ 9](0) = Real(3)/8; _points[ 9](1) = Real(1)/8; _points[ 9](2) = Real(1)/8; _weights[ 9] = -Real(256)/2835;
123  _points[10](0) = Real(3)/10; _points[10](1) = Real(3)/10; _points[10](2) = Real(3)/10; _weights[10] = Real(3125)/72576;
124  _points[11](0) = Real(1)/10; _points[11](1) = Real(1)/2; _points[11](2) = Real(3)/10; _weights[11] = Real(3125)/72576;
125  _points[12](0) = Real(1)/10; _points[12](1) = Real(1)/10; _points[12](2) = Real(7)/10; _weights[12] = Real(3125)/72576;
126  _points[13](0) = Real(1)/2; _points[13](1) = Real(1)/10; _points[13](2) = Real(3)/10; _weights[13] = Real(3125)/72576;
127  _points[14](0) = Real(1)/10; _points[14](1) = Real(7)/10; _points[14](2) = Real(1)/10; _weights[14] = Real(3125)/72576;
128  _points[15](0) = Real(1)/10; _points[15](1) = Real(1)/2; _points[15](2) = Real(1)/10; _weights[15] = Real(3125)/72576;
129  _points[16](0) = Real(1)/6; _points[16](1) = Real(1)/6; _points[16](2) = Real(1)/6; _weights[16] = Real(243)/4480;
130  _points[17](0) = Real(3)/8; _points[17](1) = Real(1)/8; _points[17](2) = Real(3)/8; _weights[17] = -Real(256)/2835;
131  _points[18](0) = Real(1)/8; _points[18](1) = Real(1)/8; _points[18](2) = Real(5)/8; _weights[18] = -Real(256)/2835;
132  _points[19](0) = Real(1)/10; _points[19](1) = Real(1)/10; _points[19](2) = Real(3)/10; _weights[19] = Real(3125)/72576;
133  _points[20](0) = Real(1)/8; _points[20](1) = Real(3)/8; _points[20](2) = Real(3)/8; _weights[20] = -Real(256)/2835;
134  _points[21](0) = Real(5)/8; _points[21](1) = Real(1)/8; _points[21](2) = Real(1)/8; _weights[21] = -Real(256)/2835;
135  _points[22](0) = Real(1)/8; _points[22](1) = Real(5)/8; _points[22](2) = Real(1)/8; _weights[22] = -Real(256)/2835;
136  _points[23](0) = Real(1)/10; _points[23](1) = Real(3)/10; _points[23](2) = Real(1)/10; _weights[23] = Real(3125)/72576;
137  _points[24](0) = Real(1)/4; _points[24](1) = Real(1)/4; _points[24](2) = Real(1)/4; _weights[24] = -Real(8)/945;
138  _points[25](0) = Real(1)/8; _points[25](1) = Real(1)/8; _points[25](2) = Real(3)/8; _weights[25] = -Real(256)/2835;
139  _points[26](0) = Real(3)/8; _points[26](1) = Real(3)/8; _points[26](2) = Real(1)/8; _weights[26] = -Real(256)/2835;
140  _points[27](0) = Real(1)/8; _points[27](1) = Real(3)/8; _points[27](2) = Real(1)/8; _weights[27] = -Real(256)/2835;
141  _points[28](0) = Real(1)/10; _points[28](1) = Real(3)/10; _points[28](2) = Real(1)/2; _weights[28] = Real(3125)/72576;
142  _points[29](0) = Real(3)/10; _points[29](1) = Real(1)/2; _points[29](2) = Real(1)/10; _weights[29] = Real(3125)/72576;
143  _points[30](0) = Real(1)/10; _points[30](1) = Real(1)/10; _points[30](2) = Real(1)/2; _weights[30] = Real(3125)/72576;
144  _points[31](0) = Real(1)/2; _points[31](1) = Real(3)/10; _points[31](2) = Real(1)/10; _weights[31] = Real(3125)/72576;
145  _points[32](0) = Real(1)/8; _points[32](1) = Real(1)/8; _points[32](2) = Real(1)/8; _weights[32] = -Real(256)/2835;
146  _points[33](0) = Real(1)/10; _points[33](1) = Real(3)/10; _points[33](2) = Real(3)/10; _weights[33] = Real(3125)/72576;
147  _points[34](0) = Real(1)/10; _points[34](1) = Real(1)/10; _points[34](2) = Real(1)/10; _weights[34] = Real(3125)/72576;
148  return;
149  }
150 
151  default:
152  {
153  // Untested above _order=23 but should work...
154  gm_rule(get_order()/2, /*dim=*/3);
155  return;
156  }
157  } // end switch (order)
158  } // end case TET4, TET10
159 
160  default:
161  libmesh_error_msg("ERROR: Unsupported element type: " << Utility::enum_to_string(_type));
162  } // end switch (_type)
163 }

References libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::QBase::allow_rules_with_negative_weights, libMesh::CONSTANT, libMesh::Utility::enum_to_string(), libMesh::FIFTH, libMesh::FIRST, libMesh::FOURTH, libMesh::QBase::get_order(), gm_rule(), libMesh::Real, libMesh::SECOND, libMesh::SEVENTH, libMesh::SIXTH, libMesh::TET10, libMesh::TET4, and libMesh::THIRD.

◆ 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(), A0::interior_assembly(), B::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]

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

◆ operator=() [2/2]

QGrundmann_Moller& libMesh::QGrundmann_Moller::operator= ( QGrundmann_Moller &&  )
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::QGrundmann_Moller::type ( ) const
overridevirtual
Returns
QGRUNDMANN_MOLLER.

Implements libMesh::QBase.

Definition at line 32 of file quadrature_gm.C.

33 {
34  return QGRUNDMANN_MOLLER;
35 }

References libMesh::QGRUNDMANN_MOLLER.

◆ 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 init_2D(), libMesh::QGauss::init_3D(), libMesh::QMonomial::init_3D(), and init_3D().


The documentation for this class was generated from the following files:
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::SIXTH
Definition: enum_order.h:47
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::FIFTH
Definition: enum_order.h:46
libMesh::SECOND
Definition: enum_order.h:43
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::QGrundmann_Moller::gm_rule
void gm_rule(unsigned int s, unsigned int dim)
This routine is called from init_2D() and init_3D().
Definition: quadrature_gm.C:52
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::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::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::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::CONSTANT
Definition: enum_order.h:41
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::FOURTH
Definition: enum_order.h:45
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::SEVENTH
Definition: enum_order.h:48
libMesh::QBase::allow_rules_with_negative_weights
bool allow_rules_with_negative_weights
Flag (default true) controlling the use of quadrature rules with negative weights.
Definition: quadrature.h:254
libMesh::QJACOBI_2_0
Definition: enum_quadrature_type.h:36
libMesh::THIRD
Definition: enum_order.h:44
libMesh::TRI3SUBDIVISION
Definition: enum_elem_type.h:69
libMesh::QGrundmann_Moller::compose_all
void compose_all(unsigned int s, unsigned int p, std::vector< std::vector< unsigned int >> &result)
Routine which generates p-compositions of a given order, s, as well as permutations thereof.
Definition: quadrature_gm.C:150
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QBase::get_order
Order get_order() const
Definition: quadrature.h:211
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::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::MeshTools::weight
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:236
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::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