$extrastylesheet
#include <quadrature_conical.h>

Public Member Functions | |
| QConical (const unsigned int _dim, const Order _order=INVALID_ORDER) | |
| ~QConical () | |
| 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 |
| virtual void | init (const ElemType type=INVALID_ELEM, unsigned int p_level=0) |
| virtual void | init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0) |
| Order | get_order () const |
| void | print_info (std::ostream &os=libMesh::out) const |
| void | scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range) |
| virtual bool | shapes_need_reinit () |
Static Public Member Functions | |
| static UniquePtr< QBase > | build (const std::string &name, const unsigned int _dim, const Order _order=INVALID_ORDER) |
| static UniquePtr< QBase > | build (const QuadratureType _qt, const unsigned int _dim, const Order _order=INVALID_ORDER) |
| static void | print_info (std::ostream &out=libMesh::out) |
| static std::string | get_info () |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
Public Attributes | |
| bool | allow_rules_with_negative_weights |
Protected Types | |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| virtual void | init_0D (const ElemType type=INVALID_ELEM, unsigned int p_level=0) |
| libmesh_error_msg ("ERROR: Seems as if this quadrature rule \nis not implemented for 2D.") | |
| void | tensor_product_hex (const QBase &q1D) |
| void | tensor_product_prism (const QBase &q1D, const QBase &q2D) |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| const unsigned int | _dim |
| const Order | _order |
| ElemType | _type |
| unsigned int | _p_level |
| std::vector< Point > | _points |
| std::vector< Real > | _weights |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
| static bool | _enable_print_counter = true |
Private Member Functions | |
| void | init_1D (const ElemType, unsigned int=0) |
| void | init_2D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0) |
| void | init_3D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0) |
| void | conical_product_tri (unsigned int p) |
| void | conical_product_tet (unsigned int p) |
| void | conical_product_pyramid (unsigned int p) |
Friends | |
| std::ostream & | operator<< (std::ostream &os, const QBase &q) |
This class implements the so-called conical product quadrature rules for Tri and Tet elements. These rules are generally non-optimal in the number of evaluation points, but have the nice property of having all positive weights and being well-defined to any order for which their underlying 1D Gauss and Jacobi quadrature rules are available.
The construction of these rules is given by e.g.
Stroud, A.H. "Approximate Calculation of Multiple Integrals.", 1972
Definition at line 43 of file quadrature_conical.h.
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
| libMesh::QConical::QConical | ( | const unsigned int | _dim, |
| const Order | _order = INVALID_ORDER |
||
| ) |
Constructor. Declares the order of the quadrature rule.
Definition at line 35 of file quadrature_conical.C.
: QBase(d,o) { }
| UniquePtr< QBase > libMesh::QBase::build | ( | const std::string & | name, |
| const unsigned int | _dim, | ||
| const Order | _order = INVALID_ORDER |
||
| ) | [static, inherited] |
Builds a specific quadrature rule, identified through the name string. An UniquePtr<QBase> is returned to prevent a memory leak. This way the user need not remember to delete the object. Enables run-time decision of the quadrature rule. The input parameter name must be mappable through the Utility::string_to_enum<>() function.
Definition at line 42 of file quadrature_build.C.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().
{
return QBase::build (Utility::string_to_enum<QuadratureType> (type),
_dim,
_order);
}
| UniquePtr< QBase > libMesh::QBase::build | ( | const QuadratureType | _qt, |
| const unsigned int | _dim, | ||
| const Order | _order = INVALID_ORDER |
||
| ) | [static, inherited] |
Builds a specific quadrature rule, identified through the QuadratureType. An UniquePtr<QBase> is returned to prevent a memory leak. This way the user need not remember to delete the object. Enables run-time decision of the quadrature rule.
Definition at line 53 of file quadrature_build.C.
References libMesh::FIRST, libMesh::FORTYTHIRD, libMesh::QBase::libmesh_error_msg(), libMesh::out, libMesh::QCLOUGH, libMesh::QCONICAL, libMesh::QGAUSS, libMesh::QGAUSS_LOBATTO, libMesh::QGRID, libMesh::QGRUNDMANN_MOLLER, libMesh::QJACOBI_1_0, libMesh::QJACOBI_2_0, libMesh::QMONOMIAL, libMesh::QSIMPSON, libMesh::QTRAP, libMesh::THIRD, and libMesh::TWENTYTHIRD.
{
switch (_qt)
{
case QCLOUGH:
{
#ifdef DEBUG
if (_order > TWENTYTHIRD)
{
libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
<< " up to TWENTYTHIRD order." << std::endl;
}
#endif
return UniquePtr<QBase>(new QClough(_dim, _order));
}
case QGAUSS:
{
#ifdef DEBUG
if (_order > FORTYTHIRD)
{
libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
<< " up to FORTYTHIRD order." << std::endl;
}
#endif
return UniquePtr<QBase>(new QGauss(_dim, _order));
}
case QJACOBI_1_0:
{
#ifdef DEBUG
if (_order > FORTYTHIRD)
{
libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
<< " up to FORTYTHIRD order." << std::endl;
}
if (_dim > 1)
{
libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
<< " in 1D only." << std::endl;
}
#endif
return UniquePtr<QBase>(new QJacobi(_dim, _order, 1, 0));
}
case QJACOBI_2_0:
{
#ifdef DEBUG
if (_order > FORTYTHIRD)
{
libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
<< " up to FORTYTHIRD order." << std::endl;
}
if (_dim > 1)
{
libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
<< " in 1D only." << std::endl;
}
#endif
return UniquePtr<QBase>(new QJacobi(_dim, _order, 2, 0));
}
case QSIMPSON:
{
#ifdef DEBUG
if (_order > THIRD)
{
libMesh::out << "WARNING: Simpson rule provides only" << std::endl
<< " THIRD order!" << std::endl;
}
#endif
return UniquePtr<QBase>(new QSimpson(_dim));
}
case QTRAP:
{
#ifdef DEBUG
if (_order > FIRST)
{
libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
<< " FIRST order!" << std::endl;
}
#endif
return UniquePtr<QBase>(new QTrap(_dim));
}
case QGRID:
return UniquePtr<QBase>(new QGrid(_dim, _order));
case QGRUNDMANN_MOLLER:
return UniquePtr<QBase>(new QGrundmann_Moller(_dim, _order));
case QMONOMIAL:
return UniquePtr<QBase>(new QMonomial(_dim, _order));
case QGAUSS_LOBATTO:
return UniquePtr<QBase>(new QGaussLobatto(_dim, _order));
case QCONICAL:
return UniquePtr<QBase>(new QConical(_dim, _order));
default:
libmesh_error_msg("ERROR: Bad qt=" << _qt);
}
libmesh_error_msg("We'll never get here!");
return UniquePtr<QBase>();
}
| void libMesh::QConical::conical_product_pyramid | ( | unsigned int | p | ) | [private] |
Implementation of conical product rule for a Pyramid in 3D of order = _order+2*p.
Definition at line 180 of file quadrature_conical.C.
References libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::get_dim(), libMesh::QBase::n_points(), libMesh::QBase::qp(), libMesh::Real, and libMesh::QBase::w().
Referenced by init_3D().
{
// Be sure the underlying rule object was built with the same dimension as the
// rule we are about to construct.
libmesh_assert_equal_to (this->get_dim(), 3);
QGauss gauss1D(1,static_cast<Order>(_order+2*p));
QJacobi jac1D(1,static_cast<Order>(_order+2*p),2,0);
// These rules should have the same number of points
libmesh_assert_equal_to (gauss1D.n_points(), jac1D.n_points());
// Save the number of points as a convenient variable
const unsigned int np = gauss1D.n_points();
// Resize the points and weights vectors
_points.resize(np * np * np);
_weights.resize(np * np * np);
// Compute the conical product
unsigned int q = 0;
for (unsigned int i=0; i<np; ++i)
for (unsigned int j=0; j<np; ++j)
for (unsigned int k=0; k<np; ++k, ++q)
{
const Real xi=gauss1D.qp(i)(0);
const Real yj=gauss1D.qp(j)(0);
const Real zk=jac1D.qp(k)(0);
_points[q](0) = (1.-zk) * xi;
_points[q](1) = (1.-zk) * yj;
_points[q](2) = zk;
_weights[q] = gauss1D.w(i) * gauss1D.w(j) * jac1D.w(k);
}
}
| void libMesh::QConical::conical_product_tet | ( | unsigned int | p | ) | [private] |
Implementation of conical product rule for a Tet in 3D of order = _order+2*p.
Definition at line 103 of file quadrature_conical.C.
References libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::get_dim(), libMesh::QBase::n_points(), libMesh::QBase::qp(), libMesh::QBase::scale(), and libMesh::QBase::w().
Referenced by init_3D().
{
// Be sure the underlying rule object was built with the same dimension as the
// rule we are about to construct.
libmesh_assert_equal_to (this->get_dim(), 3);
QGauss gauss1D(1,static_cast<Order>(_order+2*p));
QJacobi jacA1D(1,static_cast<Order>(_order+2*p),1,0);
QJacobi jacB1D(1,static_cast<Order>(_order+2*p),2,0);
// The Gauss rule needs to be scaled to [0,1]
std::pair<Real, Real> old_range(-1.0L, 1.0L);
std::pair<Real, Real> new_range( 0.0L, 1.0L);
gauss1D.scale(old_range,
new_range);
// Now construct the points and weights for the conical product rule.
// All rules should have the same number of points
libmesh_assert_equal_to (gauss1D.n_points(), jacA1D.n_points());
libmesh_assert_equal_to (jacA1D.n_points(), jacB1D.n_points());
// Save the number of points as a convenient variable
const unsigned int np = gauss1D.n_points();
// All rules should be between x=0 and x=1
libmesh_assert_greater_equal (gauss1D.qp(0)(0), 0.0);
libmesh_assert_less_equal (gauss1D.qp(np-1)(0), 1.0);
libmesh_assert_greater_equal (jacA1D.qp(0)(0), 0.0);
libmesh_assert_less_equal (jacA1D.qp(np-1)(0), 1.0);
libmesh_assert_greater_equal (jacB1D.qp(0)(0), 0.0);
libmesh_assert_less_equal (jacB1D.qp(np-1)(0), 1.0);
// Resize the points and weights vectors
_points.resize(np * np * np);
_weights.resize(np * np * np);
// Compute the conical product
unsigned int gp = 0;
for (unsigned int i=0; i<np; i++)
for (unsigned int j=0; j<np; j++)
for (unsigned int k=0; k<np; k++)
{
_points[gp](0) = jacB1D.qp(k)(0); //t[k];
_points[gp](1) = jacA1D.qp(j)(0) * (1.-jacB1D.qp(k)(0)); //s[j]*(1.-t[k]);
_points[gp](2) = gauss1D.qp(i)(0) * (1.-jacA1D.qp(j)(0)) * (1.-jacB1D.qp(k)(0)); //r[i]*(1.-s[j])*(1.-t[k]);
_weights[gp] = gauss1D.w(i) * jacA1D.w(j) * jacB1D.w(k); //A[i]*B[j]*C[k];
gp++;
}
}
| void libMesh::QConical::conical_product_tri | ( | unsigned int | p | ) | [private] |
Implementation of conical product rule for a Tri in 2D of order = _order+2*p.
Definition at line 52 of file quadrature_conical.C.
References libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::get_dim(), libMesh::QBase::n_points(), libMesh::QBase::qp(), libMesh::QBase::scale(), and libMesh::QBase::w().
Referenced by init_2D().
{
// Be sure the underlying rule object was built with the same dimension as the
// rule we are about to construct.
libmesh_assert_equal_to (this->get_dim(), 2);
QGauss gauss1D(1,static_cast<Order>(_order+2*p));
QJacobi jac1D(1,static_cast<Order>(_order+2*p),1,0);
// The Gauss rule needs to be scaled to [0,1]
std::pair<Real, Real> old_range(-1.0L, 1.0L);
std::pair<Real, Real> new_range( 0.0L, 1.0L);
gauss1D.scale(old_range,
new_range);
// Now construct the points and weights for the conical product rule.
// Both rules should have the same number of points.
libmesh_assert_equal_to (gauss1D.n_points(), jac1D.n_points());
// Save the number of points as a convenient variable
const unsigned int np = gauss1D.n_points();
// Both rules should be between x=0 and x=1
libmesh_assert_greater_equal (gauss1D.qp(0)(0), 0.0);
libmesh_assert_less_equal (gauss1D.qp(np-1)(0), 1.0);
libmesh_assert_greater_equal (jac1D.qp(0)(0), 0.0);
libmesh_assert_less_equal (jac1D.qp(np-1)(0), 1.0);
// Resize the points and weights vectors
_points.resize(np * np);
_weights.resize(np * np);
// Compute the conical product
unsigned int gp = 0;
for (unsigned int i=0; i<np; i++)
for (unsigned int j=0; j<np; j++)
{
_points[gp](0) = jac1D.qp(j)(0); //s[j];
_points[gp](1) = gauss1D.qp(i)(0) * (1.-jac1D.qp(j)(0)); //r[i]*(1.-s[j]);
_weights[gp] = gauss1D.w(i) * jac1D.w(j); //A[i]*B[j];
gp++;
}
}
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::LibMeshInit::LibMeshInit().
{
_enable_print_counter = false;
return;
}
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
{
_enable_print_counter = true;
return;
}
| unsigned int libMesh::QBase::get_dim | ( | ) | const [inline, inherited] |
Definition at line 125 of file quadrature.h.
References libMesh::QBase::_dim.
Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), conical_product_pyramid(), conical_product_tet(), and conical_product_tri().
{ return _dim; }
| ElemType libMesh::QBase::get_elem_type | ( | ) | const [inline, inherited] |
Definition at line 106 of file quadrature.h.
References libMesh::QBase::_type.
{ return _type; }
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
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().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::ostringstream oss;
oss << '\n'
<< " ---------------------------------------------------------------------------- \n"
<< "| Reference count information |\n"
<< " ---------------------------------------------------------------------------- \n";
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;
oss << "| " << name << " reference count information:\n"
<< "| Creations: " << creations << '\n'
<< "| Destructions: " << destructions << '\n';
}
oss << " ---------------------------------------------------------------------------- \n";
return oss.str();
#else
return "";
#endif
}
| Order libMesh::QBase::get_order | ( | ) | const [inline, inherited] |
Definition at line 183 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().
| unsigned int libMesh::QBase::get_p_level | ( | ) | const [inline, inherited] |
Definition at line 112 of file quadrature.h.
References libMesh::QBase::_p_level.
{ return _p_level; }
| const std::vector<Point>& libMesh::QBase::get_points | ( | ) | const [inline, inherited] |
std::vector containing the quadrature point locations on a reference object. Definition at line 131 of file quadrature.h.
References libMesh::QBase::_points.
Referenced by libMesh::QClough::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QMonomial::init_3D().
{ return _points; }
| std::vector<Point>& libMesh::QBase::get_points | ( | ) | [inline, inherited] |
std::vector containing the quadrature point locations on a reference object as a writeable reference. Definition at line 137 of file quadrature.h.
References libMesh::QBase::_points.
{ return _points; }
| const std::vector<Real>& libMesh::QBase::get_weights | ( | ) | const [inline, inherited] |
std::vector containing the quadrature weights. Definition at line 142 of file quadrature.h.
References libMesh::QBase::_weights.
Referenced by libMesh::QClough::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QMonomial::init_3D().
{ return _weights; }
| std::vector<Real>& libMesh::QBase::get_weights | ( | ) | [inline, inherited] |
std::vector containing the quadrature weights. Definition at line 147 of file quadrature.h.
References libMesh::QBase::_weights.
{ return _weights; }
| void libMesh::ReferenceCounter::increment_constructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.first++;
}
| void libMesh::ReferenceCounter::increment_destructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.second++;
}
| void libMesh::QBase::init | ( | const ElemType | type = INVALID_ELEM, |
| unsigned int | p_level = 0 |
||
| ) | [virtual, inherited] |
Initializes the data structures to contain a quadrature rule for an object of type type.
Definition at line 28 of file quadrature.C.
References libMesh::QBase::_dim, libMesh::QBase::_p_level, libMesh::QBase::_type, libMesh::QBase::init_0D(), libMesh::QBase::init_1D(), libMesh::QBase::init_2D(), libMesh::QBase::init_3D(), and libMesh::QBase::libmesh_error_msg().
Referenced by libMesh::QBase::init(), libMesh::QClough::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QJacobi::QJacobi(), libMesh::QSimpson::QSimpson(), and libMesh::QTrap::QTrap().
{
// 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_msg("Invalid dimension _dim = " << _dim);
}
}
| void libMesh::QBase::init | ( | const Elem & | elem, |
| const std::vector< Real > & | vertex_distance_func, | ||
| unsigned int | p_level = 0 |
||
| ) | [virtual, inherited] |
Initializes the data structures for a specific, potentially cut element. The array vertex_distance_func contains vertex values of a signed distance function that cuts the element. This interface is indended to be extended by derived classes that can cut the element into subelements, for example, and constuct a composite quadrature rule for the cut element.
Definition at line 72 of file quadrature.C.
References libMesh::QBase::init(), and libMesh::Elem::type().
{
// dispatch generic implementation
this->init(elem.type(), p_level);
}
| void libMesh::QBase::init_0D | ( | const ElemType | type = INVALID_ELEM, |
| unsigned int | p_level = 0 |
||
| ) | [protected, virtual, inherited] |
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate values. Generally this is just one point with weight 1.
Definition at line 82 of file quadrature.C.
References libMesh::QBase::_points, and libMesh::QBase::_weights.
Referenced by libMesh::QBase::init().
| void libMesh::QConical::init_1D | ( | const ElemType | type, |
| unsigned | p_level = 0 |
||
| ) | [inline, private, virtual] |
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 65 of file quadrature_conical.h.
{
// See about making this non-pure virtual in the base class
libmesh_not_implemented();
}
| void libMesh::QConical::init_2D | ( | const ElemType | _type = INVALID_ELEM, |
| unsigned int | p_level = 0 |
||
| ) | [private, virtual] |
The conical product rules are defined in 2D only for Tris.
Reimplemented from libMesh::QBase.
Definition at line 27 of file quadrature_conical_2D.C.
References conical_product_tri(), libMesh::QBase::libmesh_error_msg(), libMesh::TRI3, and libMesh::TRI6.
{
switch (type_in)
{
case TRI3:
case TRI6:
{
this->conical_product_tri(p);
return;
} // end case TRI3, TRI6
//---------------------------------------------
// Unsupported element type
default:
libmesh_error_msg("ERROR: Unsupported element type: " << type_in);
} // end switch (type_in)
}
| void libMesh::QConical::init_3D | ( | const ElemType | _type = INVALID_ELEM, |
| unsigned int | p_level = 0 |
||
| ) | [private, virtual] |
The conical product rules are defined in 3D only for Tets.
Reimplemented from libMesh::QBase.
Definition at line 27 of file quadrature_conical_3D.C.
References conical_product_pyramid(), conical_product_tet(), libMesh::QBase::libmesh_error_msg(), libMesh::PYRAMID13, libMesh::PYRAMID14, libMesh::PYRAMID5, libMesh::TET10, and libMesh::TET4.
{
switch (type_in)
{
case TET4:
case TET10:
{
this->conical_product_tet(p);
return;
} // end case TET4, TET10
case PYRAMID5:
case PYRAMID13:
case PYRAMID14:
{
this->conical_product_pyramid(p);
return;
} // end case PYRAMID5
//---------------------------------------------
// Unsupported element type
default:
libmesh_error_msg("ERROR: Unsupported element type: " << type_in);
} // end switch (type_in)
}
| libMesh::QBase::libmesh_error_msg | ( | "ERROR: Seems as if this quadrature rule \nis not implemented for 2D." | ) | [protected, inherited] |
Referenced by libMesh::QBase::build(), libMesh::QGauss::dunavant_rule(), libMesh::QGauss::dunavant_rule2(), libMesh::QBase::init(), libMesh::QGaussLobatto::init_1D(), libMesh::QGauss::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), init_2D(), libMesh::QGrid::init_2D(), libMesh::QGrundmann_Moller::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QClough::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::init_3D(), init_3D(), libMesh::QGrundmann_Moller::init_3D(), libMesh::QGauss::keast_rule(), libMesh::QMonomial::kim_rule(), libMesh::QMonomial::stroud_rule(), and libMesh::QJacobi::type().
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{ return _n_objects; }
| unsigned int libMesh::QBase::n_points | ( | ) | const [inline, inherited] |
Definition at line 118 of file quadrature.h.
References libMesh::QBase::_points, and libMesh::libmesh_assert().
Referenced by libMesh::ExactSolution::_compute_error(), conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::ProjectFEMSolution::operator()(), libMesh::QBase::print_info(), libMesh::QBase::tensor_product_hex(), and libMesh::QBase::tensor_product_prism().
{ libmesh_assert (!_points.empty());
return cast_int<unsigned int>(_points.size()); }
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out | ) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{
if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
}
| void libMesh::QBase::print_info | ( | std::ostream & | os = libMesh::out | ) | const [inline, inherited] |
Prints information relevant to the quadrature rule, by default to libMesh::out.
Definition at line 372 of file quadrature.h.
References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::libmesh_assert(), libMesh::QBase::n_points(), and libMesh::Real.
Referenced by libMesh::operator<<().
{
libmesh_assert(!_points.empty());
libmesh_assert(!_weights.empty());
Real summed_weights=0;
os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
{
os << " Point " << qpoint << ":\n"
<< " "
<< _points[qpoint]
<< "\n Weight:\n "
<< " w=" << _weights[qpoint] << "\n" << std::endl;
summed_weights += _weights[qpoint];
}
os << "Summed Weights: " << summed_weights << std::endl;
}
| Point libMesh::QBase::qp | ( | const unsigned int | i | ) | const [inline, inherited] |
quadrature point on the reference object. Definition at line 152 of file quadrature.h.
References libMesh::QBase::_points.
Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::QBase::tensor_product_hex(), and libMesh::QBase::tensor_product_prism().
| void libMesh::QBase::scale | ( | std::pair< Real, Real > | old_range, |
| std::pair< Real, Real > | new_range | ||
| ) | [inherited] |
Maps the points of a 1D interval quadrature rule (typically [-1,1]) to any other 1D interval (typically [0,1]) and scales the weights accordingly. The quadrature rule will be mapped from the entries of old_range to the entries of new_range.
Definition at line 93 of file quadrature.C.
References libMesh::QBase::_dim, libMesh::QBase::_points, libMesh::QBase::_weights, and libMesh::Real.
Referenced by conical_product_tet(), and conical_product_tri().
{
// Make sure we are in 1D
libmesh_assert_equal_to (_dim, 1);
Real
h_new = new_range.second - new_range.first,
h_old = old_range.second - old_range.first;
// Make sure that we have sane ranges
libmesh_assert_greater (h_new, 0.);
libmesh_assert_greater (h_old, 0.);
// Make sure there are some points
libmesh_assert_greater (_points.size(), 0);
// Compute the scale factor
Real scfact = h_new/h_old;
// We're mapping from old_range -> new_range
for (unsigned int i=0; i<_points.size(); i++)
{
_points[i](0) = new_range.first +
(_points[i](0) - old_range.first) * scfact;
// Scale the weights
_weights[i] *= scfact;
}
}
| virtual bool libMesh::QBase::shapes_need_reinit | ( | ) | [inline, virtual, inherited] |
Returns true if the shape functions need to be recalculated.
This can happen if the number of points or their positions change.
By default this will return false.
Definition at line 212 of file quadrature.h.
{ return false; }
| void libMesh::QBase::tensor_product_hex | ( | const QBase & | q1D | ) | [protected, inherited] |
Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D. Used in the init_3D routines for hexahedral element types.
Definition at line 154 of file quadrature.C.
References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::n_points(), libMesh::QBase::qp(), and libMesh::QBase::w().
Referenced by libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGrid::init_3D().
{
const unsigned int np = q1D.n_points();
_points.resize(np * np * np);
_weights.resize(np * np * np);
unsigned int q=0;
for (unsigned int k=0; k<np; k++)
for (unsigned int j=0; j<np; j++)
for (unsigned int i=0; i<np; i++)
{
_points[q](0) = q1D.qp(i)(0);
_points[q](1) = q1D.qp(j)(0);
_points[q](2) = q1D.qp(k)(0);
_weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);
q++;
}
}
| void libMesh::QBase::tensor_product_prism | ( | const QBase & | q1D, |
| const QBase & | q2D | ||
| ) | [protected, inherited] |
Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule. Used in the init_3D routines for prismatic element types.
Definition at line 181 of file quadrature.C.
References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::QBase::n_points(), libMesh::QBase::qp(), and libMesh::QBase::w().
Referenced by libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGrid::init_3D().
{
const unsigned int n_points1D = q1D.n_points();
const unsigned int n_points2D = q2D.n_points();
_points.resize (n_points1D * n_points2D);
_weights.resize (n_points1D * n_points2D);
unsigned int q=0;
for (unsigned int j=0; j<n_points1D; j++)
for (unsigned int i=0; i<n_points2D; i++)
{
_points[q](0) = q2D.qp(i)(0);
_points[q](1) = q2D.qp(i)(1);
_points[q](2) = q1D.qp(j)(0);
_weights[q] = q2D.w(i) * q1D.w(j);
q++;
}
}
| QuadratureType libMesh::QConical::type | ( | ) | const [inline, virtual] |
Implements libMesh::QBase.
Definition at line 61 of file quadrature_conical.h.
References libMesh::QCONICAL.
{ return QCONICAL; }
| Real libMesh::QBase::w | ( | const unsigned int | i | ) | const [inline, inherited] |
quadrature weight. Definition at line 158 of file quadrature.h.
References libMesh::QBase::_weights.
Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::QBase::tensor_product_hex(), and libMesh::QBase::tensor_product_prism().
| std::ostream& operator<< | ( | std::ostream & | os, |
| const QBase & | q | ||
| ) | [friend, inherited] |
Same as above, but allows you to use the stream syntax.
Definition at line 208 of file quadrature.C.
{
q.print_info(os);
return os;
}
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
const unsigned int libMesh::QBase::_dim [protected, inherited] |
The dimension
Definition at line 319 of file quadrature.h.
Referenced by libMesh::QBase::get_dim(), libMesh::QBase::init(), libMesh::QGauss::QGauss(), libMesh::QGaussLobatto::QGaussLobatto(), libMesh::QJacobi::QJacobi(), libMesh::QSimpson::QSimpson(), libMesh::QTrap::QTrap(), and libMesh::QBase::scale().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
const Order libMesh::QBase::_order [protected, inherited] |
The order of the quadrature rule.
Definition at line 324 of file quadrature.h.
Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::QBase::get_order(), libMesh::QGaussLobatto::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QGrid::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QGrid::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGrundmann_Moller::init_2D(), libMesh::QGaussLobatto::init_3D(), libMesh::QGauss::init_3D(), libMesh::QGrid::init_3D(), libMesh::QMonomial::init_3D(), and libMesh::QGrundmann_Moller::init_3D().
unsigned int libMesh::QBase::_p_level [protected, inherited] |
The p level of element for which the current values have been computed.
Definition at line 336 of file quadrature.h.
Referenced by libMesh::QBase::get_order(), libMesh::QBase::get_p_level(), and libMesh::QBase::init().
std::vector<Point> libMesh::QBase::_points [protected, inherited] |
The reference element locations of the quadrature points.
Definition at line 342 of file quadrature.h.
Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::QGauss::dunavant_rule(), libMesh::QGauss::dunavant_rule2(), libMesh::QBase::get_points(), libMesh::QGrundmann_Moller::gm_rule(), libMesh::QBase::init_0D(), libMesh::QGaussLobatto::init_1D(), libMesh::QTrap::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QGrid::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::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(), and libMesh::QMonomial::wissmann_rule().
ElemType libMesh::QBase::_type [protected, inherited] |
The type of element for which the current values have been computed.
Definition at line 330 of file quadrature.h.
Referenced by libMesh::QBase::get_elem_type(), and libMesh::QBase::init().
std::vector<Real> libMesh::QBase::_weights [protected, inherited] |
The value of the quadrature weights.
Definition at line 347 of file quadrature.h.
Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), libMesh::QGauss::dunavant_rule(), libMesh::QGauss::dunavant_rule2(), libMesh::QBase::get_weights(), libMesh::QGrundmann_Moller::gm_rule(), libMesh::QBase::init_0D(), libMesh::QGaussLobatto::init_1D(), libMesh::QTrap::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QGrid::init_1D(), libMesh::QJacobi::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QGrid::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QGrid::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::w(), and libMesh::QMonomial::wissmann_rule().
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 ONLY use (potentially) safer but more expensive rules with all positive weights.
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 229 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().