Go to the documentation of this file.
   20 #include "libmesh/elem.h" 
   21 #include "libmesh/quadrature.h" 
   22 #include "libmesh/int_range.h" 
   29   allow_rules_with_negative_weights(true),
 
   42   Real summed_weights=0;
 
   43   os << 
"N_Q_Points=" << this->
n_points() << std::endl << std::endl;
 
   46       os << 
" Point " << qpoint << 
":\n" 
   50          << 
"  w=" << 
_weights[qpoint] << 
"\n" << std::endl;
 
   54   os << 
"Summed Weights: " << summed_weights << std::endl;
 
   97       libmesh_error_msg(
"Invalid dimension _dim = " << 
_dim);
 
  104                   const std::vector<Real> & ,
 
  105                   unsigned int p_level)
 
  125   libmesh_not_implemented();
 
  132   libmesh_not_implemented();
 
  138                   std::pair<Real, Real> new_range)
 
  141   libmesh_assert_equal_to (
_dim, 1);
 
  144     h_new = new_range.second - new_range.first,
 
  145     h_old = old_range.second - old_range.first;
 
  148   libmesh_assert_greater (h_new, 0.);
 
  149   libmesh_assert_greater (h_old, 0.);
 
  152   libmesh_assert_greater (
_points.size(), 0);
 
  155   Real scfact = h_new/h_old;
 
  160       _points[i](0) = new_range.first +
 
  161         (
_points[i](0) - old_range.first) * scfact;
 
  174   const unsigned int np = q1D.
n_points();
 
  182   for (
unsigned int j=0; j<np; j++)
 
  183     for (
unsigned int i=0; i<np; i++)
 
  200   const unsigned int np = q1D.
n_points();
 
  208   for (
unsigned int k=0; k<np; k++)
 
  209     for (
unsigned int j=0; j<np; j++)
 
  210       for (
unsigned int i=0; i<np; i++)
 
  227   const unsigned int n_points1D = q1D.
n_points();
 
  228   const unsigned int n_points2D = q2D.
n_points();
 
  230   _points.resize  (n_points1D * n_points2D);
 
  231   _weights.resize (n_points1D * n_points2D);
 
  235   for (
unsigned int j=0; j<n_points1D; j++)
 
  236     for (
unsigned int i=0; i<n_points2D; i++)
 
  
The QBase class provides the basic functionality from which various quadrature rules can be derived.
 
unsigned int _p_level
The p-level of the element for which the current values have been computed.
 
void tensor_product_hex(const QBase &q1D)
Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D.
 
unsigned int n_points() const
 
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
 
The libMesh namespace provides an interface to certain functionality in the library.
 
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
 
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate val...
 
Real w(const unsigned int i) const
 
virtual void init_3D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate val...
 
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the data structures for a quadrature rule for an element of type type.
 
QBase(unsigned int dim, Order order=INVALID_ORDER)
Constructor.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
ElemType _type
The type of element for which the current values have been computed.
 
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...
 
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule.
 
std::vector< Real > _weights
The quadrature weights.
 
void tensor_product_quad(const QBase &q1D)
Constructs a 2D rule from the tensor product of q1D with itself.
 
This is the base class from which all geometric element types are derived.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate val...
 
unsigned int _dim
The spatial dimension of the quadrature rule.
 
Point qp(const unsigned int i) const
 
std::vector< Point > _points
The locations of the quadrature points in reference element space.
 
virtual ElemType type() const =0
 
void print_info(std::ostream &os=libMesh::out) const
Prints information relevant to the quadrature rule, by default to libMesh::out.
 
ElemType
Defines an enum for geometric element types.
 
virtual void init_2D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate val...