19#include <deal.II/base/mpi.templates.h>
44 template <
int dim,
int spacedim>
47 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
49 const bool check_for_distorted_cells)
51 check_for_distorted_cells)
52 , mpi_communicator(mpi_communicator)
56#ifndef DEAL_II_WITH_MPI
63 template <
int dim,
int spacedim>
66 const ::Triangulation<dim, spacedim> &other_tria)
68#ifndef DEAL_II_WITH_MPI
74 if (const ::parallel::TriangulationBase<dim, spacedim> *other_tria_x =
75 dynamic_cast<const ::parallel::TriangulationBase<dim, spacedim>
88 template <
int dim,
int spacedim>
97 number_cache.n_global_active_cells) +
102 template <
int dim,
int spacedim>
111 template <
int dim,
int spacedim>
113 : n_locally_owned_active_cells(0)
117 template <
int dim,
int spacedim>
124 template <
int dim,
int spacedim>
131 template <
int dim,
int spacedim>
138 template <
int dim,
int spacedim>
145#ifdef DEAL_II_WITH_MPI
146 template <
int dim,
int spacedim>
170 if (cell->is_ghost())
187 this->mpi_communicator);
209 cell->level_subdomain_id());
222 std::vector<MPI_Request> requests(
224 unsigned int dummy = 0;
225 unsigned int req_counter = 0;
227 for (std::set<types::subdomain_id>::iterator it =
229 it != this->number_cache.level_ghost_owners.end();
232 ierr = MPI_Isend(&dummy,
238 &requests[req_counter]);
242 for (std::set<types::subdomain_id>::iterator it =
244 it != this->number_cache.level_ghost_owners.end();
248 ierr = MPI_Recv(&dummy,
258 if (requests.size() > 0)
260 ierr = MPI_Waitall(requests.size(),
262 MPI_STATUSES_IGNORE);
282 template <
int dim,
int spacedim>
291 template <
int dim,
int spacedim>
300 std::vector<unsigned int> reference_cells_ui;
303 reference_cells_ui.push_back(
static_cast<unsigned int>(i));
311 this->reference_cells.clear();
312 for (
const auto &i : reference_cells_ui)
313 this->reference_cells.emplace_back(
318 template <
int dim,
int spacedim>
327 template <
int dim,
int spacedim>
328 const std::set<types::subdomain_id> &
336 template <
int dim,
int spacedim>
337 const std::set<types::subdomain_id> &
345 template <
int dim,
int spacedim>
346 std::map<unsigned int, std::set<::types::subdomain_id>>
355 template <
int dim,
int spacedim>
356 std::vector<types::boundary_id>
366 template <
int dim,
int spacedim>
367 std::vector<types::manifold_id>
377 template <
int dim,
int spacedim>
381#ifndef DEAL_II_WITH_MPI
387 if (pst->with_artificial_cells() ==
false)
393 std::vector<unsigned int> cell_counter(
n_subdomains + 1);
397 cell_counter[cell->subdomain_id() + 1]++;
401 cell_counter[i + 1] += cell_counter[i];
410 std::make_shared<const Utilities::MPI::Partitioner>(
417 cell->set_global_active_cell_index(
418 cell_counter[cell->subdomain_id()]++);
434 MPI_Exscan(&n_locally_owned_cells,
437 Utilities::MPI::internal::mpi_type_id(&n_locally_owned_cells),
445 if (cell->is_locally_owned())
446 cell->set_global_active_cell_index(
cell_index++);
451 GridTools::exchange_cell_data_to_ghosts<types::global_cell_index>(
453 [](
const auto &cell) {
return cell->global_active_cell_index(); },
454 [](
const auto &cell,
const auto &id) {
455 cell->set_global_active_cell_index(
id);
459 std::vector<types::global_dof_index> is_local_vector;
460 std::vector<types::global_dof_index> is_ghost_vector;
463 if (!cell->is_artificial())
465 const auto index = cell->global_active_cell_index();
470 if (cell->is_locally_owned())
471 is_local_vector.push_back(index);
473 is_ghost_vector.push_back(index);
476 std::sort(is_local_vector.begin(), is_local_vector.end());
478 is_local.add_indices(is_local_vector.begin(), is_local_vector.end());
480 std::sort(is_ghost_vector.begin(), is_ghost_vector.end());
482 is_ghost.add_indices(is_ghost_vector.begin(), is_ghost_vector.end());
485 std::make_shared<const Utilities::MPI::Partitioner>(
492 std::vector<types::global_cell_index> n_locally_owned_cells(
496 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
497 n_locally_owned_cells[cell->level()]++;
500 std::vector<types::global_cell_index>
cell_index(
503 int ierr = MPI_Exscan(n_locally_owned_cells.data(),
505 this->n_global_levels(),
506 Utilities::MPI::internal::mpi_type_id(
507 n_locally_owned_cells.data()),
509 this->mpi_communicator);
513 std::vector<types::global_cell_index> n_cells_level(
520 MPI_Bcast(n_cells_level.data(),
521 this->n_global_levels(),
522 Utilities::MPI::internal::mpi_type_id(n_cells_level.data()),
523 this->n_subdomains - 1,
524 this->mpi_communicator);
530 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
531 cell->set_global_level_cell_index(
cell_index[cell->level()]++);
540 [](
const auto &cell) {
return cell->global_level_cell_index(); },
541 [](
const auto &cell,
const auto &id) {
542 return cell->set_global_level_cell_index(
id);
551 std::vector<types::global_dof_index> is_local_vector;
552 std::vector<types::global_dof_index> is_ghost_vector;
555 if (cell->level_subdomain_id() !=
558 const auto index = cell->global_level_cell_index();
563 if (cell->level_subdomain_id() ==
565 is_local_vector.push_back(index);
567 is_ghost_vector.push_back(index);
572 std::sort(is_local_vector.begin(), is_local_vector.end());
573 is_local.add_indices(is_local_vector.begin(),
574 is_local_vector.end());
577 std::sort(is_ghost_vector.begin(), is_ghost_vector.end());
578 is_ghost.add_indices(is_ghost_vector.begin(),
579 is_ghost_vector.end());
582 std::make_shared<const Utilities::MPI::Partitioner>(
592 template <
int dim,
int spacedim>
595 const std::vector<bool> &vertex_locally_moved)
600 const std::vector<bool> locally_owned_vertices =
602 for (
unsigned int i = 0; i < locally_owned_vertices.size(); ++i)
603 Assert((vertex_locally_moved[i] ==
false) ||
604 (locally_owned_vertices[i] ==
true),
605 ExcMessage(
"The vertex_locally_moved argument must not "
606 "contain vertices that are not locally owned"));
611 for (
int d = 0;
d < spacedim; ++
d)
612 invalid_point[
d] = std::numeric_limits<double>::quiet_NaN();
614 const auto pack = [&](
const auto &cell) {
615 std::vector<Point<spacedim>>
vertices(cell->n_vertices());
617 for (
const auto v : cell->vertex_indices())
618 if (vertex_locally_moved[cell->vertex_index(v)])
627 for (
const auto v : cell->vertex_indices())
628 if (std::isnan(
vertices[v][0]) ==
false)
636 GridTools::exchange_cell_data_to_ghosts<std::vector<Point<spacedim>>>(
642 template <
int dim,
int spacedim>
643 const std::weak_ptr<const Utilities::MPI::Partitioner>
649 template <
int dim,
int spacedim>
650 const std::weak_ptr<const Utilities::MPI::Partitioner>
652 const unsigned int level)
const
662 template <
int dim,
int spacedim>
665 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
672 , cell_attached_data({0, 0, {}, {}})
673 , data_transfer(mpi_communicator)
678 template <
int dim,
int spacedim>
682 cell_attached_data = {0, 0, {}, {}};
690 template <
int dim,
int spacedim>
693 const unsigned int global_first_cell,
694 const unsigned int global_num_cells,
695 const std::string &filename)
const
698 auto tria =
const_cast<
701 if (this->cell_attached_data.n_attached_data_sets > 0)
705 tria->local_cell_relations,
706 tria->cell_attached_data.pack_callbacks_fixed,
707 tria->cell_attached_data.pack_callbacks_variable);
710 tria->data_transfer.save(global_first_cell, global_num_cells, filename);
713 tria->data_transfer.clear();
719 tria->cell_attached_data.n_attached_data_sets = 0;
720 tria->cell_attached_data.pack_callbacks_fixed.clear();
721 tria->cell_attached_data.pack_callbacks_variable.clear();
726 template <
int dim,
int spacedim>
729 const unsigned int global_first_cell,
730 const unsigned int global_num_cells,
731 const unsigned int local_num_cells,
732 const std::string &filename,
733 const unsigned int n_attached_deserialize_fixed,
734 const unsigned int n_attached_deserialize_variable)
737 if (this->cell_attached_data.n_attached_deserialize > 0)
739 this->data_transfer.
load(global_first_cell,
743 n_attached_deserialize_fixed,
744 n_attached_deserialize_variable);
746 this->data_transfer.unpack_cell_status(this->local_cell_relations);
749 for (
const auto &cell_rel : this->local_cell_relations)
755 spacedim>::CELL_PERSIST),
761 template <
int dim,
int spacedim>
766 const bool returns_variable_size_data)
772 if (returns_variable_size_data)
774 handle = 2 * cell_attached_data.pack_callbacks_variable.size();
775 cell_attached_data.pack_callbacks_variable.push_back(pack_callback);
779 handle = 2 * cell_attached_data.pack_callbacks_fixed.size() + 1;
780 cell_attached_data.pack_callbacks_fixed.push_back(pack_callback);
784 ++cell_attached_data.n_attached_data_sets;
790 template <
int dim,
int spacedim>
793 const unsigned int handle,
797 const boost::iterator_range<std::vector<char>::const_iterator> &)>
801 data_transfer.unpack_data(local_cell_relations, handle, unpack_callback);
804 --cell_attached_data.n_attached_data_sets;
805 if (cell_attached_data.n_attached_deserialize > 0)
806 --cell_attached_data.n_attached_deserialize;
814 if (cell_attached_data.n_attached_data_sets == 0 &&
815 cell_attached_data.n_attached_deserialize == 0)
818 cell_attached_data.pack_callbacks_fixed.
clear();
819 cell_attached_data.pack_callbacks_variable.clear();
822 if (this->cell_attached_data.n_attached_data_sets == 0 &&
823 this->cell_attached_data.n_attached_deserialize == 0)
826 for (
auto &cell_rel : local_cell_relations)
837 template <
int dim,
int spacedim>
840 : variable_size_data_stored(false)
841 , mpi_communicator(mpi_communicator)
846 template <
int dim,
int spacedim>
849 const std::vector<cell_relation_t> &cell_relations,
850 const std::vector<typename CellAttachedData::pack_callback_t>
851 &pack_callbacks_fixed,
852 const std::vector<typename CellAttachedData::pack_callback_t>
853 &pack_callbacks_variable)
855 Assert(src_data_fixed.size() == 0,
856 ExcMessage(
"Previously packed data has not been released yet!"));
859 const unsigned int n_callbacks_fixed = pack_callbacks_fixed.size();
860 const unsigned int n_callbacks_variable = pack_callbacks_variable.size();
864 variable_size_data_stored = (n_callbacks_variable > 0);
870 std::vector<unsigned int> cell_sizes_variable_cumulative(
871 n_callbacks_variable);
885 std::vector<std::vector<std::vector<char>>> packed_fixed_size_data(
886 cell_relations.size());
887 std::vector<std::vector<std::vector<char>>> packed_variable_size_data(
888 variable_size_data_stored ? cell_relations.size() : 0);
896 auto cell_rel_it = cell_relations.cbegin();
897 auto data_cell_fixed_it = packed_fixed_size_data.begin();
898 auto data_cell_variable_it = packed_variable_size_data.begin();
899 for (; cell_rel_it != cell_relations.cend(); ++cell_rel_it)
901 const auto &dealii_cell = cell_rel_it->first;
902 const auto &cell_status = cell_rel_it->second;
923 for (
unsigned int c = 0;
924 c < GeometryInfo<dim>::max_children_per_cell;
926 Assert(dealii_cell->child(c)->is_active(),
948 const unsigned int n_fixed_size_data_sets_on_cell =
954 ((variable_size_data_stored ? 1 : 0) + n_callbacks_fixed));
955 data_cell_fixed_it->resize(n_fixed_size_data_sets_on_cell);
958 auto data_fixed_it = data_cell_fixed_it->begin();
974 for (
auto callback_it = pack_callbacks_fixed.cbegin();
975 callback_it != pack_callbacks_fixed.cend();
976 ++callback_it, ++data_fixed_it)
978 *data_fixed_it = (*callback_it)(dealii_cell, cell_status);
985 if (variable_size_data_stored)
987 const unsigned int n_variable_size_data_sets_on_cell =
992 n_callbacks_variable);
993 data_cell_variable_it->resize(
994 n_variable_size_data_sets_on_cell);
996 auto callback_it = pack_callbacks_variable.cbegin();
997 auto data_variable_it = data_cell_variable_it->begin();
998 auto sizes_variable_it =
999 cell_sizes_variable_cumulative.begin();
1000 for (; callback_it != pack_callbacks_variable.cend();
1001 ++callback_it, ++data_variable_it, ++sizes_variable_it)
1004 (*callback_it)(dealii_cell, cell_status);
1008 *sizes_variable_it = data_variable_it->size();
1012 std::partial_sum(cell_sizes_variable_cumulative.begin(),
1013 cell_sizes_variable_cumulative.end(),
1014 cell_sizes_variable_cumulative.begin());
1020 data_fixed_it->resize(n_callbacks_variable *
1021 sizeof(
unsigned int));
1022 for (
unsigned int i = 0; i < n_callbacks_variable; ++i)
1023 std::memcpy(&(data_fixed_it->at(i *
sizeof(
unsigned int))),
1024 &(cell_sizes_variable_cumulative.at(i)),
1025 sizeof(
unsigned int));
1032 Assert(data_fixed_it == data_cell_fixed_it->end(),
1036 ++data_cell_fixed_it;
1041 if (variable_size_data_stored)
1042 ++data_cell_variable_it;
1059 std::vector<unsigned int> local_sizes_fixed(
1060 1 + n_callbacks_fixed + (variable_size_data_stored ? 1 : 0));
1061 for (
const auto &data_cell : packed_fixed_size_data)
1063 if (data_cell.size() == local_sizes_fixed.size())
1065 auto sizes_fixed_it = local_sizes_fixed.begin();
1066 auto data_fixed_it = data_cell.cbegin();
1067 for (; data_fixed_it != data_cell.cend();
1068 ++data_fixed_it, ++sizes_fixed_it)
1070 *sizes_fixed_it = data_fixed_it->size();
1078 for (
auto data_cell_fixed_it = packed_fixed_size_data.cbegin();
1079 data_cell_fixed_it != packed_fixed_size_data.cend();
1080 ++data_cell_fixed_it)
1082 Assert((data_cell_fixed_it->size() == 1) ||
1083 (data_cell_fixed_it->size() == local_sizes_fixed.size()),
1090 std::vector<unsigned int> global_sizes_fixed(local_sizes_fixed.size());
1093 global_sizes_fixed);
1097 sizes_fixed_cumulative.resize(global_sizes_fixed.size());
1098 std::partial_sum(global_sizes_fixed.begin(),
1099 global_sizes_fixed.end(),
1100 sizes_fixed_cumulative.begin());
1105 if (variable_size_data_stored)
1107 src_sizes_variable.reserve(packed_variable_size_data.size());
1108 for (
const auto &data_cell : packed_variable_size_data)
1110 int variable_data_size_on_cell = 0;
1112 for (
const auto &data : data_cell)
1113 variable_data_size_on_cell += data.size();
1115 src_sizes_variable.push_back(variable_data_size_on_cell);
1122 const unsigned int expected_size_fixed =
1123 cell_relations.size() * sizes_fixed_cumulative.back();
1124 const unsigned int expected_size_variable =
1125 std::accumulate(src_sizes_variable.begin(),
1126 src_sizes_variable.end(),
1130 src_data_fixed.reserve(expected_size_fixed);
1131 for (
const auto &data_cell_fixed : packed_fixed_size_data)
1135 for (
const auto &data_fixed : data_cell_fixed)
1136 std::move(data_fixed.begin(),
1138 std::back_inserter(src_data_fixed));
1144 if ((data_cell_fixed.size() == 1) &&
1145 (sizes_fixed_cumulative.size() > 1))
1147 const std::size_t bytes_skipped =
1148 sizes_fixed_cumulative.back() - sizes_fixed_cumulative.front();
1150 src_data_fixed.insert(src_data_fixed.end(),
1152 static_cast<char>(-1));
1158 if (variable_size_data_stored)
1160 src_data_variable.reserve(expected_size_variable);
1161 for (
const auto &data_cell : packed_variable_size_data)
1165 for (
const auto &data : data_cell)
1166 std::move(data.begin(),
1168 std::back_inserter(src_data_variable));
1174 Assert(src_data_variable.size() == expected_size_variable,
1180 template <
int dim,
int spacedim>
1183 std::vector<cell_relation_t> &cell_relations)
const
1185 Assert(sizes_fixed_cumulative.size() > 0,
1187 if (cell_relations.size() > 0)
1189 Assert(dest_data_fixed.size() > 0,
1194 const unsigned int size = sizes_fixed_cumulative.front();
1200 auto cell_rel_it = cell_relations.begin();
1201 auto dest_fixed_it = dest_data_fixed.cbegin();
1202 for (; cell_rel_it != cell_relations.end();
1203 ++cell_rel_it, dest_fixed_it += sizes_fixed_cumulative.back())
1205 cell_rel_it->second =
1209 dest_fixed_it + size,
1216 template <
int dim,
int spacedim>
1219 const std::vector<cell_relation_t> &cell_relations,
1220 const unsigned int handle,
1221 const std::function<
1222 void(
const typename ::Triangulation<dim, spacedim>::cell_iterator &,
1223 const typename ::Triangulation<dim, spacedim>::CellStatus &,
1224 const boost::iterator_range<std::vector<char>::const_iterator> &)>
1225 &unpack_callback)
const
1231 const bool callback_variable_transfer = (handle % 2 == 0);
1232 const unsigned int callback_index = handle / 2;
1238 Assert(sizes_fixed_cumulative.size() > 0,
1240 if (cell_relations.size() > 0)
1242 Assert(dest_data_fixed.size() > 0,
1246 std::vector<char>::const_iterator dest_data_it;
1247 std::vector<char>::const_iterator dest_sizes_cell_it;
1257 if (callback_variable_transfer)
1270 const unsigned int offset_variable_data_sizes =
1271 sizes_fixed_cumulative[sizes_fixed_cumulative.size() - 2];
1278 dest_sizes_cell_it = dest_data_fixed.cbegin() +
1279 offset_variable_data_sizes +
1280 callback_index *
sizeof(
unsigned int);
1283 dest_data_it = dest_data_variable.cbegin();
1290 offset = sizes_fixed_cumulative[callback_index];
1291 size = sizes_fixed_cumulative[callback_index + 1] - offset;
1292 data_increment = sizes_fixed_cumulative.back();
1298 dest_data_it = dest_data_fixed.cbegin() + offset;
1302 auto cell_rel_it = cell_relations.begin();
1303 auto dest_sizes_it = dest_sizes_variable.cbegin();
1304 for (; cell_rel_it != cell_relations.end(); ++cell_rel_it)
1306 const auto &dealii_cell = cell_rel_it->first;
1307 const auto &cell_status = cell_rel_it->second;
1309 if (callback_variable_transfer)
1313 data_increment = *dest_sizes_it;
1321 if (callback_index == 0)
1324 std::memcpy(&offset,
1325 &(*(dest_sizes_cell_it -
sizeof(
unsigned int))),
1326 sizeof(
unsigned int));
1329 &(*dest_sizes_cell_it),
1330 sizeof(
unsigned int));
1337 dest_data_it += offset;
1338 data_increment -= offset;
1342 dest_sizes_cell_it += sizes_fixed_cumulative.back();
1346 switch (cell_status)
1352 unpack_callback(dealii_cell,
1355 dest_data_it + size));
1360 unpack_callback(dealii_cell->parent(),
1363 dest_data_it + size));
1376 dest_data_it += data_increment;
1382 template <
int dim,
int spacedim>
1385 const unsigned int global_first_cell,
1386 const unsigned int global_num_cells,
1387 const std::string &filename)
const
1389#ifdef DEAL_II_WITH_MPI
1394 Assert(sizes_fixed_cumulative.size() > 0,
1403 const std::string fname_fixed = std::string(filename) +
"_fixed.data";
1406 int ierr = MPI_Info_create(&info);
1412 MPI_MODE_CREATE | MPI_MODE_WRONLY,
1417 ierr = MPI_File_set_size(fh, 0);
1423 ierr = MPI_Info_free(&info);
1432 const unsigned int *data = sizes_fixed_cumulative.data();
1434 ierr = MPI_File_write_at(fh,
1437 sizes_fixed_cumulative.size(),
1444 const MPI_Offset size_header =
1445 sizes_fixed_cumulative.size() *
sizeof(
unsigned int);
1449 const MPI_Offset my_global_file_position =
1450 size_header +
static_cast<MPI_Offset
>(global_first_cell) *
1451 sizes_fixed_cumulative.back();
1453 const char *data = src_data_fixed.data();
1455 ierr = MPI_File_write_at(fh,
1456 my_global_file_position,
1458 src_data_fixed.size(),
1463 ierr = MPI_File_close(&fh);
1470 if (variable_size_data_stored)
1472 const std::string fname_variable =
1473 std::string(filename) +
"_variable.data";
1476 int ierr = MPI_Info_create(&info);
1482 MPI_MODE_CREATE | MPI_MODE_WRONLY,
1487 ierr = MPI_File_set_size(fh, 0);
1493 ierr = MPI_Info_free(&info);
1498 const int *data = src_sizes_variable.data();
1500 MPI_File_write_at(fh,
1502 sizeof(
unsigned int),
1504 src_sizes_variable.size(),
1511 const unsigned int offset_variable =
1512 global_num_cells *
sizeof(
unsigned int);
1515 const unsigned int size_on_proc = src_data_variable.size();
1518 unsigned int prefix_sum = 0;
1527 const char *data = src_data_variable.data();
1530 ierr = MPI_File_write_at(fh,
1534 src_data_variable.size(),
1539 ierr = MPI_File_close(&fh);
1543 (void)global_first_cell;
1544 (void)global_num_cells;
1553 template <
int dim,
int spacedim>
1556 const unsigned int global_first_cell,
1557 const unsigned int global_num_cells,
1558 const unsigned int local_num_cells,
1559 const std::string &filename,
1560 const unsigned int n_attached_deserialize_fixed,
1561 const unsigned int n_attached_deserialize_variable)
1563#ifdef DEAL_II_WITH_MPI
1568 Assert(dest_data_fixed.size() == 0,
1569 ExcMessage(
"Previously loaded data has not been released yet!"));
1571 variable_size_data_stored = (n_attached_deserialize_variable > 0);
1577 const std::string fname_fixed = std::string(filename) +
"_fixed.data";
1580 int ierr = MPI_Info_create(&info);
1591 ierr = MPI_Info_free(&info);
1598 sizes_fixed_cumulative.resize(1 + n_attached_deserialize_fixed +
1599 (variable_size_data_stored ? 1 : 0));
1600 ierr = MPI_File_read_at(fh,
1602 sizes_fixed_cumulative.data(),
1603 sizes_fixed_cumulative.size(),
1609 dest_data_fixed.resize(local_num_cells * sizes_fixed_cumulative.back());
1612 const MPI_Offset size_header =
1613 sizes_fixed_cumulative.size() *
sizeof(
unsigned int);
1617 const MPI_Offset my_global_file_position =
1618 size_header +
static_cast<MPI_Offset
>(global_first_cell) *
1619 sizes_fixed_cumulative.back();
1621 ierr = MPI_File_read_at(fh,
1622 my_global_file_position,
1623 dest_data_fixed.data(),
1624 dest_data_fixed.size(),
1629 ierr = MPI_File_close(&fh);
1636 if (variable_size_data_stored)
1638 const std::string fname_variable =
1639 std::string(filename) +
"_variable.data";
1642 int ierr = MPI_Info_create(&info);
1653 ierr = MPI_Info_free(&info);
1657 dest_sizes_variable.resize(local_num_cells);
1658 ierr = MPI_File_read_at(fh,
1659 global_first_cell *
sizeof(
unsigned int),
1660 dest_sizes_variable.data(),
1661 dest_sizes_variable.size(),
1666 const unsigned int offset = global_num_cells *
sizeof(
unsigned int);
1668 const unsigned int size_on_proc =
1669 std::accumulate(dest_sizes_variable.begin(),
1670 dest_sizes_variable.end(),
1674 unsigned int prefix_sum = 0;
1683 dest_data_variable.resize(size_on_proc);
1684 ierr = MPI_File_read_at(fh,
1685 offset + prefix_sum,
1686 dest_data_variable.data(),
1687 dest_data_variable.size(),
1692 ierr = MPI_File_close(&fh);
1696 (void)global_first_cell;
1697 (void)global_num_cells;
1698 (void)local_num_cells;
1700 (void)n_attached_deserialize_fixed;
1701 (void)n_attached_deserialize_variable;
1709 template <
int dim,
int spacedim>
1713 variable_size_data_stored =
false;
1716 sizes_fixed_cumulative.clear();
1717 sizes_fixed_cumulative.shrink_to_fit();
1720 src_data_fixed.clear();
1721 src_data_fixed.shrink_to_fit();
1723 dest_data_fixed.clear();
1724 dest_data_fixed.shrink_to_fit();
1727 src_sizes_variable.clear();
1728 src_sizes_variable.shrink_to_fit();
1730 src_data_variable.clear();
1731 src_data_variable.shrink_to_fit();
1733 dest_sizes_variable.clear();
1734 dest_sizes_variable.shrink_to_fit();
1736 dest_data_variable.clear();
1737 dest_data_variable.shrink_to_fit();
1745#include "tria_base.inst"
void add_range(const size_type begin, const size_type end)
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
cell_iterator begin(const unsigned int level=0) const
virtual std::size_t memory_consumption() const
std::vector< Point< spacedim > > vertices
const bool check_for_distorted_cells
unsigned int n_active_cells() const
virtual void update_reference_cells()
std::vector< ReferenceCell > reference_cells
unsigned int n_levels() const
cell_iterator end() const
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)
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
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)
virtual void clear() override
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)
typename ::Triangulation< dim, spacedim >::CellStatus CellStatus
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
virtual void load(const std::string &filename, const bool autopartition=true)=0
DataTransfer data_transfer
const std::set< types::subdomain_id > & level_ghost_owners() const
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
void release_all_unused_memory()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors(const Triangulation< dim, spacedim > &tria)
void exchange_cell_data_to_level_ghosts(const MeshType &mesh, const std::function< std_cxx17::optional< DataType >(const typename MeshType::level_cell_iterator &)> &pack, const std::function< void(const typename MeshType::level_cell_iterator &, const DataType &)> &unpack, const std::function< bool(const typename MeshType::level_cell_iterator &)> &cell_filter=always_return< typename MeshType::level_cell_iterator, bool >{ true})
IndexSet complete_index_set(const IndexSet::size_type N)
IteratorRange< BaseIterator > make_iterator_range(const BaseIterator &begin, const typename identity< BaseIterator >::type &end)
#define DEAL_II_MPI_CONST_CAST(expr)
types::global_dof_index size_type
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
T sum(const T &t, const MPI_Comm &mpi_communicator)
std::vector< T > compute_set_union(const std::vector< T > &vec, const MPI_Comm &comm)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
T max(const T &t, const MPI_Comm &mpi_communicator)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
constexpr ::ReferenceCell make_reference_cell_from_int(const std::uint8_t kind)
const types::subdomain_id artificial_subdomain_id
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
unsigned int global_cell_index
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