16#ifndef dealii_distributed_tria_base_h
17#define dealii_distributed_tria_base_h
78 template <
int dim,
int spacedim = dim>
88 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
90 const bool check_for_distorted_cells =
false);
101 get_communicator()
const override;
119 const ::Triangulation<dim, spacedim> &old_tria)
override;
140 n_locally_owned_active_cells()
const;
148 n_global_active_cells()
const override;
154 memory_consumption()
const override;
165 n_global_levels()
const override;
174 locally_owned_subdomain()
const override;
185 const std::set<types::subdomain_id> &
186 ghost_owners()
const;
203 const std::set<types::subdomain_id> &
204 level_ghost_owners()
const;
206 const std::weak_ptr<const Utilities::MPI::Partitioner>
207 global_active_cell_index_partitioner()
const override;
209 const std::weak_ptr<const Utilities::MPI::Partitioner>
210 global_level_cell_index_partitioner(
211 const unsigned int level)
const override;
219 virtual std::vector<types::boundary_id>
220 get_boundary_ids()
const override;
228 virtual std::vector<types::manifold_id>
229 get_manifold_ids()
const override;
284 communicate_locally_moved_vertices(
285 const std::vector<bool> &vertex_locally_moved);
288 n_global_coarse_cells()
const override;
353 std::shared_ptr<const Utilities::MPI::Partitioner>
359 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
371 update_number_cache();
377 update_reference_cells()
override;
383 reset_global_cell_indices();
435 template <
int dim,
int spacedim = dim>
446 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
448 const bool check_for_distorted_cells =
false);
460 typename ::Triangulation<dim, spacedim>::cell_iterator;
463 typename ::Triangulation<dim, spacedim>::CellStatus;
481 has_hanging_nodes()
const override;
491 save(
const std::string &filename)
const = 0;
499 load(
const std::string &filename) = 0;
508 load(
const std::string &filename,
const bool autopartition) = 0;
619 register_data_attach(
622 const bool returns_variable_size_data);
673 notify_ready_to_unpack(
674 const unsigned int handle,
678 const boost::iterator_range<std::vector<char>::const_iterator> &)>
690 save_attached_data(
const unsigned int global_first_cell,
691 const unsigned int global_num_cells,
692 const std::string &filename)
const;
703 load_attached_data(
const unsigned int global_first_cell,
704 const unsigned int global_num_cells,
705 const unsigned int local_num_cells,
706 const std::string &filename,
707 const unsigned int n_attached_deserialize_fixed,
708 const unsigned int n_attached_deserialize_variable);
757 typename ::Triangulation<dim, spacedim>::cell_iterator,
758 typename ::Triangulation<dim, spacedim>::CellStatus)>;
792 pack_data(
const std::vector<cell_relation_t> &cell_relations,
793 const std::vector<typename CellAttachedData::pack_callback_t>
794 &pack_callbacks_fixed,
795 const std::vector<typename CellAttachedData::pack_callback_t>
796 &pack_callbacks_variable);
808 unpack_cell_status(std::vector<cell_relation_t> &cell_relations)
const;
824 const std::vector<cell_relation_t> &cell_relations,
825 const unsigned int handle,
826 const std::function<
void(
827 const typename ::Triangulation<dim, spacedim>::cell_iterator &,
828 const typename ::Triangulation<dim, spacedim>::CellStatus &,
829 const boost::iterator_range<std::vector<char>::const_iterator> &)>
830 &unpack_callback)
const;
847 save(
const unsigned int global_first_cell,
848 const unsigned int global_num_cells,
849 const std::string &filename)
const;
870 load(
const unsigned int global_first_cell,
871 const unsigned int global_num_cells,
872 const unsigned int local_num_cells,
873 const std::string &filename,
874 const unsigned int n_attached_deserialize_fixed,
875 const unsigned int n_attached_deserialize_variable);
std::vector< char > dest_data_fixed
std::vector< int > dest_sizes_variable
MPI_Comm mpi_communicator
std::vector< int > src_sizes_variable
std::vector< char > src_data_variable
std::vector< unsigned int > sizes_fixed_cumulative
std::vector< char > dest_data_variable
bool variable_size_data_stored
std::vector< char > src_data_fixed
virtual void load(const std::string &filename)=0
virtual void load(const std::string &filename, const bool autopartition)=0
CellAttachedData cell_attached_data
virtual void update_cell_relations()=0
typename ::Triangulation< dim, spacedim >::CellStatus CellStatus
std::vector< cell_relation_t > local_cell_relations
virtual void save(const std::string &filename) const =0
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
typename std::pair< cell_iterator, CellStatus > cell_relation_t
DataTransfer data_transfer
types::subdomain_id n_subdomains
const MPI_Comm mpi_communicator
virtual bool is_multilevel_hierarchy_constructed() const =0
types::subdomain_id my_subdomain
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
unsigned int n_attached_data_sets
std::vector< pack_callback_t > pack_callbacks_fixed
std::vector< pack_callback_t > pack_callbacks_variable
std::function< std::vector< char >(typename ::Triangulation< dim, spacedim >::cell_iterator, typename ::Triangulation< dim, spacedim >::CellStatus)> pack_callback_t
unsigned int n_attached_deserialize
unsigned int n_global_levels
unsigned int n_locally_owned_active_cells
std::set< types::subdomain_id > level_ghost_owners
types::global_cell_index n_global_active_cells
std::set< types::subdomain_id > ghost_owners
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
types::coarse_cell_id number_of_global_coarse_cells