33 template <
int dim,
int spacedim>
43 namespace fullydistributed
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;
78 this->set_mesh_smoothing(
80 typename ::Triangulation<dim, spacedim>::MeshSmoothing
>(
84 this->set_mesh_smoothing(
86 typename ::Triangulation<dim, spacedim>::MeshSmoothing
>(
89 this->set_mesh_smoothing(construction_data.smoothing);
92 this->coarse_cell_id_to_coarse_cell_index_vector.clear();
93 this->coarse_cell_index_to_coarse_cell_id_vector.clear();
96 if (construction_data.coarse_cell_vertices.empty())
99 currently_processing_create_triangulation_for_internal_usage =
true;
101 currently_processing_create_triangulation_for_internal_usage =
false;
104 auto cell = this->begin();
106 cell->set_level_subdomain_id(
111 this->coarse_cell_id_to_coarse_cell_index_vector.emplace_back(
113 this->coarse_cell_index_to_coarse_cell_id_vector.emplace_back(
119 this->coarse_cell_index_to_coarse_cell_id_vector =
120 construction_data.coarse_cell_index_to_coarse_cell_id;
123 std::map<types::coarse_cell_id, unsigned int>
124 coarse_cell_id_to_coarse_cell_index_vector;
125 for (
unsigned int i = 0;
126 i < construction_data.coarse_cell_index_to_coarse_cell_id.size();
128 coarse_cell_id_to_coarse_cell_index_vector
129 [construction_data.coarse_cell_index_to_coarse_cell_id[i]] = i;
131 for (
auto i : coarse_cell_id_to_coarse_cell_index_vector)
132 this->coarse_cell_id_to_coarse_cell_index_vector.emplace_back(i);
135 currently_processing_prepare_coarsening_and_refinement_for_internal_usage =
137 currently_processing_create_triangulation_for_internal_usage =
true;
140 currently_processing_prepare_coarsening_and_refinement_for_internal_usage =
142 currently_processing_create_triangulation_for_internal_usage =
false;
145 auto cell_infos = construction_data.cell_infos;
150 std::sort(cell_info.begin(),
157 const auto a_coarse_cell_index =
158 this->coarse_cell_id_to_coarse_cell_index(
160 const auto b_coarse_cell_index =
161 this->coarse_cell_id_to_coarse_cell_index(
169 if (a_coarse_cell_index != b_coarse_cell_index)
170 return a_coarse_cell_index < b_coarse_cell_index;
177 for (
const auto &cell : this->cell_iterators())
179 if (cell->is_active())
180 cell->set_subdomain_id(
183 cell->set_level_subdomain_id(
188 for (
unsigned int level = 0;
192 auto cell = this->begin(
level);
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);
212 this->update_number_cache();
213 this->update_cell_relations();
218 template <
int dim,
int spacedim>
221 const
std::vector<
Point<spacedim>> &vertices,
226 currently_processing_create_triangulation_for_internal_usage,
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)
280 this->partitioner(serial_tria,
282 this->mpi_communicator));
285 if (this->is_multilevel_hierarchy_constructed())
289 other_tria_ptr = &serial_tria;
295 this->mpi_communicator,
299 this->create_triangulation(construction_data);
304 template <
int dim,
int spacedim>
308 const
unsigned int)> &partitioner,
311 this->partitioner = partitioner;
312 this->settings = settings;
317 template <
int dim,
int spacedim>
323 this->partitioner_distributed = &partitioner;
324 this->settings = settings;
329 template <
int dim,
int spacedim>
334 this->signals.pre_distributed_repartition();
340 this->partitioner_distributed->partition(*
this),
345 this->coarse_cell_id_to_coarse_cell_index_vector.clear();
346 this->coarse_cell_index_to_coarse_cell_id_vector.clear();
349 this->create_triangulation(construction_data);
352 this->signals.post_distributed_repartition();
357 template <
int dim,
int spacedim>
366 template <
int dim,
int spacedim>
371 currently_processing_prepare_coarsening_and_refinement_for_internal_usage,
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 =
388 coarse_cell_id_to_coarse_cell_index_vector) +
390 coarse_cell_index_to_coarse_cell_id_vector);
396 template <
int dim,
int spacedim>
408 template <
int dim,
int spacedim>
411 coarse_cell_id_to_coarse_cell_index(
412 const
types::coarse_cell_id coarse_cell_id)
const
414 const auto coarse_cell_index = std::lower_bound(
415 coarse_cell_id_to_coarse_cell_index_vector.begin(),
416 coarse_cell_id_to_coarse_cell_index_vector.end(),
418 [](
const std::pair<types::coarse_cell_id, unsigned int> &pair,
420 if (coarse_cell_index !=
421 coarse_cell_id_to_coarse_cell_index_vector.cend())
422 return coarse_cell_index->second;
429 template <
int dim,
int spacedim>
433 const unsigned int coarse_cell_index)
const
436 coarse_cell_index_to_coarse_cell_id_vector.size());
438 const auto coarse_cell_id =
439 coarse_cell_index_to_coarse_cell_id_vector[coarse_cell_index];
441 ExcMessage(
"You are trying to access a dummy cell!"));
442 return coarse_cell_id;
446 template <
int dim,
int spacedim>
451 this->local_cell_relations.
clear();
452 this->local_cell_relations.reserve(this->n_locally_owned_active_cells());
454 for (
const auto &cell : this->active_cell_iterators())
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
469 this->cell_attached_data.n_attached_deserialize == 0,
471 "Not all SolutionTransfer objects have been deserialized after the last call to load()."));
472 Assert(this->n_cells() > 0,
473 ExcMessage(
"Can not save() an empty Triangulation."));
481 unsigned int n_locally_owned_cells = this->n_locally_owned_active_cells();
483 unsigned int global_first_cell = 0;
485 int ierr = MPI_Exscan(&n_locally_owned_cells,
490 this->mpi_communicator);
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"
505 <<
" " << this->cell_attached_data.pack_callbacks_fixed.size()
506 <<
" " << this->cell_attached_data.pack_callbacks_variable.size()
507 <<
" " << this->n_global_active_cells() << std::endl;
511 this->save_attached_data(global_first_cell,
512 this->n_global_active_cells(),
518 int ierr = MPI_Info_create(&info);
521 const std::string fname_tria = filename +
"_triangulation.data";
525 ierr = MPI_File_open(this->mpi_communicator,
527 MPI_MODE_CREATE | MPI_MODE_WRONLY,
532 ierr = MPI_File_set_size(fh, 0);
536 ierr = MPI_Barrier(this->mpi_communicator);
538 ierr = MPI_Info_free(&info);
545 this->mpi_communicator,
549 std::vector<char> buffer;
553 const std::uint64_t buffer_size = buffer.size();
555 std::uint64_t offset = 0;
563 this->mpi_communicator);
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
607 Assert(this->n_cells() == 0,
608 ExcMessage(
"load() only works if the Triangulation is empty!"));
611 unsigned int version, numcpus, attached_count_fixed,
612 attached_count_variable, n_global_active_cells;
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 >>
620 attached_count_variable >> n_global_active_cells;
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";
646 ierr = MPI_File_open(this->mpi_communicator,
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;
676 this->mpi_communicator);
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<
702 construction_data.
comm = this->mpi_communicator;
704 this->create_triangulation(construction_data);
708 unsigned int n_locally_owned_cells = this->n_locally_owned_active_cells();
710 unsigned int global_first_cell = 0;
712 int ierr = MPI_Exscan(&n_locally_owned_cells,
717 this->mpi_communicator);
720 global_first_cell *=
sizeof(
unsigned int);
722 Assert(this->n_global_active_cells() == n_global_active_cells,
723 ExcMessage(
"Number of global active cells differ!"));
727 this->cell_attached_data.n_attached_data_sets = 0;
728 this->cell_attached_data.n_attached_deserialize =
729 attached_count_fixed + attached_count_variable;
732 this->load_attached_data(global_first_cell,
733 this->n_global_active_cells(),
734 this->n_locally_owned_active_cells(),
736 attached_count_fixed,
737 attached_count_variable);
739 this->update_cell_relations();
740 this->update_periodic_face_map();
741 this->update_number_cache();
751 template <
int dim,
int spacedim>
760 for (
const auto &cell : this->active_cell_iterators())
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 =
768 this->mpi_communicator) +
771 this->number_cache.number_of_global_coarse_cells =
772 number_of_global_coarse_cells;
782#include "fully_distributed_tria.inst"
types::coarse_cell_id get_coarse_cell_id() const
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
virtual std::size_t memory_consumption() const override
virtual void update_number_cache()
virtual void clear() override
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const override
void create_triangulation(const TriangulationDescription::Description< dim, spacedim > &construction_data) override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
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)
#define DEAL_II_NOT_IMPLEMENTED()
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
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
const types::coarse_cell_id invalid_coarse_cell_id
const types::subdomain_id artificial_subdomain_id
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
std::vector< std::vector< CellData< dim > > > cell_infos