33 template <
int dim,
int spacedim>
45 template <
int dim,
int spacedim>
51 const
unsigned int n_partitions) {
54 , currently_processing_create_triangulation_for_internal_usage(
false)
55 , currently_processing_prepare_coarsening_and_refinement_for_internal_usage(
61 template <
int dim,
int spacedim>
69 Assert(construction_data.comm == this->mpi_communicator,
70 ExcMessage(
"MPI communicators do not match!"));
73 settings = construction_data.settings;
80 typename ::Triangulation<dim, spacedim>::MeshSmoothing
>(
86 typename ::Triangulation<dim, spacedim>::MeshSmoothing
>(
96 if (construction_data.coarse_cell_vertices.empty())
104 auto cell = this->
begin();
106 cell->set_level_subdomain_id(
120 construction_data.coarse_cell_index_to_coarse_cell_id;
123 std::map<types::coarse_cell_id, unsigned int>
125 for (
unsigned int i = 0;
126 i < construction_data.coarse_cell_index_to_coarse_cell_id.size();
129 [construction_data.coarse_cell_index_to_coarse_cell_id[i]] = i;
132 this->coarse_cell_id_to_coarse_cell_index_vector.emplace_back(i);
145 auto cell_infos = construction_data.cell_infos;
149 for (
auto &cell_info : cell_infos)
150 std::sort(cell_info.begin(),
154 const CellId a_id(a.id);
155 const CellId b_id(b.id);
157 const auto a_coarse_cell_index =
158 this->coarse_cell_id_to_coarse_cell_index(
159 a_id.get_coarse_cell_id());
160 const auto b_coarse_cell_index =
161 this->coarse_cell_id_to_coarse_cell_index(
162 b_id.get_coarse_cell_id());
169 if (a_coarse_cell_index != b_coarse_cell_index)
170 return a_coarse_cell_index < b_coarse_cell_index;
179 if (cell->is_active())
180 cell->set_subdomain_id(
183 cell->set_level_subdomain_id(
188 for (
unsigned int level = 0;
189 level < cell_infos.size() && !cell_infos[level].empty();
192 auto cell = this->
begin(level);
193 auto cell_info = cell_infos[level].begin();
194 for (; cell_info != cell_infos[level].end(); ++cell_info)
197 while (cell_info->id != cell->id().template to_binary<dim>())
201 if (cell->is_active())
202 cell->set_subdomain_id(cell_info->subdomain_id);
206 construct_multigrid_hierarchy)
207 cell->set_level_subdomain_id(cell_info->level_subdomain_id);
218 template <
int dim,
int spacedim>
228 "You have called the overload of\n"
230 " parallel::fullydistributed::Triangulation::"
231 "create_triangulation()\n"
233 "which takes 3 arguments. This function is not yet implemented for "
234 "this class. If you have not called this function directly, it "
235 "might have been called via a function from the GridGenerator or "
236 "GridIn namespace. To set up a fully-distributed Triangulation with "
237 "these utility functions, please start by using the same process to "
238 "set up a serial Triangulation, parallel::shared::Triangulation, or "
239 "a parallel::distributed::Triangulation. Once that is complete use "
240 "the copy_triangulation() member function to finish setting up the "
241 "original fully distributed Triangulation. Alternatively, you can "
242 "use TriangulationDescription::Utilities::"
243 "create_description_from_triangulation() or "
244 "create_description_from_triangulation_in_groups() to create the "
245 "description of the local partition, and pass that description to "
246 "parallel::fullydistributed::Triangulation::create_triangulation()."));
255 template <
int dim,
int spacedim>
264 const ::Triangulation<dim, spacedim> *other_tria_ptr = &other_tria;
273 if (
dynamic_cast<const ::parallel::TriangulationBase<dim, spacedim>
274 *
>(&other_tria) ==
nullptr)
289 other_tria_ptr = &serial_tria;
304 template <
int dim,
int spacedim>
317 template <
int dim,
int spacedim>
329 template <
int dim,
int spacedim>
334 this->
signals.pre_distributed_repartition();
352 this->
signals.post_distributed_repartition();
357 template <
int dim,
int spacedim>
366 template <
int dim,
int spacedim>
372 ExcMessage(
"No coarsening and refinement is supported!"));
374 return ::Triangulation<dim, spacedim>::
375 prepare_coarsening_and_refinement();
380 template <
int dim,
int spacedim>
384 const std::size_t mem =
396 template <
int dim,
int spacedim>
408 template <
int dim,
int spacedim>
412 const
types::coarse_cell_id coarse_cell_id)
const
414 const auto coarse_cell_index = std::lower_bound(
418 [](
const std::pair<types::coarse_cell_id, unsigned int> &pair,
420 if (coarse_cell_index !=
422 return coarse_cell_index->second;
429 template <
int dim,
int spacedim>
433 const unsigned int coarse_cell_index)
const
438 const auto coarse_cell_id =
441 ExcMessage(
"You are trying to access a dummy cell!"));
442 return coarse_cell_id;
446 template <
int dim,
int spacedim>
455 if (cell->is_locally_owned())
456 this->local_cell_relations.emplace_back(
462 template <
int dim,
int spacedim>
466#ifdef DEAL_II_WITH_MPI
471 "Not all SolutionTransfer objects have been deserialized after the last call to load()."));
473 ExcMessage(
"Can not save() an empty Triangulation."));
483 unsigned int global_first_cell = 0;
485 int ierr = MPI_Exscan(&n_locally_owned_cells,
493 global_first_cell *=
sizeof(
unsigned int);
498 std::string fname = std::string(filename) +
".info";
499 std::ofstream f(fname);
500 f <<
"version nproc n_attached_fixed_size_objs n_attached_variable_size_objs n_global_active_cells"
518 int ierr = MPI_Info_create(&info);
521 const std::string fname_tria = filename +
"_triangulation.data";
527 MPI_MODE_CREATE | MPI_MODE_WRONLY,
532 ierr = MPI_File_set_size(fh, 0);
538 ierr = MPI_Info_free(&info);
549 std::vector<char> buffer;
553 const std::uint64_t buffer_size = buffer.size();
555 std::uint64_t offset = 0;
567 ierr = MPI_File_write_at(
569 myrank *
sizeof(std::uint64_t),
577 const std::uint64_t global_position =
578 mpisize *
sizeof(std::uint64_t) + offset;
590 ierr = MPI_File_close(&fh);
602 template <
int dim,
int spacedim>
606#ifdef DEAL_II_WITH_MPI
608 ExcMessage(
"load() only works if the Triangulation is empty!"));
611 unsigned int version, numcpus, attached_count_fixed,
614 std::string fname = std::string(filename) +
".info";
615 std::ifstream f(fname);
617 std::string firstline;
618 getline(f, firstline);
619 f >> version >> numcpus >> attached_count_fixed >>
623 const auto expected_version = ::internal::
624 CellAttachedDataSerializer<dim, spacedim>::version_number;
627 ExcMessage(
"Incompatible version found in .info file."));
640 int ierr = MPI_Info_create(&info);
643 const std::string fname_tria = filename +
"_triangulation.data";
653 ierr = MPI_Info_free(&info);
657 std::uint64_t buffer_size;
659 ierr = MPI_File_read_at(
661 myrank *
sizeof(std::uint64_t),
668 std::uint64_t offset = 0;
680 const std::uint64_t global_position =
681 mpisize *
sizeof(std::uint64_t) + offset;
684 std::vector<char> buffer(buffer_size);
694 ierr = MPI_File_close(&fh);
697 auto construction_data = ::Utilities::template unpack<
710 unsigned int global_first_cell = 0;
712 int ierr = MPI_Exscan(&n_locally_owned_cells,
720 global_first_cell *=
sizeof(
unsigned int);
723 ExcMessage(
"Number of global active cells differ!"));
729 attached_count_fixed + attached_count_variable;
733 this->n_global_active_cells(),
736 attached_count_fixed,
737 attached_count_variable);
751 template <
int dim,
int spacedim>
761 if (!cell->is_artificial())
762 number_of_global_coarse_cells =
763 std::max(number_of_global_coarse_cells,
764 cell->id().get_coarse_cell_id());
766 number_of_global_coarse_cells =
772 number_of_global_coarse_cells;
782#include "distributed/fully_distributed_tria.inst"
@ limit_level_difference_at_vertices
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
std::vector< typename internal::CellAttachedDataSerializer< dim, spacedim >::cell_relation_t > local_cell_relations
cell_iterator begin(const unsigned int level=0) const
virtual void set_mesh_smoothing(const MeshSmoothing mesh_smoothing)
std::vector< Point< spacedim > > vertices
void update_periodic_face_map()
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &file_basename) const
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 &file_basename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
internal::CellAttachedData< dim, spacedim > cell_attached_data
unsigned int n_cells() const
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)
virtual std::size_t memory_consumption() const override
virtual types::global_cell_index n_global_active_cells() const override
const MPI_Comm mpi_communicator
unsigned int n_locally_owned_active_cells() const
virtual void update_number_cache()
virtual void clear() override
virtual void execute_coarsening_and_refinement() override
virtual void update_number_cache() override
bool currently_processing_create_triangulation_for_internal_usage
TriangulationDescription::Settings settings
Triangulation(const MPI_Comm mpi_communicator)
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id_vector
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const override
bool currently_processing_prepare_coarsening_and_refinement_for_internal_usage
void copy_triangulation(const ::Triangulation< dim, spacedim > &other_tria) override
void set_partitioner(const std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> &partitioner, const TriangulationDescription::Settings &settings)
ObserverPointer< const RepartitioningPolicyTools::Base< dim, spacedim > > partitioner_distributed
void update_cell_relations()
virtual bool prepare_coarsening_and_refinement() override
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const override
std::vector< std::pair< types::coarse_cell_id, unsigned int > > coarse_cell_id_to_coarse_cell_index_vector
virtual void save(const std::string &filename) const override
std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> partitioner
virtual bool is_multilevel_hierarchy_constructed() const override
void create_triangulation(const TriangulationDescription::Description< dim, spacedim > &construction_data) override
virtual void load(const std::string &filename) override
virtual std::size_t memory_consumption() const override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NOT_IMPLEMENTED()
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators() const
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Description< dim, spacedim > create_description_from_triangulation(const ::Triangulation< dim, spacedim > &tria, const MPI_Comm comm, const TriangulationDescription::Settings settings=TriangulationDescription::Settings::default_setting, const unsigned int my_rank_in=numbers::invalid_unsigned_int)
@ construct_multigrid_hierarchy
int File_write_at_c(MPI_File fh, MPI_Offset offset, const void *buf, MPI_Count count, MPI_Datatype datatype, MPI_Status *status)
int File_read_at_c(MPI_File fh, MPI_Offset offset, void *buf, MPI_Count count, MPI_Datatype datatype, MPI_Status *status)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
T max(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
std::size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
constexpr unsigned int invalid_unsigned_int
constexpr types::subdomain_id artificial_subdomain_id
constexpr types::coarse_cell_id invalid_coarse_cell_id
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
global_cell_index coarse_cell_id