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

This class implements Gauss-Lobatto quadrature for 1D elements and 2D/3D tensor product elements. More...

#include <quadrature_gauss_lobatto.h>

Inheritance diagram for libMesh::QGaussLobatto:
[legend]

Public Member Functions

 QGaussLobatto (unsigned int dim, Order order=INVALID_ORDER)
 Constructor. More...
 
 QGaussLobatto (const QGaussLobatto &)=default
 Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class. More...
 
 QGaussLobatto (QGaussLobatto &&)=default
 
QGaussLobattooperator= (const QGaussLobatto &)=default
 
QGaussLobattooperator= (QGaussLobatto &&)=default
 
virtual ~QGaussLobatto ()=default
 
virtual QuadratureType type () const override
 
virtual std::unique_ptr< QBaseclone () const override
 
ElemType get_elem_type () const
 
unsigned int get_p_level () const
 
unsigned int n_points () const
 
unsigned int size () const
 Alias for n_points() to enable use in index_range. More...
 
unsigned int get_dim () const
 
const std::vector< Point > & get_points () const
 
std::vector< Point > & get_points ()
 
const std::vector< Real > & get_weights () const
 
std::vector< Real > & get_weights ()
 
Point qp (const unsigned int i) const
 
Real w (const unsigned int i) const
 
virtual void init (const Elem &e, unsigned int p_level=invalid_uint)
 Initializes the data structures for a quadrature rule for the element e. More...
 
virtual void init (const ElemType type=INVALID_ELEM, unsigned int p_level=0, bool simple_type_only=false)
 Initializes the data structures for a quadrature rule for an element of type type. More...
 
virtual void init (const QBase &other_rule)
 Initializes the data structures for a quadrature rule based on the element, element type, and p_level settings of other_rule. More...
 
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0)
 Initializes the data structures for an element potentially "cut" by a signed distance function. More...
 
Order get_order () const
 
Order get_base_order () const
 
void print_info (std::ostream &os=libMesh::out) const
 Prints information relevant to the quadrature rule, by default to libMesh::out. More...
 
void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
 Maps the points of a 1D quadrature rule defined by "old_range" to another 1D interval defined by "new_range" and scales the weights accordingly. More...
 
virtual bool shapes_need_reinit ()
 

Static Public Member Functions

static std::unique_ptr< QBasebuild (std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
 Builds a specific quadrature rule based on the name string. More...
 
static std::unique_ptr< QBasebuild (const QuadratureType qt, const unsigned int dim, const Order order=INVALID_ORDER)
 Builds a specific quadrature rule based on the QuadratureType. More...
 
static void print_info (std::ostream &out_stream=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...
 
bool allow_nodal_pyramid_quadrature
 The flag's value defaults to false so that one does not accidentally use a nodal quadrature rule on Pyramid elements, since evaluating the inverse element Jacobian (e.g. More...
 

Protected Types

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

Protected Member Functions

virtual void init_0D ()
 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) noexcept
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name) noexcept
 Increments the destruction counter. More...
 

Protected Attributes

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

Static Protected Attributes

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

Private Member Functions

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

Detailed Description

This class implements Gauss-Lobatto quadrature for 1D elements and 2D/3D tensor product elements.

Properties of Gauss-Lobatto quadrature rules: .) Include the "end-points" of the domain (have points on edges/faces in 2D/3D). .) Rules with n points can exactly integrate polynomials of degree 2n-3.

http://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Lobatto_rules

Author
John W. Peterson
Date
2014 Implements 1D and 2/3D tensor product Gauss-Lobatto quadrature rules.

Definition at line 41 of file quadrature_gauss_lobatto.h.

Member Typedef Documentation

◆ Counts

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

Data structure to log the information.

The log is identified by the class name.

Definition at line 119 of file reference_counter.h.

Constructor & Destructor Documentation

◆ QGaussLobatto() [1/3]

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

Constructor.

Declares the order of the quadrature rule.

Definition at line 33 of file quadrature_gauss_lobatto.C.

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

34  : QBase(d,o)
35 {
36  // explicitly call the init function in 1D since the
37  // other tensor-product rules require this one.
38  // note that EDGE will not be used internally, however
39  // if we called the function with INVALID_ELEM it would try to
40  // be smart and return, thinking it had already done the work.
41  if (_dim == 1)
42  init(EDGE2);
43 }
QBase(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
Definition: quadrature.C:27
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e.
Definition: quadrature.C:65

◆ QGaussLobatto() [2/3]

libMesh::QGaussLobatto::QGaussLobatto ( const QGaussLobatto )
default

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

◆ QGaussLobatto() [3/3]

libMesh::QGaussLobatto::QGaussLobatto ( QGaussLobatto &&  )
default

◆ ~QGaussLobatto()

virtual libMesh::QGaussLobatto::~QGaussLobatto ( )
virtualdefault

Member Function Documentation

◆ build() [1/2]

std::unique_ptr< QBase > libMesh::QBase::build ( std::string_view  name,
const unsigned int  dim,
const Order  order = INVALID_ORDER 
)
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 43 of file quadrature_build.C.

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::QBase::clone(), main(), libMesh::OverlapCoupling::OverlapCoupling(), libMesh::InfFE< Dim, T_radial, T_map >::reinit(), QuadratureTest::testBuild(), QuadratureTest::testJacobi(), and QuadratureTest::testPolynomials().

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

◆ build() [2/2]

std::unique_ptr< QBase > libMesh::QBase::build ( const QuadratureType  qt,
const unsigned int  dim,
const Order  order = INVALID_ORDER 
)
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 54 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::QNODAL, libMesh::QSIMPSON, libMesh::QTRAP, libMesh::THIRD, and libMesh::TWENTYTHIRD.

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

◆ clone()

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

Reimplemented from libMesh::QBase.

Definition at line 51 of file quadrature_gauss_lobatto.C.

52 {
53  return std::make_unique<QGaussLobatto>(*this);
54 }

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 94 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ get_base_order()

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

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

Definition at line 258 of file quadrature.h.

References libMesh::QBase::_order.

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

◆ get_dim()

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

Definition at line 150 of file quadrature.h.

References libMesh::QBase::_dim.

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

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

◆ get_elem_type()

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

Definition at line 121 of file quadrature.h.

References libMesh::QBase::_type.

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

◆ get_info()

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

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & [name, cd] : _counts)
59  oss << "| " << name << " reference count information:\n"
60  << "| Creations: " << cd.first << '\n'
61  << "| Destructions: " << cd.second << '\n';
62 
63  oss << " ---------------------------------------------------------------------------- \n";
64 
65  return oss.str();
66 
67 #else
68 
69  return "";
70 
71 #endif
72 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
static Counts _counts
Actually holds the data.

◆ get_order()

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

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

Definition at line 249 of file quadrature.h.

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

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

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

◆ get_p_level()

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

Definition at line 126 of file quadrature.h.

References libMesh::QBase::_p_level.

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

◆ get_points() [1/2]

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

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

References libMesh::QBase::_points.

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

◆ get_weights() [1/2]

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

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

References libMesh::QBase::_weights.

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

◆ increment_constructor_count()

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

Increments the construction counter.

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

Definition at line 183 of file reference_counter.h.

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

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

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

◆ increment_destructor_count()

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

Increments the destruction counter.

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

Definition at line 207 of file reference_counter.h.

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

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

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

◆ init() [1/4]

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

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

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

Definition at line 65 of file quadrature.C.

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

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

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

◆ init() [2/4]

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

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

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

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

Definition at line 121 of file quadrature.C.

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

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

◆ init() [3/4]

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

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

Definition at line 178 of file quadrature.C.

References libMesh::QBase::_elem, libMesh::QBase::_p_level, libMesh::QBase::_type, and libMesh::QBase::init().

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

◆ init() [4/4]

void libMesh::QBase::init ( const Elem elem,
const std::vector< Real > &  vertex_distance_func,
unsigned int  p_level = 0 
)
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 188 of file quadrature.C.

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

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

◆ init_0D()

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

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

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

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

◆ init_1D()

void libMesh::QGaussLobatto::init_1D ( )
overrideprivatevirtual

Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate values.

The order of the rule will be defined by the implementing class. It is assumed that derived quadrature rules will at least define the init_1D function, therefore it is pure virtual.

Implements libMesh::QBase.

Definition at line 28 of file quadrature_gauss_lobatto_1D.C.

References libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::CONSTANT, libMesh::EIGHTH, libMesh::EIGHTTEENTH, libMesh::ELEVENTH, libMesh::FIFTEENTH, libMesh::FIFTH, libMesh::FIRST, libMesh::FORTIETH, libMesh::FORTYFIRST, libMesh::FORTYSECOND, libMesh::FORTYTHIRD, libMesh::FOURTEENTH, libMesh::FOURTH, libMesh::QBase::get_order(), libMesh::NINETEENTH, libMesh::NINTH, libMesh::Real, libMesh::SECOND, libMesh::SEVENTEENTH, libMesh::SEVENTH, libMesh::SIXTEENTH, libMesh::SIXTH, libMesh::TENTH, libMesh::THIRD, libMesh::THIRTEENTH, libMesh::THIRTIETH, libMesh::THIRTYEIGHTH, libMesh::THIRTYFIFTH, libMesh::THIRTYFIRST, libMesh::THIRTYFOURTH, libMesh::THIRTYNINTH, libMesh::THIRTYSECOND, libMesh::THIRTYSEVENTH, libMesh::THIRTYSIXTH, libMesh::THIRTYTHIRD, libMesh::TWELFTH, libMesh::TWENTIETH, libMesh::TWENTYEIGHTH, libMesh::TWENTYFIFTH, libMesh::TWENTYFIRST, libMesh::TWENTYFOURTH, libMesh::TWENTYNINTH, libMesh::TWENTYSECOND, libMesh::TWENTYSEVENTH, libMesh::TWENTYSIXTH, and libMesh::TWENTYTHIRD.

29 {
30  //----------------------------------------------------------------------
31  // 1D quadrature rules
32  switch(get_order())
33  {
34  // Since Gauss-Lobatto rules must include the endpoints of the
35  // domain, there is no 1-point rule. The two-point
36  // Gauss-Lobatto rule is equivalent to the trapezoidal rule.
37  case CONSTANT:
38  case FIRST:
39  {
40  _points.resize (2);
41  _weights.resize(2);
42 
43  _points[0](0) = -1;
44  _points[1] = -_points[0];
45 
46  _weights[0] = 1.;
47  _weights[1] = _weights[0];
48 
49  return;
50  }
51 
52  // The three-point Gauss-Lobatto rule is equivalent to Simpsons' rule.
53  // It can integrate cubic polynomials exactly.
54  case SECOND:
55  case THIRD:
56  {
57  _points.resize (3);
58  _weights.resize(3);
59 
60  _points[0](0) = -1;
61  _points[1] = 0;
62  _points[2] = -_points[0];
63 
64  _weights[0] = Real(1)/3;
65  _weights[1] = Real(4)/3;
66  _weights[2] = _weights[0];
67  return;
68  }
69 
70  // The four-point Gauss-Lobatto rule can integrate 2*4-3 =
71  // 5th-order polynomials exactly.
72  case FOURTH:
73  case FIFTH:
74  {
75  _points.resize (4);
76  _weights.resize(4);
77 
78  _points[ 0](0) = -1;
79  _points[ 1](0) = -std::sqrt(Real(1)/5); // ~ -0.4472135955
80  _points[ 2] = -_points[1];
81  _points[ 3] = -_points[0];
82 
83  _weights[ 0] = Real(1)/6;
84  _weights[ 1] = Real(5)/6;
85  _weights[ 2] = _weights[1];
86  _weights[ 3] = _weights[0];
87 
88  return;
89  }
90 
91  // The five-point Gauss-Lobatto rule can integrate 2*5-3 =
92  // 7th-order polynomials exactly.
93  case SIXTH:
94  case SEVENTH:
95  {
96  _points.resize (5);
97  _weights.resize(5);
98 
99  _points[ 0](0) = -1;
100  _points[ 1](0) = -std::sqrt(Real(3)/7); // ~ -0.6546536707
101  _points[ 2](0) = 0.;
102  _points[ 3] = -_points[1];
103  _points[ 4] = -_points[0];
104 
105  _weights[ 0] = 1/Real(10);
106  _weights[ 1] = Real(49)/90;
107  _weights[ 2] = Real(32)/45;
108  _weights[ 3] = _weights[1];
109  _weights[ 4] = _weights[0];
110 
111  return;
112  }
113 
114  // 2*6-3 = 9
115  case EIGHTH:
116  case NINTH:
117  {
118  _points.resize (6);
119  _weights.resize(6);
120 
121  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
122  _points[ 1](0) = -7.6505532392946469285100297395934e-01_R;
123  _points[ 2](0) = -2.8523151648064509631415099404088e-01_R;
124  _points[ 3] = -_points[2];
125  _points[ 4] = -_points[1];
126  _points[ 5] = -_points[0];
127 
128  _weights[ 0] = 6.6666666666666666666666666666667e-02_R;
129  _weights[ 1] = 3.7847495629784698031661280821202e-01_R;
130  _weights[ 2] = 5.5485837703548635301672052512131e-01_R;
131  _weights[ 3] = _weights[2];
132  _weights[ 4] = _weights[1];
133  _weights[ 5] = _weights[0];
134 
135  return;
136  }
137 
138  // 2*7-3 = 11
139  case TENTH:
140  case ELEVENTH:
141  {
142  _points.resize (7);
143  _weights.resize(7);
144 
145  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
146  _points[ 1](0) = -8.3022389627856692987203221396747e-01_R;
147  _points[ 2](0) = -4.6884879347071421380377188190877e-01_R;
148  _points[ 3](0) = 0.;
149  _points[ 4] = -_points[2];
150  _points[ 5] = -_points[1];
151  _points[ 6] = -_points[0];
152 
153  _weights[ 0] = 4.7619047619047619047619047619048e-02_R;
154  _weights[ 1] = 2.7682604736156594801070040629007e-01_R;
155  _weights[ 2] = 4.3174538120986262341787102228136e-01_R;
156  _weights[ 3] = 4.8761904761904761904761904761905e-01_R;
157  _weights[ 4] = _weights[2];
158  _weights[ 5] = _weights[1];
159  _weights[ 6] = _weights[0];
160 
161  return;
162  }
163 
164  // 2*8-3 = 13
165  case TWELFTH:
166  case THIRTEENTH:
167  {
168  _points.resize (8);
169  _weights.resize(8);
170 
171  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
172  _points[ 1](0) = -8.7174014850960661533744576122066e-01_R;
173  _points[ 2](0) = -5.9170018143314230214451073139795e-01_R;
174  _points[ 3](0) = -2.0929921790247886876865726034535e-01_R;
175  _points[ 4] = -_points[3];
176  _points[ 5] = -_points[2];
177  _points[ 6] = -_points[1];
178  _points[ 7] = -_points[0];
179 
180  _weights[ 0] = 3.5714285714285714285714285714286e-02_R;
181  _weights[ 1] = 2.1070422714350603938299206577576e-01_R;
182  _weights[ 2] = 3.4112269248350436476424067710775e-01_R;
183  _weights[ 3] = 4.1245879465870388156705297140221e-01_R;
184  _weights[ 4] = _weights[3];
185  _weights[ 5] = _weights[2];
186  _weights[ 6] = _weights[1];
187  _weights[ 7] = _weights[0];
188 
189  return;
190  }
191 
192  // 2*9-3 = 15
193  case FOURTEENTH:
194  case FIFTEENTH:
195  {
196  _points.resize (9);
197  _weights.resize(9);
198 
199  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
200  _points[ 1](0) = -8.9975799541146015731234524441834e-01_R;
201  _points[ 2](0) = -6.7718627951073775344588542709134e-01_R;
202  _points[ 3](0) = -3.6311746382617815871075206870866e-01_R;
203  _points[ 4](0) = 0.;
204  _points[ 5] = -_points[3];
205  _points[ 6] = -_points[2];
206  _points[ 7] = -_points[1];
207  _points[ 8] = -_points[0];
208 
209  _weights[ 0] = 2.7777777777777777777777777777778e-02_R;
210  _weights[ 1] = 1.6549536156080552504633972002921e-01_R;
211  _weights[ 2] = 2.7453871250016173528070561857937e-01_R;
212  _weights[ 3] = 3.4642851097304634511513153213972e-01_R;
213  _weights[ 4] = 3.7151927437641723356009070294785e-01_R;
214  _weights[ 5] = _weights[3];
215  _weights[ 6] = _weights[2];
216  _weights[ 7] = _weights[1];
217  _weights[ 8] = _weights[0];
218 
219  return;
220  }
221 
222  // 2*10-3 = 17
223  case SIXTEENTH:
224  case SEVENTEENTH:
225  {
226  _points.resize (10);
227  _weights.resize(10);
228 
229  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
230  _points[ 1](0) = -9.1953390816645881382893266082234e-01_R;
231  _points[ 2](0) = -7.3877386510550507500310617485983e-01_R;
232  _points[ 3](0) = -4.7792494981044449566117509273126e-01_R;
233  _points[ 4](0) = -1.6527895766638702462621976595817e-01_R;
234  _points[ 5] = -_points[4];
235  _points[ 6] = -_points[3];
236  _points[ 7] = -_points[2];
237  _points[ 8] = -_points[1];
238  _points[ 9] = -_points[0];
239 
240  _weights[ 0] = 2.2222222222222222222222222222222e-02_R;
241  _weights[ 1] = 1.3330599085107011112622717075539e-01_R;
242  _weights[ 2] = 2.2488934206312645211945782173105e-01_R;
243  _weights[ 3] = 2.9204268367968375787558225737444e-01_R;
244  _weights[ 4] = 3.2753976118389745665651052791689e-01_R;
245  _weights[ 5] = _weights[4];
246  _weights[ 6] = _weights[3];
247  _weights[ 7] = _weights[2];
248  _weights[ 8] = _weights[1];
249  _weights[ 9] = _weights[0];
250 
251  return;
252  }
253 
254  // 2*11-3 = 19
255  case EIGHTTEENTH:
256  case NINETEENTH:
257  {
258  _points.resize (11);
259  _weights.resize(11);
260 
261  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
262  _points[ 1](0) = -9.3400143040805913433227413609938e-01_R;
263  _points[ 2](0) = -7.8448347366314441862241781610846e-01_R;
264  _points[ 3](0) = -5.6523532699620500647096396947775e-01_R;
265  _points[ 4](0) = -2.9575813558693939143191151555906e-01_R;
266  _points[ 5](0) = 0.;
267  _points[ 6] = -_points[4];
268  _points[ 7] = -_points[3];
269  _points[ 8] = -_points[2];
270  _points[ 9] = -_points[1];
271  _points[10] = -_points[0];
272 
273  _weights[ 0] = 1.8181818181818181818181818181818e-02_R;
274  _weights[ 1] = 1.0961227326699486446140344958035e-01_R;
275  _weights[ 2] = 1.8716988178030520410814152189943e-01_R;
276  _weights[ 3] = 2.4804810426402831404008486642187e-01_R;
277  _weights[ 4] = 2.8687912477900808867922240333154e-01_R;
278  _weights[ 5] = 3.0021759545569069378593188116998e-01_R;
279  _weights[ 6] = _weights[4];
280  _weights[ 7] = _weights[3];
281  _weights[ 8] = _weights[2];
282  _weights[ 9] = _weights[1];
283  _weights[10] = _weights[0];
284 
285  return;
286  }
287 
288  // 2*12-3 = 21
289  case TWENTIETH:
290  case TWENTYFIRST:
291  {
292  _points.resize (12);
293  _weights.resize(12);
294 
295  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
296  _points[ 1](0) = -9.4489927222288222340758013830322e-01_R;
297  _points[ 2](0) = -8.1927932164400667834864158171690e-01_R;
298  _points[ 3](0) = -6.3287615303186067766240485444366e-01_R;
299  _points[ 4](0) = -3.9953094096534893226434979156697e-01_R;
300  _points[ 5](0) = -1.3655293285492755486406185573969e-01_R;
301  _points[ 6] = -_points[5];
302  _points[ 7] = -_points[4];
303  _points[ 8] = -_points[3];
304  _points[ 9] = -_points[2];
305  _points[10] = -_points[1];
306  _points[11] = -_points[0];
307 
308  _weights[ 0] = 1.5151515151515151515151515151515e-02_R;
309  _weights[ 1] = 9.1684517413196130668342594134079e-02_R;
310  _weights[ 2] = 1.5797470556437011516467106270034e-01_R;
311  _weights[ 3] = 2.1250841776102114535830207736687e-01_R;
312  _weights[ 4] = 2.5127560319920128029324441214760e-01_R;
313  _weights[ 5] = 2.7140524091069617700028833849960e-01_R;
314  _weights[ 6] = _weights[5];
315  _weights[ 7] = _weights[4];
316  _weights[ 8] = _weights[3];
317  _weights[ 9] = _weights[2];
318  _weights[10] = _weights[1];
319  _weights[11] = _weights[0];
320 
321  return;
322  }
323 
324  // 2*13-3 = 23
325  case TWENTYSECOND:
326  case TWENTYTHIRD:
327  {
328  _points.resize (13);
329  _weights.resize(13);
330 
331  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
332  _points[ 1](0) = -9.5330984664216391189690546475545e-01_R;
333  _points[ 2](0) = -8.4634756465187231686592560709875e-01_R;
334  _points[ 3](0) = -6.8618846908175742607275903956636e-01_R;
335  _points[ 4](0) = -4.8290982109133620174693723363693e-01_R;
336  _points[ 5](0) = -2.4928693010623999256867370037423e-01_R;
337  _points[ 6](0) = 0.;
338  _points[ 7] = -_points[5];
339  _points[ 8] = -_points[4];
340  _points[ 9] = -_points[3];
341  _points[10] = -_points[2];
342  _points[11] = -_points[1];
343  _points[12] = -_points[0];
344 
345  _weights[ 0] = 1.2820512820512820512820512820513e-02_R;
346  _weights[ 1] = 7.7801686746818927793588988333134e-02_R;
347  _weights[ 2] = 1.3498192668960834911991476258937e-01_R;
348  _weights[ 3] = 1.8364686520355009200749425874681e-01_R;
349  _weights[ 4] = 2.2076779356611008608553400837940e-01_R;
350  _weights[ 5] = 2.4401579030667635645857814836016e-01_R;
351  _weights[ 6] = 2.5193084933344673604413864154124e-01_R;
352  _weights[ 7] = _weights[5];
353  _weights[ 8] = _weights[4];
354  _weights[ 9] = _weights[3];
355  _weights[10] = _weights[2];
356  _weights[11] = _weights[1];
357  _weights[12] = _weights[0];
358 
359  return;
360  }
361 
362  // 2*14-3 = 25
363  case TWENTYFOURTH:
364  case TWENTYFIFTH:
365  {
366  _points.resize (14);
367  _weights.resize(14);
368 
369  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
370  _points[ 1](0) = -9.5993504526726090135510016201542e-01_R;
371  _points[ 2](0) = -8.6780105383034725100022020290826e-01_R;
372  _points[ 3](0) = -7.2886859909132614058467240052088e-01_R;
373  _points[ 4](0) = -5.5063940292864705531662270585908e-01_R;
374  _points[ 5](0) = -3.4272401334271284504390340364167e-01_R;
375  _points[ 6](0) = -1.1633186888370386765877670973616e-01_R;
376  _points[ 7] = -_points[6];
377  _points[ 8] = -_points[5];
378  _points[ 9] = -_points[4];
379  _points[10] = -_points[3];
380  _points[11] = -_points[2];
381  _points[12] = -_points[1];
382  _points[13] = -_points[0];
383 
384  _weights[ 0] = 1.0989010989010989010989010989011e-02_R;
385  _weights[ 1] = 6.6837284497681284634070660746053e-02_R;
386  _weights[ 2] = 1.1658665589871165154099667065465e-01_R;
387  _weights[ 3] = 1.6002185176295214241282099798759e-01_R;
388  _weights[ 4] = 1.9482614937341611864033177837588e-01_R;
389  _weights[ 5] = 2.1912625300977075487116252395417e-01_R;
390  _weights[ 6] = 2.3161279446845705888962835729264e-01_R;
391  _weights[ 7] = _weights[6];
392  _weights[ 8] = _weights[5];
393  _weights[ 9] = _weights[4];
394  _weights[10] = _weights[3];
395  _weights[11] = _weights[2];
396  _weights[12] = _weights[1];
397  _weights[13] = _weights[0];
398 
399  return;
400  }
401 
402  // 2*15-3 = 27
403  case TWENTYSIXTH:
404  case TWENTYSEVENTH:
405  {
406  _points.resize (15);
407  _weights.resize(15);
408 
409  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
410  _points[ 1](0) = -9.6524592650383857279585139206960e-01_R;
411  _points[ 2](0) = -8.8508204422297629882540163148223e-01_R;
412  _points[ 3](0) = -7.6351968995181520070411847597629e-01_R;
413  _points[ 4](0) = -6.0625320546984571112352993863673e-01_R;
414  _points[ 5](0) = -4.2063805471367248092189693873858e-01_R;
415  _points[ 6](0) = -2.1535395536379423822567944627292e-01_R;
416  _points[ 7](0) = 0.;
417  _points[ 8] = -_points[6];
418  _points[ 9] = -_points[5];
419  _points[10] = -_points[4];
420  _points[11] = -_points[3];
421  _points[12] = -_points[2];
422  _points[13] = -_points[1];
423  _points[14] = -_points[0];
424 
425  _weights[ 0] = 9.5238095238095238095238095238095e-03_R;
426  _weights[ 1] = 5.8029893028601249096880584025282e-02_R;
427  _weights[ 2] = 1.0166007032571806760366617078880e-01_R;
428  _weights[ 3] = 1.4051169980242810946044680564367e-01_R;
429  _weights[ 4] = 1.7278964725360094905207709940835e-01_R;
430  _weights[ 5] = 1.9698723596461335609250034650741e-01_R;
431  _weights[ 6] = 2.1197358592682092012743007697722e-01_R;
432  _weights[ 7] = 2.1704811634881564951495021425091e-01_R;
433  _weights[ 8] = _weights[6];
434  _weights[ 9] = _weights[5];
435  _weights[10] = _weights[4];
436  _weights[11] = _weights[3];
437  _weights[12] = _weights[2];
438  _weights[13] = _weights[1];
439  _weights[14] = _weights[0];
440 
441  return;
442  }
443 
444  // 2*16-3 = 29
445  case TWENTYEIGHTH:
446  case TWENTYNINTH:
447  {
448  _points.resize (16);
449  _weights.resize(16);
450 
451  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
452  _points[ 1](0) = -9.6956804627021793295224273836746e-01_R;
453  _points[ 2](0) = -8.9920053309347209299462826151985e-01_R;
454  _points[ 3](0) = -7.9200829186181506393108827096315e-01_R;
455  _points[ 4](0) = -6.5238870288249308946788321964058e-01_R;
456  _points[ 5](0) = -4.8605942188713761178189078584687e-01_R;
457  _points[ 6](0) = -2.9983046890076320809835345472230e-01_R;
458  _points[ 7](0) = -1.0132627352194944784303300504592e-01_R;
459  _points[ 8] = -_points[7];
460  _points[ 9] = -_points[6];
461  _points[10] = -_points[5];
462  _points[11] = -_points[4];
463  _points[12] = -_points[3];
464  _points[13] = -_points[2];
465  _points[14] = -_points[1];
466  _points[15] = -_points[0];
467 
468  _weights[ 0] = 8.3333333333333333333333333333333e-03_R;
469  _weights[ 1] = 5.0850361005919905403244919565455e-02_R;
470  _weights[ 2] = 8.9393697325930800991052080166084e-02_R;
471  _weights[ 3] = 1.2425538213251409834953633265731e-01_R;
472  _weights[ 4] = 1.5402698080716428081564494048499e-01_R;
473  _weights[ 5] = 1.7749191339170412530107566952836e-01_R;
474  _weights[ 6] = 1.9369002382520358431691359885352e-01_R;
475  _weights[ 7] = 2.0195830817822987148919912541094e-01_R;
476  _weights[ 8] = _weights[7];
477  _weights[ 9] = _weights[6];
478  _weights[10] = _weights[5];
479  _weights[11] = _weights[4];
480  _weights[12] = _weights[3];
481  _weights[13] = _weights[2];
482  _weights[14] = _weights[1];
483  _weights[15] = _weights[0];
484 
485  return;
486  }
487 
488  // 2*17-3 = 31
489  case THIRTIETH:
490  case THIRTYFIRST:
491  {
492  _points.resize (17);
493  _weights.resize(17);
494 
495  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
496  _points[ 1](0) = -9.7313217663141831415697950187372e-01_R;
497  _points[ 2](0) = -9.1087999591557359562380250639773e-01_R;
498  _points[ 3](0) = -8.1569625122177030710675055323753e-01_R;
499  _points[ 4](0) = -6.9102898062768470539491935737245e-01_R;
500  _points[ 5](0) = -5.4138539933010153912373340750406e-01_R;
501  _points[ 6](0) = -3.7217443356547704190723468073526e-01_R;
502  _points[ 7](0) = -1.8951197351831738830426301475311e-01_R;
503  _points[ 8](0) = 0.;
504  _points[ 9] = -_points[7];
505  _points[10] = -_points[6];
506  _points[11] = -_points[5];
507  _points[12] = -_points[4];
508  _points[13] = -_points[3];
509  _points[14] = -_points[2];
510  _points[15] = -_points[1];
511  _points[16] = -_points[0];
512 
513  _weights[ 0] = 7.3529411764705882352941176470588e-03_R;
514  _weights[ 1] = 4.4921940543254209647400954623212e-02_R;
515  _weights[ 2] = 7.9198270503687119190264429952835e-02_R;
516  _weights[ 3] = 1.1059290900702816137577270522008e-01_R;
517  _weights[ 4] = 1.3798774620192655905620157495403e-01_R;
518  _weights[ 5] = 1.6039466199762153951632836586475e-01_R;
519  _weights[ 6] = 1.7700425351565787043694574536329e-01_R;
520  _weights[ 7] = 1.8721633967761923589208848286062e-01_R;
521  _weights[ 8] = 1.9066187475346943329940724702825e-01_R;
522  _weights[ 9] = _weights[7];
523  _weights[10] = _weights[6];
524  _weights[11] = _weights[5];
525  _weights[12] = _weights[4];
526  _weights[13] = _weights[3];
527  _weights[14] = _weights[2];
528  _weights[15] = _weights[1];
529  _weights[16] = _weights[0];
530 
531  return;
532  }
533 
534  // 2*18-3 = 33
535  case THIRTYSECOND:
536  case THIRTYTHIRD:
537  {
538  _points.resize (18);
539  _weights.resize(18);
540 
541  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
542  _points[ 1](0) = -9.7610555741219854286451892434170e-01_R;
543  _points[ 2](0) = -9.2064918534753387383785462543128e-01_R;
544  _points[ 3](0) = -8.3559353521809021371364636232794e-01_R;
545  _points[ 4](0) = -7.2367932928324268130621036530207e-01_R;
546  _points[ 5](0) = -5.8850483431866176117353589319356e-01_R;
547  _points[ 6](0) = -4.3441503691212397534228713674067e-01_R;
548  _points[ 7](0) = -2.6636265287828098416766533202560e-01_R;
549  _points[ 8](0) = -8.9749093484652111022645010088562e-02_R;
550  _points[ 9] = -_points[8];
551  _points[10] = -_points[7];
552  _points[11] = -_points[6];
553  _points[12] = -_points[5];
554  _points[13] = -_points[4];
555  _points[14] = -_points[3];
556  _points[15] = -_points[2];
557  _points[16] = -_points[1];
558  _points[17] = -_points[0];
559 
560  _weights[ 0] = 6.5359477124183006535947712418301e-03_R;
561  _weights[ 1] = 3.9970628810914066137599176410101e-02_R;
562  _weights[ 2] = 7.0637166885633664999222960167786e-02_R;
563  _weights[ 3] = 9.9016271717502802394423605318672e-02_R;
564  _weights[ 4] = 1.2421053313296710026339635889675e-01_R;
565  _weights[ 5] = 1.4541196157380226798300321049443e-01_R;
566  _weights[ 6] = 1.6193951723760248926432670670023e-01_R;
567  _weights[ 7] = 1.7326210948945622601061440382668e-01_R;
568  _weights[ 8] = 1.7901586343970308229381880694353e-01_R;
569  _weights[ 9] = _weights[8];
570  _weights[10] = _weights[7];
571  _weights[11] = _weights[6];
572  _weights[12] = _weights[5];
573  _weights[13] = _weights[4];
574  _weights[14] = _weights[3];
575  _weights[15] = _weights[2];
576  _weights[16] = _weights[1];
577  _weights[17] = _weights[0];
578 
579  return;
580  }
581 
582  // 2*19-3 = 35
583  case THIRTYFOURTH:
584  case THIRTYFIFTH:
585  {
586  _points.resize (19);
587  _weights.resize(19);
588 
589  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
590  _points[ 1](0) = -9.7861176622208009515263406311022e-01_R;
591  _points[ 2](0) = -9.2890152815258624371794025879655e-01_R;
592  _points[ 3](0) = -8.5246057779664609308595597004106e-01_R;
593  _points[ 4](0) = -7.5149420255261301416363748963394e-01_R;
594  _points[ 5](0) = -6.2890813726522049776683230622873e-01_R;
595  _points[ 6](0) = -4.8822928568071350277790963762492e-01_R;
596  _points[ 7](0) = -3.3350484782449861029850010384493e-01_R;
597  _points[ 8](0) = -1.6918602340928157137515415344488e-01_R;
598  _points[ 9](0) = 0.;
599  _points[10] = -_points[8];
600  _points[11] = -_points[7];
601  _points[12] = -_points[6];
602  _points[13] = -_points[5];
603  _points[14] = -_points[4];
604  _points[15] = -_points[3];
605  _points[16] = -_points[2];
606  _points[17] = -_points[1];
607  _points[18] = -_points[0];
608 
609  _weights[ 0] = 5.8479532163742690058479532163743e-03_R;
610  _weights[ 1] = 3.5793365186176477115425569035122e-02_R;
611  _weights[ 2] = 6.3381891762629736851695690418317e-02_R;
612  _weights[ 3] = 8.9131757099207084448008790556153e-02_R;
613  _weights[ 4] = 1.1231534147730504407091001546378e-01_R;
614  _weights[ 5] = 1.3226728044875077692604673390973e-01_R;
615  _weights[ 6] = 1.4841394259593888500968064366841e-01_R;
616  _weights[ 7] = 1.6029092404406124197991096818359e-01_R;
617  _weights[ 8] = 1.6755658452714286727013727774026e-01_R;
618  _weights[ 9] = 1.7000191928482723464467271561652e-01_R;
619  _weights[10] = _weights[8];
620  _weights[11] = _weights[7];
621  _weights[12] = _weights[6];
622  _weights[13] = _weights[5];
623  _weights[14] = _weights[4];
624  _weights[15] = _weights[3];
625  _weights[16] = _weights[2];
626  _weights[17] = _weights[1];
627  _weights[18] = _weights[0];
628 
629  return;
630  }
631 
632  // 2*20-3 = 37
633  case THIRTYSIXTH:
634  case THIRTYSEVENTH:
635  {
636  _points.resize (20);
637  _weights.resize(20);
638 
639  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
640  _points[ 1](0) = -9.8074370489391417192544643858423e-01_R;
641  _points[ 2](0) = -9.3593449881266543571618158493063e-01_R;
642  _points[ 3](0) = -8.6687797808995014130984721461629e-01_R;
643  _points[ 4](0) = -7.7536826095205587041431752759469e-01_R;
644  _points[ 5](0) = -6.6377640229031128984640332297116e-01_R;
645  _points[ 6](0) = -5.3499286403188626164813596182898e-01_R;
646  _points[ 7](0) = -3.9235318371390929938647470381582e-01_R;
647  _points[ 8](0) = -2.3955170592298649518240135692709e-01_R;
648  _points[ 9](0) = -8.0545937238821837975944518159554e-02_R;
649  _points[10] = -_points[9];
650  _points[11] = -_points[8];
651  _points[12] = -_points[7];
652  _points[13] = -_points[6];
653  _points[14] = -_points[5];
654  _points[15] = -_points[4];
655  _points[16] = -_points[3];
656  _points[17] = -_points[2];
657  _points[18] = -_points[1];
658  _points[19] = -_points[0];
659 
660  _weights[ 0] = 5.2631578947368421052631578947368e-03_R;
661  _weights[ 1] = 3.2237123188488941491605028117294e-02_R;
662  _weights[ 2] = 5.7181802127566826004753627173243e-02_R;
663  _weights[ 3] = 8.0631763996119603144776846113721e-02_R;
664  _weights[ 4] = 1.0199149969945081568378120573289e-01_R;
665  _weights[ 5] = 1.2070922762867472509942970500239e-01_R;
666  _weights[ 6] = 1.3630048235872418448978079298903e-01_R;
667  _weights[ 7] = 1.4836155407091682581471301373397e-01_R;
668  _weights[ 8] = 1.5658010264747548715816989679364e-01_R;
669  _weights[ 9] = 1.6074328638784574900772672644908e-01_R;
670  _weights[10] = _weights[9];
671  _weights[11] = _weights[8];
672  _weights[12] = _weights[7];
673  _weights[13] = _weights[6];
674  _weights[14] = _weights[5];
675  _weights[15] = _weights[4];
676  _weights[16] = _weights[3];
677  _weights[17] = _weights[2];
678  _weights[18] = _weights[1];
679  _weights[19] = _weights[0];
680 
681  return;
682  }
683 
684  // 2*21-3 = 39
685  case THIRTYEIGHTH:
686  case THIRTYNINTH:
687  {
688  _points.resize (21);
689  _weights.resize(21);
690 
691  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
692  _points[ 1](0) = -9.8257229660454802823448127655541e-01_R;
693  _points[ 2](0) = -9.4197629695974553429610265066144e-01_R;
694  _points[ 3](0) = -8.7929475532359046445115359630494e-01_R;
695  _points[ 4](0) = -7.9600192607771240474431258966036e-01_R;
696  _points[ 5](0) = -6.9405102606222323262731639319467e-01_R;
697  _points[ 6](0) = -5.7583196026183068692702187033809e-01_R;
698  _points[ 7](0) = -4.4411578327900210119451634960735e-01_R;
699  _points[ 8](0) = -3.0198985650876488727535186785875e-01_R;
700  _points[ 9](0) = -1.5278551580218546600635832848567e-01_R;
701  _points[10](0) = 0.;
702  _points[11] = -_points[9];
703  _points[12] = -_points[8];
704  _points[13] = -_points[7];
705  _points[14] = -_points[6];
706  _points[15] = -_points[5];
707  _points[16] = -_points[4];
708  _points[17] = -_points[3];
709  _points[18] = -_points[2];
710  _points[19] = -_points[1];
711  _points[20] = -_points[0];
712 
713  _weights[ 0] = 4.7619047619047619047619047619048e-03_R;
714  _weights[ 1] = 2.9184840098505458609458543613171e-02_R;
715  _weights[ 2] = 5.1843169000849625072722971852830e-02_R;
716  _weights[ 3] = 7.3273918185074144252547861041894e-02_R;
717  _weights[ 4] = 9.2985467957886065301137664149214e-02_R;
718  _weights[ 5] = 1.1051708321912333526700048678439e-01_R;
719  _weights[ 6] = 1.2545812119086894801515753570800e-01_R;
720  _weights[ 7] = 1.3745846286004134358089961741515e-01_R;
721  _weights[ 8] = 1.4623686244797745926727053063439e-01_R;
722  _weights[ 9] = 1.5158757511168138445325068150529e-01_R;
723  _weights[10] = 1.5338519033217494855158440506754e-01_R;
724  _weights[11] = _weights[9];
725  _weights[12] = _weights[8];
726  _weights[13] = _weights[7];
727  _weights[14] = _weights[6];
728  _weights[15] = _weights[5];
729  _weights[16] = _weights[4];
730  _weights[17] = _weights[3];
731  _weights[18] = _weights[2];
732  _weights[19] = _weights[1];
733  _weights[20] = _weights[0];
734 
735  return;
736  }
737 
738  // 2*22-3 = 41
739  case FORTIETH:
740  case FORTYFIRST:
741  {
742  _points.resize (22);
743  _weights.resize(22);
744 
745  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
746  _points[ 1](0) = -9.8415243845764617655228962221207e-01_R;
747  _points[ 2](0) = -9.4720428399922868052421376661573e-01_R;
748  _points[ 3](0) = -8.9006229019090447052965782577909e-01_R;
749  _points[ 4](0) = -8.1394892761192113604544184805614e-01_R;
750  _points[ 5](0) = -7.2048723996120215811988189639847e-01_R;
751  _points[ 6](0) = -6.1166943828425897122621160586993e-01_R;
752  _points[ 7](0) = -4.8981487518990234980875123568327e-01_R;
753  _points[ 8](0) = -3.5752071013891953806095728024018e-01_R;
754  _points[ 9](0) = -2.1760658515928504178795509346539e-01_R;
755  _points[10](0) = -7.3054540010898334761088790464107e-02_R;
756  _points[11] = -_points[10];
757  _points[12] = -_points[9];
758  _points[13] = -_points[8];
759  _points[14] = -_points[7];
760  _points[15] = -_points[6];
761  _points[16] = -_points[5];
762  _points[17] = -_points[4];
763  _points[18] = -_points[3];
764  _points[19] = -_points[2];
765  _points[20] = -_points[1];
766  _points[21] = -_points[0];
767 
768  _weights[ 0] = 4.3290043290043290043290043290043e-03_R;
769  _weights[ 1] = 2.6545747682501757911627904520543e-02_R;
770  _weights[ 2] = 4.7214465293740752123775734864792e-02_R;
771  _weights[ 3] = 6.6865605864553076012404194157097e-02_R;
772  _weights[ 4] = 8.5090060391838447815711236095748e-02_R;
773  _weights[ 5] = 1.0150057480164767437243730374960e-01_R;
774  _weights[ 6] = 1.1574764465393906659003636772146e-01_R;
775  _weights[ 7] = 1.2752769665343027553084445930883e-01_R;
776  _weights[ 8] = 1.3658968861374142668617736220617e-01_R;
777  _weights[ 9] = 1.4274049227136140033623599356679e-01_R;
778  _weights[10] = 1.4584901944424179361642043947997e-01_R;
779  _weights[11] = _weights[10];
780  _weights[12] = _weights[9];
781  _weights[13] = _weights[8];
782  _weights[14] = _weights[7];
783  _weights[15] = _weights[6];
784  _weights[16] = _weights[5];
785  _weights[17] = _weights[4];
786  _weights[18] = _weights[3];
787  _weights[19] = _weights[2];
788  _weights[20] = _weights[1];
789  _weights[21] = _weights[0];
790 
791  return;
792  }
793 
794  // 2*23-3 = 43
795  case FORTYSECOND:
796  case FORTYTHIRD:
797  {
798  _points.resize (23);
799  _weights.resize(23);
800 
801  _points[ 0](0) = -1.0000000000000000000000000000000e+00_R;
802  _points[ 1](0) = -9.8552715587873257808146276673810e-01_R;
803  _points[ 2](0) = -9.5175795571071020413563967985143e-01_R;
804  _points[ 3](0) = -8.9945855804034501095016032034737e-01_R;
805  _points[ 4](0) = -8.2965109665128588622320061929000e-01_R;
806  _points[ 5](0) = -7.4369504117206068394516354306700e-01_R;
807  _points[ 6](0) = -6.4326364446013620847614553360277e-01_R;
808  _points[ 7](0) = -5.3031177113684416813011532015230e-01_R;
809  _points[ 8](0) = -4.0703793791447482919595048821510e-01_R;
810  _points[ 9](0) = -2.7584154894579306710687763267914e-01_R;
811  _points[10](0) = -1.3927620404066839859186261298277e-01_R;
812  _points[11](0) = 0.;
813  _points[12] = -_points[10];
814  _points[13] = -_points[9];
815  _points[14] = -_points[8];
816  _points[15] = -_points[7];
817  _points[16] = -_points[6];
818  _points[17] = -_points[5];
819  _points[18] = -_points[4];
820  _points[19] = -_points[3];
821  _points[20] = -_points[2];
822  _points[21] = -_points[1];
823  _points[22] = -_points[0];
824 
825  _weights[ 0] = 3.9525691699604743083003952569170e-03_R;
826  _weights[ 1] = 2.4248600771531736517399658937097e-02_R;
827  _weights[ 2] = 4.3175871170241834748876465612042e-02_R;
828  _weights[ 3] = 6.1252477129554206381382847440355e-02_R;
829  _weights[ 4] = 7.8135449475569989741934255347965e-02_R;
830  _weights[ 5] = 9.3497246163512341833500706906697e-02_R;
831  _weights[ 6] = 1.0703910172433651153518362791547e-01_R;
832  _weights[ 7] = 1.1849751066274913130212600472426e-01_R;
833  _weights[ 8] = 1.2764947470175887663614855305567e-01_R;
834  _weights[ 9] = 1.3431687263860381990156489770071e-01_R;
835  _weights[10] = 1.3836993638580739452350273386294e-01_R;
836  _weights[11] = 1.3972978001274736514015970647975e-01_R;
837  _weights[12] = _weights[10];
838  _weights[13] = _weights[9];
839  _weights[14] = _weights[8];
840  _weights[15] = _weights[7];
841  _weights[16] = _weights[6];
842  _weights[17] = _weights[5];
843  _weights[18] = _weights[4];
844  _weights[19] = _weights[3];
845  _weights[20] = _weights[2];
846  _weights[21] = _weights[1];
847  _weights[22] = _weights[0];
848 
849  return;
850  }
851 
852  default:
853  libmesh_error_msg("Quadrature rule " << _order << " not supported!");
854  }
855 }
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
Order get_order() const
Definition: quadrature.h:249
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ init_2D()

void libMesh::QGaussLobatto::init_2D ( )
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.

Reimplemented from libMesh::QBase.

Definition at line 28 of file quadrature_gauss_lobatto_2D.C.

References libMesh::QBase::_dim, libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::Utility::enum_to_string(), libMesh::QBase::get_order(), libMesh::QBase::get_points(), libMesh::QBase::get_weights(), libMesh::QBase::init(), libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::QUADSHELL4, libMesh::QUADSHELL8, libMesh::QUADSHELL9, and libMesh::QBase::tensor_product_quad().

29 {
30  switch (_type)
31  {
32  case QUAD4:
33  case QUADSHELL4:
34  case QUAD8:
35  case QUADSHELL8:
36  case QUAD9:
37  case QUADSHELL9:
38  {
39  // We compute the 2D quadrature rule as a tensor
40  // product of the 1D quadrature rule.
41  QGaussLobatto q1D(1, get_order());
43  return;
44  }
45 
46  // We fall back on a Gauss type rule for other types of elements,
47  // but we warn the user (once) that we are doing this instead of
48  // silently switching out quadrature rules on them.
49  default:
50  {
51  libmesh_warning("Warning: QGaussLobatto falling back on QGauss rule "
52  "for unsupported Elem type: " << Utility::enum_to_string(_type));
53 
54  QGauss gauss_rule(_dim, _order);
55  gauss_rule.init(*this);
56 
57  // Swap points and weights with the about-to-be destroyed rule.
58  _points.swap (gauss_rule.get_points() );
59  _weights.swap(gauss_rule.get_weights());
60  }
61  }
62 }
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
Order get_order() const
Definition: quadrature.h:249
std::string enum_to_string(const T e)
QGaussLobatto(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385
void tensor_product_quad(const QBase &q1D)
Constructs a 2D rule from the tensor product of q1D with itself.
Definition: quadrature.C:256

◆ init_3D()

void libMesh::QGaussLobatto::init_3D ( )
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.

Reimplemented from libMesh::QBase.

Definition at line 28 of file quadrature_gauss_lobatto_3D.C.

References libMesh::QBase::_dim, libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_type, libMesh::QBase::_weights, libMesh::Utility::enum_to_string(), libMesh::QBase::get_order(), libMesh::QBase::get_points(), libMesh::QBase::get_weights(), libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::QBase::init(), and libMesh::QBase::tensor_product_hex().

29 {
30  switch (_type)
31  {
32  case HEX8:
33  case HEX20:
34  case HEX27:
35  {
36  // We compute the 3D quadrature rule as a tensor
37  // product of the 1D quadrature rule.
38  QGaussLobatto q1D(1, get_order());
39  tensor_product_hex(q1D);
40  return;
41  }
42 
43  // We fall back on a Gauss type rule for other types of elements,
44  // but we warn the user (once) that we are doing this instead of
45  // silently switching out quadrature rules on them.
46  default:
47  {
48  libmesh_warning("Warning: QGaussLobatto falling back on QGauss rule "
49  "for unsupported Elem type: " << Utility::enum_to_string(_type));
50 
51  QGauss gauss_rule(_dim, _order);
52  gauss_rule.init(*this);
53 
54  // Swap points and weights with the about-to-be destroyed rule.
55  _points.swap (gauss_rule.get_points() );
56  _weights.swap(gauss_rule.get_weights());
57  }
58  }
59 }
ElemType _type
The type of element for which the current values have been computed.
Definition: quadrature.h:391
unsigned int _dim
The spatial dimension of the quadrature rule.
Definition: quadrature.h:379
std::vector< Point > _points
The locations of the quadrature points in reference element space.
Definition: quadrature.h:409
std::vector< Real > _weights
The quadrature weights.
Definition: quadrature.h:415
void tensor_product_hex(const QBase &q1D)
Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D.
Definition: quadrature.C:283
Order get_order() const
Definition: quadrature.h:249
std::string enum_to_string(const T e)
QGaussLobatto(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
Definition: quadrature.h:385

◆ n_objects()

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

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

Definition at line 85 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

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

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

◆ n_points()

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

Definition at line 131 of file quadrature.h.

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(), AssemblyA1::boundary_assembly(), AssemblyF0::boundary_assembly(), AssemblyF1::boundary_assembly(), AssemblyA2::boundary_assembly(), AssemblyF2::boundary_assembly(), A2::boundary_assembly(), A3::boundary_assembly(), F0::boundary_assembly(), Output0::boundary_assembly(), compute_enriched_soln(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::HDGProblem::create_identity_jacobian(), libMesh::HDGProblem::create_identity_residual(), SecondOrderScalarSystemSecondOrderTimeSolverBase::damping_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::damping_residual(), NavierSystem::element_constraint(), CoupledSystem::element_constraint(), PoissonSystem::element_postprocess(), LaplaceSystem::element_postprocess(), LaplaceQoI::element_qoi(), HeatSystem::element_qoi(), LaplaceQoI::element_qoi_derivative(), LaplaceSystem::element_qoi_derivative(), HeatSystem::element_qoi_derivative(), NavierSystem::element_time_derivative(), SolidSystem::element_time_derivative(), PoissonSystem::element_time_derivative(), LaplaceSystem::element_time_derivative(), CurlCurlSystem::element_time_derivative(), ElasticitySystem::element_time_derivative(), CoupledSystem::element_time_derivative(), HilbertSystem::element_time_derivative(), HeatSystem::element_time_derivative(), FirstOrderScalarSystemBase::element_time_derivative(), SecondOrderScalarSystemFirstOrderTimeSolverBase::element_time_derivative(), fe_assembly(), A0::interior_assembly(), B::interior_assembly(), M0::interior_assembly(), A1::interior_assembly(), AssemblyA0::interior_assembly(), AcousticsInnerProduct::interior_assembly(), AssemblyA1::interior_assembly(), A2::interior_assembly(), AssemblyA2::interior_assembly(), F0::interior_assembly(), OutputAssembly::interior_assembly(), EIM_IP_assembly::interior_assembly(), EIM_F::interior_assembly(), AssemblyEIM::interior_assembly(), InnerProductAssembly::interior_assembly(), AssemblyF0::interior_assembly(), AssemblyF1::interior_assembly(), Ex6InnerProduct::interior_assembly(), LaplaceYoung::jacobian(), NavierSystem::mass_residual(), ElasticitySystem::mass_residual(), libMesh::FEMPhysics::mass_residual(), FirstOrderScalarSystemBase::mass_residual(), SecondOrderScalarSystemSecondOrderTimeSolverBase::mass_residual(), SecondOrderScalarSystemFirstOrderTimeSolverBase::mass_residual(), libMesh::FEAbstract::n_quadrature_points(), 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::size(), libMesh::QBase::tensor_product_hex(), libMesh::QBase::tensor_product_prism(), libMesh::QBase::tensor_product_quad(), and QuadratureTest::testPolynomial().

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

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ print_info() [1/2]

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

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

Definition at line 81 of file reference_counter.C.

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

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

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

◆ print_info() [2/2]

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

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

Definition at line 43 of file quadrature.C.

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

Referenced by libMesh::operator<<().

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

◆ qp()

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

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

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

◆ scale()

void libMesh::QBase::scale ( std::pair< Real, Real old_range,
std::pair< Real, Real new_range 
)
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 222 of file quadrature.C.

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

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

◆ shapes_need_reinit()

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

286 { return false; }

◆ size()

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

Alias for n_points() to enable use in index_range.

Returns
The number of points associated with the quadrature rule.

Definition at line 142 of file quadrature.h.

References libMesh::QBase::n_points().

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

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

◆ tensor_product_hex()

void libMesh::QBase::tensor_product_hex ( const QBase q1D)
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 283 of file quadrature.C.

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

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

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

◆ tensor_product_prism()

void libMesh::QBase::tensor_product_prism ( const QBase q1D,
const QBase q2D 
)
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 310 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::QGauss::init_3D(), libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGrid::init_3D().

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

◆ tensor_product_quad()

void libMesh::QBase::tensor_product_quad ( const QBase q1D)
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 256 of file quadrature.C.

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

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

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

◆ type()

QuadratureType libMesh::QGaussLobatto::type ( ) const
overridevirtual
Returns
QGAUSS_LOBATTO.

Implements libMesh::QBase.

Definition at line 46 of file quadrature_gauss_lobatto.C.

References libMesh::QGAUSS_LOBATTO.

◆ w()

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

Definition at line 188 of file quadrature.h.

References libMesh::QBase::_weights.

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(), libMesh::QBase::tensor_product_quad(), and QuadratureTest::testPolynomial().

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

Member Data Documentation

◆ _counts

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

Actually holds the data.

Definition at line 124 of file reference_counter.h.

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

◆ _dim

unsigned int libMesh::QBase::_dim
protectedinherited

◆ _elem

const Elem* libMesh::QBase::_elem
protectedinherited

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

Definition at line 397 of file quadrature.h.

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

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 143 of file reference_counter.h.

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

◆ _mutex

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

◆ _n_objects

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

The number of objects.

Print the reference count information when the number returns to 0.

Definition at line 132 of file reference_counter.h.

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

◆ _order

Order libMesh::QBase::_order
protectedinherited

◆ _p_level

unsigned int libMesh::QBase::_p_level
protectedinherited

◆ _points

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

The locations of the quadrature points in reference element space.

Definition at line 409 of file quadrature.h.

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::QGauss::dunavant_rule(), libMesh::QGauss::dunavant_rule2(), libMesh::QBase::get_points(), libMesh::QGrundmann_Moller::gm_rule(), libMesh::QBase::init_0D(), init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QTrap::init_1D(), libMesh::QGrid::init_1D(), libMesh::QConical::init_1D(), libMesh::QNodal::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QGrundmann_Moller::init_1D(), init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QTrap::init_2D(), libMesh::QGrid::init_2D(), libMesh::QNodal::init_2D(), libMesh::QMonomial::init_2D(), init_3D(), libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGauss::init_3D(), libMesh::QGrid::init_3D(), libMesh::QNodal::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGrundmann_Moller::init_3D(), libMesh::QGauss::keast_rule(), libMesh::QMonomial::kim_rule(), libMesh::QBase::n_points(), libMesh::QBase::print_info(), libMesh::QBase::qp(), libMesh::QBase::scale(), libMesh::QMonomial::stroud_rule(), libMesh::QBase::tensor_product_hex(), libMesh::QBase::tensor_product_prism(), libMesh::QBase::tensor_product_quad(), and libMesh::QMonomial::wissmann_rule().

◆ _type

ElemType libMesh::QBase::_type
protectedinherited

◆ _weights

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

The quadrature weights.

The order of the weights matches the ordering of the _points vector.

Definition at line 415 of file quadrature.h.

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::QGauss::dunavant_rule(), libMesh::QGauss::dunavant_rule2(), libMesh::QBase::get_weights(), libMesh::QGrundmann_Moller::gm_rule(), libMesh::QBase::init_0D(), init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QTrap::init_1D(), libMesh::QGrid::init_1D(), libMesh::QConical::init_1D(), libMesh::QNodal::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QGrundmann_Moller::init_1D(), init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QTrap::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QNodal::init_2D(), libMesh::QMonomial::init_2D(), init_3D(), libMesh::QTrap::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGauss::init_3D(), libMesh::QGrid::init_3D(), libMesh::QNodal::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGrundmann_Moller::init_3D(), libMesh::QGauss::keast_rule(), libMesh::QMonomial::kim_rule(), libMesh::QBase::print_info(), libMesh::QBase::scale(), libMesh::QMonomial::stroud_rule(), libMesh::QBase::tensor_product_hex(), libMesh::QBase::tensor_product_prism(), libMesh::QBase::tensor_product_quad(), libMesh::QBase::w(), and libMesh::QMonomial::wissmann_rule().

◆ allow_nodal_pyramid_quadrature

bool libMesh::QBase::allow_nodal_pyramid_quadrature
inherited

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

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

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

Definition at line 314 of file quadrature.h.

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

◆ allow_rules_with_negative_weights

bool libMesh::QBase::allow_rules_with_negative_weights
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 301 of file quadrature.h.

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


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