19 #include "libmesh/libmesh_config.h"    20 #if defined(LIBMESH_HAVE_TRIANGLE) && defined(LIBMESH_HAVE_TETGEN)    22 #include "libmesh/fe.h"    23 #include "libmesh/quadrature_gauss.h"    24 #include "libmesh/quadrature_trap.h"    25 #include "libmesh/quadrature_simpson.h"    26 #include "libmesh/quadrature_composite.h"    27 #include "libmesh/elem.h"    28 #include "libmesh/enum_quadrature_type.h"    36 template <
class QSubCell>
    44 template <
class QSubCell>
    66 template <
class QSubCell>
    68                                  const std::vector<Real> & vertex_distance_func,
    71   libmesh_assert_equal_to (vertex_distance_func.size(), elem.
n_vertices());
    72   libmesh_assert_equal_to (_dim, elem.
dim());
    79   if (!_elem_cutter.is_cut (elem, vertex_distance_func))
    81       _q_subcell.init (elem, p_level);
    82       _points  = _q_subcell.get_points();
    83       _weights = _q_subcell.get_weights();
    98   _elem_cutter(*reference_elem, vertex_distance_func);
   107     const std::vector<Elem const *> & inside_elem (_elem_cutter.inside_elements());
   110     this->add_subelem_values(inside_elem);
   115     const std::vector<Elem const *> & outside_elem (_elem_cutter.outside_elements());
   118     this->add_subelem_values(outside_elem);
   126 template <
class QSubCell>
   130   const std::vector<Real>  & subelem_weights = _lagrange_fe->get_JxW();
   131   const std::vector<Point> & subelem_points  = _lagrange_fe->get_xyz();
   133   for (
const auto & elem : subelem)
   140           _lagrange_fe->reinit(elem);
   141           _weights.insert(_weights.end(),
   142                           subelem_weights.begin(), subelem_weights.end());
   144           _points.insert(_points.end(),
   145                          subelem_points.begin(), subelem_points.end());
   151           for (
auto n : elem->node_index_range())
   154           libmesh_error_msg(
"Tetgen may have created a 0-volume cell during Cutcell integration.");
   169 #endif // LIBMESH_HAVE_TRIANGLE && LIBMESH_HAVE_TETGEN class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Order
defines an enum for polynomial orders. 
QSubCell _q_subcell
Subcell quadrature object. 
This is the base class from which all geometric element types are derived. 
QComposite(unsigned int dim, Order order=INVALID_ORDER)
Constructor. 
QuadratureType
Defines an enum for currently available quadrature rules. 
virtual void init(const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0) override
Overrides the base class init() function, and uses the ElemCutter to subdivide the element into "insi...
The libMesh namespace provides an interface to certain functionality in the library. 
ElemMappingType mapping_type() const
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary. 
const Elem * reference_elem() const
This class implements generic composite quadrature rules. 
std::unique_ptr< FEBase > _lagrange_fe
Lagrange FE to use for subcell mapping. 
virtual unsigned int n_vertices() const =0
virtual unsigned short dim() const =0
virtual QuadratureType type() const override
This class forms the foundation from which generic finite elements may be derived. 
void add_subelem_values(const std::vector< Elem const *> &subelem)
Helper function called from init() to collect all the points and weights of the subelement quadrature...