#include <agglomeration_handler.h>
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_iterator > | polytope_iterators () const |
void | distribute_agglomerated_dofs (const FiniteElement< dim > &fe_space) |
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) |
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) |
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) |
![]() | |
EnableObserverPointer () | |
EnableObserverPointer (const EnableObserverPointer &) | |
EnableObserverPointer (EnableObserverPointer &&) noexcept | |
virtual | ~EnableObserverPointer () |
EnableObserverPointer & | operator= (const EnableObserverPointer &) |
EnableObserverPointer & | operator= (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 () |
Friends | |
template<int, int> | |
class | AgglomerationIterator |
template<int, int> | |
class | AgglomerationAccessor |
class | internal::AgglomerationHandlerImplementation< dim, spacedim > |
Additional Inherited Members | |
![]() | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
#TODO: Documentation.
Definition at line 172 of file agglomeration_handler.h.
using AgglomerationHandler< dim, spacedim >::agglomeration_iterator = AgglomerationIterator<dim, spacedim> |
Definition at line 175 of file agglomeration_handler.h.
using AgglomerationHandler< dim, spacedim >::AgglomerationContainer |
Definition at line 177 of file agglomeration_handler.h.
|
private |
Definition at line 652 of file agglomeration_handler.h.
enum AgglomerationHandler::CellAgglomerationType |
Enumerator | |
---|---|
master | |
slave |
Definition at line 181 of file agglomeration_handler.h.
|
explicit |
Definition at line 21 of file agglomeration_handler.cc.
References cached_tria, communicator, get_mapping(), get_triangulation(), hybrid_mesh, initialize_agglomeration_data(), and n_agglomerations.
|
default |
|
inline |
Definition at line 194 of file agglomeration_handler.h.
References tria_listener.
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 438 of file agglomeration_handler.cc.
References bboxes, is_master_cell(), master2polygon, and no_values.
Referenced by reinit().
|
inlineprivate |
Returns true if the two given cells are agglomerated together.
Definition at line 1082 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().
|
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 1015 of file agglomeration_handler.h.
References is_slave_cell(), master2polygon, and polytope_cache.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin | ( | ) |
Iterator to the first polytope.
Definition at line 1129 of file agglomeration_handler.h.
References master_cells_container, and n_agglomerations.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin | ( | ) | const |
Iterator to the first polytope.
Definition at line 1118 of file agglomeration_handler.h.
References master_cells_container, and n_agglomerations.
Referenced by polytope_iterators().
|
inline |
Definition at line 927 of file agglomeration_handler.h.
References master2polygon.
Referenced by dealii::PolyUtils::fill_interpolation_matrix(), and dealii::PolyUtils::internal::interpolate_to_fine_grid().
void AgglomerationHandler< dim, spacedim >::connect_hierarchy | ( | const CellsAgglomerator< dim, RtreeType > & | agglomerator | ) |
Definition at line 1184 of file agglomeration_handler.h.
References dealii::CellsAgglomerator< dim, RtreeType >::extraction_level, parent_child_info, dealii::CellsAgglomerator< dim, RtreeType >::parent_node_to_children_nodes, and present_extraction_level.
|
inlineprivate |
Reinitialize the agglomeration data.
Definition at line 549 of file agglomeration_handler.h.
References cached_tria, connect_to_tria_signals(), initialize_agglomeration_data(), tria, and tria_listener.
Referenced by connect_to_tria_signals(), and initialize_agglomeration_data().
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 670 of file agglomeration_handler.cc.
References AffineConstraints< typename number >::add_entries_local_to_global(), agglo_dh, communicator, DoFTools::extract_locally_relevant_dofs(), AgglomerationAccessor< dim, spacedim >::get_dof_indices(), DoFTools::make_sparsity_pattern(), n_agglomerations, AgglomerationAccessor< dim, spacedim >::neighbor(), polytope_iterators(), IndexSet::size(), AgglomerationIterator< dim, spacedim >::state(), and IteratorState::valid.
|
private |
Construct bounding boxes for an agglomeration described by a sequence of cells. This fills also the euler vector
Definition at line 292 of file agglomeration_handler.cc.
Referenced by define_agglomerate().
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
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, slave, and tria.
Referenced by define_agglomerate_with_check(), and dealii::PolyUtils::partition_locally_owned_regions().
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 110 of file agglomeration_handler.cc.
References Utils::compute_connected_components(), Utils::create_graph_from_agglomerate(), and define_agglomerate().
void AgglomerationHandler< dim, spacedim >::distribute_agglomerated_dofs | ( | const FiniteElement< dim > & | fe_space | ) |
Distribute degrees of freedom on a grid where some cells have been agglomerated.
Definition at line 235 of file agglomeration_handler.cc.
References FiniteElement< int dim, int spacedim >::degree, and fe.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end | ( | ) |
Iterator to one past the last polygonal element.
Definition at line 1151 of file agglomeration_handler.h.
References master_cells_container, and n_agglomerations.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end | ( | ) | const |
Iterator to one past the last polygonal element.
Definition at line 1140 of file agglomeration_handler.h.
References master_cells_container, and n_agglomerations.
Referenced by polytope_iterators().
void AgglomerationHandler< dim, spacedim >::exchange_interface_values | ( | ) |
Definition at line 347 of file agglomeration_handler.cc.
References fe, local_jxws, local_normals, local_qpoints, polytope_iterators(), and reinit().
|
inline |
TODO: remove this in favour of the accessor version.
master_cell |
Definition at line 955 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().
|
inline |
Definition at line 901 of file agglomeration_handler.h.
References box_mapping.
|
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 969 of file agglomeration_handler.h.
References agglo_dh.
Referenced by dealii::PolyUtils::compute_quality_metrics().
|
inline |
Definition at line 883 of file agglomeration_handler.h.
References fe.
Referenced by Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
inline |
Definition at line 937 of file agglomeration_handler.h.
References get_interface(), and polytope_cache.
Referenced by get_interface().
|
inline |
Definition at line 918 of file agglomeration_handler.h.
References bboxes.
Referenced by dealii::PolyUtils::compute_quality_metrics(), Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
inline |
Definition at line 892 of file agglomeration_handler.h.
References mapping.
Referenced by AgglomerationHandler(), Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
inlineprivate |
Definition at line 1068 of file agglomeration_handler.h.
References master_slave_relationships.
Referenced by are_cells_agglomerated().
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.
|
inline |
Definition at line 946 of file agglomeration_handler.h.
References master_slave_relationships.
|
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 979 of file agglomeration_handler.h.
References master2slaves.
Referenced by dealii::PolyUtils::fill_interpolation_matrix(), get_agglomerate(), and dealii::PolyUtils::internal::interpolate_to_fine_grid().
|
inline |
Definition at line 910 of file agglomeration_handler.h.
References tria.
Referenced by AgglomerationHandler(), dealii::PolyUtils::compute_global_error(), dealii::PolyUtils::compute_quality_metrics(), Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
private |
Initialize connectivity informations
Definition at line 199 of file agglomeration_handler.cc.
References agglo_dh, bboxes, cached_tria, communicator, connect_to_tria_signals(), mapping, master_slave_relationships, n_agglomerations, polytope_cache, and tria.
Referenced by AgglomerationHandler(), and connect_to_tria_signals().
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 148 of file agglomeration_handler.cc.
References agglomeration_face_flags, agglomeration_face_quad, agglomeration_flags, agglomeration_quad, dummy_fe, internal_agglomeration_face_flags, mapping, no_face_values, and no_values.
|
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 482 of file agglomeration_handler.cc.
References agglo_dh.
|
inline |
Helper function to determine whether or not a cell is a master or a slave
Definition at line 990 of file agglomeration_handler.h.
References master_slave_relationships.
Referenced by agglomerated_quadrature(), dealii::PolyUtils::fill_interpolation_matrix(), get_agglomerate(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), reinit(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
inlineprivate |
Helper function to determine whether or not a cell is a slave cell, without any information about his parents.
Definition at line 1005 of file agglomeration_handler.h.
References master_slave_relationships.
Referenced by at_boundary(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
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 1058 of file agglomeration_handler.h.
References master_slave_relationships_iterators.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::last | ( | ) |
Iterator to the last polygonal element.
Definition at line 1162 of file agglomeration_handler.h.
References master_cells_container, and n_agglomerations.
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 178 of file agglomeration_handler.cc.
References are_cells_agglomerated().
|
inline |
Returns the number of agglomerate cells in the grid.
Definition at line 1098 of file agglomeration_handler.h.
References n_agglomerations.
Referenced by dealii::PolyUtils::partition_locally_owned_regions().
|
inlinenoexcept |
Definition at line 1039 of file agglomeration_handler.h.
References agglo_dh.
Referenced by Utils::fill_injection_matrix().
|
inlinenoexcept |
Definition at line 1030 of file agglomeration_handler.h.
References fe.
Referenced by dealii::PolyUtils::compute_global_error(), dealii::PolyUtils::interpolate_to_fine_grid(), and setup_ghost_polytopes().
|
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 1048 of file agglomeration_handler.h.
References polygon_boundary.
IteratorRange< typename AgglomerationHandler< dim, spacedim >::agglomeration_iterator > AgglomerationHandler< dim, spacedim >::polytope_iterators | ( | ) | const |
Returns an IteratorRange that makes up all the polygonal elements in the mesh.
Definition at line 1174 of file agglomeration_handler.h.
References begin(), and end().
Referenced by dealii::PolyUtils::assemble_dg_matrix(), dealii::PolyUtils::compute_global_error(), dealii::PolyUtils::compute_quality_metrics(), create_agglomeration_sparsity_pattern(), exchange_interface_values(), dealii::PolyUtils::export_polygon_to_csv_file(), Utils::fill_injection_matrix(), dealii::PolyUtils::interpolate_to_fine_grid(), and setup_ghost_polytopes().
|
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 1107 of file agglomeration_handler.h.
References agglo_dh, and master_cells_container.
Referenced by Utils::fill_injection_matrix().
|
inline |
Display the indices of the vector identifying which cell is agglomerated with which master.
Definition at line 363 of file agglomeration_handler.h.
References master_slave_relationships, and tria.
const FEValues< dim, spacedim > & AgglomerationHandler< dim, spacedim >::reinit | ( | const AgglomerationIterator< dim, spacedim > & | polytope | ) | const |
Construct a finite element space on the agglomeration.
Definition at line 500 of file agglomeration_handler.cc.
References agglo_dh, agglomerated_quadrature(), agglomerated_scratch, agglomeration_flags, AgglomerationIterator, AgglomerationAccessor< dim, spacedim >::as_dof_handler_iterator(), box_mapping, fe_collection, and AgglomerationAccessor< dim, spacedim >::get_agglomerate().
Referenced by dealii::PolyUtils::assemble_dg_matrix(), dealii::PolyUtils::compute_global_error(), dealii::PolyUtils::compute_quality_metrics(), and exchange_interface_values().
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 544 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().
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 |
Return a pair of FEValuesBase object reinited from the two sides of the agglomeration.
Definition at line 565 of file agglomeration_handler.cc.
References AgglomerationIterator.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
private |
Helper function to call reinit on a master cell.
Definition at line 530 of file agglomeration_handler.cc.
|
private |
Initialize all the necessary connectivity information for an agglomeration. #TODO: loop over polytopes, avoid using explicitly master cells.
Definition at line 311 of file agglomeration_handler.cc.
References agglo_dh, communicator, Utilities::MPI::job_supports_mpi(), local_bdary_info, local_ghosted_master_id, local_n_faces, master2polygon, master_cells_container, number_of_agglomerated_faces, recv_bdary_info, recv_ghosted_master_id, recv_n_faces, dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity(), and Utilities::MPI::some_to_some().
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 741 of file agglomeration_handler.cc.
References agglo_dh, bboxes, communicator, fe, local_cell_ids_neigh_cell, local_ghost_dofs, local_ghosted_bbox, n_dofs_per_cell(), polytope_iterators(), recv_cell_ids_neigh_cell, recv_ghost_dofs, recv_ghosted_bbox, Utilities::MPI::some_to_some(), and tria.
|
private |
References update_agglomerate().
Referenced by update_agglomerate().
Definition at line 241 of file agglomeration_handler.h.
References AgglomerationAccessor, and numbers::invalid_subdomain_id.
Referenced by AgglomerationAccessor.
Definition at line 238 of file agglomeration_handler.h.
References AgglomerationIterator.
Referenced by AgglomerationIterator, reinit(), reinit(), and reinit_interface().
|
friend |
Definition at line 859 of file agglomeration_handler.h.
References master_cells_container.
DoFHandler<dim, spacedim> AgglomerationHandler< dim, spacedim >::agglo_dh |
DoFHandler for the agglomerated space
Definition at line 468 of file agglomeration_handler.h.
Referenced by create_agglomeration_sparsity_pattern(), define_agglomerate(), Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), get_dof_handler(), initialize_agglomeration_data(), initialize_hp_structure(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), n_dofs(), polytope_to_dh_iterator(), reinit(), reinit(), setup_connectivity_of_agglomeration(), and setup_ghost_polytopes().
|
mutableprivate |
Definition at line 799 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 805 of file agglomeration_handler.h.
Referenced by reinit().
|
mutableprivate |
Definition at line 802 of file agglomeration_handler.h.
|
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 795 of file agglomeration_handler.h.
Referenced by reinit().
|
private |
Definition at line 815 of file agglomeration_handler.h.
Referenced by initialize_fe_values(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
private |
Definition at line 823 of file agglomeration_handler.h.
Referenced by initialize_fe_values(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
private |
Definition at line 809 of file agglomeration_handler.h.
Referenced by initialize_fe_values(), and reinit().
|
private |
Definition at line 821 of file agglomeration_handler.h.
Referenced by initialize_fe_values().
|
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
Definition at line 673 of file agglomeration_handler.h.
Referenced by agglomerated_quadrature(), get_local_bboxes(), initialize_agglomeration_data(), dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master(), and setup_ghost_polytopes().
std::unique_ptr<MappingBox<dim> > AgglomerationHandler< dim, spacedim >::box_mapping |
Definition at line 475 of file agglomeration_handler.h.
Referenced by get_agglomeration_mapping(), reinit(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
private |
Definition at line 768 of file agglomeration_handler.h.
Referenced by AgglomerationHandler(), connect_to_tria_signals(), and initialize_agglomeration_data().
|
private |
Definition at line 770 of file agglomeration_handler.h.
Referenced by AgglomerationHandler(), create_agglomeration_sparsity_pattern(), initialize_agglomeration_data(), setup_connectivity_of_agglomeration(), and setup_ghost_polytopes().
|
private |
Dummy FE_Nothing
Definition at line 843 of file agglomeration_handler.h.
Referenced by initialize_fe_values().
|
private |
Eulerian vector describing the new cells obtained by the bounding boxes
Definition at line 781 of file agglomeration_handler.h.
|
private |
Definition at line 774 of file agglomeration_handler.h.
Referenced by distribute_agglomerated_dofs(), exchange_interface_values(), get_fe(), n_dofs_per_cell(), dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master(), and setup_ghost_polytopes().
|
private |
Definition at line 776 of file agglomeration_handler.h.
Referenced by reinit().
|
private |
Bool that keeps track whether the mesh is composed also by standard deal.II cells as (trivial) polytopes.
Definition at line 869 of file agglomeration_handler.h.
Referenced by AgglomerationHandler(), and define_agglomerate().
|
private |
Definition at line 817 of file agglomeration_handler.h.
Referenced by initialize_fe_values().
|
private |
Definition at line 811 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 709 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 688 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes().
|
mutableprivate |
Definition at line 725 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes().
|
mutableprivate |
Definition at line 717 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes().
|
mutableprivate |
Definition at line 700 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 758 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 741 of file agglomeration_handler.h.
Referenced by exchange_interface_values().
|
mutableprivate |
Definition at line 680 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 747 of file agglomeration_handler.h.
Referenced by exchange_interface_values().
|
mutableprivate |
Definition at line 735 of file agglomeration_handler.h.
Referenced by exchange_interface_values().
|
mutableprivate |
Definition at line 753 of file agglomeration_handler.h.
|
private |
Definition at line 766 of file agglomeration_handler.h.
Referenced by get_mapping(), initialize_agglomeration_data(), and initialize_fe_values().
|
private |
Definition at line 832 of file agglomeration_handler.h.
Referenced by agglomerated_quadrature(), at_boundary(), cell_to_polytope_index(), define_agglomerate(), dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master(), setup_connectivity_of_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Definition at line 829 of file agglomeration_handler.h.
Referenced by define_agglomerate(), and get_slaves_of_idx().
|
private |
A contiguous container for all of the master cells.
Definition at line 859 of file agglomeration_handler.h.
Referenced by begin(), begin(), define_agglomerate(), end(), end(), internal::AgglomerationHandlerImplementation< dim, spacedim >, last(), polytope_to_dh_iterator(), and setup_connectivity_of_agglomeration().
|
private |
Definition at line 836 of file agglomeration_handler.h.
|
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 642 of file agglomeration_handler.h.
Referenced by define_agglomerate(), get_master_idx_of_cell(), get_relationships(), initialize_agglomeration_data(), is_master_cell(), is_slave_cell(), print_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Same as the one above, but storing cell iterators rather than indices.
Definition at line 650 of file agglomeration_handler.h.
Referenced by define_agglomerate(), is_slave_cell_of(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Record the number of agglomerations on the grid.
Definition at line 633 of file agglomeration_handler.h.
Referenced by AgglomerationHandler(), begin(), begin(), create_agglomeration_sparsity_pattern(), define_agglomerate(), end(), end(), initialize_agglomeration_data(), last(), and n_agglomerates().
|
private |
Dummy FEFaceValues, needed for face quadratures.
Definition at line 853 of file agglomeration_handler.h.
Referenced by initialize_fe_values(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master().
|
private |
Dummy FEValues, needed for cell quadratures.
Definition at line 848 of file agglomeration_handler.h.
Referenced by agglomerated_quadrature(), and initialize_fe_values().
|
mutableprivate |
Definition at line 654 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
DoFHandler<dim, spacedim> AgglomerationHandler< dim, spacedim >::output_dh |
DoFHandler for the finest space: classical deal.II space
Definition at line 473 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
private |
Definition at line 873 of file agglomeration_handler.h.
Referenced by connect_hierarchy().
|
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 664 of file agglomeration_handler.h.
Referenced by polytope_boundary(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Definition at line 863 of file agglomeration_handler.h.
Referenced by at_boundary(), get_interface(), initialize_agglomeration_data(), dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::reinit_master(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Definition at line 875 of file agglomeration_handler.h.
Referenced by connect_hierarchy().
|
mutableprivate |
Definition at line 713 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration().
|
mutableprivate |
Definition at line 691 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 729 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes().
|
mutableprivate |
Definition at line 720 of file agglomeration_handler.h.
Referenced by setup_ghost_polytopes().
|
mutableprivate |
Definition at line 704 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration().
|
mutable |
Definition at line 515 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
mutable |
Definition at line 500 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
mutableprivate |
Definition at line 683 of file agglomeration_handler.h.
Referenced by setup_connectivity_of_agglomeration().
|
mutable |
Definition at line 505 of file agglomeration_handler.h.
|
mutable |
Definition at line 495 of file agglomeration_handler.h.
|
mutable |
Definition at line 510 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
mutableprivate |
Use this in reinit(cell) for (non-agglomerated, standard) cells, and return the result of scratch.reinit(cell) for cells
Definition at line 788 of file agglomeration_handler.h.
|
private |
Definition at line 764 of file agglomeration_handler.h.
Referenced by connect_to_tria_signals(), define_agglomerate(), get_triangulation(), initialize_agglomeration_data(), print_agglomeration(), setup_ghost_polytopes(), and dealii::internal::AgglomerationHandlerImplementation< dim, spacedim >::setup_master_neighbor_connectivity().
|
private |
Definition at line 807 of file agglomeration_handler.h.
Referenced by connect_to_tria_signals(), and ~AgglomerationHandler().