libMesh::QGrundmann_Moller Class Referencefinal

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

#include <quadrature_gm.h>

Inheritance diagram for libMesh::QGrundmann_Moller:

Public Member Functions

 QGrundmann_Moller (unsigned int dim, Order order=INVALID_ORDER)
 
 QGrundmann_Moller (const QGrundmann_Moller &)=default
 
 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)
 
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0)
 
Order get_order () const
 
void print_info (std::ostream &os=libMesh::out) const
 
void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
 
virtual bool shapes_need_reinit ()
 

Static Public Member Functions

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

Public Attributes

bool allow_rules_with_negative_weights
 

Protected Types

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

Protected Member Functions

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

Protected Attributes

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

Static Protected Attributes

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

Private Member Functions

virtual void init_1D (const ElemType, unsigned int=0) override
 
virtual void init_3D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
 
virtual void init_2D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0) override
 
void gm_rule (unsigned int s, unsigned int dim)
 
void compose_all (unsigned int s, unsigned int p, std::vector< std::vector< unsigned int >> &result)
 

Detailed Description

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

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

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  {}
QBase(unsigned int dim, Order order=INVALID_ORDER)
Definition: quadrature.h:364

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

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

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

43 {
44  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
45  _dim,
46  _order);
47 }
unsigned int _dim
Definition: quadrature.h:325
virtual QuadratureType type() const =0
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)

◆ build() [2/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 51 of file quadrature_build.C.

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

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

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

Referenced by gm_rule().

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 }

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

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

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

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ 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 117 of file quadrature.h.

References libMesh::QBase::_type.

117 { return _type; }
ElemType _type
Definition: quadrature.h:337

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & 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 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ get_order()

Order libMesh::QBase::get_order ( ) const
inlineinherited
Returns
The order of the quadrature rule.

Definition at line 203 of file quadrature.h.

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

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

203 { return static_cast<Order>(_order + _p_level); }
unsigned int _p_level
Definition: quadrature.h:343

◆ get_p_level()

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

Definition at line 122 of file quadrature.h.

References libMesh::QBase::_p_level.

122 { return _p_level; }
unsigned int _p_level
Definition: quadrature.h:343

◆ get_points() [1/2]

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

Definition at line 142 of file quadrature.h.

References libMesh::QBase::_points.

Referenced by libMesh::QClough::init_1D(), libMesh::QConical::init_1D(), libMesh::QMonomial::init_1D(), init_1D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QMonomial::init_3D().

142 { return _points; }
std::vector< Point > _points
Definition: quadrature.h:349

◆ get_points() [2/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 148 of file quadrature.h.

References libMesh::QBase::_points.

148 { return _points; }
std::vector< Point > _points
Definition: quadrature.h:349

◆ get_weights() [1/2]

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

Definition at line 154 of file quadrature.h.

References libMesh::QBase::_weights.

Referenced by libMesh::QClough::init_1D(), libMesh::QConical::init_1D(), libMesh::QMonomial::init_1D(), init_1D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QMonomial::init_3D().

154 { return _weights; }
std::vector< Real > _weights
Definition: quadrature.h:355

◆ get_weights() [2/2]

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

Definition at line 160 of file quadrature.h.

References libMesh::QBase::_weights.

160 { return _weights; }
std::vector< Real > _weights
Definition: quadrature.h:355

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

References libMesh::QBase::_points, libMesh::QBase::_weights, compose_all(), std::max(), libMesh::Real, and libMesh::MeshTools::weight().

Referenced by init_2D(), and init_3D().

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 (std::size_t p=0; p<permutations.size(); ++p)
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 (std::size_t j=0; j<permutations.size(); ++j)
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 }
std::vector< Point > _points
Definition: quadrature.h:349
long double max(long double a, double b)
std::vector< Real > _weights
Definition: quadrature.h:355
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:233
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void compose_all(unsigned int s, unsigned int p, std::vector< std::vector< unsigned int >> &result)

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

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

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

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 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

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

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

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

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 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init() [1/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 28 of file quadrature.C.

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::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGaussLobatto::init_2D(), libMesh::QGrid::init_2D(), libMesh::QTrap::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGauss::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::QMonomial::init_3D(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QJacobi::QJacobi(), libMesh::QSimpson::QSimpson(), and libMesh::QTrap::QTrap().

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

◆ init() [2/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 72 of file quadrature.C.

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

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

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

Definition at line 82 of file quadrature.C.

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

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

84 {
85  _points.resize(1);
86  _weights.resize(1);
87  _points[0] = Point(0.);
88  _weights[0] = 1.0;
89 }
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355

◆ init_1D()

void libMesh::QGrundmann_Moller::init_1D ( const ElemType  ,
unsigned int  p = 0 
)
overrideprivatevirtual

In 1D, simply use a Gauss rule.

Implements libMesh::QBase.

Definition at line 38 of file quadrature_gm.C.

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

40 {
41  QGauss gauss1D(1, static_cast<Order>(_order+2*p));
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 }
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355

◆ init_2D()

void libMesh::QGrundmann_Moller::init_2D ( const ElemType  _type = INVALID_ELEM,
unsigned int  p_level = 0 
)
overrideprivatevirtual

Initialize a 2D GM rule. Only makes sense for Tris.

Reimplemented from libMesh::QBase.

Definition at line 28 of file quadrature_gm_2D.C.

References libMesh::QBase::_order, libMesh::QBase::allow_rules_with_negative_weights, gm_rule(), libMesh::TRI3, and libMesh::TRI6.

30 {
31  // Nearly all GM rules contain negative weights, so if you are not
32  // allowing rules with negative weights, we cannot continue!
34  libmesh_error_msg("You requested a Grundmann-Moller rule but\n" \
35  << "are not allowing rules with negative weights!\n" \
36  << "Either select a different quadrature class or\n" \
37  << "set allow_rules_with_negative_weights==true.");
38 
39  switch (type_in)
40  {
41  case TRI3:
42  case TRI6:
43  {
44  switch(_order + 2*p)
45  {
46 
47  default:
48  {
49  // Untested above _order=23 but should work...
50  gm_rule((_order + 2*p)/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: " << type_in);
58  } // end switch (type_in)
59 }
bool allow_rules_with_negative_weights
Definition: quadrature.h:246
void gm_rule(unsigned int s, unsigned int dim)
Definition: quadrature_gm.C:52

◆ init_3D()

void libMesh::QGrundmann_Moller::init_3D ( const ElemType  _type = INVALID_ELEM,
unsigned int  p_level = 0 
)
overrideprivatevirtual

Initialize a 3D GM rule. Only makes sense for Tets.

Reimplemented from libMesh::QBase.

Definition at line 28 of file quadrature_gm_3D.C.

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

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

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

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ n_points()

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

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

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

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

◆ print_info() [2/2]

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

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

Definition at line 378 of file quadrature.h.

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

Referenced by libMesh::operator<<().

379 {
380  libmesh_assert(!_points.empty());
381  libmesh_assert(!_weights.empty());
382 
383  Real summed_weights=0;
384  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
385  for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
386  {
387  os << " Point " << qpoint << ":\n"
388  << " "
389  << _points[qpoint]
390  << "\n Weight:\n "
391  << " w=" << _weights[qpoint] << "\n" << std::endl;
392 
393  summed_weights += _weights[qpoint];
394  }
395  os << "Summed Weights: " << summed_weights << std::endl;
396 }
std::vector< Point > _points
Definition: quadrature.h:349
std::vector< Real > _weights
Definition: quadrature.h:355
unsigned int n_points() const
Definition: quadrature.h:127
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ qp()

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

Definition at line 165 of file quadrature.h.

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().

166  {
167  libmesh_assert_less (i, _points.size());
168  return _points[i];
169  }
std::vector< Point > _points
Definition: quadrature.h:349

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

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

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

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

◆ 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 231 of file quadrature.h.

231 { 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 154 of file quadrature.C.

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().

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

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

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().

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

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

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().

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

◆ type()

QuadratureType libMesh::QGrundmann_Moller::type ( ) const
overridevirtual
Returns
QGRUNDMANN_MOLLER.

Implements libMesh::QBase.

Definition at line 31 of file quadrature_gm.C.

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

◆ _p_level

unsigned int libMesh::QBase::_p_level
protectedinherited

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

Definition at line 343 of file quadrature.h.

Referenced by libMesh::QBase::get_order(), libMesh::QBase::get_p_level(), and libMesh::QBase::init().

◆ _points

◆ _type

ElemType libMesh::QBase::_type
protectedinherited

The type of element for which the current values have been computed.

Definition at line 337 of file quadrature.h.

Referenced by libMesh::QBase::get_elem_type(), and libMesh::QBase::init().

◆ _weights

◆ 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 246 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: