34template <
int structdim>
38 , manifold_id(
numbers::flat_manifold_id)
43template <
int structdim>
50 for (
unsigned int i = 0; i <
vertices.size(); ++i)
90 template <
int dim,
int spacedim>
91 struct DescriptionTemp
97 template <
class Archive>
99 serialize(Archive &ar,
const unsigned int )
112 const std::vector<unsigned int> & relevant_processes,
113 const std::vector<DescriptionTemp<dim, spacedim>> &description_temp,
115 const bool vertices_have_unique_ids)
117 const auto create_request = [&](
const unsigned int other_rank) {
118 const auto ptr = std::find(relevant_processes.begin(),
119 relevant_processes.end(),
124 const auto other_rank_index =
125 std::distance(relevant_processes.begin(), ptr);
127 return description_temp[other_rank_index];
130 const auto process_request =
131 [&](
const unsigned int,
132 const DescriptionTemp<dim, spacedim> &request) {
133 this->
merge(request, vertices_have_unique_ids);
137 DescriptionTemp<dim, spacedim>>(relevant_processes,
149 merge(
const DescriptionTemp<dim, spacedim> &other,
150 const bool vertices_have_unique_ids)
153 std::max(other.cell_infos.size(), this->cell_infos.size()));
155 if (vertices_have_unique_ids ==
false)
158 std::map<Point<spacedim>,
161 map_point_to_local_vertex_index(
167 map_point_to_local_vertex_index[coarse_cell_vertices[i]
171 std::map<unsigned int, unsigned int>
172 map_old_to_new_local_vertex_index;
177 if (map_point_to_local_vertex_index.find(p.
second) ==
178 map_point_to_local_vertex_index.end())
181 map_point_to_local_vertex_index[p.second] =
182 map_old_to_new_local_vertex_index[p.first] = counter++;
185 map_old_to_new_local_vertex_index[p.first] =
186 map_point_to_local_vertex_index[p.second];
189 auto other_coarse_cells_copy = other.coarse_cells;
191 for (
auto &cell : other_coarse_cells_copy)
193 v = map_old_to_new_local_vertex_index[v];
196 other_coarse_cells_copy.begin(),
197 other_coarse_cells_copy.end());
202 other.coarse_cells.begin(),
203 other.coarse_cells.end());
206 other.coarse_cell_vertices.begin(),
207 other.coarse_cell_vertices.end());
212 other.coarse_cell_index_to_coarse_cell_id.begin(),
213 other.coarse_cell_index_to_coarse_cell_id.end());
215 for (
unsigned int i = 0; i < this->
cell_infos.size(); ++i)
216 this->cell_infos[i].
insert(this->cell_infos[i].
end(),
217 other.cell_infos[i].begin(),
218 other.cell_infos[i].end());
234 for (
unsigned int i = 0; i < this->
coarse_cells.size(); ++i)
235 temp.emplace_back(this->coarse_cell_index_to_coarse_cell_id[i],
236 this->coarse_cells[i],
239 std::sort(temp.begin(),
241 [](
const auto &a,
const auto &b) {
242 return std::get<0>(a) < std::get<0>(b);
244 temp.erase(std::unique(temp.begin(),
246 [](
const auto &a,
const auto &b) {
247 return std::get<0>(a) == std::get<0>(b);
250 std::sort(temp.begin(),
252 [](
const auto &a,
const auto &b) {
253 return std::get<2>(a) < std::get<2>(b);
259 for (
unsigned int i = 0; i < temp.size(); ++i)
262 std::get<0>(temp[i]);
270 this->coarse_cell_vertices.end(),
271 [](
const auto &a,
const auto &b) {
272 return a.first < b.first;
276 this->coarse_cell_vertices.end(),
277 [](
const auto &a,
const auto &b) {
278 if (a.first == b.first)
280 Assert(a.second.distance(b.second) < 10e-8,
290 for (
unsigned int i = 0; i < this->
cell_infos.size(); ++i)
292 if (this->cell_infos[i].size() == 0)
295 std::sort(this->cell_infos[i].
begin(),
296 this->cell_infos[i].
end(),
297 [](
const auto &a,
const auto &b) {
301 std::vector<CellData<dim>> temp;
302 temp.push_back(this->cell_infos[i][0]);
304 for (
unsigned int j = 1; j < this->
cell_infos[i].size(); ++j)
305 if (temp.back().id == cell_infos[i][j].id)
307 temp.back().subdomain_id =
309 this->cell_infos[i][j].subdomain_id);
310 temp.back().level_subdomain_id =
311 std::min(temp.back().level_subdomain_id,
312 this->cell_infos[i][j].level_subdomain_id);
316 temp.push_back(this->cell_infos[i][j]);
328 Description<dim, spacedim>
334 Description<dim, spacedim> description;
337 description.comm =
comm;
342 description.smoothing = mesh_smoothing;
344 std::map<unsigned int, unsigned int> map;
348 description.coarse_cell_vertices.push_back(
356 for (unsigned
int v = 0; v < cell.vertices.size(); ++v)
357 cell.vertices[v] = map[cell.vertices[v]];
359 description.coarse_cell_index_to_coarse_cell_id =
368 std::vector<std::pair<unsigned int, Point<spacedim>>>
379 template <
int dim,
int spacedim>
381 mark_cell_and_its_parents(
383 std::vector<std::vector<bool>> & cell_marked)
385 cell_marked[cell->level()][cell->index()] =
true;
386 if (cell->level() != 0)
387 mark_cell_and_its_parents(cell->parent(), cell_marked);
395 template <
int dim,
int spacedim>
396 class CreateDescriptionFromTriangulationHelper
399 CreateDescriptionFromTriangulationHelper(
400 const ::Triangulation<dim, spacedim> &tria,
402 const typename ::Triangulation<dim, spacedim>::cell_iterator
403 &)> & subdomain_id_function,
405 const typename ::Triangulation<dim, spacedim>::cell_iterator
406 &)> & level_subdomain_id_function,
421 (
tria.get_mesh_smoothing() &
423 spacedim>::limit_level_difference_at_vertices),
425 "Source triangulation has to be set up with "
426 "limit_level_difference_at_vertices if the construction of the "
427 "multigrid hierarchy is requested!"));
430 tria, coinciding_vertex_groups, vertex_to_coinciding_vertex_group);
433 template <
typename T>
435 create_description_for_rank(
const unsigned int my_rank)
const
439 set_additional_data(construction_data);
444 add_vertices_of_cell_to_vertices_owned_by_locally_owned_cells =
445 [
this](
const auto & cell,
446 std::vector<bool> &vertices_owned_by_locally_owned_cells) {
450 vertices_owned_by_locally_owned_cells[cell->vertex_index(
452 const auto coinciding_vertex_group =
454 cell->vertex_index(v));
455 if (coinciding_vertex_group !=
458 coinciding_vertex_group->
second))
459 vertices_owned_by_locally_owned_cells[co_vertex] = true;
465 std::vector<std::vector<bool>> cell_marked(
tria.
n_levels());
475 std::vector<bool> vertices_owned_by_locally_owned_cells_on_level(
477 for (
const auto &cell :
tria.cell_iterators_on_level(
level))
480 add_vertices_of_cell_to_vertices_owned_by_locally_owned_cells(
481 cell, vertices_owned_by_locally_owned_cells_on_level);
483 for (
const auto &cell :
tria.active_cell_iterators())
485 add_vertices_of_cell_to_vertices_owned_by_locally_owned_cells(
486 cell, vertices_owned_by_locally_owned_cells_on_level);
491 const auto is_locally_relevant_on_level = [&](
const auto &cell) {
493 if (vertices_owned_by_locally_owned_cells_on_level
494 [cell->vertex_index(v)])
500 for (
const auto &cell :
tria.cell_iterators_on_level(
level))
501 if (is_locally_relevant_on_level(cell))
502 mark_cell_and_its_parents(cell, cell_marked);
511 for (
const auto &cell :
tria.cell_iterators_on_level(0))
513 if (!cell_marked[cell->level()][cell->index()])
518 cell_data.material_id = cell->material_id();
519 cell_data.manifold_id = cell->manifold_id();
521 cell_data.
vertices[v] = cell->vertex_index(v);
522 construction_data.coarse_cells.push_back(cell_data);
526 vertices_locally_relevant[cell->vertex_index(v)] = true;
529 construction_data.coarse_cell_index_to_coarse_cell_id.push_back(
530 cell->id().get_coarse_cell_id());
533 add_vertices(construction_data, vertices_locally_relevant);
538 construction_data.cell_infos.resize(
542 std::vector<bool> vertices_owned_by_locally_owned_active_cells(
544 for (
const auto &cell :
tria.active_cell_iterators())
546 add_vertices_of_cell_to_vertices_owned_by_locally_owned_cells(
547 cell, vertices_owned_by_locally_owned_active_cells);
551 const auto is_locally_relevant_on_active_level =
552 [&](
const auto &cell) {
553 if (cell->is_active())
555 if (vertices_owned_by_locally_owned_active_cells
556 [cell->vertex_index(v)])
561 for (
unsigned int level = 0;
566 std::vector<bool> vertices_owned_by_locally_owned_cells_on_level(
568 for (
const auto &cell :
tria.cell_iterators_on_level(
level))
571 (cell->is_active() &&
573 add_vertices_of_cell_to_vertices_owned_by_locally_owned_cells(
574 cell, vertices_owned_by_locally_owned_cells_on_level);
578 const auto is_locally_relevant_on_level = [&](
const auto &cell) {
580 if (vertices_owned_by_locally_owned_cells_on_level
581 [cell->vertex_index(v)])
586 auto &level_cell_infos = construction_data.cell_infos[
level];
587 for (
const auto &cell :
tria.cell_iterators_on_level(
level))
590 if (!cell_marked[cell->level()][cell->index()])
596 cell_info.id = cell->id().template to_binary<dim>();
599 for (
const auto f : cell->face_indices())
602 cell->face(f)->boundary_id();
604 cell_info.boundary_ids.emplace_back(f, boundary_ind);
614 for (
const auto line : cell->line_indices())
615 cell_info.manifold_line_ids[line] =
620 for (
const auto f : cell->face_indices())
621 cell_info.manifold_quad_ids[f] =
627 cell_info.level_subdomain_id =
630 if (is_locally_relevant_on_active_level(cell))
634 cell_info.level_subdomain_id =
637 else if (is_locally_relevant_on_level(cell))
639 cell_info.level_subdomain_id =
647 level_cell_infos.emplace_back(cell_info);
651 return construction_data;
656 set_additional_data(Description<dim, spacedim> &construction_data)
const
658 construction_data.comm =
comm;
660 construction_data.settings =
settings;
664 set_additional_data(DescriptionTemp<dim, spacedim> &)
const
670 add_vertices(Description<dim, spacedim> &construction_data,
671 const std::vector<bool> &vertices_locally_relevant)
const
673 std::vector<unsigned int> vertices_locally_relevant_indices(
674 vertices_locally_relevant.size());
677 unsigned int vertex_counter = 0;
678 for (
unsigned int i = 0; i < vertices_locally_relevant.size(); ++i)
679 if (vertices_locally_relevant[i])
681 construction_data.coarse_cell_vertices.push_back(
683 vertices_locally_relevant_indices[i] = vertex_counter++;
688 for (unsigned
int v = 0; v < cell.vertices.size(); ++v)
690 vertices_locally_relevant_indices[cell.vertices[v]];
694 add_vertices(DescriptionTemp<dim, spacedim> &construction_data,
695 const std::vector<bool> &vertices_locally_relevant)
const
697 for (
unsigned int i = 0; i < vertices_locally_relevant.size(); ++i)
698 if (vertices_locally_relevant[i])
699 construction_data.coarse_cell_vertices.emplace_back(
704 const ::Triangulation<dim, spacedim> &
tria;
706 const typename ::Triangulation<dim, spacedim>::cell_iterator &)>
709 const typename ::Triangulation<dim, spacedim>::cell_iterator &)>
716 std::map<unsigned int, std::vector<unsigned int>>
724 template <
int dim,
int spacedim>
725 Description<dim, spacedim>
727 const ::Triangulation<dim, spacedim> &
tria,
730 const unsigned int my_rank_in)
732 if (
const auto tria_pdt =
dynamic_cast<
735 ExcMessage(
"MPI communicators do not match."));
739 unsigned int my_rank = my_rank_in;
742 ExcMessage(
"Rank has to be smaller than available processes."));
744 if (
auto tria_pdt =
dynamic_cast<
750 "For creation from a parallel::distributed::Triangulation, "
751 "my_rank has to equal global rank."));
755 else if (
auto tria_serial =
756 dynamic_cast<const ::Triangulation<dim, spacedim> *
>(
765 ExcMessage(
"This type of triangulation is not supported!"));
769 return cell->subdomain_id();
773 return cell->level_subdomain_id();
776 return CreateDescriptionFromTriangulationHelper<dim, spacedim>(
782 .template create_description_for_rank<Description<dim, spacedim>>(
788 template <
int dim,
int spacedim>
792 & serial_grid_generator,
795 const unsigned int)> &serial_grid_partitioner,
797 const int group_size,
801#ifndef DEAL_II_WITH_MPI
802 (void)serial_grid_generator;
803 (void)serial_grid_partitioner;
811 const unsigned int my_rank =
813 const unsigned int group_root = (my_rank / group_size) * group_size;
819 if (my_rank == group_root)
826 typename ::Triangulation<dim, spacedim>::MeshSmoothing
>(
829 spacedim>::limit_level_difference_at_vertices) :
831 serial_grid_generator(
tria);
834 serial_grid_partitioner(
tria,
comm, group_size);
841 const unsigned int end_group =
849 for (
unsigned int other_rank = group_root + 1; other_rank < end_group;
853 const auto construction_data =
859 std::vector<char> buffer;
863 const auto ierr = MPI_Send(buffer.data(),
885 auto ierr = MPI_Probe(group_root, mpi_tag,
comm, &status);
889 MPI_Get_count(&status, MPI_CHAR, &len);
891 std::vector<char> buf(len);
892 ierr = MPI_Recv(buf.data(),
903 auto construction_data =
904 ::Utilities::template unpack<Description<dim, spacedim>>(
909 construction_data.comm =
comm;
911 return construction_data;
918 template <
int dim,
int spacedim>
926 (partition.size() > 0) &&
935 "Source triangulation has to be set up with "
936 "limit_level_difference_at_vertices if the construction of the "
937 "multigrid hierarchy is requested!"));
939 std::vector<LinearAlgebra::distributed::Vector<double>> partitions_mg;
943 const auto tria_parallel =
949 partition.update_ghost_values();
954 partitions_mg[l].reinit(
955 tria_parallel->global_level_cell_index_partitioner(l).lock());
961 if (cell->is_locally_owned_on_level() ==
false)
964 if (cell->is_active())
965 partitions_mg[
level][cell->global_level_cell_index()] =
966 partition[cell->global_active_cell_index()];
968 partitions_mg[
level][cell->global_level_cell_index()] =
969 partitions_mg[
level + 1]
970 [cell->child(0)->global_level_cell_index()];
973 partitions_mg[
level].update_ghost_values();
985 template <
int dim,
int spacedim>
994#ifdef DEAL_II_WITH_MPI
999 if (partition.size() == 0)
1007 partition.update_ghost_values();
1009 for (
const auto &partition : partitions_mg)
1010 partition.update_ghost_values();
1013 const std::vector<unsigned int> relevant_processes = [&]() {
1014 std::set<unsigned int> relevant_processes;
1016 for (
unsigned int i = 0; i < partition.locally_owned_size(); ++i)
1017 relevant_processes.insert(
1018 static_cast<unsigned int>(partition.local_element(i)));
1020 for (
const auto &partition : partitions_mg)
1021 for (
unsigned int i = 0; i < partition.locally_owned_size(); ++i)
1022 relevant_processes.insert(
1023 static_cast<unsigned int>(partition.local_element(i)));
1025 return std::vector<unsigned int>(relevant_processes.begin(),
1026 relevant_processes.end());
1039 if ((cell->is_active() && (cell->is_artificial() ==
false)))
1040 return static_cast<unsigned int>(
1041 partition[cell->global_active_cell_index()]);
1049 return static_cast<unsigned int>(
1050 partitions_mg[cell->level()][cell->global_level_cell_index()]);
1055 CreateDescriptionFromTriangulationHelper<dim, spacedim> helper(
1064 std::vector<DescriptionTemp<dim, spacedim>> descriptions_per_rank;
1065 descriptions_per_rank.reserve(relevant_processes.size());
1067 for (
const auto rank : relevant_processes)
1068 descriptions_per_rank.emplace_back(
1069 helper.template create_description_for_rank<
1070 DescriptionTemp<dim, spacedim>>(rank));
1074 DescriptionTemp<dim, spacedim> description_merged;
1075 description_merged.collect(
1077 descriptions_per_rank,
1078 partition.get_mpi_communicator(),
1084 description_merged.reduce();
1087 return description_merged.convert(partition.get_mpi_communicator(),
1098#include "tria_description.inst"
virtual MPI_Comm get_communicator() const
virtual const MeshSmoothing & get_mesh_smoothing() const
const std::vector< Point< spacedim > > & get_vertices() const
unsigned int n_levels() const
unsigned int n_raw_cells(const unsigned int level) const
virtual unsigned int n_global_levels() const
Triangulation< dim, spacedim > & get_triangulation()
unsigned int n_vertices() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
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)
Description< dim, spacedim > create_description_from_triangulation_in_groups(const std::function< void(::Triangulation< dim, spacedim > &)> &serial_grid_generator, const std::function< void(::Triangulation< dim, spacedim > &, const MPI_Comm, const unsigned int)> &serial_grid_partitioner, const MPI_Comm comm, const int group_size=1, const typename Triangulation< dim, spacedim >::MeshSmoothing smoothing=::Triangulation< dim, spacedim >::none, const TriangulationDescription::Settings setting=TriangulationDescription::Settings::default_setting)
@ construct_multigrid_hierarchy
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::vector< unsigned int > selector(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
DEAL_II_HOST constexpr TableIndices< 2 > merge(const TableIndices< 2 > &previous_indices, const unsigned int new_index, const unsigned int position)
const types::boundary_id internal_face_boundary_id
const types::subdomain_id artificial_subdomain_id
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
global_cell_index coarse_cell_id
unsigned int subdomain_id
std::vector< unsigned int > vertices
bool operator==(const CellData< structdim > &other) const
types::manifold_id manifold_id
types::material_id material_id
types::boundary_id boundary_id
CellData(const unsigned int n_vertices=ReferenceCells::get_hypercube< structdim >().n_vertices())
std::vector< CellData< 2 > > boundary_quads
bool check_consistency(const unsigned int dim) const
std::vector< CellData< 1 > > boundary_lines
std::vector< std::pair< unsigned int, Point< spacedim > > > coarse_cell_vertices
const TriangulationDescription::Settings settings
const bool construct_multigrid
std::map< unsigned int, unsigned int > vertex_to_coinciding_vertex_group
const std::function< types::subdomain_id(const typename ::Triangulation< dim, spacedim >::cell_iterator &)> subdomain_id_function
const std::function< types::subdomain_id(const typename ::Triangulation< dim, spacedim >::cell_iterator &)> level_subdomain_id_function
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id
std::vector< std::vector< CellData< dim > > > cell_infos
std::map< unsigned int, std::vector< unsigned int > > coinciding_vertex_groups
const ::Triangulation< dim, spacedim > & tria
std::vector<::CellData< dim > > coarse_cells