#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) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | serialize (Archive &ar, const unsigned int version) |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") 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 Public Member Functions inherited from Subscriptor | |
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 171 of file agglomeration_handler.h.
using AgglomerationHandler< dim, spacedim >::agglomeration_iterator = AgglomerationIterator<dim, spacedim> |
Definition at line 174 of file agglomeration_handler.h.
using AgglomerationHandler< dim, spacedim >::AgglomerationContainer |
Definition at line 176 of file agglomeration_handler.h.
|
private |
Definition at line 651 of file agglomeration_handler.h.
enum AgglomerationHandler::CellAgglomerationType |
Enumerator | |
---|---|
master | |
slave |
Definition at line 180 of file agglomeration_handler.h.
|
explicit |
Definition at line 21 of file agglomeration_handler.cc.
References Assert, AgglomerationHandler< dim, spacedim >::cached_tria, AgglomerationHandler< dim, spacedim >::hybrid_mesh, AgglomerationHandler< dim, spacedim >::initialize_agglomeration_data(), and AgglomerationHandler< dim, spacedim >::n_agglomerations.
|
default |
|
inline |
Definition at line 193 of file agglomeration_handler.h.
References AgglomerationHandler< dim, spacedim >::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.
|
inlineprivate |
Returns true if the two given cells are agglomerated together.
Definition at line 1081 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::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 1014 of file agglomeration_handler.h.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin | ( | ) |
Iterator to the first polytope.
Definition at line 1128 of file agglomeration_handler.h.
References Assert.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::begin | ( | ) | const |
Iterator to the first polytope.
Definition at line 1117 of file agglomeration_handler.h.
References Assert.
|
inline |
Definition at line 926 of file agglomeration_handler.h.
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 1183 of file agglomeration_handler.h.
References dealii::CellsAgglomerator< dim, RtreeType >::extraction_level, and dealii::CellsAgglomerator< dim, RtreeType >::parent_node_to_children_nodes.
|
inlineprivate |
Reinitialize the agglomeration data.
Definition at line 548 of file agglomeration_handler.h.
References AgglomerationHandler< dim, spacedim >::cached_tria, AgglomerationHandler< dim, spacedim >::initialize_agglomeration_data(), AgglomerationHandler< dim, spacedim >::tria, and AgglomerationHandler< dim, spacedim >::tria_listener.
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(), Assert, AssertThrow, DoFTools::extract_locally_relevant_dofs(), DoFTools::make_sparsity_pattern(), 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.
References Assert.
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 Assert, and numbers::invalid_subdomain_id.
Referenced by 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 Assert, Utils::compute_connected_components(), and Utils::create_graph_from_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 AssertThrow, and FiniteElement< int dim, int spacedim >::degree.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end | ( | ) |
Iterator to one past the last polygonal element.
Definition at line 1150 of file agglomeration_handler.h.
References Assert.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::end | ( | ) | const |
Iterator to one past the last polygonal element.
Definition at line 1139 of file agglomeration_handler.h.
References Assert.
void AgglomerationHandler< dim, spacedim >::exchange_interface_values | ( | ) |
Definition at line 347 of file agglomeration_handler.cc.
|
inline |
TODO: remove this in favour of the accessor version.
master_cell |
Definition at line 954 of file agglomeration_handler.h.
References Assert.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
inline |
Definition at line 900 of file agglomeration_handler.h.
|
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 968 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::compute_quality_metrics().
|
inline |
Definition at line 882 of file agglomeration_handler.h.
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 936 of file agglomeration_handler.h.
|
inline |
Definition at line 917 of file agglomeration_handler.h.
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 891 of file agglomeration_handler.h.
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().
|
inlineprivate |
Definition at line 1067 of file agglomeration_handler.h.
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 945 of file agglomeration_handler.h.
|
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 978 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::fill_interpolation_matrix(), and dealii::PolyUtils::internal::interpolate_to_fine_grid().
|
inline |
Definition at line 909 of file agglomeration_handler.h.
Referenced by 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.
Referenced by AgglomerationHandler< dim, spacedim >::AgglomerationHandler(), and AgglomerationHandler< dim, spacedim >::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 update_JxW_values, update_normal_vectors, and update_quadrature_points.
|
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 Assert.
|
inline |
Helper function to determine whether or not a cell is a master or a slave
Definition at line 989 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::fill_interpolation_matrix(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), and dealii::internal::AgglomerationHandlerImplementation< int, int >::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 1004 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::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 1057 of file agglomeration_handler.h.
AgglomerationIterator< dim, spacedim > AgglomerationHandler< dim, spacedim >::last | ( | ) |
Iterator to the last polygonal element.
Definition at line 1161 of file agglomeration_handler.h.
References Assert.
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.
|
inline |
Returns the number of agglomerate cells in the grid.
Definition at line 1097 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::partition_locally_owned_regions().
|
inlinenoexcept |
Definition at line 1038 of file agglomeration_handler.h.
Referenced by Utils::fill_injection_matrix().
|
inlinenoexcept |
Definition at line 1029 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::compute_global_error(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
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 1047 of file agglomeration_handler.h.
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 1173 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix(), dealii::PolyUtils::compute_global_error(), dealii::PolyUtils::compute_quality_metrics(), dealii::PolyUtils::export_polygon_to_csv_file(), Utils::fill_injection_matrix(), and dealii::PolyUtils::interpolate_to_fine_grid().
|
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 1106 of file agglomeration_handler.h.
Referenced by Utils::fill_injection_matrix().
|
inline |
Display the indices of the vector identifying which cell is agglomerated with which master.
Definition at line 362 of file agglomeration_handler.h.
References AgglomerationHandler< dim, spacedim >::master_slave_relationships, and AgglomerationHandler< dim, spacedim >::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.
Referenced by dealii::PolyUtils::assemble_dg_matrix(), dealii::PolyUtils::compute_global_error(), and dealii::PolyUtils::compute_quality_metrics().
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 Assert, and dealii::internal::AgglomerationHandlerImplementation< int, int >::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.
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 Assert, Utilities::MPI::job_supports_mpi(), dealii::internal::AgglomerationHandlerImplementation< int, int >::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 Assert, and Utilities::MPI::some_to_some().
|
private |
|
friend |
Definition at line 240 of file agglomeration_handler.h.
|
friend |
Definition at line 237 of file agglomeration_handler.h.
|
friend |
Definition at line 858 of file agglomeration_handler.h.
DoFHandler<dim, spacedim> AgglomerationHandler< dim, spacedim >::agglo_dh |
DoFHandler for the agglomerated space
Definition at line 467 of file agglomeration_handler.h.
Referenced by Utils::fill_injection_matrix(), dealii::PolyUtils::fill_interpolation_matrix(), and dealii::PolyUtils::internal::interpolate_to_fine_grid().
|
mutableprivate |
Definition at line 798 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 804 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 801 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 794 of file agglomeration_handler.h.
|
private |
Definition at line 814 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
|
private |
Definition at line 822 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
|
private |
Definition at line 808 of file agglomeration_handler.h.
|
private |
Definition at line 820 of file agglomeration_handler.h.
|
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 672 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
std::unique_ptr<MappingBox<dim> > AgglomerationHandler< dim, spacedim >::box_mapping |
Definition at line 474 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
|
private |
Definition at line 767 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::AgglomerationHandler(), and AgglomerationHandler< dim, spacedim >::connect_to_tria_signals().
|
private |
Definition at line 769 of file agglomeration_handler.h.
|
private |
Dummy FE_Nothing
Definition at line 842 of file agglomeration_handler.h.
|
private |
Eulerian vector describing the new cells obtained by the bounding boxes
Definition at line 780 of file agglomeration_handler.h.
|
private |
Definition at line 773 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
|
private |
Definition at line 775 of file agglomeration_handler.h.
|
private |
Bool that keeps track whether the mesh is composed also by standard deal.II cells as (trivial) polytopes.
Definition at line 868 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::AgglomerationHandler().
|
private |
Definition at line 816 of file agglomeration_handler.h.
|
private |
Definition at line 810 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 708 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 687 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 724 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 716 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 699 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 757 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 740 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 679 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 746 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 734 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 752 of file agglomeration_handler.h.
|
private |
Definition at line 765 of file agglomeration_handler.h.
|
private |
|
private |
Definition at line 828 of file agglomeration_handler.h.
|
private |
A contiguous container for all of the master cells.
Definition at line 858 of file agglomeration_handler.h.
|
private |
Definition at line 835 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 641 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::print_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
private |
Same as the one above, but storing cell iterators rather than indices.
Definition at line 649 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
private |
Record the number of agglomerations on the grid.
Definition at line 632 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::AgglomerationHandler().
|
private |
Dummy FEFaceValues, needed for face quadratures.
Definition at line 852 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::reinit_master().
|
private |
Dummy FEValues, needed for cell quadratures.
Definition at line 847 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 653 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
DoFHandler<dim, spacedim> AgglomerationHandler< dim, spacedim >::output_dh |
DoFHandler for the finest space: classical deal.II space
Definition at line 472 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 872 of file agglomeration_handler.h.
|
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 663 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
private |
|
private |
Definition at line 874 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 712 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 690 of file agglomeration_handler.h.
Referenced by dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
mutableprivate |
Definition at line 728 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 719 of file agglomeration_handler.h.
|
mutableprivate |
Definition at line 703 of file agglomeration_handler.h.
|
mutable |
Definition at line 514 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
mutable |
Definition at line 499 of file agglomeration_handler.h.
Referenced by dealii::PolyUtils::assemble_dg_matrix().
|
mutableprivate |
Definition at line 682 of file agglomeration_handler.h.
|
mutable |
Definition at line 504 of file agglomeration_handler.h.
|
mutable |
Definition at line 494 of file agglomeration_handler.h.
|
mutable |
Definition at line 509 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 787 of file agglomeration_handler.h.
|
private |
Definition at line 763 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::connect_to_tria_signals(), AgglomerationHandler< dim, spacedim >::print_agglomeration(), and dealii::internal::AgglomerationHandlerImplementation< int, int >::setup_master_neighbor_connectivity().
|
private |
Definition at line 806 of file agglomeration_handler.h.
Referenced by AgglomerationHandler< dim, spacedim >::connect_to_tria_signals(), and AgglomerationHandler< dim, spacedim >::~AgglomerationHandler().