20 #include "libmesh/quadrature_clough.h"    21 #include "libmesh/quadrature_gauss.h"    22 #include "libmesh/quadrature_gm.h"    23 #include "libmesh/quadrature_grid.h"    24 #include "libmesh/quadrature_jacobi.h"    25 #include "libmesh/quadrature_monomial.h"    26 #include "libmesh/quadrature_simpson.h"    27 #include "libmesh/quadrature_trap.h"    28 #include "libmesh/quadrature_gauss_lobatto.h"    29 #include "libmesh/quadrature_conical.h"    30 #include "libmesh/quadrature_nodal.h"    31 #include "libmesh/string_to_enum.h"    32 #include "libmesh/enum_quadrature_type.h"    44                                      const unsigned int _dim,
    55                                     const unsigned int _dim,
    66             libMesh::out << 
"WARNING: Clough quadrature implemented" << std::endl
    67                          << 
" up to TWENTYTHIRD order." << std::endl;
    71         return std::make_unique<QClough>(
_dim, 
_order);
    80             libMesh::out << 
"WARNING: Gauss quadrature implemented" << std::endl
    81                          << 
" up to FORTYTHIRD order." << std::endl;
    94             libMesh::out << 
"WARNING: Jacobi(1,0) quadrature implemented" << std::endl
    95                          << 
" up to FORTYTHIRD order." << std::endl;
   100             libMesh::out << 
"WARNING: Jacobi(1,0) quadrature implemented" << std::endl
   101                          << 
" in 1D only." << std::endl;
   105         return std::make_unique<QJacobi>(
_dim, 
_order, 1, 0);
   114             libMesh::out << 
"WARNING: Jacobi(2,0) quadrature implemented" << std::endl
   115                          << 
" up to FORTYTHIRD order." << std::endl;
   120             libMesh::out << 
"WARNING: Jacobi(2,0) quadrature implemented" << std::endl
   121                          << 
" in 1D only." << std::endl;
   125         return std::make_unique<QJacobi>(
_dim, 
_order, 2, 0);
   134             libMesh::out << 
"WARNING: Simpson rule provides only" << std::endl
   135                          << 
" THIRD order!" << std::endl;
   139         return std::make_unique<QSimpson>(
_dim);
   148             libMesh::out << 
"WARNING: Trapezoidal rule provides only" << std::endl
   149                          << 
" FIRST order!" << std::endl;
   153         return std::make_unique<QTrap>(
_dim);
   160       return std::make_unique<QGrundmann_Moller>(
_dim, 
_order);
   163       return std::make_unique<QMonomial>(
_dim, 
_order);
   166       return std::make_unique<QGaussLobatto>(
_dim, 
_order);
   169       return std::make_unique<QConical>(
_dim, 
_order);
   172       return std::make_unique<QNodal>(
_dim, 
_order);
   175       libmesh_error_msg(
"ERROR: Bad qt=" << _qt);
 
Order
defines an enum for polynomial orders. 
QuadratureType
Defines an enum for currently available quadrature rules. 
unsigned int _dim
The spatial dimension of the quadrature rule. 
The libMesh namespace provides an interface to certain functionality in the library. 
virtual QuadratureType type() const =0
Order _order
The polynomial order which the quadrature rule is capable of integrating exactly. ...
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.