libMesh::QComposite< QSubCell > Class Template Referencefinal

A quadrature rule for subdivided elements. More...

#include <quadrature_composite.h>

Inheritance diagram for libMesh::QComposite< QSubCell >:

Public Member Functions

 QComposite (unsigned int dim, Order order=INVALID_ORDER)
 
 QComposite (const QComposite &)=delete
 
QCompositeoperator= (const QComposite &)=delete
 
 QComposite (QComposite &&)=default
 
QCompositeoperator= (QComposite &&)=default
 
virtual ~QComposite ()=default
 
virtual QuadratureType type () const override
 
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0) override
 

Private Member Functions

void add_subelem_values (const std::vector< Elem const *> &subelem)
 

Private Attributes

QSubCell _q_subcell
 
ElemCutter _elem_cutter
 
std::unique_ptr< FEBase_lagrange_fe
 

Detailed Description

template<class QSubCell>
class libMesh::QComposite< QSubCell >

A quadrature rule for subdivided elements.

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.

Author
Benjamin Kirk
Date
2013

Definition at line 52 of file quadrature_composite.h.

Constructor & Destructor Documentation

◆ QComposite() [1/3]

template<class QSubCell >
libMesh::QComposite< QSubCell >::QComposite ( unsigned int  dim,
Order  order = INVALID_ORDER 
)

Constructor. Declares the order of the quadrature rule.

Definition at line 45 of file quadrature_composite.C.

References libMesh::QComposite< QSubCell >::_lagrange_fe, libMesh::QComposite< QSubCell >::_q_subcell, libMesh::EDGE2, and libMesh::TriangleWrapper::init().

46  :
47  QSubCell(d,o), // explicitly call base class constructor
48  _q_subcell(d,o),
50 {
51  // explicitly call the init function in 1D since the
52  // other tensor-product rules require this one.
53  // note that EDGE will not be used internally, however
54  // if we called the function with INVALID_ELEM it would try to
55  // be smart and return, thinking it had already done the work.
56  if (_dim == 1)
58 
59  libmesh_assert (_lagrange_fe.get() != nullptr);
60 
61  _lagrange_fe->attach_quadrature_rule (&_q_subcell);
62 }
void init(triangulateio &t)
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
std::unique_ptr< FEBase > _lagrange_fe

◆ QComposite() [2/3]

template<class QSubCell >
libMesh::QComposite< QSubCell >::QComposite ( const QComposite< QSubCell > &  )
delete

This class contains a unique_ptr member, so it can't be default copy constructed or assigned.

◆ QComposite() [3/3]

template<class QSubCell >
libMesh::QComposite< QSubCell >::QComposite ( QComposite< QSubCell > &&  )
default

Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class.

◆ ~QComposite()

template<class QSubCell >
virtual libMesh::QComposite< QSubCell >::~QComposite ( )
virtualdefault

Member Function Documentation

◆ add_subelem_values()

template<class QSubCell >
void libMesh::QComposite< QSubCell >::add_subelem_values ( const std::vector< Elem const *> &  subelem)
private

Helper function called from init() to collect all the points and weights of the subelement quadrature rules.

Definition at line 123 of file quadrature_composite.C.

References libMesh::err.

125 {
126  const std::vector<Real> & subelem_weights = _lagrange_fe->get_JxW();
127  const std::vector<Point> & subelem_points = _lagrange_fe->get_xyz();
128 
129  for (const auto & elem : subelem)
130  {
131  // tetgen seems to create 0-volume cells on occasion, but we *should*
132  // be catching that appropriately now inside the ElemCutter class.
133  // Just in case trap here, describe the error, and abort.
134  libmesh_try
135  {
136  _lagrange_fe->reinit(elem);
137  _weights.insert(_weights.end(),
138  subelem_weights.begin(), subelem_weights.end());
139 
140  _points.insert(_points.end(),
141  subelem_points.begin(), subelem_points.end());
142  }
143  libmesh_catch (...)
144  {
145  libMesh::err << "ERROR: found a bad cut cell!\n";
146 
147  for (unsigned int n=0; n<elem->n_nodes(); n++)
148  libMesh::err << elem->point(n) << std::endl;
149 
150  libmesh_error_msg("Tetgen may have created a 0-volume cell during Cutcell integration.");
151  }
152  }
153 }
OStreamProxy err(std::cerr)
std::unique_ptr< FEBase > _lagrange_fe

◆ init()

template<class QSubCell >
void libMesh::QComposite< QSubCell >::init ( const Elem elem,
const std::vector< Real > &  vertex_distance_func,
unsigned int  p_level = 0 
)
overridevirtual

Overrides the base class init() function, and uses the ElemCutter to subdivide the element into "inside" and "outside" subelements.

Definition at line 67 of file quadrature_composite.C.

References libMesh::Elem::dim(), libMesh::Elem::n_vertices(), libMesh::Elem::reference_elem(), and libMesh::Elem::type().

70 {
71  libmesh_assert_equal_to (vertex_distance_func.size(), elem.n_vertices());
72  libmesh_assert_equal_to (_dim, elem.dim());
73 
74  // if we are not cut, revert to simple base class init() method.
75  if (!_elem_cutter.is_cut (elem, vertex_distance_func))
76  {
77  _q_subcell.init (elem.type(), p_level);
78  _points = _q_subcell.get_points();
79  _weights = _q_subcell.get_weights();
80 
81  //this->print_info();
82  return;
83  }
84 
85  // Get a pointer to the element's reference element. We want to
86  // perform cutting on the reference element such that the quadrature
87  // point locations of the subelements live in the reference
88  // coordinate system, thereby eliminating the need for inverse
89  // mapping.
90  const Elem * reference_elem = elem.reference_elem();
91 
92  libmesh_assert (reference_elem != nullptr);
93 
94  _elem_cutter(*reference_elem, vertex_distance_func);
95  //_elem_cutter(elem, vertex_distance_func);
96 
97  // clear our state & accumulate points from subelements
98  _points.clear();
99  _weights.clear();
100 
101  // inside subelem
102  {
103  const std::vector<Elem const *> & inside_elem (_elem_cutter.inside_elements());
104  std::cout << inside_elem.size() << " elements inside\n";
105 
106  this->add_subelem_values(inside_elem);
107  }
108 
109  // outside subelem
110  {
111  const std::vector<Elem const *> & outside_elem (_elem_cutter.outside_elements());
112  std::cout << outside_elem.size() << " elements outside\n";
113 
114  this->add_subelem_values(outside_elem);
115  }
116 
117  this->print_info();
118 }
bool is_cut(const Elem &elem, const std::vector< Real > &vertex_distance_func) const
Definition: elem_cutter.C:88
const std::vector< Elem const * > & outside_elements() const
Definition: elem_cutter.h:124
const std::vector< Elem const * > & inside_elements() const
Definition: elem_cutter.h:116
void add_subelem_values(const std::vector< Elem const *> &subelem)

◆ operator=() [1/2]

template<class QSubCell >
QComposite& libMesh::QComposite< QSubCell >::operator= ( const QComposite< QSubCell > &  )
delete

◆ operator=() [2/2]

template<class QSubCell >
QComposite& libMesh::QComposite< QSubCell >::operator= ( QComposite< QSubCell > &&  )
default

◆ type()

template<class QSubCell >
QuadratureType libMesh::QComposite< QSubCell >::type ( ) const
overridevirtual
Returns
QCOMPOSITE.

Definition at line 37 of file quadrature_composite.C.

References libMesh::QCOMPOSITE.

38 {
39  return QCOMPOSITE;
40 }

Member Data Documentation

◆ _elem_cutter

template<class QSubCell >
ElemCutter libMesh::QComposite< QSubCell >::_elem_cutter
private

ElemCutter object.

Definition at line 114 of file quadrature_composite.h.

◆ _lagrange_fe

template<class QSubCell >
std::unique_ptr<FEBase> libMesh::QComposite< QSubCell >::_lagrange_fe
private

Lagrange FE to use for subcell mapping.

Definition at line 119 of file quadrature_composite.h.

Referenced by libMesh::QComposite< QSubCell >::QComposite().

◆ _q_subcell

template<class QSubCell >
QSubCell libMesh::QComposite< QSubCell >::_q_subcell
private

Subcell quadrature object.

Definition at line 109 of file quadrature_composite.h.

Referenced by libMesh::QComposite< QSubCell >::QComposite().


The documentation for this class was generated from the following files: