16#ifndef dealii_distributed_tria_base_h
17#define dealii_distributed_tria_base_h
77 template <
int dim,
int spacedim = dim>
86 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
117 const ::Triangulation<dim, spacedim> &old_tria)
override;
183 const std::set<types::subdomain_id> &
201 const std::set<types::subdomain_id> &
208 const std::weak_ptr<const Utilities::MPI::Partitioner>
215 const std::weak_ptr<const Utilities::MPI::Partitioner>
227 std::set<::types::subdomain_id>>
236 virtual std::vector<types::boundary_id>
245 virtual std::vector<types::manifold_id>
302 const std::vector<bool> &vertex_locally_moved);
358 std::shared_ptr<const Utilities::MPI::Partitioner>
364 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
395 template <
int dim,
int spacedim = dim>
445 template <
int dim,
int spacedim = dim>
455 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
469 typename ::Triangulation<dim, spacedim>::cell_iterator;
472 typename ::Triangulation<dim, spacedim>::CellStatus;
482 save(
const std::string &filename)
const = 0;
490 load(
const std::string &filename,
const bool autopartition =
true) = 0;
604 const bool returns_variable_size_data);
656 const unsigned int handle,
660 const boost::iterator_range<std::vector<char>::const_iterator> &)>
672 const unsigned int global_num_cells,
673 const std::string &filename)
const;
684 const unsigned int global_num_cells,
685 const unsigned int local_num_cells,
686 const std::string &filename,
687 const unsigned int n_attached_deserialize_fixed,
688 const unsigned int n_attached_deserialize_variable);
737 typename ::Triangulation<dim, spacedim>::cell_iterator,
738 typename ::Triangulation<dim, spacedim>::CellStatus)>;
772 pack_data(
const std::vector<cell_relation_t> &cell_relations,
773 const std::vector<typename CellAttachedData::pack_callback_t>
774 &pack_callbacks_fixed,
775 const std::vector<typename CellAttachedData::pack_callback_t>
776 &pack_callbacks_variable);
804 const std::vector<cell_relation_t> &cell_relations,
805 const unsigned int handle,
806 const std::function<
void(
807 const typename ::Triangulation<dim, spacedim>::cell_iterator &,
808 const typename ::Triangulation<dim, spacedim>::CellStatus &,
809 const boost::iterator_range<std::vector<char>::const_iterator> &)>
810 &unpack_callback)
const;
827 save(
const unsigned int global_first_cell,
828 const unsigned int global_num_cells,
829 const std::string &filename)
const;
850 load(
const unsigned int global_first_cell,
851 const unsigned int global_num_cells,
852 const unsigned int local_num_cells,
853 const std::string &filename,
854 const unsigned int n_attached_deserialize_fixed,
855 const unsigned int n_attached_deserialize_variable);
const bool check_for_distorted_cells
MeshSmoothing smooth_grid
void pack_data(const std::vector< cell_relation_t > &cell_relations, const std::vector< typename CellAttachedData::pack_callback_t > &pack_callbacks_fixed, const std::vector< typename CellAttachedData::pack_callback_t > &pack_callbacks_variable)
std::vector< char > dest_data_fixed
void save(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) const
void unpack_cell_status(std::vector< cell_relation_t > &cell_relations) const
void unpack_data(const std::vector< cell_relation_t > &cell_relations, const unsigned int handle, const std::function< void(const typename ::Triangulation< dim, spacedim >::cell_iterator &, const typename ::Triangulation< dim, spacedim >::CellStatus &, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback) const
DataTransfer(const MPI_Comm &mpi_communicator)
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
void load(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
bool variable_size_data_stored
std::vector< char > src_data_fixed
virtual void clear() override
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) const
CellAttachedData cell_attached_data
DistributedTriangulationBase(const MPI_Comm &mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false)
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback)
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
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const CellStatus)> &pack_callback, const bool returns_variable_size_data)
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
virtual void load(const std::string &filename, const bool autopartition=true)=0
typename std::pair< cell_iterator, CellStatus > cell_relation_t
void load_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
DataTransfer data_transfer
const std::set< types::subdomain_id > & level_ghost_owners() const
virtual std::size_t memory_consumption() const override
virtual types::global_cell_index n_global_active_cells() const override
void update_reference_cells() override
const std::set< types::subdomain_id > & ghost_owners() const
types::subdomain_id n_subdomains
const std::weak_ptr< const Utilities::MPI::Partitioner > global_level_cell_index_partitioner(const unsigned int level) const
types::subdomain_id locally_owned_subdomain() const override
const MPI_Comm mpi_communicator
TriangulationBase(const MPI_Comm &mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false)
virtual unsigned int n_global_levels() const override
virtual std::vector< types::boundary_id > get_boundary_ids() const override
unsigned int n_locally_owned_active_cells() const
virtual std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors() const
virtual bool is_multilevel_hierarchy_constructed() const =0
virtual void update_number_cache()
const std::weak_ptr< const Utilities::MPI::Partitioner > global_active_cell_index_partitioner() const
void communicate_locally_moved_vertices(const std::vector< bool > &vertex_locally_moved)
virtual MPI_Comm get_communicator() const override
void reset_global_cell_indices()
types::subdomain_id my_subdomain
virtual std::vector< types::manifold_id > get_manifold_ids() const override
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria) override
virtual ~TriangulationBase() override
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#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