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

#include <agglomeration_handler.h>

Inheritance diagram for AgglomerationHandler< dim, spacedim >:
EnableObserverPointer

Public Types

enum  CellAgglomerationType { master = 0 , slave = 1 }
using agglomeration_iterator = AgglomerationIterator<dim, spacedim>
using AgglomerationContainer

Public Member Functions

 AgglomerationHandler (const GridTools::Cache< dim, spacedim > &cached_tria)
 AgglomerationHandler ()=default
 ~AgglomerationHandler ()
agglomeration_iterator begin () const
agglomeration_iterator begin ()
agglomeration_iterator end () const
agglomeration_iterator end ()
agglomeration_iterator last ()
IteratorRange< agglomeration_iteratorpolytope_iterators () const
void distribute_agglomerated_dofs (const FiniteElement< dim > &fe_space)
void distribute_agglomerated_dofs (const hp::FECollection< dim, spacedim > &fe_collection_in)
void initialize_fe_values (const Quadrature< dim > &cell_quadrature=QGauss< dim >(1), const UpdateFlags &flags=UpdateFlags::update_default, const Quadrature< dim - 1 > &face_quadrature=QGauss< dim - 1 >(1), const UpdateFlags &face_flags=UpdateFlags::update_default)
void initialize_fe_values (const hp::QCollection< dim > &cell_qcollection=hp::QCollection< dim >(QGauss< dim >(1)), const UpdateFlags &flags=UpdateFlags::update_default, const hp::QCollection< dim - 1 > &face_qcollection=hp::QCollection< dim - 1 >(QGauss< dim - 1 >(1)), const UpdateFlags &face_flags=UpdateFlags::update_default)
template<typename SparsityPatternType, typename Number = double>
void create_agglomeration_sparsity_pattern (SparsityPatternType &sparsity_pattern, const AffineConstraints< Number > &constraints=AffineConstraints< Number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
agglomeration_iterator define_agglomerate (const AgglomerationContainer &cells)
agglomeration_iterator define_agglomerate (const AgglomerationContainer &cells, const unsigned int fecollection_size)
void define_agglomerate_with_check (const AgglomerationContainer &cells)
const Triangulation< dim, spacedim > & get_triangulation () const
const FiniteElement< dim, spacedim > & get_fe () const
const Mapping< dim > & get_mapping () const
const MappingBox< dim > & get_agglomeration_mapping () const
const std::vector< BoundingBox< dim > > & get_local_bboxes () const
double get_mesh_size () const
types::global_cell_index cell_to_polytope_index (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
decltype(auto) get_interface () const
template<typename CellIterator>
bool is_master_cell (const CellIterator &cell) const
const std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > & get_slaves_of_idx (types::global_cell_index idx) const
const LinearAlgebra::distributed::Vector< float > & get_relationships () const
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > get_agglomerate (const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
template<class StreamType>
void print_agglomeration (StreamType &out)
const DoFHandler< dim, spacedim > & get_dof_handler () const
unsigned int n_agglomerates () const
unsigned int n_agglomerated_faces_per_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
const FEValues< dim, spacedim > & reinit (const AgglomerationIterator< dim, spacedim > &polytope) const
const FEValuesBase< dim, spacedim > & reinit (const AgglomerationIterator< dim, spacedim > &polytope, const unsigned int face_index) const
std::pair< const FEValuesBase< dim, spacedim > &, const FEValuesBase< dim, spacedim > & > reinit_interface (const AgglomerationIterator< dim, spacedim > &polytope_in, const AgglomerationIterator< dim, spacedim > &neigh_polytope, const unsigned int local_in, const unsigned int local_outside) const
Quadrature< dim > agglomerated_quadrature (const AgglomerationContainer &cells, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
bool at_boundary (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int f) const
unsigned int n_dofs_per_cell () const noexcept
types::global_dof_index n_dofs () const noexcept
const std::vector< typename Triangulation< dim >::active_face_iterator > & polytope_boundary (const typename Triangulation< dim >::active_cell_iterator &cell)
void setup_ghost_polytopes ()
void exchange_interface_values ()
const DoFHandler< dim, spacedim >::active_cell_iterator polytope_to_dh_iterator (const types::global_cell_index polytope_index) const
template<typename RtreeType>
void connect_hierarchy (const CellsAgglomerator< dim, RtreeType > &agglomerator)
const hp::FECollection< dim, spacedim > & get_fe_collection () const
bool used_fe_collection () const
Public Member Functions inherited from EnableObserverPointer
 EnableObserverPointer ()
 EnableObserverPointer (const EnableObserverPointer &)
 EnableObserverPointer (EnableObserverPointer &&) noexcept
virtual ~EnableObserverPointer ()
EnableObserverPointeroperator= (const EnableObserverPointer &)
EnableObserverPointeroperator= (EnableObserverPointer &&) noexcept
void serialize (Archive &ar, const unsigned int version)
unsigned int n_subscriptions () const
void list_subscribers (StreamType &stream) const
void list_subscribers () const
unsigned int n_subscriptions () const
void list_subscribers (StreamType &stream) const
void list_subscribers () const

Public Attributes

DoFHandler< dim, spacedim > agglo_dh
DoFHandler< dim, spacedim > output_dh
std::unique_ptr< MappingBox< dim > > box_mapping
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Point< spacedim > > > > recv_qpoints
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< double > > > recv_jxws
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Tensor< 1, spacedim > > > > recv_normals
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< double > > > > recv_values
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< std::vector< Tensor< 1, spacedim > > > > > recv_gradients

Private Types

using ScratchData = MeshWorker::ScratchData<dim, spacedim>

Private Member Functions

void initialize_agglomeration_data (const std::unique_ptr< GridTools::Cache< dim, spacedim > > &cache_tria)
void update_agglomerate (AgglomerationContainer &polytope, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell)
void connect_to_tria_signals ()
Triangulation< dim, spacedim >::active_cell_iterator & is_slave_cell_of (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
void create_bounding_box (const AgglomerationContainer &polytope)
types::global_cell_index get_master_idx_of_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
bool are_cells_agglomerated (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const typename Triangulation< dim, spacedim >::active_cell_iterator &other_cell) const
void initialize_hp_structure ()
const FEValuesBase< dim, spacedim > & reinit_master (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_number, std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > &agglo_isv_ptr) const
template<typename CellIterator>
bool is_slave_cell (const CellIterator &cell) const
void setup_connectivity_of_agglomeration ()

Private Attributes

unsigned int n_agglomerations
LinearAlgebra::distributed::Vector< float > master_slave_relationships
std::map< types::global_cell_index, typename Triangulation< dim, spacedim >::active_cell_iterator > master_slave_relationships_iterators
std::vector< types::global_cell_indexnumber_of_agglomerated_faces
std::map< const typename Triangulation< dim, spacedim >::active_cell_iterator, std::vector< typename Triangulation< dim >::active_face_iterator > > polygon_boundary
LinearAlgebra::distributed::Vector< double > euler_vector
std::unique_ptr< ScratchDatastandard_scratch
std::unique_ptr< ScratchDataagglomerated_scratch
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_neigh
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_bdary
boost::signals2::connection tria_listener
UpdateFlags agglomeration_flags
const UpdateFlags internal_agglomeration_flags
UpdateFlags agglomeration_face_flags
const UpdateFlags internal_agglomeration_face_flags
Quadrature< dim > agglomeration_quad
Quadrature< dim - 1 > agglomeration_face_quad
std::unordered_map< types::global_cell_index, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > > master2slaves
std::map< types::global_cell_index, types::global_cell_indexmaster2polygon
std::vector< typename Triangulation< dim >::active_cell_iterator > master_disconnected
FE_Nothing< dim, spacedim > dummy_fe
std::unique_ptr< FEValues< dim, spacedim > > no_values
std::unique_ptr< FEFaceValues< dim, spacedim > > no_face_values
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > master_cells_container
internal::PolytopeCache< dim, spacedim > polytope_cache
bool hybrid_mesh
std::map< std::pair< types::global_cell_index, types::global_cell_index >, std::vector< types::global_cell_index > > parent_child_info
unsigned int present_extraction_level
bool is_hp_collection = false
std::unique_ptr< hp::FECollection< dim, spacedim > > hp_fe_collection
hp::QCollection< dim > agglomeration_quad_collection
hp::QCollection< dim - 1 > agglomeration_face_quad_collection
hp::MappingCollection< dim > mapping_collection
hp::FECollection< dim, spacedim > dummy_fe_collection
std::unique_ptr< hp::FEValues< dim, spacedim > > hp_no_values
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > hp_no_face_values

Friends

template<int, int>
class AgglomerationIterator
template<int, int>
class AgglomerationAccessor
class internal::AgglomerationHandlerImplementation< dim, spacedim >

Additional Inherited Members

Static Public Member Functions inherited from EnableObserverPointer
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)

Detailed Description

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

#TODO: Documentation.

Definition at line 172 of file agglomeration_handler.h.

Member Typedef Documentation

◆ agglomeration_iterator

template<int dim, int spacedim = dim>
using AgglomerationHandler< dim, spacedim >::agglomeration_iterator = AgglomerationIterator<dim, spacedim>

Definition at line 175 of file agglomeration_handler.h.

◆ AgglomerationContainer

template<int dim, int spacedim = dim>
using AgglomerationHandler< dim, spacedim >::AgglomerationContainer
Initial value:
typename AgglomerationAccessor< dim, spacedim >::AgglomerationContainer AgglomerationContainer

Definition at line 177 of file agglomeration_handler.h.

◆ ScratchData

template<int dim, int spacedim = dim>
using AgglomerationHandler< dim, spacedim >::ScratchData = MeshWorker::ScratchData<dim, spacedim>
private

Definition at line 698 of file agglomeration_handler.h.

Member Enumeration Documentation

◆ CellAgglomerationType

template<int dim, int spacedim = dim>
enum AgglomerationHandler::CellAgglomerationType
Enumerator
master 
slave 

Definition at line 181 of file agglomeration_handler.h.

Constructor & Destructor Documentation

◆ AgglomerationHandler() [1/2]

template<int dim, int spacedim>
AgglomerationHandler< dim, spacedim >::AgglomerationHandler ( const GridTools::Cache< dim, spacedim > & cached_tria)
explicit

◆ AgglomerationHandler() [2/2]

template<int dim, int spacedim = dim>
AgglomerationHandler< dim, spacedim >::AgglomerationHandler ( )
default

◆ ~AgglomerationHandler()

template<int dim, int spacedim = dim>
AgglomerationHandler< dim, spacedim >::~AgglomerationHandler ( )
inline

Definition at line 194 of file agglomeration_handler.h.

References tria_listener.

Member Function Documentation

◆ agglomerated_quadrature()

template<int dim, int spacedim = dim>
Quadrature< dim > AgglomerationHandler< dim, spacedim >::agglomerated_quadrature ( const AgglomerationContainer & cells,
const typename Triangulation< dim, spacedim >::active_cell_iterator & master_cell ) const

Return the agglomerated quadrature for the given agglomeration. This amounts to loop over all cells in an agglomeration and collecting together all the rules.

Definition at line 624 of file agglomeration_handler.cc.

References is_hp_collection, is_master_cell(), and no_values.

Referenced by reinit().

◆ are_cells_agglomerated()

template<int dim, int spacedim>
bool AgglomerationHandler< dim, spacedim >::are_cells_agglomerated ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell,
const typename Triangulation< dim, spacedim >::active_cell_iterator & other_cell ) const
inlineprivate

Returns true if the two given cells are agglomerated together.

Definition at line 1151 of file agglomeration_handler.h.

References get_master_idx_of_cell().

Referenced by n_agglomerated_faces_per_cell(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ at_boundary()

template<int dim, int spacedim>
bool AgglomerationHandler< dim, spacedim >::at_boundary ( const typename DoFHandler< dim, spacedim >::active_cell_iterator & cell,
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 1084 of file agglomeration_handler.h.

References is_slave_cell(), master2polygon, and polytope_cache.

◆ begin() [1/2]

template<int dim, int spacedim>
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin ( )

Iterator to the first polytope.

Definition at line 1198 of file agglomeration_handler.h.

References master_cells_container, and n_agglomerations.

◆ begin() [2/2]

template<int dim, int spacedim>
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin ( ) const

Iterator to the first polytope.

Definition at line 1187 of file agglomeration_handler.h.

References master_cells_container, and n_agglomerations.

Referenced by polytope_iterators().

◆ cell_to_polytope_index()

template<int dim, int spacedim>
types::global_cell_index AgglomerationHandler< dim, spacedim >::cell_to_polytope_index ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell) const
inline

◆ connect_hierarchy()

template<int dim, int spacedim>
template<typename RtreeType>
void AgglomerationHandler< dim, spacedim >::connect_hierarchy ( const CellsAgglomerator< dim, RtreeType > & agglomerator)

◆ connect_to_tria_signals()

template<int dim, int spacedim = dim>
void AgglomerationHandler< dim, spacedim >::connect_to_tria_signals ( )
inlineprivate

Reinitialize the agglomeration data.

Definition at line 595 of file agglomeration_handler.h.

References connect_to_tria_signals(), initialize_agglomeration_data(), and tria_listener.

Referenced by connect_to_tria_signals().

◆ create_agglomeration_sparsity_pattern()

template<int dim, int spacedim>
template<typename SparsityPatternType, typename Number>
template void AgglomerationHandler< dim, spacedim >::create_agglomeration_sparsity_pattern ( SparsityPatternType & sparsity_pattern,
const AffineConstraints< Number > & constraints = AffineConstraints< Number >(),
const bool keep_constrained_dofs = true,
const types::subdomain_id subdomain_id = numbers::invalid_subdomain_id )

Given a Triangulation with some agglomerated cells, create the sparsity pattern corresponding to a Discontinuous Galerkin discretization where the agglomerated cells are seen as one unique cell, with only the DoFs associated to the master cell of the agglomeration.

Definition at line 913 of file agglomeration_handler.cc.

References AffineConstraints< typename number >::add_entries_local_to_global(), agglo_dh, DoFTools::extract_locally_relevant_dofs(), is_hp_collection, DoFTools::make_sparsity_pattern(), n_agglomerations, polytope_iterators(), IndexSet::size(), and IteratorState::valid.

◆ create_bounding_box()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::create_bounding_box ( const AgglomerationContainer & polytope)
private

Construct bounding boxes for an agglomeration described by a sequence of cells. This fills also the euler vector

Definition at line 478 of file agglomeration_handler.cc.

References n_agglomerations.

Referenced by define_agglomerate(), and define_agglomerate().

◆ define_agglomerate() [1/2]

template<int dim, int spacedim>
AgglomerationHandler< dim, spacedim >::agglomeration_iterator AgglomerationHandler< dim, spacedim >::define_agglomerate ( const AgglomerationContainer & cells)

Store internally that the given cells are agglomerated. The convenction we take is the following: -1: cell is a master cell

Note
cells are assumed to be adjacent one to each other, and no check about this is done.

Definition at line 46 of file agglomeration_handler.cc.

References agglo_dh, create_bounding_box(), hybrid_mesh, numbers::invalid_subdomain_id, master, master2polygon, master2slaves, master_cells_container, master_slave_relationships, master_slave_relationships_iterators, n_agglomerations, and slave.

Referenced by dealii::PolyUtils::partition_locally_owned_regions().

◆ define_agglomerate() [2/2]

template<int dim, int spacedim>
AgglomerationHandler< dim, spacedim >::agglomeration_iterator AgglomerationHandler< dim, spacedim >::define_agglomerate ( const AgglomerationContainer & cells,
const unsigned int fecollection_size )

Overload for hp::FECollection.

The parameter fecollection_size provides the number of finite elements in the collection, allowing Polydeal to insert an empty element for slave cells internally.

When fecollection_size equals 1, this function behaves identically to define_agglomerate(const AgglomerationContainer &cells).

Definition at line 108 of file agglomeration_handler.cc.

References agglo_dh, create_bounding_box(), hybrid_mesh, numbers::invalid_subdomain_id, master, master2polygon, master2slaves, master_cells_container, master_slave_relationships, master_slave_relationships_iterators, and n_agglomerations.

◆ define_agglomerate_with_check()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::define_agglomerate_with_check ( const AgglomerationContainer & cells)

Same as above, but checking that every vector of cells is connected. If not, each connected component is agglomerated by calling the define_agglomerate() function defined above.

Definition at line 174 of file agglomeration_handler.cc.

References Utils::compute_connected_components(), and Utils::create_graph_from_agglomerate().

◆ distribute_agglomerated_dofs() [1/2]

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::distribute_agglomerated_dofs ( const FiniteElement< dim > & fe_space)

◆ distribute_agglomerated_dofs() [2/2]

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::distribute_agglomerated_dofs ( const hp::FECollection< dim, spacedim > & fe_collection_in)

◆ end() [1/2]

template<int dim, int spacedim>
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end ( )

Iterator to one past the last polygonal element.

Definition at line 1220 of file agglomeration_handler.h.

References master_cells_container, and n_agglomerations.

◆ end() [2/2]

template<int dim, int spacedim>
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end ( ) const

Iterator to one past the last polygonal element.

Definition at line 1209 of file agglomeration_handler.h.

References master_cells_container, and n_agglomerations.

Referenced by polytope_iterators().

◆ exchange_interface_values()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::exchange_interface_values ( )

Definition at line 533 of file agglomeration_handler.cc.

References polytope_iterators(), and reinit().

Referenced by distribute_agglomerated_dofs().

◆ get_agglomerate()

template<int dim, int spacedim>
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > AgglomerationHandler< dim, spacedim >::get_agglomerate ( const typename Triangulation< dim, spacedim >::active_cell_iterator & master_cell) const
inline

TODO: remove this in favour of the accessor version.

Parameters
master_cell
Returns
std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator>

Definition at line 1024 of file agglomeration_handler.h.

References get_slaves_of_idx(), and is_master_cell().

Referenced by dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ get_agglomeration_mapping()

template<int dim, int spacedim>
const MappingBox< dim > & AgglomerationHandler< dim, spacedim >::get_agglomeration_mapping ( ) const
inline

Definition at line 970 of file agglomeration_handler.h.

References box_mapping.

◆ get_dof_handler()

template<int dim, int spacedim>
const DoFHandler< dim, spacedim > & AgglomerationHandler< dim, spacedim >::get_dof_handler ( ) const
inline

Return a constant reference to the DoFHandler underlying the agglomeration. It knows which cell have been agglomerated, and which FE spaces are present on each cell of the triangulation.

Definition at line 1038 of file agglomeration_handler.h.

References agglo_dh.

Referenced by dealii::PolyUtils::compute_quality_metrics().

◆ get_fe()

template<int dim, int spacedim>
const FiniteElement< dim, spacedim > & AgglomerationHandler< dim, spacedim >::get_fe ( ) const
inline

◆ get_fe_collection()

template<int dim, int spacedim>
const hp::FECollection< dim, spacedim > & AgglomerationHandler< dim, spacedim >::get_fe_collection ( ) const
inline

Return the finite element collection passed to distribute_agglomerated_dofs().

Definition at line 1262 of file agglomeration_handler.h.

References hp_fe_collection.

Referenced by dealii::PolyUtils::interpolate_to_fine_grid().

◆ get_interface()

template<int dim, int spacedim>
decltype(auto) AgglomerationHandler< dim, spacedim >::get_interface ( ) const
inline

Definition at line 1006 of file agglomeration_handler.h.

References get_interface(), and polytope_cache.

Referenced by get_interface().

◆ get_local_bboxes()

◆ get_mapping()

◆ get_master_idx_of_cell()

template<int dim, int spacedim>
types::global_cell_index AgglomerationHandler< dim, spacedim >::get_master_idx_of_cell ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell) const
inlineprivate

Definition at line 1137 of file agglomeration_handler.h.

References master_slave_relationships.

Referenced by are_cells_agglomerated().

◆ get_mesh_size()

template<int dim, int spacedim = dim>
double AgglomerationHandler< dim, spacedim >::get_mesh_size ( ) const

Return the mesh size of the polytopal mesh. It simply takes the maximum diameter over all the polytopes.

◆ get_relationships()

template<int dim, int spacedim>
const LinearAlgebra::distributed::Vector< float > & AgglomerationHandler< dim, spacedim >::get_relationships ( ) const
inline

Definition at line 1015 of file agglomeration_handler.h.

References master_slave_relationships.

◆ get_slaves_of_idx()

template<int dim, int spacedim>
const std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > & AgglomerationHandler< dim, spacedim >::get_slaves_of_idx ( types::global_cell_index idx) const
inline

Find (if any) the cells that have the given master index. Note that idx is as it can be equal to -1 (meaning that the cell is a master one).

Definition at line 1048 of file agglomeration_handler.h.

References master2slaves.

Referenced by dealii::PolyUtils::fill_interpolation_matrix(), get_agglomerate(), and dealii::PolyUtils::internal::interpolate_to_fine_grid().

◆ get_triangulation()

◆ initialize_agglomeration_data()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::initialize_agglomeration_data ( const std::unique_ptr< GridTools::Cache< dim, spacedim > > & cache_tria)
private

Initialize connectivity informations

Definition at line 292 of file agglomeration_handler.cc.

References agglo_dh.

Referenced by AgglomerationHandler(), and connect_to_tria_signals().

◆ initialize_fe_values() [1/2]

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::initialize_fe_values ( const hp::QCollection< dim > & cell_qcollection = hp::QCollection<dim>(QGauss<dim>(1)),
const UpdateFlags & flags = UpdateFlags::update_default,
const hp::QCollection< dim - 1 > & face_qcollection = hp::QCollection<dim - 1>(QGauss<dim - 1>(1)),
const UpdateFlags & face_flags = UpdateFlags::update_default )

Overload for hp::FECollection.

Definition at line 240 of file agglomeration_handler.cc.

References agglomeration_flags, and agglomeration_quad_collection.

◆ initialize_fe_values() [2/2]

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::initialize_fe_values ( const Quadrature< dim > & cell_quadrature = QGauss<dim>(1),
const UpdateFlags & flags = UpdateFlags::update_default,
const Quadrature< dim - 1 > & face_quadrature = QGauss<dim - 1>(1),
const UpdateFlags & face_flags = UpdateFlags::update_default )

Set the degree of the quadrature formula to be used and the proper flags for the FEValues object on the agglomerated cell.

Definition at line 212 of file agglomeration_handler.cc.

References agglomeration_face_flags, agglomeration_face_quad, agglomeration_flags, agglomeration_quad, dummy_fe, internal_agglomeration_face_flags, no_face_values, and no_values.

◆ initialize_hp_structure()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::initialize_hp_structure ( )
private

Assign a finite element index on each cell of a triangulation, depending if it is a master cell, a slave cell, or a standard deal.II cell. A user doesn't need to know the internals of this, the only thing that is relevant is that after the call to the present function, DoFs are distributed in a different way if a cell is a master, slave, or standard cell.

Definition at line 713 of file agglomeration_handler.cc.

References agglo_dh, and n_agglomerations.

Referenced by distribute_agglomerated_dofs().

◆ is_master_cell()

template<int dim, int spacedim>
template<typename CellIterator>
bool AgglomerationHandler< dim, spacedim >::is_master_cell ( const CellIterator & cell) const
inline

◆ is_slave_cell()

template<int dim, int spacedim>
template<typename CellIterator>
bool AgglomerationHandler< dim, spacedim >::is_slave_cell ( const CellIterator & cell) const
inlineprivate

Helper function to determine whether or not a cell is a slave cell, without any information about his parents.

Definition at line 1074 of file agglomeration_handler.h.

References master_slave_relationships.

Referenced by at_boundary(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ is_slave_cell_of()

template<int dim, int spacedim>
Triangulation< dim, spacedim >::active_cell_iterator & AgglomerationHandler< dim, spacedim >::is_slave_cell_of ( const typename Triangulation< dim, spacedim >::active_cell_iterator & cell)
inlineprivate

Helper function to determine whether or not a cell is a slave cell. Instead of returning a boolean, it gives the index of the master cell. If it's a master cell, then the it returns -1, by construction.

Definition at line 1127 of file agglomeration_handler.h.

References master_slave_relationships_iterators.

◆ last()

template<int dim, int spacedim>
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::last ( )

Iterator to the last polygonal element.

Definition at line 1231 of file agglomeration_handler.h.

References master_cells_container, and n_agglomerations.

◆ n_agglomerated_faces_per_cell()

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

Return the number of agglomerated faces for a generic deal.II cell.

Definition at line 271 of file agglomeration_handler.cc.

References are_cells_agglomerated().

◆ n_agglomerates()

template<int dim, int spacedim>
unsigned int AgglomerationHandler< dim, spacedim >::n_agglomerates ( ) const
inline

Returns the number of agglomerate cells in the grid.

Definition at line 1167 of file agglomeration_handler.h.

References n_agglomerations.

Referenced by dealii::PolyUtils::partition_locally_owned_regions().

◆ n_dofs()

template<int dim, int spacedim>
types::global_dof_index AgglomerationHandler< dim, spacedim >::n_dofs ( ) const
inlinenoexcept

Definition at line 1108 of file agglomeration_handler.h.

References agglo_dh.

Referenced by Utils::fill_injection_matrix().

◆ n_dofs_per_cell()

template<int dim, int spacedim>
unsigned int AgglomerationHandler< dim, spacedim >::n_dofs_per_cell ( ) const
inlinenoexcept

◆ polytope_boundary()

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

Return the collection of vertices describing the boundary of the polytope associated to the master cell cell. The return type is meant to describe a sequence of edges (in 2D) or faces (in 3D).

Definition at line 1117 of file agglomeration_handler.h.

References polygon_boundary.

◆ polytope_iterators()

◆ polytope_to_dh_iterator()

template<int dim, int spacedim>
const DoFHandler< dim, spacedim >::active_cell_iterator AgglomerationHandler< dim, spacedim >::polytope_to_dh_iterator ( const types::global_cell_index polytope_index) const
inline

Given the index of a polytopic element, return a DoFHandler iterator for which DoFs associated to that polytope can be queried.

Definition at line 1176 of file agglomeration_handler.h.

References agglo_dh, and master_cells_container.

Referenced by Utils::fill_injection_matrix().

◆ print_agglomeration()

template<int dim, int spacedim = dim>
template<class StreamType>
void AgglomerationHandler< dim, spacedim >::print_agglomeration ( StreamType & out)
inline

Display the indices of the vector identifying which cell is agglomerated with which master.

Definition at line 396 of file agglomeration_handler.h.

References master_slave_relationships.

◆ reinit() [1/2]

◆ reinit() [2/2]

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & AgglomerationHandler< dim, spacedim >::reinit ( const AgglomerationIterator< dim, spacedim > & polytope,
const unsigned int face_index ) const

For a given polytope and face index, initialize shape functions, normals and quadratures rules to integrate there.

Definition at line 787 of file agglomeration_handler.cc.

References agglo_dh, agglomerated_isv_bdary, AgglomerationIterator, AgglomerationAccessor< dim, spacedim >::as_dof_handler_iterator(), is_master_cell(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().

◆ reinit_interface()

template<int dim, int spacedim>
std::pair< const FEValuesBase< dim, spacedim > &, const FEValuesBase< dim, spacedim > & > AgglomerationHandler< dim, spacedim >::reinit_interface ( const AgglomerationIterator< dim, spacedim > & polytope_in,
const AgglomerationIterator< dim, spacedim > & neigh_polytope,
const unsigned int local_in,
const unsigned int local_outside ) const

◆ reinit_master()

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & AgglomerationHandler< dim, spacedim >::reinit_master ( const typename DoFHandler< dim, spacedim >::active_cell_iterator & cell,
const unsigned int face_number,
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > & agglo_isv_ptr ) const
private

Helper function to call reinit on a master cell.

Definition at line 773 of file agglomeration_handler.cc.

References dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().

◆ setup_connectivity_of_agglomeration()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::setup_connectivity_of_agglomeration ( )
private

Initialize all the necessary connectivity information for an agglomeration. #TODO: loop over polytopes, avoid using explicitly master cells.

Definition at line 497 of file agglomeration_handler.cc.

References agglo_dh, Utilities::MPI::job_supports_mpi(), master2polygon, master_cells_container, number_of_agglomerated_faces, dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity(), and Utilities::MPI::some_to_some().

Referenced by distribute_agglomerated_dofs().

◆ setup_ghost_polytopes()

template<int dim, int spacedim>
void AgglomerationHandler< dim, spacedim >::setup_ghost_polytopes ( )

This function stores the information needed to identify which polytopes are ghosted w.r.t the local partition. The issue this function addresses is due to the fact that the layer of ghost cells is made by just one layer of deal.II cells. Therefore, the neighboring polytopes will always be made by some ghost cells and artificial ones. This implies that we need to communicate the missing information from the neighboring rank.

Definition at line 1028 of file agglomeration_handler.cc.

References agglo_dh, n_dofs_per_cell(), polytope_iterators(), and Utilities::MPI::some_to_some().

Referenced by distribute_agglomerated_dofs().

◆ update_agglomerate()

template<int dim, int spacedim = dim>
void AgglomerationHandler< dim, spacedim >::update_agglomerate ( AgglomerationContainer & polytope,
const typename Triangulation< dim, spacedim >::active_cell_iterator & master_cell )
private

References update_agglomerate().

Referenced by update_agglomerate().

◆ used_fe_collection()

template<int dim, int spacedim>
bool AgglomerationHandler< dim, spacedim >::used_fe_collection ( ) const
inline

Return whether a hp::FECollection is being used.

Definition at line 1269 of file agglomeration_handler.h.

References is_hp_collection.

Referenced by dealii::PolyUtils::interpolate_to_fine_grid().

◆ AgglomerationAccessor

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

Definition at line 241 of file agglomeration_handler.h.

References AgglomerationAccessor, and numbers::invalid_subdomain_id.

Referenced by AgglomerationAccessor.

◆ AgglomerationIterator

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

Definition at line 238 of file agglomeration_handler.h.

References AgglomerationIterator.

Referenced by AgglomerationIterator, reinit(), reinit(), and reinit_interface().

◆ internal::AgglomerationHandlerImplementation< dim, spacedim >

template<int dim, int spacedim = dim>
friend class internal::AgglomerationHandlerImplementation< dim, spacedim >
friend

Definition at line 905 of file agglomeration_handler.h.

References master_cells_container.

Member Data Documentation

◆ agglo_dh

◆ agglomerated_isv

template<int dim, int spacedim = dim>
std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim> > AgglomerationHandler< dim, spacedim >::agglomerated_isv
mutableprivate

Definition at line 845 of file agglomeration_handler.h.

Referenced by reinit_interface().

◆ agglomerated_isv_bdary

template<int dim, int spacedim = dim>
std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim> > AgglomerationHandler< dim, spacedim >::agglomerated_isv_bdary
mutableprivate

Definition at line 851 of file agglomeration_handler.h.

Referenced by reinit().

◆ agglomerated_isv_neigh

template<int dim, int spacedim = dim>
std::unique_ptr<NonMatching::FEImmersedSurfaceValues<spacedim> > AgglomerationHandler< dim, spacedim >::agglomerated_isv_neigh
mutableprivate

Definition at line 848 of file agglomeration_handler.h.

Referenced by reinit_interface().

◆ agglomerated_scratch

template<int dim, int spacedim = dim>
std::unique_ptr<ScratchData> AgglomerationHandler< dim, spacedim >::agglomerated_scratch
mutableprivate

Fill this up in reinit(cell), for agglomerated cells, using the custom quadrature, and return the result of scratch.reinit(cell);

Definition at line 841 of file agglomeration_handler.h.

Referenced by reinit().

◆ agglomeration_face_flags

template<int dim, int spacedim = dim>
UpdateFlags AgglomerationHandler< dim, spacedim >::agglomeration_face_flags
private

◆ agglomeration_face_quad

template<int dim, int spacedim = dim>
Quadrature<dim - 1> AgglomerationHandler< dim, spacedim >::agglomeration_face_quad
private

◆ agglomeration_face_quad_collection

template<int dim, int spacedim = dim>
hp::QCollection<dim - 1> AgglomerationHandler< dim, spacedim >::agglomeration_face_quad_collection
private

◆ agglomeration_flags

template<int dim, int spacedim = dim>
UpdateFlags AgglomerationHandler< dim, spacedim >::agglomeration_flags
private

Definition at line 855 of file agglomeration_handler.h.

Referenced by initialize_fe_values(), initialize_fe_values(), and reinit().

◆ agglomeration_quad

template<int dim, int spacedim = dim>
Quadrature<dim> AgglomerationHandler< dim, spacedim >::agglomeration_quad
private

Definition at line 867 of file agglomeration_handler.h.

Referenced by initialize_fe_values().

◆ agglomeration_quad_collection

template<int dim, int spacedim = dim>
hp::QCollection<dim> AgglomerationHandler< dim, spacedim >::agglomeration_quad_collection
private

Definition at line 930 of file agglomeration_handler.h.

Referenced by initialize_fe_values().

◆ box_mapping

◆ dummy_fe

template<int dim, int spacedim = dim>
FE_Nothing<dim, spacedim> AgglomerationHandler< dim, spacedim >::dummy_fe
private

Dummy FE_Nothing

Definition at line 889 of file agglomeration_handler.h.

Referenced by initialize_fe_values().

◆ dummy_fe_collection

template<int dim, int spacedim = dim>
hp::FECollection<dim, spacedim> AgglomerationHandler< dim, spacedim >::dummy_fe_collection
private

Definition at line 936 of file agglomeration_handler.h.

◆ euler_vector

template<int dim, int spacedim = dim>
LinearAlgebra::distributed::Vector<double> AgglomerationHandler< dim, spacedim >::euler_vector
private

Vector of BoundingBoxes s.t. bboxes[idx] equals BBOx associated to the agglomeration with master cell indexed by ìdx`. Othwerwise default BBox is empty

*/ std::vector<BoundingBox<spacedim>> bboxes;

//////////////////////////////////////////////////////

n_faces mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>> local_n_faces;

mutable std::map<types::subdomain_id, std::map<CellId, unsigned int>> recv_n_faces;

CellId (including slaves) mutable std::map<types::subdomain_id, std::map<CellId, CellId>> local_cell_ids_neigh_cell;

mutable std::map<types::subdomain_id, std::map<CellId, CellId>> recv_cell_ids_neigh_cell;

send to neighborign rank the information that

  • current polytope id
  • face f has the following neighboring id. mutable std::map<types::subdomain_id, std::map<CellId, std::map<unsigned int, CellId>>> local_ghosted_master_id;

    mutable std::map<types::subdomain_id, std::map<CellId, std::map<unsigned int, CellId>>> recv_ghosted_master_id;

CellIds from neighboring rank mutable std::map<types::subdomain_id, std::map<CellId, std::map<unsigned int, bool>>> local_bdary_info;

mutable std::map<types::subdomain_id, std::map<CellId, std::map<unsigned int, bool>>> recv_bdary_info;

Exchange neighboring bounding boxes mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>> local_ghosted_bbox;

mutable std::map<types::subdomain_id, std::map<CellId, BoundingBox<dim>>> recv_ghosted_bbox;

Exchange DoF indices with ghosted polytopes mutable std::map<types::subdomain_id, std::map<CellId, std::vector<types::global_dof_index>>> local_ghost_dofs;

mutable std::map<types::subdomain_id, std::map<CellId, std::vector<types::global_dof_index>>> recv_ghost_dofs;

Exchange qpoints mutable std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim>>>> local_qpoints;

Exchange jxws mutable std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<double>>> local_jxws;

Exchange normals mutable std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim>>>> local_normals;

Exchange values mutable std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double>>>> local_values;

mutable std::map<types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<Tensor<1, spacedim>>>>> local_gradients;

//////////////////////////////////////////////////////

ObserverPointer<const Triangulation<dim, spacedim>> tria;

ObserverPointer<const Mapping<dim, spacedim>> mapping;

std::unique_ptr<GridTools::Cache<dim, spacedim>> cached_tria;

const MPI_Comm communicator;

The FiniteElement space we have on each cell. Currently supported types are FE_DGQ and FE_DGP elements. std::unique_ptr<FiniteElement<dim>> fe;

hp::FECollection<dim, spacedim> fe_collection;

/** Eulerian vector describing the new cells obtained by the bounding boxes

Definition at line 827 of file agglomeration_handler.h.

◆ hp_fe_collection

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FECollection<dim, spacedim> > AgglomerationHandler< dim, spacedim >::hp_fe_collection
private

Definition at line 926 of file agglomeration_handler.h.

Referenced by distribute_agglomerated_dofs(), and get_fe_collection().

◆ hp_no_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FEFaceValues<dim, spacedim> > AgglomerationHandler< dim, spacedim >::hp_no_face_values
private

◆ hp_no_values

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FEValues<dim, spacedim> > AgglomerationHandler< dim, spacedim >::hp_no_values
private

Definition at line 943 of file agglomeration_handler.h.

◆ hybrid_mesh

template<int dim, int spacedim = dim>
bool AgglomerationHandler< dim, spacedim >::hybrid_mesh
private

Bool that keeps track whether the mesh is composed also by standard deal.II cells as (trivial) polytopes.

Definition at line 915 of file agglomeration_handler.h.

Referenced by AgglomerationHandler(), define_agglomerate(), define_agglomerate(), distribute_agglomerated_dofs(), and distribute_agglomerated_dofs().

◆ internal_agglomeration_face_flags

template<int dim, int spacedim = dim>
const UpdateFlags AgglomerationHandler< dim, spacedim >::internal_agglomeration_face_flags
private
Initial value:

Definition at line 863 of file agglomeration_handler.h.

Referenced by initialize_fe_values().

◆ internal_agglomeration_flags

template<int dim, int spacedim = dim>
const UpdateFlags AgglomerationHandler< dim, spacedim >::internal_agglomeration_flags
private

◆ is_hp_collection

◆ mapping_collection

template<int dim, int spacedim = dim>
hp::MappingCollection<dim> AgglomerationHandler< dim, spacedim >::mapping_collection
private

Definition at line 934 of file agglomeration_handler.h.

◆ master2polygon

◆ master2slaves

template<int dim, int spacedim = dim>
std::unordered_map< types::global_cell_index, std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator> > AgglomerationHandler< dim, spacedim >::master2slaves
private

◆ master_cells_container

template<int dim, int spacedim = dim>
std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator> AgglomerationHandler< dim, spacedim >::master_cells_container
private

◆ master_disconnected

template<int dim, int spacedim = dim>
std::vector<typename Triangulation<dim>::active_cell_iterator> AgglomerationHandler< dim, spacedim >::master_disconnected
private

Definition at line 882 of file agglomeration_handler.h.

◆ master_slave_relationships

template<int dim, int spacedim = dim>
LinearAlgebra::distributed::Vector<float> AgglomerationHandler< dim, spacedim >::master_slave_relationships
private

Vector of indices such that v[cell->active_cell_index()] returns { -1 if cell is a master cell { cell_master->active_cell_index(), i.e. the index of the master cell if cell is a slave cell.

Definition at line 688 of file agglomeration_handler.h.

Referenced by define_agglomerate(), define_agglomerate(), get_master_idx_of_cell(), get_relationships(), is_master_cell(), is_slave_cell(), print_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ master_slave_relationships_iterators

template<int dim, int spacedim = dim>
std::map<types::global_cell_index, typename Triangulation<dim, spacedim>::active_cell_iterator> AgglomerationHandler< dim, spacedim >::master_slave_relationships_iterators
private

Same as the one above, but storing cell iterators rather than indices.

Definition at line 696 of file agglomeration_handler.h.

Referenced by define_agglomerate(), define_agglomerate(), is_slave_cell_of(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ n_agglomerations

template<int dim, int spacedim = dim>
unsigned int AgglomerationHandler< dim, spacedim >::n_agglomerations
private

◆ no_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEFaceValues<dim, spacedim> > AgglomerationHandler< dim, spacedim >::no_face_values
private

◆ no_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEValues<dim, spacedim> > AgglomerationHandler< dim, spacedim >::no_values
private

Dummy FEValues, needed for cell quadratures.

Definition at line 894 of file agglomeration_handler.h.

Referenced by agglomerated_quadrature(), and initialize_fe_values().

◆ number_of_agglomerated_faces

template<int dim, int spacedim = dim>
std::vector<types::global_cell_index> AgglomerationHandler< dim, spacedim >::number_of_agglomerated_faces
mutableprivate

◆ output_dh

template<int dim, int spacedim = dim>
DoFHandler<dim, spacedim> AgglomerationHandler< dim, spacedim >::output_dh

◆ parent_child_info

template<int dim, int spacedim = dim>
std::map<std::pair<types::global_cell_index, types::global_cell_index>, std::vector<types::global_cell_index> > AgglomerationHandler< dim, spacedim >::parent_child_info
private

Definition at line 919 of file agglomeration_handler.h.

Referenced by connect_hierarchy().

◆ polygon_boundary

template<int dim, int spacedim = dim>
std::map< const typename Triangulation<dim, spacedim>::active_cell_iterator, std::vector<typename Triangulation<dim>::active_face_iterator> > AgglomerationHandler< dim, spacedim >::polygon_boundary
mutableprivate

Associate a master cell (hence, a given polytope) to its boundary faces. The boundary is described through a vector of face iterators.

Definition at line 710 of file agglomeration_handler.h.

Referenced by polytope_boundary(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().

◆ polytope_cache

◆ present_extraction_level

template<int dim, int spacedim = dim>
unsigned int AgglomerationHandler< dim, spacedim >::present_extraction_level
private

Definition at line 921 of file agglomeration_handler.h.

Referenced by connect_hierarchy().

◆ recv_gradients

template<int dim, int spacedim = dim>
std::map<types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<Tensor<1, spacedim> > > > > AgglomerationHandler< dim, spacedim >::recv_gradients
mutable

Definition at line 548 of file agglomeration_handler.h.

Referenced by dealii::PolyUtils::assemble_dg_matrix().

◆ recv_jxws

template<int dim, int spacedim = dim>
std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<double> > > AgglomerationHandler< dim, spacedim >::recv_jxws
mutable

◆ recv_normals

template<int dim, int spacedim = dim>
std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<Tensor<1, spacedim> > > > AgglomerationHandler< dim, spacedim >::recv_normals
mutable

Definition at line 538 of file agglomeration_handler.h.

Referenced by reinit_interface().

◆ recv_qpoints

template<int dim, int spacedim = dim>
std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<Point<spacedim> > > > AgglomerationHandler< dim, spacedim >::recv_qpoints
mutable

Definition at line 528 of file agglomeration_handler.h.

Referenced by reinit_interface().

◆ recv_values

template<int dim, int spacedim = dim>
std::map< types::subdomain_id, std::map<std::pair<CellId, unsigned int>, std::vector<std::vector<double> > > > AgglomerationHandler< dim, spacedim >::recv_values
mutable

Definition at line 543 of file agglomeration_handler.h.

Referenced by dealii::PolyUtils::assemble_dg_matrix().

◆ standard_scratch

template<int dim, int spacedim = dim>
std::unique_ptr<ScratchData> AgglomerationHandler< dim, spacedim >::standard_scratch
mutableprivate

Use this in reinit(cell) for (non-agglomerated, standard) cells, and return the result of scratch.reinit(cell) for cells

Definition at line 834 of file agglomeration_handler.h.

Referenced by distribute_agglomerated_dofs().

◆ tria_listener

template<int dim, int spacedim = dim>
boost::signals2::connection AgglomerationHandler< dim, spacedim >::tria_listener
private

Definition at line 853 of file agglomeration_handler.h.

Referenced by connect_to_tria_signals(), and ~AgglomerationHandler().


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