PolyDEAL
 
Loading...
Searching...
No Matches
AgglomerationAccessor< dim, spacedim > Class Template Reference

#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 AgglomerationContainerget_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
 

Detailed Description

template<int dim, int spacedim = dim>
class AgglomerationAccessor< dim, spacedim >

Accessor class used by AgglomerationIterator to access agglomeration data.

Definition at line 42 of file agglomeration_accessor.h.

Member Typedef Documentation

◆ AgglomerationContainer

template<int dim, int spacedim = dim>
using AgglomerationAccessor< dim, spacedim >::AgglomerationContainer
Initial value:
std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>

Type for storing the polygons in an agglomerate.

Definition at line 48 of file agglomeration_accessor.h.

Constructor & Destructor Documentation

◆ AgglomerationAccessor() [1/3]

template<int dim, int spacedim>
AgglomerationAccessor< dim, spacedim >::AgglomerationAccessor ( )
inlineprivate

Private default constructor. This is not supposed to be used and hence will throw.

Definition at line 462 of file agglomeration_accessor.h.

◆ AgglomerationAccessor() [2/3]

template<int dim, int spacedim>
AgglomerationAccessor< dim, spacedim >::AgglomerationAccessor ( const typename Triangulation< dim, spacedim >::active_cell_iterator & master_cell,
const AgglomerationHandler< dim, spacedim > * ah )
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 Triangulation< int dim, int spacedim >::end(), numbers::invalid_subdomain_id, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ AgglomerationAccessor() [3/3]

template<int dim, int spacedim>
AgglomerationAccessor< dim, spacedim >::AgglomerationAccessor ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell,
const CellId & cell_id,
const AgglomerationHandler< dim, spacedim > * ah )
inlineprivate

Same as above, but needed when the argument cells is a ghost cell.

Definition at line 492 of file agglomeration_accessor.h.

References Assert, numbers::invalid_unsigned_int, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ ~AgglomerationAccessor()

template<int dim, int spacedim = dim>
AgglomerationAccessor< dim, spacedim >::~AgglomerationAccessor ( )
privatedefault

Default destructor.

Member Function Documentation

◆ as_dof_handler_iterator()

template<int dim, int spacedim>
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 AgglomerationIterator< dim, spacedim >::master_cell().

◆ at_boundary()

template<int dim, int spacedim>
bool AgglomerationAccessor< dim, spacedim >::at_boundary ( const unsigned int f) const
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 Assert, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ children()

template<int dim, int spacedim>
const std::vector< types::global_cell_index > & AgglomerationAccessor< dim, spacedim >::children ( ) const
inline

Returns a vector of indices identifying the children polytopes.

Definition at line 779 of file agglomeration_accessor.h.

References Assert.

◆ diameter()

template<int dim, int spacedim>
double AgglomerationAccessor< dim, spacedim >::diameter ( ) const
inline

Return the diameter of the present polytopal element.

Definition at line 560 of file agglomeration_accessor.h.

References Assert, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ get_agglomerate()

template<int dim, int spacedim>
AgglomerationAccessor< dim, spacedim >::AgglomerationContainer AgglomerationAccessor< dim, spacedim >::get_agglomerate ( ) const
inline

Returns the deal.II cells that build the agglomerate.

Definition at line 540 of file agglomeration_accessor.h.

References AgglomerationIterator< dim, spacedim >::master_cell().

◆ get_bounding_box()

template<int dim, int spacedim>
const BoundingBox< dim > & AgglomerationAccessor< dim, spacedim >::get_bounding_box ( ) const
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.

◆ get_dof_indices()

template<int dim, int spacedim>
void AgglomerationAccessor< dim, spacedim >::get_dof_indices ( std::vector< types::global_dof_index > & dof_indices) const
inline

Get the DoFs indices associated to the agglomerate.

Definition at line 512 of file agglomeration_accessor.h.

References Assert, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ get_slaves()

template<int dim, int spacedim>
const AgglomerationAccessor< dim, spacedim >::AgglomerationContainer & AgglomerationAccessor< dim, spacedim >::get_slaves ( ) const
inlineprivate

Returns the slaves of the present agglomeration.

Definition at line 682 of file agglomeration_accessor.h.

References AgglomerationIterator< dim, spacedim >::master_cell().

◆ id()

template<int dim, int spacedim>
CellId AgglomerationAccessor< dim, spacedim >::id ( ) const
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.

◆ index()

template<int dim, int spacedim>
types::global_cell_index AgglomerationAccessor< dim, spacedim >::index ( ) const
inline

Return the index of the present polytope.

Definition at line 661 of file agglomeration_accessor.h.

◆ is_locally_owned()

template<int dim, int spacedim>
bool AgglomerationAccessor< dim, spacedim >::is_locally_owned ( ) const
inline

◆ n_agglomerated_faces()

template<int dim, int spacedim>
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.

◆ n_agglomerated_faces_per_cell()

template<int dim, int spacedim>
unsigned int AgglomerationAccessor< dim, spacedim >::n_agglomerated_faces_per_cell ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell) const
private

Definition at line 281 of file agglomeration_accessor.h.

◆ n_background_cells()

template<int dim, int spacedim>
unsigned int AgglomerationAccessor< dim, spacedim >::n_background_cells ( ) const
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 AssertThrow.

◆ n_faces()

template<int dim, int spacedim>
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 Assert, and AgglomerationIterator< dim, spacedim >::master_cell().

◆ neighbor()

template<int dim, int spacedim>
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 AgglomerationIterator< dim, spacedim >::master_cell(), and second.

◆ neighbor_of_agglomerated_neighbor()

template<int dim, int spacedim>
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 AssertThrow, numbers::invalid_unsigned_int, and IteratorState::valid.

◆ next()

template<int dim, int spacedim>
void AgglomerationAccessor< dim, spacedim >::next ( )
inlineprivate

Move to the next cell in the polytopal mesh.

Definition at line 614 of file agglomeration_accessor.h.

References AgglomerationIterator< dim, spacedim >::master_cell().

◆ operator!=()

template<int dim, int spacedim>
bool AgglomerationAccessor< dim, spacedim >::operator!= ( const AgglomerationAccessor< dim, spacedim > & other) const
inlineprivate

Compare for inequality.

Definition at line 651 of file agglomeration_accessor.h.

◆ operator==()

template<int dim, int spacedim>
bool AgglomerationAccessor< dim, spacedim >::operator== ( const AgglomerationAccessor< dim, spacedim > & other) const
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< dim, spacedim >::present_index.

◆ polytope_boundary()

template<int dim, int spacedim>
const std::vector< typename Triangulation< dim >::active_face_iterator > & AgglomerationAccessor< dim, spacedim >::polytope_boundary ( ) const
inline

Return a vector of face iterators describing the boundary of agglomerate.

Definition at line 551 of file agglomeration_accessor.h.

References AgglomerationIterator< dim, spacedim >::master_cell().

◆ prev()

template<int dim, int spacedim>
void AgglomerationAccessor< dim, spacedim >::prev ( )
inlineprivate

Move to the previous cell in the polytopal mesh.

Definition at line 632 of file agglomeration_accessor.h.

References AgglomerationIterator< dim, spacedim >::master_cell().

◆ subdomain_id()

template<int dim, int spacedim>
types::subdomain_id AgglomerationAccessor< dim, spacedim >::subdomain_id ( ) const
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.

◆ volume()

template<int dim, int spacedim>
double AgglomerationAccessor< dim, spacedim >::volume ( ) const
inline

Return the volume of a polytope.

Definition at line 595 of file agglomeration_accessor.h.

References Assert, and AgglomerationIterator< dim, spacedim >::master_cell().

Friends And Related Symbol Documentation

◆ AgglomerationIterator

template<int dim, int spacedim = dim>
template<int , int >
friend class AgglomerationIterator
friend

Definition at line 274 of file agglomeration_accessor.h.

Member Data Documentation

◆ handler

template<int dim, int spacedim = dim>
AgglomerationHandler<dim, spacedim>* AgglomerationAccessor< dim, spacedim >::handler
private

A pointer to the Handler.

Definition at line 235 of file agglomeration_accessor.h.

◆ master_cell

template<int dim, int spacedim = dim>
Triangulation<dim,spacedim>::active_cell_iterator AgglomerationAccessor< dim, spacedim >::master_cell
private

The unique deal.II cell associated to the present polytope.

Definition at line 215 of file agglomeration_accessor.h.

◆ present_id

template<int dim, int spacedim = dim>
CellId AgglomerationAccessor< dim, spacedim >::present_id
private

The index of the present polytope.

Definition at line 225 of file agglomeration_accessor.h.

◆ present_index

template<int dim, int spacedim = dim>
types::global_cell_index AgglomerationAccessor< dim, spacedim >::present_index
private

The index of the present polytope.

Definition at line 220 of file agglomeration_accessor.h.

Referenced by AgglomerationAccessor< dim, spacedim >::operator==().

◆ present_subdomain_id

template<int dim, int spacedim = dim>
types::subdomain_id AgglomerationAccessor< dim, spacedim >::present_subdomain_id
private

The rank owning of the present polytope.

Definition at line 230 of file agglomeration_accessor.h.


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