$extrastylesheet
#include <quadrature_jacobi.h>

Public Member Functions | |
| QJacobi (const unsigned int _dim, const Order _order=INVALID_ORDER, const unsigned int a=1, const unsigned int b=0) | |
| ~QJacobi () | |
| 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) |
| virtual void | init_2D (const ElemType, unsigned int=0) |
| libmesh_error_msg ("ERROR: Seems as if this quadrature rule \nis not implemented for 2D.") | |
| virtual void | init_3D (const ElemType, unsigned int=0) |
| 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 _type=INVALID_ELEM, unsigned int p_level=0) |
Private Attributes | |
| const unsigned int | _alpha |
| const unsigned int | _beta |
Friends | |
| std::ostream & | operator<< (std::ostream &os, const QBase &q) |
This class implements two (for now) Jacobi-Gauss quadrature rules. These rules have the same order of accuracy as the normal Gauss quadrature rules, but instead of integrating a weight function w(x)=1, they integrate w(x)=(1-x)^alpha * (1+x)^beta. The reason that they are useful is that they allow you to implement conical product rules for triangles and tetrahedra. Although these product rules are non-optimal (use more points than necessary) they are automatically constructable for high orders of accuracy where other formulae may not exist.
There is not much sense in using this class directly, since it only provides 1D rules, weighted, as described before. Still, this class is particularly helpful: check QGauss for triangles and tetrahedra, with orders beyond THIRTIETH.
Definition at line 52 of file quadrature_jacobi.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::QJacobi::QJacobi | ( | const unsigned int | _dim, |
| const Order | _order = INVALID_ORDER, |
||
| const unsigned int | a = 1, |
||
| const unsigned int | b = 0 |
||
| ) | [inline] |
Constructor. Currently, only one-dimensional rules provided.
Definition at line 91 of file quadrature_jacobi.h.
References libMesh::QBase::_dim, libMesh::EDGE2, and libMesh::QBase::init().
| libMesh::QJacobi::~QJacobi | ( | ) | [inline] |
| 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::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(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), and libMesh::QConical::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(), 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::QJacobi::init_1D | ( | const ElemType | type = INVALID_ELEM, |
| unsigned int | p_level = 0 |
||
| ) | [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 28 of file quadrature_jacobi_1D.C.
References _alpha, _beta, libMesh::QBase::_order, libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::CONSTANT, libMesh::EIGHTH, libMesh::EIGHTTEENTH, libMesh::ELEVENTH, libMesh::err, libMesh::FIFTEENTH, libMesh::FIFTH, libMesh::FIRST, libMesh::FORTIETH, libMesh::FORTYFIRST, libMesh::FORTYSECOND, libMesh::FORTYTHIRD, libMesh::FOURTEENTH, libMesh::FOURTH, libMesh::QBase::libmesh_error_msg(), libMesh::NINETEENTH, libMesh::NINTH, 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.
{
//----------------------------------------------------------------------
// 1D quadrature rules
// Note that we need to check _alpha and _beta
// to determine which quadrature rule to implement.
// These rules have been pre-scaled based on their
// intended use. The weights of the (alpha=1,beta=0)
// rule sum to 0.5. The weights of the (alpha=2,beta=0)
// rule sum to 0.333333333333.
if ((_alpha == 1) && (_beta == 0))
{
switch(_order + 2*p)
{
case CONSTANT:
case FIRST:
{
_points.resize (1);
_weights.resize(1);
_points[0](0) = 1./3.;
_weights[0] = 0.5;
return;
}
case SECOND:
case THIRD:
{
_points.resize (2);
_weights.resize(2);
_points[ 0](0)=1.5505102572168219018027159252941e-01L;
_points[ 1](0)=6.4494897427831780981972840747059e-01L;
_weights[ 0]=3.1804138174397716939436900207516e-01L;
_weights[ 1]=1.8195861825602283060563099792484e-01L;
return;
}
case FOURTH:
case FIFTH:
{
_points.resize (3);
_weights.resize(3);
_points[ 0](0)=8.8587959512703947395546143769456e-02L;
_points[ 1](0)=4.0946686444073471086492625206883e-01L;
_points[ 2](0)=7.8765946176084705602524188987600e-01L;
_weights[ 0]=2.0093191373895963077219813326460e-01L;
_weights[ 1]=2.2924110635958624669392059455632e-01L;
_weights[ 2]=6.9826979901454122533881272179077e-02L;
return;
}
case SIXTH:
case SEVENTH:
{
_points.resize (4);
_weights.resize(4);
_points[ 0](0)=5.7104196114517682193121192554116e-02L;
_points[ 1](0)=2.7684301363812382768004599768563e-01L;
_points[ 2](0)=5.8359043236891682005669766866292e-01L;
_points[ 3](0)=8.6024013565621944784791291887512e-01L;
_weights[ 0]=1.3550691343148811620826417407794e-01L;
_weights[ 1]=2.0346456801027136079140447593585e-01L;
_weights[ 2]=1.2984754760823244082645620288963e-01L;
_weights[ 3]=3.1180970950008082173875147096575e-02L;
return;
}
case EIGHTH:
case NINTH:
{
_points.resize (5);
_weights.resize(5);
_points[ 0](0)=3.9809857051468742340806690093333e-02L;
_points[ 1](0)=1.9801341787360817253579213679530e-01L;
_points[ 2](0)=4.3797481024738614400501252000523e-01L;
_points[ 3](0)=6.9546427335363609451461482372117e-01L;
_points[ 4](0)=9.0146491420117357387650110211225e-01L;
_weights[ 0]=9.6781590226651679274360971636151e-02L;
_weights[ 1]=1.6717463809436956549167562309770e-01L;
_weights[ 2]=1.4638698708466980869803786935596e-01L;
_weights[ 3]=7.3908870072616670350633219341704e-02L;
_weights[ 4]=1.5747914521692276185292316568488e-02L;
return;
}
case TENTH:
case ELEVENTH:
{
_points.resize (6);
_weights.resize(6);
_points[ 0](0)=2.9316427159784891972050276913165e-02L;
_points[ 1](0)=1.4807859966848429184997685249598e-01L;
_points[ 2](0)=3.3698469028115429909705297208078e-01L;
_points[ 3](0)=5.5867151877155013208139334180552e-01L;
_points[ 4](0)=7.6923386203005450091688336011565e-01L;
_points[ 5](0)=9.2694567131974111485187396581968e-01L;
_weights[ 0]=7.2310330725508683655454326124832e-02L;
_weights[ 1]=1.3554249723151861684069039663805e-01L;
_weights[ 2]=1.4079255378819892811907683907093e-01L;
_weights[ 3]=9.8661150890655264120584510548360e-02L;
_weights[ 4]=4.3955165550508975508176624305430e-02L;
_weights[ 5]=8.7383018136095317560173033123989e-03L;
return;
}
case TWELFTH:
case THIRTEENTH:
{
_points.resize (7);
_weights.resize(7);
_points[ 0](0)=2.2479386438712498108825499570845e-02L;
_points[ 1](0)=1.1467905316090423190964023706751e-01L;
_points[ 2](0)=2.6578982278458946847678939169339e-01L;
_points[ 3](0)=4.5284637366944461699855144234765e-01L;
_points[ 4](0)=6.4737528288683036262609222982876e-01L;
_points[ 5](0)=8.1975930826310763501242005719082e-01L;
_points[ 6](0)=9.4373743946307785353434780896769e-01L;
_weights[ 0]=5.5967363423491051151562113429926e-02L;
_weights[ 1]=1.1050925819087460242112073636961e-01L;
_weights[ 2]=1.2739089729958832691873448577496e-01L;
_weights[ 3]=1.0712506569587366999091225299986e-01L;
_weights[ 4]=6.6384696465491469983647988763640e-02L;
_weights[ 5]=2.7408356721873475444331048358226e-02L;
_weights[ 6]=5.2143622028074040896913743037830e-03L;
return;
}
case FOURTEENTH:
case FIFTEENTH:
{
_points.resize (8);
_weights.resize(8);
_points[ 0](0)=1.7779915147363451813205101037679e-02L;
_points[ 1](0)=9.1323607899793956003741458074541e-02L;
_points[ 2](0)=2.1430847939563075835754126758167e-01L;
_points[ 3](0)=3.7193216458327230243085396048263e-01L;
_points[ 4](0)=5.4518668480342664903227222995321e-01L;
_points[ 5](0)=7.1317524285556948105131376025091e-01L;
_points[ 6](0)=8.5563374295785442851478147977179e-01L;
_points[ 7](0)=9.5536604471003014926687897814169e-01L;
_weights[ 0]=4.4550804361555931326215619533963e-02L;
_weights[ 1]=9.1119023636373626345722461782880e-02L;
_weights[ 2]=1.1250579947088736617177209860428e-01L;
_weights[ 3]=1.0604735943593001070453109641116e-01L;
_weights[ 4]=7.9199599492319160120408189359892e-02L;
_weights[ 5]=4.5439319504698898083055421095694e-02L;
_weights[ 6]=1.7842902655986208393552797201042e-02L;
_weights[ 7]=3.2951914422487988547423160110891e-03L;
return;
}
case SIXTEENTH:
case SEVENTEENTH:
{
_points.resize (9);
_weights.resize(9);
_points[ 0](0)=1.4412409648876548632826740810813e-02L;
_points[ 1](0)=7.4387389709196044635918185955999e-02L;
_points[ 2](0)=1.7611665616299528186317574607257e-01L;
_points[ 3](0)=3.0966757992763781705962046722948e-01L;
_points[ 4](0)=4.6197040108101093488314310868056e-01L;
_points[ 5](0)=6.1811723469529402463922975160350e-01L;
_points[ 6](0)=7.6282301518503961468269330714673e-01L;
_points[ 7](0)=8.8192102121000129980771488800526e-01L;
_points[ 8](0)=9.6374218711679053905883569923192e-01L;
_weights[ 0]=3.6278003523279871459649597941160e-02L;
_weights[ 1]=7.6074255109308162580079303753964e-02L;
_weights[ 2]=9.8533742172345705160173020369274e-02L;
_weights[ 3]=1.0030880919336828965415022473259e-01L;
_weights[ 4]=8.4358321844920349394250019958587e-02L;
_weights[ 5]=5.8401195295165110565731522901640e-02L;
_weights[ 6]=3.1804821491054001261690106935870e-02L;
_weights[ 7]=1.2060004284785379051727152286900e-02L;
_weights[ 8]=2.1808470857731308725490511200170e-03L;
return;
}
case EIGHTTEENTH:
case NINETEENTH:
{
_points.resize (10);
_weights.resize(10);
_points[ 0](0)=1.1917613432415596909745586958986e-02L;
_points[ 1](0)=6.1732071877148125522629367194250e-02L;
_points[ 2](0)=1.4711144964307024042759943557982e-01L;
_points[ 3](0)=2.6115967600845624026605165805485e-01L;
_points[ 4](0)=3.9463984688578684296195210507761e-01L;
_points[ 5](0)=5.3673876571566063287309517771189e-01L;
_points[ 6](0)=6.7594446167666510735715050893476e-01L;
_points[ 7](0)=8.0097892103689884513794630161724e-01L;
_points[ 8](0)=9.0171098779014677034879897840977e-01L;
_points[ 9](0)=9.6997096783851350295693564236559e-01L;
_weights[ 0]=3.0099508024037023300654755523879e-02L;
_weights[ 1]=6.4287154509082572874826538628865e-02L;
_weights[ 2]=8.6211300289176036428333360711914e-02L;
_weights[ 3]=9.2696893677723415844756991929711e-02L;
_weights[ 4]=8.4557109690827332430645582165501e-02L;
_weights[ 5]=6.6053075563350379166074502504878e-02L;
_weights[ 6]=4.3401906407150625922456117586384e-02L;
_weights[ 7]=2.2774591453263032586487649426060e-02L;
_weights[ 8]=8.4193197829831874175163655187011e-03L;
_weights[ 9]=1.4991406024063940282481360041064e-03L;
return;
}
case TWENTIETH:
case TWENTYFIRST:
{
_points.resize (11);
_weights.resize(11);
_points[ 0](0)=1.0018280461680405843024729867805e-02L;
_points[ 1](0)=5.2035451127180552583542695772964e-02L;
_points[ 2](0)=1.2461922514444307373529958726391e-01L;
_points[ 3](0)=2.2284060704383785550783145345736e-01L;
_points[ 4](0)=3.4000815791466518823360523389702e-01L;
_points[ 5](0)=4.6813761308958404208311038075775e-01L;
_points[ 6](0)=5.9849727976713918322772071367304e-01L;
_points[ 7](0)=7.2220328489096792556332130780461e-01L;
_points[ 8](0)=8.3082489962281857403056654390566e-01L;
_points[ 9](0)=9.1695838655259485329313462701803e-01L;
_points[10](0)=9.7472637960247965024666881353837e-01L;
_weights[ 0]=2.5367340688141426184944658517383e-02L;
_weights[ 1]=5.4938091132871498258746681665651e-02L;
_weights[ 2]=7.5620048057186998321377054841660e-02L;
_weights[ 3]=8.4659422884017621219239388551910e-02L;
_weights[ 4]=8.1879102988063470389586851943461e-02L;
_weights[ 5]=6.9531875158183049799050889633934e-02L;
_weights[ 6]=5.1591360672297580193086641744608e-02L;
_weights[ 7]=3.2641546713833346017953865389764e-02L;
_weights[ 8]=1.6663623451680692415747124573826e-02L;
_weights[ 9]=6.0439209604797757927329621765769e-03L;
_weights[10]=1.0636672932445414075338809612253e-03L;
return;
}
case TWENTYSECOND:
case TWENTYTHIRD:
{
_points.resize (12);
_weights.resize(12);
_points[ 0](0)=8.5390549884274193686644608778398e-03L;
_points[ 1](0)=4.4446463155407723025466798785553e-02L;
_points[ 2](0)=1.0685449088347665763410677043250e-01L;
_points[ 3](0)=1.9215105452985404099105725622861e-01L;
_points[ 4](0)=2.9538088426258022162291683437601e-01L;
_points[ 5](0)=4.1054508120145768248903435055933e-01L;
_points[ 6](0)=5.3095084931281767062893024289682e-01L;
_points[ 7](0)=6.4960065027725499276629172334298e-01L;
_points[ 8](0)=7.5959888952522705374260257404342e-01L;
_points[ 9](0)=8.5455254376493588079021191640562e-01L;
_points[10](0)=9.2894210126441101784881015513432e-01L;
_points[11](0)=9.7843793683414963909190691691700e-01L;
_weights[ 0]=2.1664860886871767324583472077976e-02L;
_weights[ 1]=4.7427851980377470134309182197087e-02L;
_weights[ 2]=6.6606750626673198130834485351975e-02L;
_weights[ 3]=7.6896602680041006281204469112373e-02L;
_weights[ 4]=7.7696316815531026628765631099760e-02L;
_weights[ 5]=7.0109339997630094495921805029671e-02L;
_weights[ 6]=5.6613843713667631943701105156121e-02L;
_weights[ 7]=4.0451653706913616760646424250646e-02L;
_weights[ 8]=2.4876780409265432012508172290016e-02L;
_weights[ 9]=1.2436009166422265607483293056674e-02L;
_weights[10]=4.4448077956684623227961458658283e-03L;
_weights[11]=7.7518222093802835724581451187414e-04L;
return;
}
case TWENTYFOURTH:
case TWENTYFIFTH:
{
_points.resize (13);
_weights.resize(13);
_points[ 0](0)=7.3646510260893216506914984138223e-03L;
_points[ 1](0)=3.8398138739678350376832524864142e-02L;
_points[ 2](0)=9.2595224699002635282891375438616e-02L;
_points[ 3](0)=1.6725101139155773149552247875936e-01L;
_points[ 4](0)=2.5862354070576251658979073282269e-01L;
_points[ 5](0)=3.6213139728223880040868129822732e-01L;
_points[ 6](0)=4.7258438600411772517509464267860e-01L;
_points[ 7](0)=5.8444396402134045550422084781091e-01L;
_points[ 8](0)=6.9210100171960165689704195168767e-01L;
_points[ 9](0)=7.9015702827343748555286333249938e-01L;
_points[10](0)=8.7369482130668941936771456713174e-01L;
_points[11](0)=9.3852445910073101239763338676530e-01L;
_points[12](0)=9.8138963498901214856028062215971e-01L;
_weights[ 0]=1.8714731585615175193374408664994e-02L;
_weights[ 1]=4.1320289419839298476611706952624e-02L;
_weights[ 2]=5.8953937956647082536161503261961e-02L;
_weights[ 3]=6.9713468282144398251416629950427e-02L;
_weights[ 4]=7.2849696515818223079384033843973e-02L;
_weights[ 5]=6.8815529487565417217837889436694e-02L;
_weights[ 6]=5.9120481314092327419290155953032e-02L;
_weights[ 7]=4.5995780676052393313259784530213e-02L;
_weights[ 8]=3.1936787162203148842381371158461e-02L;
_weights[ 9]=1.9213957164719820645895497863743e-02L;
_weights[10]=9.4489244795927503661800865492413e-03L;
_weights[11]=3.3383592507737718066104206420419e-03L;
_weights[12]=5.7805670493619285159651119259460e-04L;
return;
}
case TWENTYSIXTH:
case TWENTYSEVENTH:
{
_points.resize (14);
_weights.resize(14);
_points[ 0](0)=6.4167607928184568108204640945072e-03L;
_points[ 1](0)=3.3501404532013140035963928570572e-02L;
_points[ 2](0)=8.0985499681955184392451307739982e-02L;
_points[ 3](0)=1.4680486768121372992366016526080e-01L;
_points[ 4](0)=2.2808427064925799153464409861994e-01L;
_points[ 5](0)=3.2127174398893617440234039741304e-01L;
_points[ 6](0)=4.2229465730757025784040898751810e-01L;
_points[ 7](0)=5.2673786133987303205372694481288e-01L;
_points[ 8](0)=6.3003668837040395788448059413143e-01L;
_points[ 9](0)=7.2767645288926468543602672799033e-01L;
_points[10](0)=8.1538973944347464197307421861206e-01L;
_points[11](0)=8.8934280881951553969087166094649e-01L;
_points[12](0)=9.4630270006027538353340594342472e-01L;
_points[13](0)=9.8377523409860023828122800914100e-01L;
_weights[ 0]=1.6326754322596839819449429350767e-02L;
_weights[ 1]=3.6296078055443693187652118821852e-02L;
_weights[ 2]=5.2445953501916745591081806162678e-02L;
_weights[ 3]=6.3213025574329203442529115357932e-02L;
_weights[ 4]=6.7831659604182837102347339834146e-02L;
_weights[ 5]=6.6392350626085425132821358651277e-02L;
_weights[ 6]=5.9785559022531584222614059168587e-02L;
_weights[ 7]=4.9519454523203004580936179335753e-02L;
_weights[ 8]=3.7443096458985826655938301826886e-02L;
_weights[ 9]=2.5423465407148762115069898757459e-02L;
_weights[10]=1.5032303743928643503408714751315e-02L;
_weights[11]=7.2964841332126214434099222674776e-03L;
_weights[12]=2.5541013176876539479794996099427e-03L;
_weights[13]=4.3971370874715925476225610392609e-04L;
return;
}
case TWENTYEIGHTH:
case TWENTYNINTH:
{
_points.resize (15);
_weights.resize(15);
_points[ 0](0)=5.6406889725117097575064665794442e-03L;
_points[ 1](0)=2.9482298647942486631640984718004e-02L;
_points[ 2](0)=7.1412953115158840054783901746735e-02L;
_points[ 3](0)=1.2983102555359105597125453809429e-01L;
_points[ 4](0)=2.0249275505010404756821998832023e-01L;
_points[ 5](0)=2.8660608625752705998402604548440e-01L;
_points[ 6](0)=3.7893868864697803810652900914277e-01L;
_points[ 7](0)=4.7594230846323486057258899201585e-01L;
_points[ 8](0)=5.7388916090668588067145385523753e-01L;
_points[ 9](0)=6.6901519502995987653869563905314e-01L;
_points[10](0)=7.5766473903134271199707427774029e-01L;
_points[11](0)=8.3643096060561012397555908565418e-01L;
_points[12](0)=9.0228670067937802663123850466201e-01L;
_points[13](0)=9.5270040990583331433804420752426e-01L;
_points[14](0)=9.8573054526317422526590063305913e-01L;
_weights[ 0]=1.4367069621712523019113300206269e-02L;
_weights[ 1]=3.2119044874965915878547666076448e-02L;
_weights[ 2]=4.6891639944059925679345570546797e-02L;
_weights[ 3]=5.7398825506941289374618916116414e-02L;
_weights[ 4]=6.2918106648596909592665349230551e-02L;
_weights[ 5]=6.3343955161511392519350738069297e-02L;
_weights[ 6]=5.9174192059712931565021168664588e-02L;
_weights[ 7]=5.1412362311886468530400169371359e-02L;
_weights[ 8]=4.1400735184422578228813150524935e-02L;
_weights[ 9]=3.0609769586468540804291808789747e-02L;
_weights[10]=2.0416309272619835691190995270654e-02L;
_weights[11]=1.1904190355097578843030186960250e-02L;
_weights[12]=5.7172223947211761716246502974268e-03L;
_weights[13]=1.9862346931450380410338051145760e-03L;
_weights[14]=3.4034238413789606095252476068833e-04L;
return;
}
case THIRTIETH:
case THIRTYFIRST:
{
_points.resize (16);
_weights.resize(16);
_points[ 0](0)=4.9972996637719210426511496015207e-03L;
_points[ 1](0)=2.6143513681394051978549323380543e-02L;
_points[ 2](0)=6.3430945583836745492528222770701e-02L;
_points[ 3](0)=1.1559843757981594017892362824740e-01L;
_points[ 4](0)=1.8087055965788366624757339125570e-01L;
_points[ 5](0)=2.5702480784517085992294423955927e-01L;
_points[ 6](0)=3.4146792754782306862589929793869e-01L;
_points[ 7](0)=4.3132434359562537508359226958274e-01L;
_points[ 8](0)=5.2353411602516832051747016507566e-01L;
_points[ 9](0)=6.1495715187648951748964083381908e-01L;
_points[10](0)=7.0248013792504104445480697553762e-01L;
_points[11](0)=7.8312255396486808385218178256205e-01L;
_points[12](0)=8.5413814777521063366511368381932e-01L;
_points[13](0)=9.1310837653714127462958554670327e-01L;
_points[14](0)=9.5802441769047545454495934455429e-01L;
_points[15](0)=9.8735302062604161803115590316790e-01L;
_weights[ 0]=1.2739310657334184747549845879578e-02L;
_weights[ 1]=2.8612452590143149357558068987418e-02L;
_weights[ 2]=4.2129431157428284757515740541878e-02L;
_weights[ 3]=5.2228594403959147363509678893029e-02L;
_weights[ 4]=5.8254240209855835379301779422393e-02L;
_weights[ 5]=6.0000854888113507055385942002166e-02L;
_weights[ 6]=5.7718879609267823617681995467969e-02L;
_weights[ 7]=5.2064137843806583754973352189273e-02L;
_weights[ 8]=4.3997382189243923173481417304490e-02L;
_weights[ 9]=3.4647816579285799103084877873731e-02L;
_weights[10]=2.5159233106470040585653957496918e-02L;
_weights[11]=1.6539584355635793501225242083105e-02L;
_weights[12]=9.5341678846659366549669638207733e-03L;
_weights[13]=4.5392323861386350682211111786128e-03L;
_weights[14]=1.5671884735958197030724006964051e-03L;
_weights[15]=2.6749366505553617681762616226073e-04L;
return;
}
case THIRTYSECOND:
case THIRTYTHIRD:
{
_points.resize (17);
_weights.resize(17);
_points[ 0](0)=4.4579935677873917545062291115464e-03L;
_points[ 1](0)=2.3340094123774273967021528207683e-02L;
_points[ 2](0)=5.6707968769078239532713331556521e-02L;
_points[ 3](0)=1.0355543293519706341349122696131e-01L;
_points[ 4](0)=1.6246000342813654998162954391410e-01L;
_points[ 5](0)=2.3163212577717218294679722573836e-01L;
_points[ 6](0)=3.0897011775219590331752426128626e-01L;
_points[ 7](0)=3.9212413472209565803361946333190e-01L;
_points[ 8](0)=4.7856759685307958878371046285937e-01L;
_points[ 9](0)=5.6567396688663433190330697600965e-01L;
_points[10](0)=6.5079655845338226473148740732315e-01L;
_points[11](0)=7.3134895273405138213302149556298e-01L;
_points[12](0)=8.0488357852496649199917770139703e-01L;
_points[13](0)=8.6916605956741361920444649801116e-01L;
_points[14](0)=9.2224303459229840848710732823038e-01L;
_points[15](0)=9.6250119782335004125046077538851e-01L;
_points[16](0)=9.8871404063224375141712140225296e-01L;
_weights[ 0]=1.1372703113196894148880104403419e-02L;
_weights[ 1]=2.5642427920205179185750205247311e-02L;
_weights[ 2]=3.8025728825742400201474917006665e-02L;
_weights[ 3]=4.7641869374618441920575209265709e-02L;
_weights[ 4]=5.3907931667527480191616020246527e-02L;
_weights[ 5]=5.6573611382471512986991034182850e-02L;
_weights[ 6]=5.5734917923560549577129019490753e-02L;
_weights[ 7]=5.1809741929296333255214704762820e-02L;
_weights[ 8]=4.5477790399142277497438924164211e-02L;
_weights[ 9]=3.7592325487112238031431998524319e-02L;
_weights[10]=2.9074524087705050563187355467924e-02L;
_weights[11]=2.0803278010295362060012866375773e-02L;
_weights[12]=1.3513697452300780909638958310491e-02L;
_weights[13]=7.7164222148059386958884227660543e-03L;
_weights[14]=3.6472304834214064450344565237838e-03L;
_weights[15]=1.2526834393912207489106525201067e-03L;
_weights[16]=2.1311628920693358082515074128342e-04L;
return;
}
case THIRTYFOURTH:
case THIRTYFIFTH:
{
_points.resize (18);
_weights.resize(18);
_points[ 0](0)=4.0014793838673867970427199102708e-03L;
_points[ 1](0)=2.0963648393766477995135600986561e-02L;
_points[ 2](0)=5.0994041587890553966774921299853e-02L;
_points[ 3](0)=9.3280395928544950777609058826066e-02L;
_points[ 4](0)=1.4667010367761795715160882973461e-01L;
_points[ 5](0)=2.0970703958884296672188171145203e-01L;
_points[ 6](0)=2.8067179006017169910681926203455e-01L;
_points[ 7](0)=3.5762864997663462939832714298390e-01L;
_points[ 8](0)=4.3847844922413811765185426038741e-01L;
_points[ 9](0)=5.2101582087077614731981479158547e-01L;
_points[10](0)=6.0298936055983207213748070418747e-01L;
_points[11](0)=6.8216303911365156501882207998153e-01L;
_points[12](0)=7.5637719340615062679816733590978e-01L;
_points[13](0)=8.2360742977486956587394885462062e-01L;
_points[14](0)=8.8201982536207914673124120886230e-01L;
_points[15](0)=9.3002088969969321648195850472769e-01L;
_points[16](0)=9.6630075194563254241150674649068e-01L;
_points[17](0)=9.8986684820259713441676302277596e-01L;
_weights[ 0]=1.0214339435040633872267566728543e-02L;
_weights[ 1]=2.3106350163377066931289331193145e-02L;
_weights[ 2]=3.4471118921922024965038695158493e-02L;
_weights[ 3]=4.3573907074643088262181650935306e-02L;
_weights[ 4]=4.9902226917702713891102223856106e-02L;
_weights[ 5]=5.3192334631835660026223777017865e-02L;
_weights[ 6]=5.3445287802814027020361560753181e-02L;
_weights[ 7]=5.0916928856374093181562996165432e-02L;
_weights[ 8]=4.6082400990241280762086995392822e-02L;
_weights[ 9]=3.9579162288607878726422969126448e-02L;
_weights[10]=3.2134806127628549240737903166505e-02L;
_weights[11]=2.4487590001291716597869859208141e-02L;
_weights[12]=1.7308300480899289815470186505791e-02L;
_weights[13]=1.1131870097409362452892993478578e-02L;
_weights[14]=6.3060388823805337332782426625400e-03L;
_weights[15]=2.9624430644328375036535062298419e-03L;
_weights[16]=1.0130251793010471252471709101709e-03L;
_weights[17]=1.7186908409819589231237151109171e-04L;
return;
}
case THIRTYSIXTH:
case THIRTYSEVENTH:
{
_points.resize (19);
_weights.resize(19);
_points[ 0](0)=3.6116428185568930344532101262520e-03L;
_points[ 1](0)=1.8931837031588217250264524316949e-02L;
_points[ 2](0)=4.6097933048431083573924770835459e-02L;
_points[ 3](0)=8.4447222784209820599847900451170e-02L;
_points[ 4](0)=1.3303618855809877542408742397522e-01L;
_points[ 5](0)=1.9066859490476335108747466488719e-01L;
_points[ 6](0)=2.5592540348210296025204739729528e-01L;
_points[ 7](0)=3.2719980076381171572842996527566e-01L;
_points[ 8](0)=4.0273678616507709512304701110488e-01L;
_points[ 9](0)=4.8067639357855190263322135248271e-01L;
_points[10](0)=5.5909949264903198784692797587908e-01L;
_points[11](0)=6.3607504487927089469109105709816e-01L;
_points[12](0)=7.0970765165390239884374678311164e-01L;
_points[13](0)=7.7818422297676161420929964704208e-01L;
_points[14](0)=8.3981861570871107783863629442929e-01L;
_points[15](0)=8.9309313498184497075098115456491e-01L;
_points[16](0)=9.3669584807436508054390550208195e-01L;
_points[17](0)=9.6955263708022088547192878623227e-01L;
_points[18](0)=9.9085180527095568535309483522012e-01L;
_weights[ 0]=9.2240377250330766917668248692203e-03L;
_weights[ 1]=2.0924651751979890481445228332988e-02L;
_weights[ 2]=3.1376275285732898854356814457609e-02L;
_weights[ 3]=3.9962370048991741940293076442317e-02L;
_weights[ 4]=4.6234928168578411765380785132191e-02L;
_weights[ 5]=4.9934138240419379553284105838773e-02L;
_weights[ 6]=5.1004186626719396295755229568486e-02L;
_weights[ 7]=4.9591046468706580141316052533532e-02L;
_weights[ 8]=4.6021985716138412908362881060678e-02L;
_weights[ 9]=4.0768870871836685572301729841664e-02L;
_weights[10]=3.4398957951907391429559705435925e-02L;
_weights[11]=2.7518058144347978141481142479516e-02L;
_weights[12]=2.0711672562757804879837159838775e-02L;
_weights[13]=1.4489847112109724657001233005036e-02L;
_weights[14]=9.2410919293165977626976877702861e-03L;
_weights[15]=5.1997813282783796868055465441831e-03L;
_weights[16]=2.4300917120432807075681660438661e-03L;
_weights[17]=8.2788076092398219334600301771448e-04L;
_weights[18]=1.4012759417838633744062778723898e-04L;
return;
}
case THIRTYEIGHTH:
case THIRTYNINTH:
{
_points.resize (20);
_weights.resize(20);
_points[ 0](0)=3.2761066690500989752646761128711e-03L;
_points[ 1](0)=1.7181218145255713738247704631582e-02L;
_points[ 2](0)=4.1871431117765195244234634143709e-02L;
_points[ 3](0)=7.6800837089621970481870523019947e-02L;
_points[ 4](0)=1.2118986732367606561768615229454e-01L;
_points[ 5](0)=1.7404711263554216522415793419966e-01L;
_points[ 6](0)=2.3419188631359219089445571283230e-01L;
_points[ 7](0)=3.0028067683595016436625841871717e-01L;
_points[ 8](0)=3.7083718058439689644532207953771e-01L;
_points[ 9](0)=4.4428528696301116212882391783598e-01L;
_points[10](0)=5.1898428870357372961841176083024e-01L;
_points[11](0)=5.9326553348350648796739308654437e-01L;
_points[12](0)=6.6546969890551364795372958813385e-01L;
_points[13](0)=7.3398385823565944973202719329584e-01L;
_points[14](0)=7.9727750833659468770218885137527e-01L;
_points[15](0)=8.5393675303589045708900439161802e-01L;
_points[16](0)=9.0269587179345369227634694169967e-01L;
_points[17](0)=9.4246554236318633869089187742766e-01L;
_points[18](0)=9.7235694664743691118470187166465e-01L;
_points[19](0)=9.9169995579293273076654365969470e-01L;
_weights[ 0]=8.3708441762725569742347797346057e-03L;
_weights[ 1]=1.9034969463701666952317752319364e-02L;
_weights[ 2]=2.8668237379130596303519483134423e-02L;
_weights[ 3]=3.6750239836647855880170014918540e-02L;
_weights[ 4]=4.2890295837978204214146023790200e-02L;
_weights[ 5]=4.6841691862557827847492882468251e-02L;
_weights[ 6]=4.8516105059854090810181874268668e-02L;
_weights[ 7]=4.7985461609943912515256158210751e-02L;
_weights[ 8]=4.5470588941353475125260875269943e-02L;
_weights[ 9]=4.1317706816665737825832945172262e-02L;
_weights[10]=3.5964919153862747854063187859799e-02L;
_weights[11]=2.9901741842629013049082489189589e-02L;
_weights[12]=2.3625296173548892813066360441583e-02L;
_weights[13]=1.7597066879206016267686128655074e-02L;
_weights[14]=1.2204046126680104278798703245791e-02L;
_weights[15]=7.7276679300514628735780671415898e-03L;
_weights[16]=4.3232191126901531337390520742896e-03L;
_weights[17]=2.0114576139640780069071394684122e-03L;
_weights[18]=6.8306227608592198451771987553801e-04L;
_weights[19]=1.1538190717568529014836276132844e-04L;
return;
}
case FORTIETH:
case FORTYFIRST:
{
_points.resize (21);
_weights.resize(21);
_points[ 0](0)=2.9852372832130627418887840757498e-03L;
_points[ 1](0)=1.5662280557573547332895610655678e-02L;
_points[ 2](0)=3.8198288245073553344832819395658e-02L;
_points[ 3](0)=7.0139619062325798139842656763758e-02L;
_points[ 4](0)=1.1083667332036350223762530758616e-01L;
_points[ 5](0)=1.5946112944406525026381380836555e-01L;
_points[ 6](0)=2.1502318535802028508091216764592e-01L;
_points[ 7](0)=2.7639177909378317240968620976306e-01L;
_points[ 8](0)=3.4231763295873079372348132371246e-01L;
_points[ 9](0)=4.1145869156190519892697533612651e-01L;
_points[10](0)=4.8240744461243888053531852185958e-01L;
_points[11](0)=5.5371958072829359140036561243486e-01L;
_points[12](0)=6.2394338972128194749630581474409e-01L;
_points[13](0)=6.9164931501296225721944667896925e-01L;
_points[14](0)=7.5545905444968698626584997281486e-01L;
_points[15](0)=8.1407361651269741678894410957934e-01L;
_points[16](0)=8.6629975897095659707636330435792e-01L;
_points[17](0)=9.1107426580261577477966643276619e-01L;
_points[18](0)=9.4748554440639734668609897166113e-01L;
_points[19](0)=9.7479197566036565307915209128675e-01L;
_points[20](0)=9.9243549072562147749379027938900e-01L;
_weights[ 0]=7.6306059065231890951709176867159e-03L;
_weights[ 1]=1.7387926076195012545310092581185e-02L;
_weights[ 2]=2.6287260083203306468723755830943e-02L;
_weights[ 3]=3.3886790344403175515158905218748e-02L;
_weights[ 4]=3.9845977571398966444570950002831e-02L;
_weights[ 5]=4.3935570635150377464859835729052e-02L;
_weights[ 6]=4.6050270723257502074077985626262e-02L;
_weights[ 7]=4.6212552611022412660350368837561e-02L;
_weights[ 8]=4.4566874686154606866984823191992e-02L;
_weights[ 9]=4.1364760396997617354403702904192e-02L;
_weights[10]=3.6942012648499385551903205764287e-02L;
_weights[11]=3.1689953050127340376763333466278e-02L;
_weights[12]=2.6023046795895555117823162899995e-02L;
_weights[13]=2.0345549389821508292035088652126e-02L;
_weights[14]=1.5019872613192784763520392446339e-02L;
_weights[15]=1.0339209757970312663832726410241e-02L;
_weights[16]=6.5065974236456053354067512991831e-03L;
_weights[17]=3.6220528697449398652362737797211e-03L;
_weights[18]=1.6787560215669155362763365937815e-03L;
_weights[19]=5.6849852669976417861179123569953e-04L;
_weights[20]=9.5861868529721828979599842865244e-05L;
return;
}
case FORTYSECOND:
case FORTYTHIRD:
{
_points.resize (22);
_weights.resize(22);
_points[ 0](0)=2.7314460088842296840089704503200e-03L;
_points[ 1](0)=1.4335933483699706185308039285053e-02L;
_points[ 2](0)=3.4986328351216449387218814071085e-02L;
_points[ 3](0)=6.4302641132158882531275071395686e-02L;
_points[ 4](0)=1.0173934345216392401258226593794e-01L;
_points[ 5](0)=1.4659920015766541692656584266952e-01L;
_points[ 6](0)=1.9804660411818261484660504028001e-01L;
_points[ 7](0)=2.5512320697168696957682873507625e-01L;
_points[ 8](0)=3.1676578871044727183421117602379e-01L;
_points[ 9](0)=3.8182606925675020305095127628112e-01L;
_points[10](0)=4.4909210098681733251387195042393e-01L;
_points[11](0)=5.1731084595099112005136286684174e-01L;
_points[12](0)=5.8521151795129912975524533467408e-01L;
_points[13](0)=6.5152925489581447403580682231185e-01L;
_points[14](0)=7.1502868047476240261215783933800e-01L;
_points[15](0)=7.7452691606482263532020273841099e-01L;
_points[16](0)=8.2891561362172651235400549854277e-01L;
_points[17](0)=8.7718159746542002099358323827803e-01L;
_points[18](0)=9.1842572498984926327926994425432e-01L;
_points[19](0)=9.5187959322971991452332541323047e-01L;
_points[20](0)=9.7691966323712886717816657113704e-01L;
_points[21](0)=9.9307748504434821490300210664154e-01L;
_weights[ 0]=6.9842561136978136538216549700282e-03L;
_weights[ 1]=1.5944053454778496864565929126309e-02L;
_weights[ 2]=2.4184236671890499758870393233781e-02L;
_weights[ 3]=3.1327596999375227892408958559682e-02L;
_weights[ 4]=3.7077054093202864483217335909271e-02L;
_weights[ 5]=4.1222404906411452577079119693783e-02L;
_weights[ 6]=4.3651203786303462995293138094828e-02L;
_weights[ 7]=4.4353370900419806188626977124445e-02L;
_weights[ 8]=4.3418750145417827908334525766667e-02L;
_weights[ 9]=4.1027795781940901638265438618895e-02L;
_weights[10]=3.7436116466930019239662981188891e-02L;
_weights[11]=3.2954058785385534039887856528676e-02L;
_weights[12]=2.7922870823442631249079662456458e-02L;
_weights[13]=2.2689228598160995380212860776143e-02L;
_weights[14]=1.7580017560152292742735027157072e-02L;
_weights[15]=1.2879230281365339696348205046681e-02L;
_weights[16]=8.8086724401109070881737532508054e-03L;
_weights[17]=5.5138751021276410109163746634529e-03L;
_weights[18]=3.0562141131504037344643025112362e-03L;
_weights[19]=1.4117668531394242829582650189559e-03L;
_weights[20]=4.7692790565900675904839748933142e-04L;
_weights[21]=8.0298216937450816028842814607537e-05L;
return;
}
default:
libmesh_error_msg("Quadrature rule " << _order << " not supported!");
} // end switch(_order + 2*p)
} // end if ((_alpha == 1) && (_beta == 0))
else if ((_alpha == 2) && (_beta == 0))
{
switch(_order + 2*p)
{
case CONSTANT:
case FIRST:
{
_points.resize (1);
_weights.resize(1);
_points[0](0) = 0.25;
_weights[0] = 1./3.;
return;
}
case SECOND:
case THIRD:
{
_points.resize (2);
_weights.resize(2);
_points[ 0](0)=1.2251482265544137786674043037115e-01L;
_points[ 1](0)=5.4415184401122528879992623629551e-01L;
_weights[ 0]=2.3254745125350790274997694884235e-01L;
_weights[ 1]=1.0078588207982543058335638449099e-01L;
return;
}
case FOURTH:
case FIFTH:
{
_points.resize (3);
_weights.resize(3);
_points[ 0](0)=7.2994024073149732155837979012003e-02L;
_points[ 1](0)=3.4700376603835188472176354340395e-01L;
_points[ 2](0)=7.0500220988849838312239847758405e-01L;
_weights[ 0]=1.5713636106488661331834482221327e-01L;
_weights[ 1]=1.4624626925986602200351202036424e-01L;
_weights[ 2]=2.9950703008580698011476490755827e-02L;
return;
}
case SIXTH:
case SEVENTH:
{
_points.resize (4);
_weights.resize(4);
_points[ 0](0)=4.8500549446997329297067257098986e-02L;
_points[ 1](0)=2.3860073755186230505898141272470e-01L;
_points[ 2](0)=5.1704729510436750234057336938307e-01L;
_points[ 3](0)=7.9585141789677286330337796079324e-01L;
_weights[ 0]=1.1088841561127798368323131746895e-01L;
_weights[ 1]=1.4345878979921420904832801427594e-01L;
_weights[ 2]=6.8633887172923075317376345041811e-02L;
_weights[ 3]=1.0352240749918065284397656546639e-02L;
return;
}
case EIGHTH:
case NINTH:
{
_points.resize (5);
_weights.resize(5);
_points[ 0](0)=3.4578939918215091524457428631527e-02L;
_points[ 1](0)=1.7348032077169572310459241798618e-01L;
_points[ 2](0)=3.8988638706551932824089541038499e-01L;
_points[ 3](0)=6.3433347263088677234716388892062e-01L;
_points[ 4](0)=8.5105421294701641811622418741001e-01L;
_weights[ 0]=8.1764784285770917904880732922352e-02L;
_weights[ 1]=1.2619896189991148802883293516467e-01L;
_weights[ 2]=8.9200161221590000186254493070384e-02L;
_weights[ 3]=3.2055600722961919254748930556633e-02L;
_weights[ 4]=4.1138252030990079586162416192983e-03L;
return;
}
case TENTH:
case ELEVENTH:
{
_points.resize (6);
_weights.resize(6);
_points[ 0](0)=2.5904555093667192754643606997235e-02L;
_points[ 1](0)=1.3156394165798513398691085074097e-01L;
_points[ 2](0)=3.0243691802289123274990557791855e-01L;
_points[ 3](0)=5.0903641316475208401103990516772e-01L;
_points[ 4](0)=7.1568112731171391876766262459361e-01L;
_points[ 5](0)=8.8680561617756186630126600601049e-01L;
_weights[ 0]=6.2538702726580937878526556468332e-02L;
_weights[ 1]=1.0737649973678063260575568234795e-01L;
_weights[ 2]=9.4577186748541203568292051720052e-02L;
_weights[ 3]=5.1289571129616210220129325076919e-02L;
_weights[ 4]=1.5720297184945051327851262130020e-02L;
_weights[ 5]=1.8310758068692977327784555900562e-03L;
return;
}
case TWELFTH:
case THIRTEENTH:
{
_points.resize (7);
_weights.resize(7);
_points[ 0](0)=2.0132773773400507230501687117472e-02L;
_points[ 1](0)=1.0308902914804901475222678600595e-01L;
_points[ 2](0)=2.4055412604805753665369914150100e-01L;
_points[ 3](0)=4.1400214459705974641828724885057e-01L;
_points[ 4](0)=6.0002151327899293019396877288994e-01L;
_points[ 5](0)=7.7351724659143750111199899642603e-01L;
_points[ 6](0)=9.1118316656300276363931736720905e-01L;
_weights[ 0]=4.9276501776438120823429129084905e-02L;
_weights[ 1]=9.0698824612686144536967479972374e-02L;
_weights[ 2]=9.1733803279795254111416595896867e-02L;
_weights[ 3]=6.3146378708891297709533940481199e-02L;
_weights[ 4]=2.9422211289528614840745351985011e-02L;
_weights[ 5]=8.1629256323046945072318090687213e-03L;
_weights[ 6]=8.9268803368920680400902684425561e-04L;
return;
}
case FOURTEENTH:
case FIFTEENTH:
{
_points.resize (8);
_weights.resize(8);
_points[ 0](0)=1.6097759551921033532013550096618e-02L;
_points[ 1](0)=8.2900617485651102700366353380388e-02L;
_points[ 2](0)=1.9547516848873991732426660953061e-01L;
_points[ 3](0)=3.4165199147720222027296226216240e-01L;
_points[ 4](0)=5.0559707818448917194006186500608e-01L;
_points[ 5](0)=6.6955227182436145183011451055437e-01L;
_points[ 6](0)=8.1577170358328376075475178697587e-01L;
_points[ 7](0)=9.2850896495990689720101861784923e-01L;
_weights[ 0]=3.9778957806690559830106455926398e-02L;
_weights[ 1]=7.6818093267222624804506735068290e-02L;
_weights[ 2]=8.5284769171938765209441226265186e-02L;
_weights[ 3]=6.8447183421653272247027200610852e-02L;
_weights[ 4]=4.0814426388544023186461596152148e-02L;
_weights[ 5]=1.7246863780234982804243655774862e-02L;
_weights[ 6]=4.4745217130144079957565210032836e-03L;
_weights[ 7]=4.6851778403469725578994253231334e-04L;
return;
}
case SIXTEENTH:
case SEVENTEENTH:
{
_points.resize (9);
_weights.resize(9);
_points[ 0](0)=1.3165885597114490545190537817972e-02L;
_points[ 1](0)=6.8084529593767587476505856986822e-02L;
_points[ 2](0)=1.6175951676407464206981091232919e-01L;
_points[ 3](0)=2.8589108833922039772798956691264e-01L;
_points[ 4](0)=4.2945364538781279250924800228659e-01L;
_points[ 5](0)=5.7969405635116312626577241331199e-01L;
_points[ 6](0)=7.2326857174033543181796015820021e-01L;
_points[ 7](0)=8.4743684201323732017318042501962e-01L;
_points[ 8](0)=9.4124586421327421141434212713498e-01L;
_weights[ 0]=3.2760145111039721461819540084274e-02L;
_weights[ 1]=6.5489537033396631781939107203912e-02L;
_weights[ 2]=7.7673569160555504078595170478444e-02L;
_weights[ 3]=6.9284395689804709215727785348892e-02L;
_weights[ 4]=4.8540627864506507959528062424865e-02L;
_weights[ 5]=2.6343280902549835737974696002673e-02L;
_weights[ 6]=1.0406116579352237052426337235938e-02L;
_weights[ 7]=2.5743986456087231023317159056037e-03L;
_weights[ 8]=2.6126234651946294299091864873302e-04L;
return;
}
case EIGHTTEENTH:
case NINETEENTH:
{
_points.resize (10);
_weights.resize(10);
_points[ 0](0)=1.0968452456174134250281032990432e-02L;
_points[ 1](0)=5.6898150533657920711434859379121e-02L;
_points[ 2](0)=1.3595023405022895426576342938519e-01L;
_points[ 3](0)=2.4228119613252356111139455536509e-01L;
_points[ 4](0)=3.6800785044933771760476956625017e-01L;
_points[ 5](0)=5.0380712641487390700553064205723e-01L;
_points[ 6](0)=6.3960948865470970932853650616696e-01L;
_points[ 7](0)=7.6534767954810789620953312724460e-01L;
_points[ 8](0)=8.7171007457440848761186134200542e-01L;
_points[ 9](0)=9.5087429264052316644634948461033e-01L;
_weights[ 0]=2.7434088710097496960451648584485e-02L;
_weights[ 1]=5.6272936402808126282415832152709e-02L;
_weights[ 2]=7.0069507708666334086372455308183e-02L;
_weights[ 3]=6.7452219381437824970525794313691e-02L;
_weights[ 4]=5.2883788766964158209884838942421e-02L;
_weights[ 5]=3.3854565016814116044101071309373e-02L;
_weights[ 6]=1.7197575046553008079137502323978e-02L;
_weights[ 7]=6.4698890685585292297306065609276e-03L;
_weights[ 8]=1.5455231947365786285853268147186e-03L;
_weights[ 9]=1.5324003669716084212825702284675e-04L;
return;
}
case TWENTIETH:
case TWENTYFIRST:
{
_points.resize (11);
_weights.resize(11);
_points[ 0](0)=9.2789738313488342029466103013448e-03L;
_points[ 1](0)=4.8249692094286258286483723731132e-02L;
_points[ 2](0)=1.1578862662939522596030591825565e-01L;
_points[ 3](0)=2.0766834159705988244389542530053e-01L;
_points[ 4](0)=3.1811795190623396637744646657914e-01L;
_points[ 5](0)=4.4019839985886246708268132217632e-01L;
_points[ 6](0)=5.6623983915456831356751247193488e-01L;
_points[ 7](0)=6.8832423986295668569553397580809e-01L;
_points[ 8](0)=7.9878435859091464197460872800735e-01L;
_points[ 9](0)=8.9069109935439218915081289990608e-01L;
_points[10](0)=9.5832514378664820192443912466615e-01L;
_weights[ 0]=2.3300850051388377415061404350855e-02L;
_weights[ 1]=4.8745861030664491864450329028080e-02L;
_weights[ 2]=6.2979543000367394338333998811129e-02L;
_weights[ 3]=6.4183559349755944992539371306173e-02L;
_weights[ 4]=5.4632222321044926070678715192708e-02L;
_weights[ 5]=3.9290217798435612246535835921329e-02L;
_weights[ 6]=2.3589663538806928649316661943808e-02L;
_weights[ 7]=1.1414590918097618249921262806781e-02L;
_weights[ 8]=4.1400042632190223547676874553245e-03L;
_weights[ 9]=9.6302057554323908320037373956771e-04L;
_weights[10]=9.3800486009778068527692777578977e-05L;
return;
}
case TWENTYSECOND:
case TWENTYTHIRD:
{
_points.resize (12);
_weights.resize(12);
_points[ 0](0)=7.9520457026384368762190933400404e-03L;
_points[ 1](0)=4.1427810454294590790770165194795e-02L;
_points[ 2](0)=9.9757625542614157066081703480300e-02L;
_points[ 3](0)=1.7981078905242324954164334430721e-01L;
_points[ 4](0)=2.7727345779931511981753891831052e-01L;
_points[ 5](0)=3.8689200999768979166563470380544e-01L;
_points[ 6](0)=5.0275736044903220327930927843912e-01L;
_points[ 7](0)=6.1862386345845604629354886189114e-01L;
_points[ 8](0)=7.2824645295307268642888033006114e-01L;
_points[ 9](0)=8.2571851421478690176408336295370e-01L;
_points[10](0)=9.0579507354454372680179647039974e-01L;
_points[11](0)=9.6420653529267155121295530627839e-01L;
_weights[ 0]=2.0031112584748752943807724216722e-02L;
_weights[ 1]=4.2556615367405637772968490900760e-02L;
_weights[ 2]=5.6584184749474069646478036154467e-02L;
_weights[ 3]=6.0250607487604518170935343500906e-02L;
_weights[ 4]=5.4573941164231582825558709780290e-02L;
_weights[ 5]=4.2764763570188902317265936785879e-02L;
_weights[ 6]=2.8908051838727360494403492180206e-02L;
_weights[ 7]=1.6547343644609412101387213540038e-02L;
_weights[ 8]=7.7163670238864266187877366246530e-03L;
_weights[ 9]=2.7208421441038393545415361663513e-03L;
_weights[10]=6.1995339854451942336673117135933e-04L;
_weights[11]=5.9550359808311663832382311702501e-05L;
return;
}
case TWENTYFOURTH:
case TWENTYFIFTH:
{
_points.resize (13);
_weights.resize(13);
_points[ 0](0)=6.8908313099587313966815341534070e-03L;
_points[ 1](0)=3.5953362700170095501132811969521e-02L;
_points[ 2](0)=8.6811780157585433284280307690601e-02L;
_points[ 3](0)=1.5709997419650927473052206555148e-01L;
_points[ 4](0)=2.4353289958712844441903410356775e-01L;
_points[ 5](0)=3.4206951264607532347450799656627e-01L;
_points[ 6](0)=4.4810263890570985362873430879398e-01L;
_points[ 7](0)=5.5667462262945501069694664236208e-01L;
_points[ 8](0)=6.6270931356038854232146315212326e-01L;
_points[ 9](0)=7.6124976630827697227665952639343e-01L;
_points[10](0)=8.4769104394723239091031404895142e-01L;
_points[11](0)=9.1799964553720667953362597028835e-01L;
_points[12](0)=9.6892889422858896211181181730274e-01L;
_weights[ 0]=1.7401228610061845488081141733526e-02L;
_weights[ 1]=3.7426708168844039953666298186474e-02L;
_weights[ 2]=5.0909830712695019092906847187640e-02L;
_weights[ 3]=5.6111537130827396185039872203667e-02L;
_weights[ 4]=5.3340021605871745725401666422880e-02L;
_weights[ 5]=4.4624600755614954430566337477243e-02L;
_weights[ 6]=3.2940451678432790083648445397381e-02L;
_weights[ 7]=2.1254996714756613031113108888658e-02L;
_weights[ 8]=1.1728259593101880662970381149369e-02L;
_weights[ 9]=5.3135693388326384480559144318474e-03L;
_weights[10]=1.8323911392393785873485644191745e-03L;
_weights[11]=4.1072089179088541678428389490865e-04L;
_weights[12]=3.9016993264146227750471940564101e-05L;
return;
}
case TWENTYSIXTH:
case TWENTYSEVENTH:
{
_points.resize (14);
_weights.resize(14);
_points[ 0](0)=6.0288088710666481700887477121600e-03L;
_points[ 1](0)=3.1494259818657590707306712766317e-02L;
_points[ 2](0)=7.6213622474853867586406786340023e-02L;
_points[ 3](0)=1.3836652595947406914513833295825e-01L;
_points[ 4](0)=2.1540975706889557490069553208515e-01L;
_points[ 5](0)=3.0418952814165630963236635490440e-01L;
_points[ 6](0)=4.0107137869715435758118603205199e-01L;
_points[ 7](0)=5.0208913974997841006932539804960e-01L;
_points[ 8](0)=6.0310739672759678368189556654036e-01L;
_points[ 9](0)=6.9999091716990649741073666393849e-01L;
_points[10](0)=7.8877423843456071534030488559891e-01L;
_points[11](0)=8.6582490073159695367114083518102e-01L;
_points[12](0)=9.2799573874675232498377812203940e-01L;
_points[13](0)=9.7277712074118323045296336316727e-01L;
_weights[ 0]=1.5255206119448295373511341972378e-02L;
_weights[ 1]=3.3139595469960559875930225836517e-02L;
_weights[ 2]=4.5914772114671178911403799944693e-02L;
_weights[ 3]=5.2025865410701701832827358849471e-02L;
_weights[ 4]=5.1389924039471503273808065266413e-02L;
_weights[ 5]=4.5253432254460881465122373566022e-02L;
_weights[ 6]=3.5739010292514612984162525810346e-02L;
_weights[ 7]=2.5216183536839055034306220202219e-02L;
_weights[ 8]=1.5694332285072957184930434798487e-02L;
_weights[ 9]=8.4129767371364622096077152778570e-03L;
_weights[10]=3.7248756752828588215538124784652e-03L;
_weights[11]=1.2617838951574880413283334197260e-03L;
_weights[12]=2.7909934094687847008945560658272e-04L;
_weights[13]=2.6276161668899854751670304157309e-05L;
return;
}
case TWENTYEIGHTH:
case TWENTYNINTH:
{
_points.resize (15);
_weights.resize(15);
_points[ 0](0)=5.3190520028437405909978402739774e-03L;
_points[ 1](0)=2.7814561918260051571980709238074e-02L;
_points[ 2](0)=6.7431865838998613384568577571818e-02L;
_points[ 3](0)=1.2274828621549783933958732401303e-01L;
_points[ 4](0)=1.9176570627765362550282684505042e-01L;
_points[ 5](0)=2.7198996592948620305633800793026e-01L;
_points[ 6](0)=3.6052169358705764815778444077840e-01L;
_points[ 7](0)=4.5416123575447709916009624996188e-01L;
_points[ 8](0)=5.4952435559355131365468962196031e-01L;
_points[ 9](0)=6.4316460176878493853949400884187e-01L;
_points[10](0)=7.3169797305252320629933948467521e-01L;
_points[11](0)=8.1192547727500966668969312120195e-01L;
_points[12](0)=8.8094952871531360238808055493652e-01L;
_points[13](0)=9.3628182173551874210251937228715e-01L;
_points[14](0)=9.7595387433502370956200384127912e-01L;
_weights[ 0]=1.3481662725134304090265733542410e-02L;
_weights[ 1]=2.9527650756798205089388266870014e-02L;
_weights[ 2]=4.1531307651604265373799456795289e-02L;
_weights[ 3]=4.8132057668529815888404876719886e-02L;
_weights[ 4]=4.9041267894409180565196427508457e-02L;
_weights[ 5]=4.4991555820615933365730152800755e-02L;
_weights[ 6]=3.7473875792101303095319296790068e-02L;
_weights[ 7]=2.8326553174241247644068974123766e-02L;
_weights[ 8]=1.9293416585274520485988288429744e-02L;
_weights[ 9]=1.1668582484776638563534331197891e-02L;
_weights[10]=6.1110501102060490412153142922240e-03L;
_weights[11]=2.6556629857743191811282114211333e-03L;
_weights[12]=8.8657363523409667398227745346013e-04L;
_weights[13]=1.9398726306394554934013378631178e-04L;
_weights[14]=1.8128785569508725971591601924058e-05L;
return;
}
case THIRTIETH:
case THIRTYFIRST:
{
_points.resize (16);
_weights.resize(16);
_points[ 0](0)=4.7276871229345922061506104450645e-03L;
_points[ 1](0)=2.4742967619434984512603456269536e-02L;
_points[ 2](0)=6.0076437716637992275242255328112e-02L;
_points[ 3](0)=1.0960060925877060286789758461950e-01L;
_points[ 4](0)=1.7172475948300210790307417791270e-01L;
_points[ 5](0)=2.4445243174451137730207656160275e-01L;
_points[ 6](0)=3.2544621002412961482683685190181e-01L;
_points[ 7](0)=4.1210296411685959691728421362251e-01L;
_points[ 8](0)=5.0163755615318759370013043235084e-01L;
_points[ 9](0)=5.9117238531650308897076051877634e-01L;
_points[10](0)=6.7782991741456081130251277264148e-01L;
_points[11](0)=7.5882525794914141335427940535114e-01L;
_points[12](0)=8.3155588302646195188920635482406e-01L;
_points[13](0)=8.9368597454880401881486969876816e-01L;
_points[14](0)=9.4322428571224559352545066239296e-01L;
_points[15](0)=9.7860643749869701257280091378128e-01L;
_weights[ 0]=1.1999405624691969586506695404143e-02L;
_weights[ 1]=2.6460842517339527040437264564557e-02L;
_weights[ 2]=3.7685824681918572833315081444311e-02L;
_weights[ 3]=4.4496573472543214976151049587710e-02L;
_weights[ 4]=4.6506239397362658797012521100333e-02L;
_weights[ 5]=4.4112657675899371113911374089320e-02L;
_weights[ 6]=3.8348046951185955673582690910292e-02L;
_weights[ 7]=3.0612037917837609541948834407964e-02L;
_weights[ 8]=2.2357185649529930223334734263644e-02L;
_weights[ 9]=1.4803776657751738340239250997704e-02L;
_weights[10]=8.7475454956120699745801506391282e-03L;
_weights[11]=4.4948488514798588832066397608164e-03L;
_weights[12]=1.9235319260667742273824398535627e-03L;
_weights[13]=6.3445999956883686801084896849641e-04L;
_weights[14]=1.3757810397555175246959527906570e-04L;
_weights[15]=1.2778410569693501244162062287056e-05L;
return;
}
case THIRTYSECOND:
case THIRTYTHIRD:
{
_points.resize (17);
_weights.resize(17);
_points[ 0](0)=4.2297654864907678749919877813007e-03L;
_points[ 1](0)=2.2152705311830232294209002101573e-02L;
_points[ 2](0)=5.3856014827184829566760322013645e-02L;
_points[ 3](0)=9.8434936830541362494246924237098e-02L;
_points[ 4](0)=1.5460792826588562987838903880788e-01L;
_points[ 5](0)=2.2075922519215138391363787768336e-01L;
_points[ 6](0)=2.9498586629802348865100598497056e-01L;
_points[ 7](0)=3.7515254164081788200834771392226e-01L;
_points[ 8](0)=4.5895305335301908085933451474375e-01L;
_points[ 9](0)=5.4397667809047570974515217435547e-01L;
_points[10](0)=6.2777753854681427256560538094776e-01L;
_points[11](0)=7.0794500138508556794362242607309e-01L;
_points[12](0)=7.8217310308122112810113854015572e-01L;
_points[13](0)=8.4832708570995060757906353137179e-01L;
_points[14](0)=9.0450542331432870037906100932836e-01L;
_points[15](0)=9.4909701659654457024654762835657e-01L;
_points[16](0)=9.8084389384741256367666372092760e-01L;
_weights[ 0]=1.0748151030768621662834062204936e-02L;
_weights[ 1]=2.3837762018703244620506074665265e-02L;
_weights[ 2]=3.4307889795671880028436869439002e-02L;
_weights[ 3]=4.1143975847068207121110777965721e-02L;
_weights[ 4]=4.3922618772592577808611328095136e-02L;
_weights[ 5]=4.2825403927934553301439133345126e-02L;
_weights[ 6]=3.8555241947102353628223582107170e-02L;
_weights[ 7]=3.2163570112378639217164962127759e-02L;
_weights[ 8]=2.4829256111976082802381065145755e-02L;
_weights[ 9]=1.7638789331050110518481458865839e-02L;
_weights[10]=1.1413642238636687747573403275566e-02L;
_weights[11]=6.6164197603617147286557116199452e-03L;
_weights[12]=3.3465172918198153715474254116460e-03L;
_weights[13]=1.4138837589561322316624033902623e-03L;
_weights[14]=4.6167405877711308077043506706165e-04L;
_weights[15]=9.9356525335473959962406835119113e-05L;
_weights[16]=9.1808042001255039722337720228143e-06L;
return;
}
case THIRTYFOURTH:
case THIRTYFIFTH:
{
_points.resize (18);
_weights.resize(18);
_points[ 0](0)=3.8065822475018763836485343548908e-03L;
_points[ 1](0)=1.9948351047343018666301273550578e-02L;
_points[ 2](0)=4.8549645304224630965525818224597e-02L;
_points[ 3](0)=8.8876259116497388048181392194957e-02L;
_points[ 4](0)=1.3988457083516926493003835526765e-01L;
_points[ 5](0)=2.0025369592363033820671255853211e-01L;
_points[ 6](0)=2.6842018278660208442014712803618e-01L;
_points[ 7](0)=3.4261859780200830045140484446480e-01L;
_points[ 8](0)=4.2092727587031795580391095608392e-01L;
_points[ 9](0)=5.0131810275027358836467598537906e-01L;
_points[10](0)=5.8170905267091792822883758529331e-01L;
_points[11](0)=6.6001812722531554655212340270381e-01L;
_points[12](0)=7.3421730813501325031921736140253e-01L;
_points[13](0)=8.0238514983641747020667379058517e-01L;
_points[14](0)=8.6275672070005501336081662021610e-01L;
_points[15](0)=9.1376986411079293616612318979159e-01L;
_points[16](0)=9.5410789551731665881890309384613e-01L;
_points[17](0)=9.8274840759428696063307389954630e-01L;
_weights[ 0]=9.6823988037660231774964807691529e-03L;
_weights[ 1]=2.1578772679193827027157750409882e-02L;
_weights[ 2]=3.1333907671729513327950397913143e-02L;
_weights[ 3]=3.8075134253447687665457823100366e-02L;
_weights[ 4]=4.1377182287519565332651504100178e-02L;
_weights[ 5]=4.1283665809957749414286828772884e-02L;
_weights[ 6]=3.8262794672349299757849656256339e-02L;
_weights[ 7]=3.3095976202072497795672769135623e-02L;
_weights[ 8]=2.6723424368960119486097972811646e-02L;
_weights[ 9]=2.0078641696369720179443819032712e-02L;
_weights[10]=1.3943897653647906425703531920973e-02L;
_weights[11]=8.8522782250949774435502470347624e-03L;
_weights[12]=5.0502487883546813166841235600508e-03L;
_weights[13]=2.5207216384362558022745702364513e-03L;
_weights[14]=1.0535502839571923712409403476880e-03L;
_weights[15]=3.4109120368765498289598525047321e-04L;
_weights[16]=7.2937008664118705328162954331465e-05L;
_weights[17]=6.7100861245431215907697266774769e-06L;
return;
}
case THIRTYSIXTH:
case THIRTYSEVENTH:
{
_points.resize (19);
_weights.resize(19);
_points[ 0](0)=3.4438904038624902806940509331722e-03L;
_points[ 1](0)=1.8056978337900562863030093954343e-02L;
_points[ 2](0)=4.3987395090842734455198832773659e-02L;
_points[ 3](0)=8.0633276366661415469721776663016e-02L;
_points[ 4](0)=1.2713640986156026953002535083969e-01L;
_points[ 5](0)=1.8240698366909185313069695387128e-01L;
_points[ 6](0)=2.4514956909516857055249217701821e-01L;
_points[ 7](0)=3.1389356741821822415240859865307e-01L;
_points[ 8](0)=3.8702770276933220341545290751683e-01L;
_points[ 9](0)=4.6283779855202027820220487054694e-01L;
_points[10](0)=5.3954696279208808766943066332998e-01L;
_points[11](0)=6.1535724480454716546694016358337e-01L;
_points[12](0)=6.8849179091409498918986551203802e-01L;
_points[13](0)=7.5723651917329793549736932794667e-01L;
_points[14](0)=8.1998035706323733537954006841029e-01L;
_points[15](0)=8.7525316298382668006594890887609e-01L;
_points[16](0)=9.2176068172914004889937421632274e-01L;
_points[17](0)=9.5841690242340714510347790926694e-01L;
_points[18](0)=9.8438280655170201067612761745569e-01L;
_weights[ 0]=8.7672970568792323661393261305840e-03L;
_weights[ 1]=1.9620856340232434512901283777202e-02L;
_weights[ 2]=2.8708149778704086122791006941433e-02L;
_weights[ 3]=3.5278117199476000434599584597235e-02L;
_weights[ 4]=3.8922386222261737550985848811514e-02L;
_weights[ 5]=3.9598248214880677015066061570123e-02L;
_weights[ 6]=3.7607500520347543950912173611332e-02L;
_weights[ 7]=3.3525068070610283960905700567460e-02L;
_weights[ 8]=2.8091607936102907649141378579235e-02L;
_weights[ 9]=2.2090563149823990465333672755377e-02L;
_weights[10]=1.6231815321996325630177612472423e-02L;
_weights[11]=1.1061440629037007897211784796273e-02L;
_weights[12]=6.9108110779321981973514541466546e-03L;
_weights[13]=3.8897848125530902102984958782958e-03L;
_weights[14]=1.9197884726923132048372710905667e-03L;
_weights[15]=7.9504432544855937325482850317179e-04L;
_weights[16]=2.5553085889097125647555166013802e-04L;
_weights[17]=5.4342643132004416773622274012570e-05L;
_weights[18]=4.9807023319691181766751703037610e-06L;
return;
}
case THIRTYEIGHTH:
case THIRTYNINTH:
{
_points.resize (20);
_weights.resize(20);
_points[ 0](0)=3.1306837407435895621501417495463e-03L;
_points[ 1](0)=1.6422088133987832902541614256258e-02L;
_points[ 2](0)=4.0036900461906780996320149828793e-02L;
_points[ 3](0)=7.3477191785285283057566545285874e-02L;
_points[ 4](0)=1.1603090765548680494669743156748e-01L;
_points[ 5](0)=1.6679125348100943001298070145849e-01L;
_points[ 6](0)=2.2467641938430503840938196483210e-01L;
_points[ 7](0)=2.8845271218127238477194070665240e-01L;
_points[ 8](0)=3.5676087010245310564705506663326e-01L;
_points[ 9](0)=4.2814504093325668461542503477245e-01L;
_points[10](0)=5.0108381548284855137062465413693e-01L;
_points[11](0)=5.7402265818282687172751219665371e-01L;
_points[12](0)=6.4540704577736756661732083476397e-01L;
_points[13](0)=7.1371561101528768074208343535393e-01L;
_points[14](0)=7.7749259203548675355112594583658e-01L;
_points[15](0)=8.3537891463156172517151839471341e-01L;
_points[16](0)=8.8614130291454714368728293036495e-01L;
_points[17](0)=9.2869901491875054534509783716727e-01L;
_points[18](0)=9.6214871215049818970916080642292e-01L;
_points[19](0)=9.8579578884064184668002313135922e-01L;
_weights[ 0]=7.9757927362766516852273101574255e-03L;
_weights[ 1]=1.7913752325738559101773958866323e-02L;
_weights[ 2]=2.6382564429125584794513408474046e-02L;
_weights[ 3]=3.2734644317573148047479656498988e-02L;
_weights[ 4]=3.6587916553227064986609295576925e-02L;
_weights[ 5]=3.7847390376922340715354694449373e-02L;
_weights[ 6]=3.6697420472018132713254249570653e-02L;
_weights[ 7]=3.3556431453475493479002714822739e-02L;
_weights[ 8]=2.9002404745589285169916357261403e-02L;
_weights[ 9]=2.3682290414246053980866450266495e-02L;
_weights[10]=1.8220503418268691551843315309759e-02L;
_weights[11]=1.3140919256661651742108133720781e-02L;
_weights[12]=8.8135560782433553289252731679767e-03L;
_weights[13]=5.4320908061230274500120111426487e-03L;
_weights[14]=3.0224891709073793749184859780611e-03L;
_weights[15]=1.4774444419714858062838576456553e-03L;
_weights[16]=6.0704570426725825415272373108771e-04L;
_weights[17]=1.9388830961751181076983851525948e-04L;
_weights[18]=4.1039102087320205549074040898061e-05L;
_weights[19]=3.7492209933371347725241368372638e-06L;
return;
}
case FORTIETH:
case FORTYFIRST:
{
_points.resize (21);
_weights.resize(21);
_points[ 0](0)=2.8583506000918943222556277168267e-03L;
_points[ 1](0)=1.4999364972075522340357135916023e-02L;
_points[ 2](0)=3.6593863281140122626545055282186e-02L;
_points[ 3](0)=6.7226441667938638526563947674675e-02L;
_points[ 4](0)=1.0630147973769311755044868847393e-01L;
_points[ 5](0)=1.5305857495571921224970195300963e-01L;
_points[ 6](0)=2.0658770608253356867698033417515e-01L;
_points[ 7](0)=2.6584701401081590981089215154372e-01L;
_points[ 8](0)=3.2968309899010272373734914205452e-01L;
_points[ 9](0)=3.9685347698561900950394342909483e-01L;
_points[10](0)=4.6605076663517464847296649003812e-01L;
_points[11](0)=5.3592813869303216093268951839176e-01L;
_points[12](0)=6.0512553397610192844810985230774e-01L;
_points[13](0)=6.7229614098692961628997049948121e-01L;
_points[14](0)=7.3613262059728734846410297670221e-01L;
_points[15](0)=7.9539257374613146232498650332774e-01L;
_points[16](0)=8.4892277353768883984401746895126e-01L;
_points[17](0)=8.9568174292546676052820257429765e-01L;
_points[18](0)=9.3476043640122994798279429359660e-01L;
_points[19](0)=9.6540160736574245809294907904429e-01L;
_points[20](0)=9.8702556657875783654690055164718e-01L;
_weights[ 0]=7.2866317313744308995629793584932e-03L;
_weights[ 1]=1.6417061144242818824757454238052e-02L;
_weights[ 2]=2.4316074956305588348553255406534e-02L;
_weights[ 3]=3.0423854654878694196479542521305e-02L;
_weights[ 4]=3.4388553550627421941340376509459e-02L;
_weights[ 5]=3.6085277384893044986268336555960e-02L;
_weights[ 6]=3.5615944774693119613624648820995e-02L;
_weights[ 7]=3.3281046762833346236473134492600e-02L;
_weights[ 8]=2.9528072702082493907799384765866e-02L;
_weights[ 9]=2.4885083875255978309567929915161e-02L;
_weights[10]=1.9889780066352815019424534213776e-02L;
_weights[11]=1.5024528041917190276613254752701e-02L;
_weights[12]=1.0666270365508680718786852036243e-02L;
_weights[13]=7.0573209743621604363614669672455e-03L;
_weights[14]=4.2993100642134152210363702100934e-03L;
_weights[15]=2.3686215346115318974092271035074e-03L;
_weights[16]=1.1482414207164748272703814817286e-03L;
_weights[17]=4.6857361253397420043294612589344e-04L;
_weights[18]=1.4884999898655056859287865070751e-04L;
_weights[19]=3.1377233396792732027953172274463e-05L;
_weights[20]=2.8584835468101709504260347382715e-06L;
return;
}
case FORTYSECOND:
case FORTYTHIRD:
{
_points.resize (22);
_weights.resize(22);
_points[ 0](0)=2.6200747203711594248593674552649e-03L;
_points[ 1](0)=1.3753657987385804025756750709059e-02L;
_points[ 2](0)=3.3575226848026692714308191640158e-02L;
_points[ 3](0)=6.1735561623069601690826574395308e-02L;
_points[ 4](0)=9.7732696675368642458479884996723e-02L;
_points[ 5](0)=1.4092439775695158519450619775816e-01L;
_points[ 6](0)=1.9053995518421662849913969152881e-01L;
_points[ 7](0)=2.4569399962450332222446519900435e-01L;
_points[ 8](0)=3.0540231844259902315260853518244e-01L;
_points[ 9](0)=3.6859942495594128250321664201293e-01L;
_points[10](0)=4.3415757484015513750820121247106e-01L;
_points[11](0)=5.0090689264793505011775015892772e-01L;
_points[12](0)=5.6765625027260978702177266771788e-01L;
_points[13](0)=6.3321452557266284279166256182373e-01L;
_points[14](0)=6.9641186298595172661126707689713e-01L;
_points[15](0)=7.5612055917893487122928793527786e-01L;
_points[16](0)=8.1127520683922050787811864631196e-01L;
_points[17](0)=8.6089175276286621730656128908876e-01L;
_points[18](0)=9.0408517849775272066170136431884e-01L;
_points[19](0)=9.4008566921081846961603077466108e-01L;
_points[20](0)=9.6825388381656365288164210087944e-01L;
_points[21](0)=9.8810245999087788318348935085438e-01L;
_weights[ 0]=6.6829315239216368796332393406135e-03L;
_weights[ 1]=1.5098062421623821319493004238262e-02L;
_weights[ 2]=2.2473713734526648003280648297671e-02L;
_weights[ 3]=2.8324458025288763364233486167273e-02L;
_weights[ 4]=3.2329473952869436965326587507885e-02L;
_weights[ 5]=3.4348599664417779736425965942178e-02L;
_weights[ 6]=3.4426309828130701949676126316990e-02L;
_weights[ 7]=3.2774586900497126738748991342107e-02L;
_weights[ 8]=2.9737336628427308246400548097974e-02L;
_weights[ 9]=2.5741821234796918381507016798653e-02L;
_weights[10]=2.1244167813524304464862872202199e-02L;
_weights[11]=1.6676482379699419616981906191727e-02L;
_weights[12]=1.2402453510959574984251948178874e-02L;
_weights[13]=8.6866984738749141556196547655377e-03L;
_weights[14]=5.6807683849250404394970034550940e-03L;
_weights[15]=3.4260685209698087638377353482425e-03L;
_weights[16]=1.8713879532585689543042377871244e-03L;
_weights[17]=9.0066625279251987400261357834378e-04L;
_weights[18]=3.6536328467284676654467859219800e-04L;
_weights[19]=1.1551470917602499323205530447797e-04L;
_weights[20]=2.4263177919766775694554800111547e-05L;
_weights[21]=2.2049570604019597784590797971169e-06L;
return;
}
default:
libmesh_error_msg("Quadrature rule " << _order << " not supported!");
} // end switch(_order + 2*p)
} // end else if ((_alpha == 2) && (_beta == 0))
else
{
libMesh::err << "Unsupported combination of (alpha,beta) = ("
<< _alpha
<< ","
<< _beta
<< ") requested in Jacobi-Gauss quadrature rule."
<< std::endl;
}
}
| virtual void libMesh::QBase::init_2D | ( | const ElemType | , |
| unsigned int | = 0 |
||
| ) | [inline, protected, virtual, inherited] |
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 redefined, gives an error (when DEBUG defined) when called.
Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QGrid, libMesh::QConical, libMesh::QSimpson, libMesh::QGauss, libMesh::QClough, libMesh::QTrap, and libMesh::QGaussLobatto.
Definition at line 260 of file quadrature.h.
Referenced by libMesh::QBase::init().
{}
| virtual void libMesh::QBase::init_3D | ( | const ElemType | , |
| unsigned int | = 0 |
||
| ) | [inline, protected, virtual, inherited] |
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 redefined, gives an error (when DEBUG defined) when called.
Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QConical, libMesh::QGrid, libMesh::QSimpson, libMesh::QGauss, libMesh::QClough, libMesh::QTrap, and libMesh::QGaussLobatto.
Definition at line 278 of file quadrature.h.
Referenced by libMesh::QBase::init().
{}
| 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(), init_1D(), libMesh::QGaussLobatto::init_2D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QConical::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(), libMesh::QConical::init_3D(), libMesh::QGrundmann_Moller::init_3D(), libMesh::QGauss::keast_rule(), libMesh::QMonomial::kim_rule(), libMesh::QMonomial::stroud_rule(), and 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(), libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::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 libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::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 libMesh::QConical::conical_product_tet(), and libMesh::QConical::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::QJacobi::type | ( | ) | const [inline, virtual] |
QuadratureType, either QJACOBI_1_0 or QJACOBI_2_0. Implements libMesh::QBase.
Definition at line 103 of file quadrature_jacobi.h.
References _alpha, _beta, libMesh::QBase::libmesh_error_msg(), libMesh::QJACOBI_1_0, and libMesh::QJACOBI_2_0.
{
if ((_alpha == 1) && (_beta == 0))
return QJACOBI_1_0;
else if ((_alpha == 2) && (_beta == 0))
return QJACOBI_2_0;
else
libmesh_error_msg("Invalid Jacobi quadrature rule: alpha = " << _alpha << ", beta = " << _beta);
}
| 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 libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::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;
}
const unsigned int libMesh::QJacobi::_alpha [private] |
Definition at line 77 of file quadrature_jacobi.h.
const unsigned int libMesh::QJacobi::_beta [private] |
Definition at line 78 of file quadrature_jacobi.h.
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(), 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 libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::QBase::get_order(), libMesh::QGaussLobatto::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QGrid::init_1D(), 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 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(), libMesh::QGaussLobatto::init_1D(), libMesh::QTrap::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QGrid::init_1D(), 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 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(), libMesh::QGaussLobatto::init_1D(), libMesh::QTrap::init_1D(), libMesh::QClough::init_1D(), libMesh::QGauss::init_1D(), libMesh::QSimpson::init_1D(), libMesh::QGrid::init_1D(), 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().