#include <agglomeration_accessor.h>
Public Types | |
using | AgglomerationContainer |
Public Member Functions | |
void | get_dof_indices (std::vector< types::global_dof_index > &) const |
unsigned int | n_faces () const |
unsigned int | n_agglomerated_faces () const |
const AgglomerationIterator< dim, spacedim > | neighbor (const unsigned int f) const |
unsigned int | neighbor_of_agglomerated_neighbor (const unsigned int f) const |
bool | at_boundary (const unsigned int f) const |
const std::vector< typename Triangulation< dim >::active_face_iterator > & | polytope_boundary () const |
double | volume () const |
double | diameter () const |
AgglomerationContainer | get_agglomerate () const |
const BoundingBox< dim > & | get_bounding_box () const |
types::global_cell_index | index () const |
DoFHandler< dim, spacedim >::active_cell_iterator | as_dof_handler_iterator (const DoFHandler< dim, spacedim > &dof_handler) const |
unsigned int | n_background_cells () const |
bool | is_locally_owned () const |
CellId | id () const |
types::subdomain_id | subdomain_id () const |
const std::vector< types::global_cell_index > & | children () const |
Private Member Functions | |
AgglomerationAccessor () | |
AgglomerationAccessor (const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell, const AgglomerationHandler< dim, spacedim > *ah) | |
AgglomerationAccessor (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const CellId &cell_id, const AgglomerationHandler< dim, spacedim > *ah) | |
~AgglomerationAccessor ()=default | |
bool | operator== (const AgglomerationAccessor< dim, spacedim > &other) const |
bool | operator!= (const AgglomerationAccessor< dim, spacedim > &other) const |
void | next () |
void | prev () |
const AgglomerationContainer & | get_slaves () const |
unsigned int | n_agglomerated_faces_per_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const |
Private Attributes | |
Triangulation< dim, spacedim >::active_cell_iterator | master_cell |
types::global_cell_index | present_index |
CellId | present_id |
types::subdomain_id | present_subdomain_id |
AgglomerationHandler< dim, spacedim > * | handler |
Friends | |
template<int, int> | |
class | AgglomerationIterator |
Accessor class used by AgglomerationIterator to access agglomeration data.
Definition at line 42 of file agglomeration_accessor.h.
using AgglomerationAccessor< dim, spacedim >::AgglomerationContainer |
Type for storing the polygons in an agglomerate.
Definition at line 48 of file agglomeration_accessor.h.
|
inlineprivate |
Private default constructor. This is not supposed to be used and hence will throw.
Definition at line 462 of file agglomeration_accessor.h.
Referenced by operator!=(), and operator==().
|
inlineprivate |
Private constructor for an agglomerate. This is meant to be invoked by the AgglomerationIterator class. It takes as input the master cell of the agglomerate and a pointer to the handler.
Definition at line 468 of file agglomeration_accessor.h.
References handler, numbers::invalid_subdomain_id, master_cell, present_id, present_index, and present_subdomain_id.
|
inlineprivate |
Same as above, but needed when the argument cells
is a ghost cell.
Definition at line 492 of file agglomeration_accessor.h.
References handler, numbers::invalid_unsigned_int, master_cell, present_id, present_index, and present_subdomain_id.
|
privatedefault |
Default destructor.
DoFHandler< dim, spacedim >::active_cell_iterator AgglomerationAccessor< dim, spacedim >::as_dof_handler_iterator | ( | const DoFHandler< dim, spacedim > & | dof_handler | ) | const |
Returns an active cell iterator for the dof_handler, matching the polytope referenced by the input iterator. The type of the returned object is a DoFHandler::active_cell_iterator which can be used to initialize FiniteElement data.
Definition at line 670 of file agglomeration_accessor.h.
References master_cell.
Referenced by AgglomerationHandler< dim, dim >::create_bounding_box(), AgglomerationHandler< dim, spacedim >::reinit(), and AgglomerationHandler< dim, spacedim >::reinit().
|
inline |
This function generalizes the behaviour of cell->face(f)->at_boundary() in the case where f is an index out of the range [0,..., n_faces). In practice, if you call this function with a standard deal.II cell, you have precisely the same result as calling cell->face(f)->at_boundary(). Otherwise, if the cell is a master one, you have a boolean returning true is that face for the agglomeration is on the boundary or not.
Definition at line 714 of file agglomeration_accessor.h.
References handler, master_cell, and present_id.
Referenced by neighbor(), and neighbor_of_agglomerated_neighbor().
|
inline |
Returns a vector of indices identifying the children polytopes.
Definition at line 779 of file agglomeration_accessor.h.
References handler, and present_index.
|
inline |
Return the diameter of the present polytopal element.
Definition at line 560 of file agglomeration_accessor.h.
References handler, master_cell, and present_index.
|
inline |
Returns the deal.II cells that build the agglomerate.
Definition at line 540 of file agglomeration_accessor.h.
References get_slaves(), and master_cell.
Referenced by n_agglomerated_faces(), n_background_cells(), and AgglomerationHandler< dim, spacedim >::reinit().
|
inline |
Return the BoundingBox which bounds the present polytope. In case the present polytope is not locally owned, it returns the BoundingBox of that ghosted polytope.
Definition at line 583 of file agglomeration_accessor.h.
References handler, is_locally_owned(), present_id, present_index, and present_subdomain_id.
|
inline |
Get the DoFs indices associated to the agglomerate.
Definition at line 512 of file agglomeration_accessor.h.
References handler, is_locally_owned(), master_cell, present_id, and present_subdomain_id.
Referenced by AgglomerationHandler< dim, spacedim >::create_agglomeration_sparsity_pattern().
|
inlineprivate |
Returns the slaves of the present agglomeration.
Definition at line 682 of file agglomeration_accessor.h.
References handler, and master_cell.
Referenced by get_agglomerate().
|
inline |
The polytopal analogue of CellAccessor::id(). It provides a way to uniquely identify cells in a parallel Triangulation such as a parallel::distributed::Triangulation.
Definition at line 763 of file agglomeration_accessor.h.
References present_id.
Referenced by dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
inline |
Return the index of the present polytope.
Definition at line 661 of file agglomeration_accessor.h.
References present_index.
|
inline |
Definition at line 754 of file agglomeration_accessor.h.
References master_cell.
Referenced by get_bounding_box(), and get_dof_indices().
unsigned int AgglomerationAccessor< dim, spacedim >::n_agglomerated_faces | ( | ) | const |
Return the number of deal.II faces that is building a polygon face.
Definition at line 701 of file agglomeration_accessor.h.
References get_agglomerate(), and n_agglomerated_faces_per_cell().
|
private |
Definition at line 281 of file agglomeration_accessor.h.
References handler.
Referenced by n_agglomerated_faces().
|
inline |
Returns the number of classical deal.II cells that are building the present polygon.
Definition at line 691 of file agglomeration_accessor.h.
References get_agglomerate().
unsigned int AgglomerationAccessor< dim, spacedim >::n_faces | ( | ) | const |
Return, for a cell, the number of faces. In case the cell is a standard cell, then the number of faces is the classical one. If it's a master cell, then it returns the number of faces of the agglomeration identified by the master cell itself.
Definition at line 302 of file agglomeration_accessor.h.
References handler, master_cell, and present_index.
const AgglomerationIterator< dim, spacedim > AgglomerationAccessor< dim, spacedim >::neighbor | ( | const unsigned int | f | ) | const |
Return the agglomerate which shares face f.
Definition at line 313 of file agglomeration_accessor.h.
References at_boundary(), handler, master_cell, and present_id.
Referenced by AgglomerationHandler< dim, spacedim >::create_agglomeration_sparsity_pattern(), neighbor_of_agglomerated_neighbor(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
unsigned int AgglomerationAccessor< dim, spacedim >::neighbor_of_agglomerated_neighbor | ( | const unsigned int | f | ) | const |
Return the present index (seen from the neighboring agglomerate) of the present face f.
Definition at line 404 of file agglomeration_accessor.h.
References at_boundary(), handler, numbers::invalid_unsigned_int, neighbor(), present_id, and IteratorState::valid.
|
inlineprivate |
Move to the next cell in the polytopal mesh.
Definition at line 614 of file agglomeration_accessor.h.
References handler, master_cell, present_id, present_index, and present_subdomain_id.
|
inlineprivate |
Compare for inequality.
Definition at line 651 of file agglomeration_accessor.h.
References AgglomerationAccessor().
|
inlineprivate |
Comparison operator for Accessor. Two accessors are equal if they refer to the same polytopal element.
Definition at line 643 of file agglomeration_accessor.h.
References AgglomerationAccessor(), and present_index.
|
inline |
Return a vector of face iterators describing the boundary of agglomerate.
Definition at line 551 of file agglomeration_accessor.h.
References handler, and master_cell.
|
inlineprivate |
Move to the previous cell in the polytopal mesh.
Definition at line 632 of file agglomeration_accessor.h.
References handler, master_cell, present_id, and present_index.
|
inline |
The polytopal analogue of CellAccessor::subdomain_id(). In case of a serial Triangulation, it returns the numbers::invalid_subdomain_id.
Definition at line 772 of file agglomeration_accessor.h.
References present_subdomain_id.
|
inline |
Return the volume of a polytope.
Definition at line 595 of file agglomeration_accessor.h.
References handler, master_cell, and present_index.
Definition at line 274 of file agglomeration_accessor.h.
|
private |
A pointer to the Handler.
Definition at line 235 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), at_boundary(), children(), diameter(), get_bounding_box(), get_dof_indices(), get_slaves(), n_agglomerated_faces_per_cell(), n_faces(), neighbor(), neighbor_of_agglomerated_neighbor(), next(), polytope_boundary(), prev(), and volume().
|
private |
The unique deal.II cell associated to the present polytope.
Definition at line 215 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), as_dof_handler_iterator(), at_boundary(), diameter(), get_agglomerate(), get_dof_indices(), get_slaves(), is_locally_owned(), n_faces(), neighbor(), next(), polytope_boundary(), prev(), and volume().
|
private |
The index of the present polytope.
Definition at line 225 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), at_boundary(), get_bounding_box(), get_dof_indices(), id(), neighbor(), neighbor_of_agglomerated_neighbor(), next(), and prev().
|
private |
The index of the present polytope.
Definition at line 220 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), children(), diameter(), get_bounding_box(), index(), n_faces(), next(), operator==(), prev(), and volume().
|
private |
The rank owning of the present polytope.
Definition at line 230 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), get_bounding_box(), get_dof_indices(), next(), and subdomain_id().