20 #include "libmesh/quadrature_conical.h"    21 #include "libmesh/quadrature_gauss.h"    22 #include "libmesh/quadrature_jacobi.h"    23 #include "libmesh/enum_quadrature_type.h"    40   return std::make_unique<QConical>(*this);
    62   libmesh_assert_equal_to (this->
get_dim(), 2);
    68   std::pair<Real, Real> old_range(-1, 1);
    69   std::pair<Real, Real> new_range( 0, 1);
    70   gauss1D.
scale(old_range,
    79   const unsigned int np = gauss1D.
n_points();
    82   libmesh_assert_greater_equal (gauss1D.
qp(0)(0), 0.0);
    83   libmesh_assert_less_equal (gauss1D.
qp(np-1)(0), 1.0);
    84   libmesh_assert_greater_equal (jac1D.
qp(0)(0), 0.0);
    85   libmesh_assert_less_equal (jac1D.
qp(np-1)(0), 1.0);
    93   for (
unsigned int i=0; i<np; i++)
    94     for (
unsigned int j=0; j<np; j++)
    97         _points[gp](1) = gauss1D.
qp(i)(0) * (1.-jac1D.
qp(j)(0)); 
   113   libmesh_assert_equal_to (this->
get_dim(), 3);
   120   std::pair<Real, Real> old_range(-1, 1);
   121   std::pair<Real, Real> new_range( 0, 1);
   122   gauss1D.
scale(old_range,
   132   const unsigned int np = gauss1D.
n_points();
   135   libmesh_assert_greater_equal (gauss1D.
qp(0)(0), 0.0);
   136   libmesh_assert_less_equal (gauss1D.
qp(np-1)(0), 1.0);
   137   libmesh_assert_greater_equal (jacA1D.
qp(0)(0), 0.0);
   138   libmesh_assert_less_equal (jacA1D.
qp(np-1)(0), 1.0);
   139   libmesh_assert_greater_equal (jacB1D.
qp(0)(0), 0.0);
   140   libmesh_assert_less_equal (jacB1D.
qp(np-1)(0), 1.0);
   148   for (
unsigned int i=0; i<np; i++)
   149     for (
unsigned int j=0; j<np; j++)
   150       for (
unsigned int k=0; k<np; k++)
   153           _points[gp](1) = jacA1D.
qp(j)(0)  * (1.-jacB1D.
qp(k)(0));                         
   154           _points[gp](2) = gauss1D.
qp(i)(0) * (1.-jacA1D.
qp(j)(0)) * (1.-jacB1D.
qp(k)(0)); 
   155           _weights[gp]   = gauss1D.
w(i)     * jacA1D.
w(j)          * jacB1D.
w(k);          
   190   libmesh_assert_equal_to (this->
get_dim(), 3);
   199   const unsigned int np = gauss1D.
n_points();
   207   for (
unsigned int i=0; i<np; ++i)
   208     for (
unsigned int j=0; j<np; ++j)
   209       for (
unsigned int k=0; k<np; ++k, ++q)
   211           const Real xi=gauss1D.
qp(i)(0);
   212           const Real yj=gauss1D.
qp(j)(0);
   213           const Real zk=jac1D.
qp(k)(0);
   218           _weights[q]   = gauss1D.
w(i) * gauss1D.
w(j) * jac1D.
w(k);
 void conical_product_tri()
Implementation of conical product rule for a Tri in 2D of order get_order(). 
Point qp(const unsigned int i) const
const std::vector< Real > & get_weights() const
QuadratureType
Defines an enum for currently available quadrature rules. 
void conical_product_tet()
Implementation of conical product rule for a Tet in 3D of order get_order(). 
The libMesh namespace provides an interface to certain functionality in the library. 
std::vector< Point > _points
The locations of the quadrature points in reference element space. 
std::vector< Real > _weights
The quadrature weights. 
void conical_product_pyramid()
Implementation of conical product rule for a Pyramid in 3D of order get_order(). 
unsigned int get_dim() const
unsigned int n_points() const
virtual std::unique_ptr< QBase > clone() const override
virtual void init_1D() override
In 1D, use a Gauss rule. 
virtual QuadratureType type() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Point > & get_points() const
This class implements two (for now) Jacobi-Gauss quadrature rules. 
This class implements specific orders of Gauss quadrature. 
virtual void init(const Elem &e, unsigned int p_level=invalid_uint)
Initializes the data structures for a quadrature rule for the element e. 
void scale(std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
Maps the points of a 1D quadrature rule defined by "old_range" to another 1D interval defined by "new...
Real w(const unsigned int i) const