$extrastylesheet
#include <quadrature_composite.h>
Public Member Functions | |
| QComposite (const unsigned int _dim, const Order _order=INVALID_ORDER) | |
| ~QComposite () | |
| QuadratureType | type () const |
| virtual void | init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0) |
Private Member Functions | |
| void | add_subelem_values (const std::vector< Elem const * > &subelem) |
Private Attributes | |
| QSubCell | _q_subcell |
| ElemCutter | _elem_cutter |
| UniquePtr< FEBase > | _lagrange_fe |
This class implements generic composite quadrature rules. Composite quadrature rules are constructed from any of the supported rules by breaking an element into subelements and applying the base rule on each subelement. This class uses the ElemCutter, which is only available if libmesh is configured with --disable-strict-lgpl.
Definition at line 49 of file quadrature_composite.h.
| libMesh::QComposite< QSubCell >::QComposite | ( | const unsigned int | _dim, |
| const Order | _order = INVALID_ORDER |
||
| ) |
Constructor. Declares the order of the quadrature rule.
Definition at line 36 of file quadrature_composite.C.
References libMesh::QComposite< QSubCell >::_lagrange_fe, libMesh::QComposite< QSubCell >::_q_subcell, libMesh::EDGE2, libMesh::QComposite< QSubCell >::init(), and libMesh::libmesh_assert().
: QSubCell(d,o), // explicitly call base class constructor _q_subcell(d,o), _lagrange_fe(FEBase::build (d, FEType (FIRST, LAGRANGE))) { // explicitly call the init function in 1D since the // other tensor-product rules require this one. // note that EDGE will not be used internally, however // if we called the function with INVALID_ELEM it would try to // be smart and return, thinking it had already done the work. if (_dim == 1) QSubCell::init(EDGE2); libmesh_assert (_lagrange_fe.get() != NULL); _lagrange_fe->attach_quadrature_rule (&_q_subcell); }
| libMesh::QComposite< QSubCell >::~QComposite | ( | ) |
| void libMesh::QComposite< QSubCell >::add_subelem_values | ( | const std::vector< Elem const * > & | subelem | ) | [private] |
Definition at line 120 of file quadrature_composite.C.
References libMesh::err.
{
const std::vector<Real> &subelem_weights = _lagrange_fe->get_JxW();
const std::vector<Point> &subelem_points = _lagrange_fe->get_xyz();
for (std::vector<Elem const*>::const_iterator it = subelem.begin();
it!=subelem.end(); ++it)
{
// tetgen seems to create 0-volume cells on occasion, but we *should*
// be catching that appropriately now inside the ElemCutter class.
// Just in case trap here, describe the error, and abort.
#ifdef LIBMESH_ENABLE_EXCEPTIONS
try
{
#endif
_lagrange_fe->reinit(*it);
_weights.insert(_weights.end(),
subelem_weights.begin(), subelem_weights.end());
_points.insert(_points.end(),
subelem_points.begin(), subelem_points.end());
#ifdef LIBMESH_ENABLE_EXCEPTIONS
}
catch (...)
{
libMesh::err << "ERROR: found a bad cut cell!\n";
for (unsigned int n=0; n<(*it)->n_nodes(); n++)
libMesh::err << (*it)->point(n) << std::endl;
libmesh_error_msg("Tetgen may have created a 0-volume cell during Cutcell integration.");
}
#endif
}
}
| void libMesh::QComposite< QSubCell >::init | ( | const Elem & | elem, |
| const std::vector< Real > & | vertex_distance_func, | ||
| unsigned int | p_level = 0 |
||
| ) | [virtual] |
Initializes the data structures for a specific, potentially cut element. The array vertex_distance_func contains vertex values of a signed distance function that cuts the element. This interface is indended to be extended by derived classes that can cut the element into subelements, for example, and constuct a composite quadrature rule for the cut element.
Definition at line 64 of file quadrature_composite.C.
References libMesh::Elem::dim(), libMesh::libmesh_assert(), libMesh::Elem::n_vertices(), libMesh::Elem::reference_elem(), and libMesh::Elem::type().
Referenced by libMesh::QComposite< QSubCell >::QComposite().
{
libmesh_assert_equal_to (vertex_distance_func.size(), elem.n_vertices());
libmesh_assert_equal_to (_dim, elem.dim());
// if we are not cut, revert to simple base class init() method.
if (!_elem_cutter.is_cut (elem, vertex_distance_func))
{
_q_subcell.init (elem.type(), p_level);
_points = _q_subcell.get_points();
_weights = _q_subcell.get_weights();
//this->print_info();
return;
}
// Get a pointer to the element's reference element. We want to
// perform cutting on the reference element such that the quadrature
// point locations of the subelements live in the reference
// coordinate system, thereby eliminating the need for inverse
// mapping.
const Elem *reference_elem = elem.reference_elem();
libmesh_assert (reference_elem != NULL);
_elem_cutter(*reference_elem, vertex_distance_func);
//_elem_cutter(elem, vertex_distance_func);
// clear our state & accumulate points from subelements
_points.clear();
_weights.clear();
// inside subelem
{
const std::vector<Elem const*> &inside_elem (_elem_cutter.inside_elements());
std::cout << inside_elem.size() << " elements inside\n";
this->add_subelem_values(inside_elem);
}
// outside subelem
{
const std::vector<Elem const*> &outside_elem (_elem_cutter.outside_elements());
std::cout << outside_elem.size() << " elements outside\n";
this->add_subelem_values(outside_elem);
}
this->print_info();
}
| QuadratureType libMesh::QComposite< QSubCell >::type | ( | ) | const [inline] |
QCOMPOSITE Definition at line 75 of file quadrature_composite.h.
References libMesh::QCOMPOSITE.
{ return QCOMPOSITE; }
ElemCutter libMesh::QComposite< QSubCell >::_elem_cutter [private] |
ElemCutter object.
Definition at line 104 of file quadrature_composite.h.
UniquePtr<FEBase> libMesh::QComposite< QSubCell >::_lagrange_fe [private] |
Lagrange FE to use for subcell mapping.
Definition at line 109 of file quadrature_composite.h.
Referenced by libMesh::QComposite< QSubCell >::QComposite().
QSubCell libMesh::QComposite< QSubCell >::_q_subcell [private] |
Subcell quadrature object.
Definition at line 99 of file quadrature_composite.h.
Referenced by libMesh::QComposite< QSubCell >::QComposite().