#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 |
const FiniteElement< dim, spacedim > & | get_fe () const |
void | set_active_fe_index (const types::fe_index index) const |
types::fe_index | active_fe_index () 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 486 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 492 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 516 of file agglomeration_accessor.h.
References handler, numbers::invalid_unsigned_int, master_cell, present_id, present_index, and present_subdomain_id.
|
privatedefault |
Default destructor.
|
inline |
Returns the index of the active finite element. When using hp::FECollection, this function retrieves the index of the finite element assigned to the current polytope.
Definition at line 835 of file agglomeration_accessor.h.
References handler, and master_cell.
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 694 of file agglomeration_accessor.h.
References master_cell.
Referenced by AgglomerationHandler< dim, spacedim >::reinit(), AgglomerationHandler< dim, spacedim >::reinit(), and AgglomerationHandler< dim, spacedim >::reinit_interface().
|
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 738 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 803 of file agglomeration_accessor.h.
References handler, and present_index.
|
inline |
Return the diameter of the present polytopal element.
Definition at line 584 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 564 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 607 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 536 of file agglomeration_accessor.h.
References handler, is_locally_owned(), master_cell, present_id, and present_subdomain_id.
|
inline |
Returns the FiniteElement object used by the current polytope. This function should only be called after the corresponding agglomeration handler has invoked distribute_agglomerated_dofs().
Definition at line 812 of file agglomeration_accessor.h.
References DoFHandler< int dim, int spacedim >::get_fe(), handler, and master_cell.
Referenced by AgglomerationHandler< dim, spacedim >::reinit().
|
inlineprivate |
Returns the slaves of the present agglomeration.
Definition at line 706 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 787 of file agglomeration_accessor.h.
References present_id.
Referenced by AgglomerationHandler< dim, spacedim >::reinit_interface().
|
inline |
Return the index of the present polytope.
Definition at line 685 of file agglomeration_accessor.h.
References present_index.
Referenced by set_active_fe_index().
|
inline |
Definition at line 778 of file agglomeration_accessor.h.
References master_cell.
Referenced by get_bounding_box(), get_dof_indices(), and AgglomerationHandler< dim, spacedim >::reinit_interface().
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 725 of file agglomeration_accessor.h.
References get_agglomerate(), and n_agglomerated_faces_per_cell().
|
private |
Definition at line 305 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 715 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 326 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 337 of file agglomeration_accessor.h.
References at_boundary(), handler, master_cell, and present_id.
Referenced by 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 428 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 638 of file agglomeration_accessor.h.
References handler, master_cell, present_id, present_index, and present_subdomain_id.
|
inlineprivate |
Compare for inequality.
Definition at line 675 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 667 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 575 of file agglomeration_accessor.h.
References handler, and master_cell.
|
inlineprivate |
Move to the previous cell in the polytopal mesh.
Definition at line 656 of file agglomeration_accessor.h.
References handler, master_cell, present_id, and present_index.
|
inline |
Sets the active finite element index. This function should be called when using hp::FECollection to specify which finite element in the collection is assigned to the current polytope.
Definition at line 822 of file agglomeration_accessor.h.
References handler, index(), and master_cell.
|
inline |
The polytopal analogue of CellAccessor::subdomain_id(). In case of a serial Triangulation, it returns the numbers::invalid_subdomain_id.
Definition at line 796 of file agglomeration_accessor.h.
References present_subdomain_id.
Referenced by AgglomerationHandler< dim, spacedim >::reinit_interface().
|
inline |
Return the volume of a polytope.
Definition at line 619 of file agglomeration_accessor.h.
References handler, master_cell, and present_index.
Definition at line 298 of file agglomeration_accessor.h.
|
private |
A pointer to the Handler.
Definition at line 259 of file agglomeration_accessor.h.
Referenced by active_fe_index(), AgglomerationAccessor(), AgglomerationAccessor(), at_boundary(), children(), diameter(), get_bounding_box(), get_dof_indices(), get_fe(), get_slaves(), n_agglomerated_faces_per_cell(), n_faces(), neighbor(), neighbor_of_agglomerated_neighbor(), next(), polytope_boundary(), prev(), set_active_fe_index(), and volume().
|
private |
The unique deal.II cell associated to the present polytope.
Definition at line 239 of file agglomeration_accessor.h.
Referenced by active_fe_index(), AgglomerationAccessor(), AgglomerationAccessor(), as_dof_handler_iterator(), at_boundary(), diameter(), get_agglomerate(), get_dof_indices(), get_fe(), get_slaves(), is_locally_owned(), n_faces(), neighbor(), next(), polytope_boundary(), prev(), set_active_fe_index(), and volume().
|
private |
The index of the present polytope.
Definition at line 249 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 244 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 254 of file agglomeration_accessor.h.
Referenced by AgglomerationAccessor(), AgglomerationAccessor(), get_bounding_box(), get_dof_indices(), next(), and subdomain_id().