20 #ifndef LIBMESH_QUADRATURE_H 21 #define LIBMESH_QUADRATURE_H 31 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS 70 QBase (
unsigned int dim,
83 virtual ~QBase() =
default;
99 static std::unique_ptr<QBase>
build (
const std::string &
name,
100 const unsigned int dim,
111 const unsigned int dim,
129 libmesh_assert (!
_points.empty());
130 return cast_int<unsigned int>(
_points.size());
167 libmesh_assert_less (i,
_points.size());
176 libmesh_assert_less (i,
_weights.size());
185 unsigned int p_level=0);
196 virtual void init (
const Elem & elem,
197 const std::vector<Real> & vertex_distance_func,
198 unsigned int p_level=0);
216 void scale(std::pair<Real, Real> old_range,
217 std::pair<Real, Real> new_range);
257 unsigned int p_level=0);
267 unsigned int p_level=0) = 0;
281 libmesh_error_msg(
"ERROR: Seems as if this quadrature rule \nis not implemented for 2D.");
297 libmesh_error_msg(
"ERROR: Seems as if this quadrature rule \nis not implemented for 3D.");
366 allow_rules_with_negative_weights(true),
380 libmesh_assert(!
_points.empty());
383 Real summed_weights=0;
384 os <<
"N_Q_Points=" << this->
n_points() << std::endl << std::endl;
385 for (
unsigned int qpoint=0; qpoint<this->
n_points(); qpoint++)
387 os <<
" Point " << qpoint <<
":\n" 391 <<
" w=" <<
_weights[qpoint] <<
"\n" << std::endl;
395 os <<
"Summed Weights: " << summed_weights << std::endl;
400 #endif // LIBMESH_QUADRATURE_H std::string name(const ElemQuality q)
virtual void init(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
virtual bool shapes_need_reinit()
QBase & operator=(const QBase &)=default
std::vector< Point > & get_points()
bool allow_rules_with_negative_weights
Point qp(const unsigned int i) const
ElemType get_elem_type() const
The base class for all geometric element types.
const std::vector< Real > & get_weights() const
QBase(unsigned int dim, Order order=INVALID_ORDER)
friend std::ostream & operator<<(std::ostream &os, const QBase &q)
std::vector< Point > _points
std::vector< Real > _weights
virtual QuadratureType type() const =0
void tensor_product_prism(const QBase &q1D, const QBase &q2D)
virtual void init_2D(const ElemType, unsigned int=0)
virtual void init_1D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
void tensor_product_hex(const QBase &q1D)
unsigned int get_p_level() const
virtual void init_0D(const ElemType type=INVALID_ELEM, unsigned int p_level=0)
unsigned int get_dim() const
unsigned int n_points() const
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
void print_info(std::ostream &os=libMesh::out) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< Point > & get_points() const
void scale(std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
void tensor_product_quad(const QBase &q1D)
Real w(const unsigned int i) const
OStreamProxy out(std::cout)
A geometric point in (x,y,z) space.
virtual void init_3D(const ElemType, unsigned int=0)
Base class for all quadrature families and orders.
std::vector< Real > & get_weights()