#include <quadrature_gm.h>
QGrundmann_Moller (const unsigned int _dim, const Order _order=INVALID_ORDER)
~QGrundmann_Moller ()
QuadratureType type () const
ElemType get_elem_type () const
unsigned int get_p_level () const
unsigned int n_points () const
unsigned int get_dim () const
const std::vector< Point > & get_points () const
std::vector< Point > & get_points ()
const std::vector< Real > & get_weights () const
std::vector< Real > & get_weights ()
Point qp (const unsigned int i) const
Real w (const unsigned int i) const
void init (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)
Order get_order () const
void print_info (std::ostream &os=std::cout) const
void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
static AutoPtr< QBase > build (const std::string &name, const unsigned int _dim, const Order _order=INVALID_ORDER)
static AutoPtr< QBase > build (const QuadratureType _qt, const unsigned int _dim, const Order _order=INVALID_ORDER)
static void print_info ()
static std::string get_info ()
static unsigned int n_objects ()
bool allow_rules_with_negative_weights
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
virtual void init_0D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)
virtual void init_2D (const ElemType, unsigned int=0)
void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)
std::cerr<< 'ERROR: Seems as if this quadrature rule'<< std::endl<< ' is not implemented for 2D.'<< std::endl;libmesh_error();}#endif virtual void init_3D(const ElemType, unsigned int=0)#ifndef DEBUG{}#else{std::cerr<< 'ERROR: Seems as if this quadrature rule'<< std::endl<< ' is not implemented for 3D.'<< std::endl;libmesh_error();}#endif void tensor_product_quad(const QBase &q1D);void tensor_product_hex(const QBase &q1D);void tensor_product_prism(const QBase &q1D, const QBase &q2D);const unsigned int _dim;const Order _order;ElemType _type;unsigned int _p_level;std::vector< Point > _points
std::vector< Real > _weights
static Counts _counts
static Threads::atomic< unsigned int > _n_objects
static Threads::spin_mutex _mutex
void init_1D (const ElemType, unsigned int=0)
void init_3D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)
void gm_rule (unsigned int s)
void compose_all (unsigned int s, unsigned int p, std::vector< std::vector< unsigned int > > &result)
std::ostream & operator<< (std::ostream &os, const QBase &q)
This class implements the Grundmann-Moller quadrature rules for tetrahedra. The GM rules are well-defined for simplices of arbitrary dimension and to any order, but the rules by Dunavant for two-dimensional simplices are in general superior. This is primarily due to the fact that the GM rules contain a significant proportion of negative weights, making them susceptible to round-off error at high-order.
The GM rules are interesting in 3D because they overlap with the conical product rules at higher order while having significantly fewer evaluation points, making them potentially much more efficient. The table below gives a comparison between the number of points in a conical product (CP) rule and the GM rule of equivalent order. The GM rules are defined to be exact for polynomials of degree d=2*s+1, s=0,1,2,3,... The table also gives the percentage of each GM rule's weights which are negative. Although the percentage of negative weights does not grow particularly quickly, the amplification factor (a measure of the effect of round-off) defined as
amp. factor = (1/V) * |w_i|,
where V is the volume of the reference element, does grow quickly. (A rule with all positive has has an amplification factor of 1.0 by definition.)
s | d | N. CP | N. GM | % neg wts | amp. factor ----------------------------------------------------------------- 0 | 1 | | 1 | | 1 | 2-3 | | 5 | | 2 | 4-5 | | 15 | | 3 | 6-7 | | 35 | 31.43 | 11.94 4 | 8-9 | 5^3=125 | 70 | 34.29 | 25.35 5 | 10-11 | 6^3=216 | 126 | 36.51 | 54.14 6 | 12-13 | 7^3=343 | 210 | 38.10 | 116.30 7 | 14-15 | 8^3=512 | 330 | 39.39 | 251.10 8 | 16-17 | 9^3=729 | 495 | 40.40 | 544.68 9 | 18-19 | 10^3=1,000 | 715 | 41.26 | 1186.16 10 | 20-21 | 11^3=1,331 | 1,001 | 41.96 | 2591.97 11 | 22-23 | 12^3=1,728 | 1,365 | 42.56 | 5680.75 ... 16 | 32-33 | 17^3=4,913 | 4,845 | 17 | 34-35 | 18^3=5,832 | 5,985 | <= Cross-over point, CP has fewer points for d >= 34 18 | 36-37 | 19^3=6,859 | 7,315 | ... 21 | 42-43 | 22^3=10,648 | 12,650 |
Reference: Axel Grundmann and Michael M"{o}ller, 'Invariant Integration Formulas for the N-Simplex
by Combinatorial Methods,' SIAM Journal on Numerical Analysis, Volume 15, Number 2, April 1978, pages 282-290.
Reference LGPL Fortran90 code by John Burkardt can be found here: http://people.scs.fsu.edu/~burkardt/f_src/gm_rules/gm_rules.html
Author:
Definition at line 96 of file quadrature_gm.h.
Definition at line 105 of file reference_counter.h.
Definition at line 31 of file quadrature_gm.C.
: QBase(d,o)
{
}
Definition at line 38 of file quadrature_gm.C.
{
}
Definition at line 36 of file quadrature_build.C.
References Utility::string_to_enum< QuadratureType >().
Referenced by InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().
{
return QBase::build (Utility::string_to_enum<QuadratureType> (type),
_dim,
_order);
}
Definition at line 47 of file quadrature_build.C.
References libMeshEnums::FIRST, libMeshEnums::FORTYTHIRD, libMeshEnums::QCLOUGH, libMeshEnums::QGAUSS, libMeshEnums::QJACOBI_1_0, libMeshEnums::QJACOBI_2_0, libMeshEnums::QSIMPSON, libMeshEnums::QTRAP, libMeshEnums::THIRD, and libMeshEnums::TWENTYTHIRD.
{
switch (_qt)
{
case QCLOUGH:
{
#ifdef DEBUG
if (_order > TWENTYTHIRD)
{
std::cout << 'WARNING: Clough quadrature implemented' << std::endl
<< ' up to TWENTYTHIRD order.' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QClough(_dim, _order));
return ap;
}
case QGAUSS:
{
#ifdef DEBUG
if (_order > FORTYTHIRD)
{
std::cout << 'WARNING: Gauss quadrature implemented' << std::endl
<< ' up to FORTYTHIRD order.' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QGauss(_dim, _order));
return ap;
}
case QJACOBI_1_0:
{
#ifdef DEBUG
if (_order > TWENTYTHIRD)
{
std::cout << 'WARNING: Jacobi(1,0) quadrature implemented' << std::endl
<< ' up to TWENTYTHIRD order.' << std::endl;
}
if (_dim > 1)
{
std::cout << 'WARNING: Jacobi(1,0) quadrature implemented' << std::endl
<< ' in 1D only.' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QJacobi(_dim, _order, 1, 0));
return ap;
}
case QJACOBI_2_0:
{
#ifdef DEBUG
if (_order > TWENTYTHIRD)
{
std::cout << 'WARNING: Jacobi(2,0) quadrature implemented' << std::endl
<< ' up to TWENTYTHIRD order.' << std::endl;
}
if (_dim > 1)
{
std::cout << 'WARNING: Jacobi(2,0) quadrature implemented' << std::endl
<< ' in 1D only.' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QJacobi(_dim, _order, 2, 0));
return ap;
}
case QSIMPSON:
{
#ifdef DEBUG
if (_order > THIRD)
{
std::cout << 'WARNING: Simpson rule provides only' << std::endl
<< ' THIRD order!' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QSimpson(_dim));
return ap;
}
case QTRAP:
{
#ifdef DEBUG
if (_order > FIRST)
{
std::cout << 'WARNING: Trapezoidal rule provides only' << std::endl
<< ' FIRST order!' << std::endl;
}
#endif
AutoPtr<QBase> ap(new QTrap(_dim));
return ap;
}
default:
{
std::cerr << 'ERROR: Bad qt=' << _qt << std::endl;
libmesh_error();
}
}
libmesh_error();
AutoPtr<QBase> ap(NULL);
return ap;
}
Definition at line 143 of file quadrature_gm.C.
Referenced by gm_rule().
{
// Clear out results remaining from previous calls
result.clear();
// Allocate storage for a workspace. The workspace will periodically
// be copied into the result container.
std::vector<unsigned int> workspace(p);
// The first result is always (s,0,...,0)
workspace[0] = s;
result.push_back(workspace);
// the value of the first non-zero entry
unsigned int head_value=s;
// When head_index=-1, it refers to 'off the front' of the array. Therefore,
// this needs to be a regular int rather than unsigned. I initially tried to
// do this with head_index unsigned and an else statement below, but then there
// is the special case: (1,0,...,0) which does not work correctly.
int head_index = -1;
// At the end, all the entries will be in the final slot of workspace
while (workspace.back() != s)
{
// Uncomment for debugging
//std::cout << 'previous head_value=' << head_value << ' -> ';
// If the previous head value is still larger than 1, reset the index
// to 'off the front' of the array
if (head_value > 1)
head_index = -1;
// Either move the index onto the front of the array or on to
// the next value.
head_index++;
// Get current value of the head entry
head_value = workspace[head_index];
// Uncomment for debugging
//std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(std::cout, ' '));
//std::cout << ', head_index=' << head_index;
//std::cout << ', head_value=' << head_value << ' -> ';
// Put a zero into the head_index of the array. If head_index==0,
// this will be overwritten in the next line with head_value-1.
workspace[head_index] = 0;
// The initial entry gets the current head value, minus 1.
// If head_value > 1, the next loop iteration will start back
// at workspace[0] again.
libmesh_assert (head_value > 0);
workspace[0] = head_value - 1;
// Increment the head+1 value
workspace[head_index+1] += 1;
// Save this composition in the results
result.push_back(workspace);
// Uncomment for debugging
//std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(std::cout, ' '));
//std::cout<<';
}
}
Definition at line 121 of file quadrature.h.
Referenced by InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), QConical::conical_product_pyramid(), QConical::conical_product_tet(), and QConical::conical_product_tri().
{ return _dim; }
Definition at line 103 of file quadrature.h.
{ return _type; }
Definition at line 45 of file reference_counter.C.
References ReferenceCounter::_counts, and Quality::name().
Referenced by ReferenceCounter::print_info().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::ostringstream out;
out << '
<< ' ----------------------------------------------------------------------------
<< '| Reference count information |
<< ' ---------------------------------------------------------------------------- ;
for (Counts::iterator it = _counts.begin();
it != _counts.end(); ++it)
{
const std::string name(it->first);
const unsigned int creations = it->second.first;
const unsigned int destructions = it->second.second;
out << '| ' << name << ' reference count information:
<< '| Creations: ' << creations << '
<< '| Destructions: ' << destructions << ';
}
out << ' ---------------------------------------------------------------------------- ;
return out.str();
#else
return '';
#endif
}
Definition at line 167 of file quadrature.h.
Referenced by InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().
{ return static_cast<Order>(_order + _p_level); }
Definition at line 109 of file quadrature.h.
{ return _p_level; }
Definition at line 127 of file quadrature.h.
References QBase::_points.
Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QMonomial::init_2D(), QGauss::init_2D(), QClough::init_2D(), QMonomial::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), InfFE< Dim, T_radial, T_map >::reinit(), and REINIT_ERROR().
{ return _points; }
Definition at line 133 of file quadrature.h.
References QBase::_points.
{ return _points; }
Definition at line 138 of file quadrature.h.
References QBase::_weights.
Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QMonomial::init_2D(), QGauss::init_2D(), QClough::init_2D(), QMonomial::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), and REINIT_ERROR().
{ return _weights; }
Definition at line 143 of file quadrature.h.
References QBase::_weights.
{ return _weights; }
Definition at line 46 of file quadrature_gm.C.
References QBase::_points, QBase::_weights, compose_all(), std::max(), and MeshTools::weight().
Referenced by init_3D().
{
// A GM rule of index s can integrate polynomials of degree 2*s+1 exactly
const unsigned int degree = 2*s+1;
// Here we are considering only tetrahedra rules, so dim==3
const unsigned int dim = 3;
// The number of points for rule of index s is
// (dim+1+s)! / (dim+1)! / s!
// In 3D, this is = 1/24 * P_{i=1}^4 (s+i)
const unsigned int n_pts = (s+4)*(s+3)*(s+2)*(s+1) / 24;
//std::cout << 'n_pts=' << n_pts << std::endl;
// Allocate space for points and weights
_points.resize(n_pts);
_weights.resize(n_pts);
// (-1)^i -> This one flips sign at each iteration of the i-loop below.
int one_pm=1;
// Where we store all the integer point compositions/permutations
std::vector<std::vector<unsigned int> > permutations;
// Index into the vector where we should start adding the next round of points/weights
unsigned int offset=0;
// Implement the GM formula 4.1 on page 286 of the paper
for (unsigned int i=0; i<=s; ++i)
{
// Get all the ordered compositions (and their permutations)
// of |beta| = s-i into dim+1=4 parts
compose_all(s-i, dim+1, permutations);
//std::cout << 'n. permutations=' << permutations.size() << std::endl;
for (unsigned int p=0; p<permutations.size(); ++p)
{
// We use the first dim=3 entries of each permutation to
// construct an integration point.
for (unsigned int j=0; j<3; ++j)
_points[offset+p](j) =
static_cast<Real>(2.*permutations[p][j] + 1.) /
static_cast<Real>( degree + dim - 2.*i );
}
// Compute the weight for this i, being careful to avoid overflow.
// This technique is borrowed from Burkardt's code as well.
// Use once for each of the points obtained from the permutations array.
Real weight = one_pm;
// This for loop needs to run for dim, degree, or dim+degree-i iterations,
// whichever is largest.
const unsigned int weight_loop_index =
std::max(dim, std::max(degree, degree+dim-i));
for (unsigned int j=1; j<=weight_loop_index; ++j)
{
if (j <= degree) // Accumulate (d+n-2i)^d term
weight *= static_cast<Real>(degree+dim-2*i);
if (j <= 2*s) // Accumulate 2^{-2s}
weight *= 0.5;
if (j <= i) // Accumulate (i!)^{-1}
weight /= static_cast<Real>(j);
if (j <= degree+dim-i) // Accumulate ( (d+n-i)! )^{-1}
weight /= static_cast<Real>(j);
}
// This is the weight for each of the points computed previously
for (unsigned int j=0; j<permutations.size(); ++j)
_weights[offset+j] = weight;
// Change sign for next iteration
one_pm = -one_pm;
// Update offset for the next set of points
offset += permutations.size();
}
}
Definition at line 149 of file reference_counter.h.
References ReferenceCounter::_counts, Quality::name(), and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.first++;
}
Definition at line 167 of file reference_counter.h.
References ReferenceCounter::_counts, Quality::name(), and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.second++;
}
Definition at line 26 of file quadrature.C.
References QBase::init_0D(), QBase::init_1D(), and QBase::init_2D().
Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), QGauss::QGauss(), QJacobi::QJacobi(), QSimpson::QSimpson(), QTrap::QTrap(), InfFE< Dim, T_radial, T_map >::reinit(), and REINIT_ERROR().
{
// check to see if we have already
// done the work for this quadrature rule
if (t == _type && p == _p_level)
return;
else
{
_type = t;
_p_level = p;
}
switch(_dim)
{
case 0:
this->init_0D(_type,_p_level);
return;
case 1:
this->init_1D(_type,_p_level);
return;
case 2:
this->init_2D(_type,_p_level);
return;
case 3:
this->init_3D(_type,_p_level);
return;
default:
libmesh_error();
}
}
Definition at line 70 of file quadrature.C.
References QBase::_points, and QBase::_weights.
Referenced by QBase::init().
{
_points.resize(1);
_weights.resize(1);
_points[0] = Point(0.);
_weights[0] = 1.0;
}
Implements QBase.
Definition at line 119 of file quadrature_gm.h.
{
// See about making this non-pure virtual in the base class
libmesh_error();
}
Reimplemented in QClough, QConical, QGauss, QGrid, QMonomial, QSimpson, and QTrap.
Definition at line 234 of file quadrature.h.
Referenced by QBase::init().
{}
Definition at line 27 of file quadrature_gm_3D.C.
References QBase::allow_rules_with_negative_weights, gm_rule(), libMeshEnums::TET10, and libMeshEnums::TET4.
{
// Nearly all GM rules contain negative weights, so if you are not
// allowing rules with negative weights, we cannot continue!
if (!allow_rules_with_negative_weights)
{
std::cerr << 'You requested a Grundmann-Moller rule but
<< 'are not allowing rules with negative weights!
<< 'Either select a different quadrature class or
<< 'set allow_rules_with_negative_weights==true.'
<< std::endl;
libmesh_error();
}
switch (_type)
{
case TET4:
case TET10:
{
// Untested above _order=23 but should work...
gm_rule( (_order + 2*p)/2 );
return;
} // end case TET4, TET10
//---------------------------------------------
// Unsupported element type
default:
{
std::cerr << 'ERROR: Unsupported element type: ' << _type << std::endl;
libmesh_error();
}
} // end switch (_type)
// We must have returned or errored-out by this point. If not,
// throw an error now.
libmesh_error();
return;
}
Definition at line 76 of file reference_counter.h.
References ReferenceCounter::_n_objects.
Referenced by System::read_serialized_blocked_dof_objects(), and System::write_serialized_blocked_dof_objects().
{ return _n_objects; }
Definition at line 115 of file quadrature.h.
References QBase::_points.
Referenced by FEBase::coarsened_dof_values(), QConical::conical_product_pyramid(), QConical::conical_product_tet(), QConical::conical_product_tri(), FEMSystem::eulerian_residual(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), FEMSystem::mass_residual(), and QBase::print_info().
{ libmesh_assert (!_points.empty()); return _points.size(); }
Definition at line 83 of file reference_counter.C.
References ReferenceCounter::get_info().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::cout << ReferenceCounter::get_info();
#endif
}
Definition at line 350 of file quadrature.h.
References QBase::_points, QBase::_weights, QBase::n_points(), and QBase::qp().
Referenced by operator<<().
{
libmesh_assert(!_points.empty());
libmesh_assert(!_weights.empty());
os << 'N_Q_Points=' << this->n_points() << std::endl << std::endl;
for (unsigned int qp=0; qp<this->n_points(); qp++)
{
os << ' Point ' << qp << ':
<< ' '
<< _points[qp]
<< ' Weight:'
<< ' w=' << _weights[qp] << ' << std::endl;
}
}
Definition at line 148 of file quadrature.h.
References QBase::_points.
Referenced by QConical::conical_product_pyramid(), QConical::conical_product_tet(), QConical::conical_product_tri(), and QBase::print_info().
{ libmesh_assert (i < _points.size()); return _points[i]; }
Definition at line 81 of file quadrature.C.
References QBase::_points, and QBase::_weights.
Referenced by QConical::conical_product_tet(), and QConical::conical_product_tri().
{
// Make sure we are in 1D
libmesh_assert(_dim == 1);
// Make sure that we have sane ranges
libmesh_assert(new_range.second > new_range.first);
libmesh_assert(old_range.second > old_range.first);
// Make sure there are some points
libmesh_assert(_points.size() > 0);
// We're mapping from old_range -> new_range
for (unsigned int i=0; i<_points.size(); i++)
{
_points[i](0) =
(_points[i](0) - old_range.first) *
(new_range.second - new_range.first) /
(old_range.second - old_range.first) +
new_range.first;
}
// Compute the scale factor and scale the weights
const Real scfact = (new_range.second - new_range.first) /
(old_range.second - old_range.first);
for (unsigned int i=0; i<_points.size(); i++)
_weights[i] *= scfact;
}
Implements QBase.
Definition at line 114 of file quadrature_gm.h.
References libMeshEnums::QGRUNDMANN_MOLLER.
{ return QGRUNDMANN_MOLLER; }
Definition at line 154 of file quadrature.h.
References QBase::_weights.
Referenced by QConical::conical_product_pyramid(), QConical::conical_product_tet(), QConical::conical_product_tri(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), and QMonomial::stroud_rule().
{ libmesh_assert (i < _weights.size()); return _weights[i]; }
Definition at line 196 of file quadrature.C.
{
q.print_info(os);
return os;
}
Definition at line 110 of file reference_counter.h.
Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().
Definition at line 123 of file reference_counter.h.
Definition at line 118 of file reference_counter.h.
Referenced by ReferenceCounter::n_objects(), ReferenceCounter::ReferenceCounter(), and ReferenceCounter::~ReferenceCounter().
Definition at line 320 of file quadrature.h.
Referenced by QConical::conical_product_pyramid(), QConical::conical_product_tet(), QConical::conical_product_tri(), QGauss::dunavant_rule(), QGauss::dunavant_rule2(), QBase::get_points(), gm_rule(), QBase::init_0D(), QTrap::init_1D(), QSimpson::init_1D(), QJacobi::init_1D(), QGrid::init_1D(), QGauss::init_1D(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), QBase::n_points(), QBase::print_info(), QBase::qp(), QBase::scale(), QMonomial::stroud_rule(), and QMonomial::wissmann_rule().
Definition at line 325 of file quadrature.h.
Referenced by QConical::conical_product_pyramid(), QConical::conical_product_tet(), QConical::conical_product_tri(), QGauss::dunavant_rule(), QGauss::dunavant_rule2(), QBase::get_weights(), gm_rule(), QBase::init_0D(), QTrap::init_1D(), QSimpson::init_1D(), QJacobi::init_1D(), QGrid::init_1D(), QGauss::init_1D(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), QBase::print_info(), QBase::scale(), QMonomial::stroud_rule(), QBase::w(), and QMonomial::wissmann_rule().
Negative weights typically appear in Gaussian quadrature rules over three-dimensional elements. 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 203 of file quadrature.h.
Referenced by QMonomial::init_3D(), init_3D(), and QGauss::init_3D().
Generated automatically by Doxygen for libMesh from the source code.