19#include <deal.II/base/mpi.templates.h>
46 template <
int dim,
int spacedim>
52 const
bool check_for_distorted_cells)
54 check_for_distorted_cells)
55 , mpi_communicator(mpi_communicator)
56 , my_subdomain(
Utilities::MPI::this_mpi_process(this->mpi_communicator))
57 , n_subdomains(
Utilities::MPI::n_mpi_processes(this->mpi_communicator))
59#ifndef DEAL_II_WITH_MPI
66 template <
int dim,
int spacedim>
71#ifndef DEAL_II_WITH_MPI
77 if (const ::parallel::TriangulationBase<dim, spacedim> *other_tria_x =
78 dynamic_cast<const ::parallel::TriangulationBase<dim, spacedim>
91 template <
int dim,
int spacedim>
100 number_cache.n_global_active_cells) +
107 template <
int dim,
int spacedim>
119 template <
int dim,
int spacedim>
122 : n_locally_owned_active_cells(0)
123 , number_of_global_coarse_cells(0)
129 template <
int dim,
int spacedim>
139 template <
int dim,
int spacedim>
148 template <
int dim,
int spacedim>
153 return number_cache.n_global_active_cells;
158 template <
int dim,
int spacedim>
162 return mpi_communicator;
167#ifdef DEAL_II_WITH_MPI
168 template <
int dim,
int spacedim>
173 number_cache.level_ghost_owners.clear();
174 number_cache.n_locally_owned_active_cells = 0;
176 if (this->n_levels() == 0)
183 number_cache.n_global_active_cells = 0;
184 number_cache.n_global_levels = 0;
191 for (
const auto &cell : this->active_cell_iterators())
192 if (cell->is_ghost())
193 number_cache.ghost_owners.insert(cell->subdomain_id());
195 Assert(number_cache.ghost_owners.size() <
200 if (this->n_levels() > 0)
201 number_cache.n_locally_owned_active_cells = std::count_if(
202 this->begin_active(),
205 [](
const auto &i) {
return i.is_locally_owned(); });
207 number_cache.n_locally_owned_active_cells = 0;
211 number_cache.n_global_active_cells =
213 number_cache.n_locally_owned_active_cells),
214 this->mpi_communicator);
216 number_cache.n_global_levels =
221 if (this->is_multilevel_hierarchy_constructed() ==
true)
223 number_cache.level_ghost_owners.clear();
226 if (this->n_levels() == 0)
230 for (
const auto &cell : this->cell_iterators())
231 if (cell->level_subdomain_id() !=
numbers::artificial_subdomain_id &&
232 cell->level_subdomain_id() != this->locally_owned_subdomain())
233 this->number_cache.level_ghost_owners.insert(
234 cell->level_subdomain_id());
240 int ierr = MPI_Barrier(this->mpi_communicator);
247 std::vector<MPI_Request> requests(
248 this->number_cache.level_ghost_owners.size());
249 unsigned int dummy = 0;
250 unsigned int req_counter = 0;
252 for (
const auto &it : this->number_cache.level_ghost_owners)
254 ierr = MPI_Isend(&dummy,
259 this->mpi_communicator,
260 &requests[req_counter]);
265 for (
const auto &it : this->number_cache.level_ghost_owners)
268 ierr = MPI_Recv(&dummy,
273 this->mpi_communicator,
278 if (requests.size() > 0)
280 ierr = MPI_Waitall(requests.size(),
282 MPI_STATUSES_IGNORE);
286 ierr = MPI_Barrier(this->mpi_communicator);
291 Assert(this->number_cache.level_ghost_owners.size() <
296 this->number_cache.number_of_global_coarse_cells = this->n_cells(0);
299 this->reset_global_cell_indices();
304 template <
int dim,
int spacedim>
313 template <
int dim,
int spacedim>
322 std::vector<unsigned int> reference_cells_ui;
324 reference_cells_ui.reserve(this->reference_cells.size());
325 for (
const auto &i : this->reference_cells)
326 reference_cells_ui.push_back(
static_cast<unsigned int>(i));
331 this->mpi_communicator);
334 this->reference_cells.clear();
335 for (
const auto &i : reference_cells_ui)
336 this->reference_cells.emplace_back(
342 template <
int dim,
int spacedim>
352 template <
int dim,
int spacedim>
354 const std::set<types::subdomain_id>
357 return number_cache.ghost_owners;
362 template <
int dim,
int spacedim>
364 const std::set<types::subdomain_id>
367 return number_cache.level_ghost_owners;
372 template <
int dim,
int spacedim>
375 get_boundary_ids()
const
379 this->mpi_communicator);
384 template <
int dim,
int spacedim>
387 get_manifold_ids()
const
391 this->mpi_communicator);
396 template <
int dim,
int spacedim>
400#ifndef DEAL_II_WITH_MPI
406 if (pst->with_artificial_cells() ==
false)
412 std::vector<unsigned int> cell_counter(n_subdomains + 1);
415 for (
const auto &cell : this->active_cell_iterators())
416 cell_counter[cell->subdomain_id() + 1]++;
419 for (
unsigned int i = 0; i < n_subdomains; ++i)
420 cell_counter[i + 1] += cell_counter[i];
425 IndexSet is_local(this->n_active_cells());
426 is_local.
add_range(cell_counter[my_subdomain],
427 cell_counter[my_subdomain + 1]);
428 number_cache.active_cell_index_partitioner =
429 std::make_shared<const Utilities::MPI::Partitioner>(
432 this->mpi_communicator);
435 for (
const auto &cell : this->active_cell_iterators())
436 cell->set_global_active_cell_index(
437 cell_counter[cell->subdomain_id()]++);
439 Assert(this->is_multilevel_hierarchy_constructed() ==
false,
447 this->n_locally_owned_active_cells();
452 const int ierr = MPI_Exscan(
453 &n_locally_owned_cells,
458 this->mpi_communicator);
463 std::pair<types::global_cell_index, types::global_cell_index> my_range;
466 for (
const auto &cell : this->active_cell_iterators())
467 if (cell->is_locally_owned())
468 cell->set_global_active_cell_index(
cell_index++);
475 std::vector<types::global_dof_index> is_ghost_vector;
476 GridTools::exchange_cell_data_to_ghosts<types::global_cell_index>(
478 [](
const auto &cell) {
return cell->global_active_cell_index(); },
479 [&is_ghost_vector](
const auto &cell,
const auto &id) {
480 cell->set_global_active_cell_index(
id);
481 is_ghost_vector.push_back(
id);
485 IndexSet is_local(this->n_global_active_cells());
486 is_local.add_range(my_range.first, my_range.second);
488 std::sort(is_ghost_vector.begin(), is_ghost_vector.end());
489 IndexSet is_ghost(this->n_global_active_cells());
490 is_ghost.add_indices(is_ghost_vector.begin(), is_ghost_vector.end());
492 number_cache.active_cell_index_partitioner =
493 std::make_shared<const Utilities::MPI::Partitioner>(
494 is_local, is_ghost, this->mpi_communicator);
497 if (this->is_multilevel_hierarchy_constructed() ==
true)
500 std::vector<types::global_cell_index> n_cells_level(
501 this->n_global_levels(), 0);
503 for (
auto cell : this->cell_iterators())
504 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
505 n_cells_level[cell->
level()]++;
508 std::vector<types::global_cell_index>
cell_index(
509 this->n_global_levels(), 0);
511 int ierr = MPI_Exscan(
512 n_cells_level.data(),
514 this->n_global_levels(),
517 this->mpi_communicator);
522 this->mpi_communicator,
528 std::pair<types::global_cell_index, types::global_cell_index>>
529 my_ranges(this->n_global_levels());
530 for (
unsigned int l = 0; l < this->n_global_levels(); ++l)
533 for (
auto cell : this->cell_iterators())
534 if (cell->level_subdomain_id() == this->locally_owned_subdomain())
539 for (
unsigned int l = 0; l < this->n_global_levels(); ++l)
543 std::vector<std::vector<types::global_dof_index>> is_ghost_vectors(
544 this->n_global_levels());
549 [](
const auto &cell) {
return cell->global_level_cell_index(); },
550 [&is_ghost_vectors](
const auto &cell,
const auto &id) {
551 cell->set_global_level_cell_index(
id);
552 is_ghost_vectors[cell->level()].push_back(
id);
555 number_cache.level_cell_index_partitioners.resize(
556 this->n_global_levels());
559 for (
unsigned int l = 0;
l < this->n_global_levels(); ++
l)
561 IndexSet is_local(n_cells_level[l]);
562 is_local.add_range(my_ranges[l].
first, my_ranges[l].
second);
564 IndexSet is_ghost(n_cells_level[l]);
565 std::sort(is_ghost_vectors[l].
begin(), is_ghost_vectors[l].
end());
566 is_ghost.add_indices(is_ghost_vectors[l].
begin(),
567 is_ghost_vectors[l].
end());
569 number_cache.level_cell_index_partitioners[
l] =
570 std::make_shared<const Utilities::MPI::Partitioner>(
571 is_local, is_ghost, this->mpi_communicator);
580 template <
int dim,
int spacedim>
583 const
std::vector<
bool> &vertex_locally_moved)
588 const std::vector<bool> locally_owned_vertices =
590 for (
unsigned int i = 0; i < locally_owned_vertices.size(); ++i)
591 Assert((vertex_locally_moved[i] ==
false) ||
592 (locally_owned_vertices[i] ==
true),
593 ExcMessage(
"The vertex_locally_moved argument must not "
594 "contain vertices that are not locally owned"));
599 for (
unsigned int d = 0; d < spacedim; ++d)
600 invalid_point[d] = std::numeric_limits<double>::quiet_NaN();
602 const auto pack = [&](
const auto &cell) {
603 std::vector<Point<spacedim>>
vertices(cell->n_vertices());
605 for (
const auto v : cell->vertex_indices())
606 if (vertex_locally_moved[cell->vertex_index(v)])
614 const auto unpack = [&](
const auto &cell,
const auto &
vertices) {
615 for (
const auto v : cell->vertex_indices())
620 if (this->is_multilevel_hierarchy_constructed())
627 GridTools::exchange_cell_data_to_ghosts<std::vector<Point<spacedim>>>(
635 template <
int dim,
int spacedim>
639 spacedim>::global_active_cell_index_partitioner()
const
641 return number_cache.active_cell_index_partitioner;
646 template <
int dim,
int spacedim>
650 spacedim>::global_level_cell_index_partitioner(
const unsigned int level)
656 return number_cache.level_cell_index_partitioners[
level];
661 template <
int dim,
int spacedim>
666 return number_cache.number_of_global_coarse_cells;
671 template <
int dim,
int spacedim>
677 const
bool check_for_distorted_cells)
681 check_for_distorted_cells)
682 , cell_attached_data({0, 0, {}, {}})
683 , data_transfer(mpi_communicator)
688 template <
int dim,
int spacedim>
692 cell_attached_data = {0, 0, {}, {}};
693 data_transfer.
clear();
700 template <
int dim,
int spacedim>
703 const
unsigned int global_first_cell,
704 const
unsigned int global_num_cells,
705 const
std::
string &filename)
const
708 auto tria =
const_cast<
711 if (this->cell_attached_data.n_attached_data_sets > 0)
714 tria->data_transfer.pack_data(
715 tria->local_cell_relations,
716 tria->cell_attached_data.pack_callbacks_fixed,
717 tria->cell_attached_data.pack_callbacks_variable);
720 tria->data_transfer.
save(global_first_cell, global_num_cells, filename);
729 tria->cell_attached_data.n_attached_data_sets = 0;
730 tria->cell_attached_data.pack_callbacks_fixed.
clear();
731 tria->cell_attached_data.pack_callbacks_variable.
clear();
737 template <
int dim,
int spacedim>
741 if (this->n_global_levels() <= 1)
750 const bool have_coarser_cell =
751 std::any_of(this->begin_active(this->n_global_levels() - 2),
752 this->end_active(this->n_global_levels() - 2),
759 this->mpi_communicator) != 0;
764 template <
int dim,
int spacedim>
767 const
unsigned int global_first_cell,
768 const
unsigned int global_num_cells,
769 const
unsigned int local_num_cells,
770 const
std::
string &filename,
771 const
unsigned int n_attached_deserialize_fixed,
772 const
unsigned int n_attached_deserialize_variable)
775 if (this->cell_attached_data.n_attached_deserialize > 0)
777 this->data_transfer.
load(global_first_cell,
781 n_attached_deserialize_fixed,
782 n_attached_deserialize_variable);
784 this->data_transfer.unpack_cell_status(this->local_cell_relations);
787 for (
const auto &cell_rel : this->local_cell_relations)
793 spacedim>::CELL_PERSIST),
801 template <
int dim,
int spacedim>
804 register_data_attach(
807 const
bool returns_variable_size_data)
813 if (returns_variable_size_data)
815 handle = 2 * cell_attached_data.pack_callbacks_variable.size();
816 cell_attached_data.pack_callbacks_variable.push_back(pack_callback);
820 handle = 2 * cell_attached_data.pack_callbacks_fixed.size() + 1;
821 cell_attached_data.pack_callbacks_fixed.push_back(pack_callback);
825 ++cell_attached_data.n_attached_data_sets;
832 template <
int dim,
int spacedim>
835 const
unsigned int handle,
839 const
boost::iterator_range<
std::vector<
char>::const_iterator> &)>
843 data_transfer.unpack_data(local_cell_relations, handle, unpack_callback);
846 --cell_attached_data.n_attached_data_sets;
847 if (cell_attached_data.n_attached_deserialize > 0)
848 --cell_attached_data.n_attached_deserialize;
857 if (cell_attached_data.n_attached_data_sets == 0 &&
858 cell_attached_data.n_attached_deserialize == 0)
861 cell_attached_data.pack_callbacks_fixed.
clear();
862 cell_attached_data.pack_callbacks_variable.clear();
863 data_transfer.clear();
866 for (
auto &cell_rel : local_cell_relations)
877 template <
int dim,
int spacedim>
881 : variable_size_data_stored(false)
882 , mpi_communicator(mpi_communicator)
887 template <
int dim,
int spacedim>
892 &pack_callbacks_fixed,
894 &pack_callbacks_variable)
896 Assert(src_data_fixed.size() == 0,
897 ExcMessage(
"Previously packed data has not been released yet!"));
900 const unsigned int n_callbacks_fixed = pack_callbacks_fixed.size();
901 const unsigned int n_callbacks_variable = pack_callbacks_variable.size();
905 variable_size_data_stored = (n_callbacks_variable > 0);
911 std::vector<unsigned int> cell_sizes_variable_cumulative(
912 n_callbacks_variable);
926 std::vector<std::vector<std::vector<char>>> packed_fixed_size_data(
927 cell_relations.size());
928 std::vector<std::vector<std::vector<char>>> packed_variable_size_data(
929 variable_size_data_stored ? cell_relations.size() : 0);
937 auto cell_rel_it = cell_relations.cbegin();
938 auto data_cell_fixed_it = packed_fixed_size_data.begin();
939 auto data_cell_variable_it = packed_variable_size_data.begin();
940 for (; cell_rel_it != cell_relations.cend(); ++cell_rel_it)
942 const auto &dealii_cell = cell_rel_it->first;
943 const auto &cell_status = cell_rel_it->second;
964 for (
unsigned int c = 0;
965 c < GeometryInfo<dim>::max_children_per_cell;
967 Assert(dealii_cell->child(c)->is_active(),
989 const unsigned int n_fixed_size_data_sets_on_cell =
993 spacedim>::CELL_INVALID) ?
995 ((variable_size_data_stored ? 1 : 0) + n_callbacks_fixed));
996 data_cell_fixed_it->resize(n_fixed_size_data_sets_on_cell);
999 auto data_fixed_it = data_cell_fixed_it->begin();
1012 spacedim>::CELL_INVALID)
1015 for (
auto callback_it = pack_callbacks_fixed.cbegin();
1016 callback_it != pack_callbacks_fixed.cend();
1017 ++callback_it, ++data_fixed_it)
1019 *data_fixed_it = (*callback_it)(dealii_cell, cell_status);
1026 if (variable_size_data_stored)
1028 const unsigned int n_variable_size_data_sets_on_cell =
1033 n_callbacks_variable);
1034 data_cell_variable_it->resize(
1035 n_variable_size_data_sets_on_cell);
1037 auto callback_it = pack_callbacks_variable.cbegin();
1038 auto data_variable_it = data_cell_variable_it->begin();
1039 auto sizes_variable_it =
1040 cell_sizes_variable_cumulative.begin();
1041 for (; callback_it != pack_callbacks_variable.cend();
1042 ++callback_it, ++data_variable_it, ++sizes_variable_it)
1045 (*callback_it)(dealii_cell, cell_status);
1049 *sizes_variable_it = data_variable_it->size();
1053 std::partial_sum(cell_sizes_variable_cumulative.begin(),
1054 cell_sizes_variable_cumulative.end(),
1055 cell_sizes_variable_cumulative.begin());
1062 data_fixed_it->resize(n_callbacks_variable *
1063 sizeof(
unsigned int));
1064 for (
unsigned int i = 0; i < n_callbacks_variable; ++i)
1065 std::memcpy(&(data_fixed_it->at(i *
sizeof(
unsigned int))),
1066 &(cell_sizes_variable_cumulative.at(i)),
1067 sizeof(
unsigned int));
1074 Assert(data_fixed_it == data_cell_fixed_it->end(),
1078 ++data_cell_fixed_it;
1083 if (variable_size_data_stored)
1084 ++data_cell_variable_it;
1101 std::vector<unsigned int> local_sizes_fixed(
1102 1 + n_callbacks_fixed + (variable_size_data_stored ? 1 : 0));
1103 for (
const auto &data_cell : packed_fixed_size_data)
1105 if (data_cell.size() == local_sizes_fixed.size())
1107 auto sizes_fixed_it = local_sizes_fixed.begin();
1108 auto data_fixed_it = data_cell.cbegin();
1109 for (; data_fixed_it != data_cell.cend();
1110 ++data_fixed_it, ++sizes_fixed_it)
1112 *sizes_fixed_it = data_fixed_it->size();
1120 for (
auto data_cell_fixed_it = packed_fixed_size_data.cbegin();
1121 data_cell_fixed_it != packed_fixed_size_data.cend();
1122 ++data_cell_fixed_it)
1124 Assert((data_cell_fixed_it->size() == 1) ||
1125 (data_cell_fixed_it->size() == local_sizes_fixed.size()),
1132 std::vector<unsigned int> global_sizes_fixed(local_sizes_fixed.size());
1134 this->mpi_communicator,
1135 global_sizes_fixed);
1139 sizes_fixed_cumulative.resize(global_sizes_fixed.size());
1140 std::partial_sum(global_sizes_fixed.begin(),
1141 global_sizes_fixed.end(),
1142 sizes_fixed_cumulative.begin());
1147 if (variable_size_data_stored)
1149 src_sizes_variable.reserve(packed_variable_size_data.size());
1150 for (
const auto &data_cell : packed_variable_size_data)
1152 int variable_data_size_on_cell = 0;
1154 for (
const auto &data : data_cell)
1155 variable_data_size_on_cell += data.size();
1157 src_sizes_variable.push_back(variable_data_size_on_cell);
1164 const unsigned int expected_size_fixed =
1165 cell_relations.size() * sizes_fixed_cumulative.back();
1166 const unsigned int expected_size_variable =
1167 std::accumulate(src_sizes_variable.begin(),
1168 src_sizes_variable.end(),
1169 std::vector<int>::size_type(0));
1173 src_data_fixed.reserve(expected_size_fixed);
1174 for (
const auto &data_cell_fixed : packed_fixed_size_data)
1178 for (
const auto &data_fixed : data_cell_fixed)
1179 std::move(data_fixed.begin(),
1181 std::back_inserter(src_data_fixed));
1187 if ((data_cell_fixed.size() == 1) &&
1188 (sizes_fixed_cumulative.size() > 1))
1190 const std::size_t bytes_skipped =
1191 sizes_fixed_cumulative.back() - sizes_fixed_cumulative.front();
1193 src_data_fixed.insert(src_data_fixed.end(),
1195 static_cast<char>(-1));
1201 if (variable_size_data_stored)
1203 src_data_variable.reserve(expected_size_variable);
1204 for (
const auto &data_cell : packed_variable_size_data)
1208 for (
const auto &data : data_cell)
1209 std::move(data.begin(),
1211 std::back_inserter(src_data_variable));
1217 Assert(src_data_variable.size() == expected_size_variable,
1223 template <
int dim,
int spacedim>
1228 Assert(sizes_fixed_cumulative.size() > 0,
1230 if (cell_relations.size() > 0)
1232 Assert(dest_data_fixed.size() > 0,
1237 const unsigned int size = sizes_fixed_cumulative.front();
1243 auto cell_rel_it = cell_relations.begin();
1244 auto dest_fixed_it = dest_data_fixed.cbegin();
1245 for (; cell_rel_it != cell_relations.end();
1246 ++cell_rel_it, dest_fixed_it += sizes_fixed_cumulative.back())
1248 cell_rel_it->second =
1251 spacedim>::CellStatus>(dest_fixed_it,
1252 dest_fixed_it + size,
1259 template <
int dim,
int spacedim>
1263 const
unsigned int handle,
1264 const
std::function<
1267 const
boost::iterator_range<
std::vector<
char>::const_iterator> &)>
1268 &unpack_callback)
const
1274 const bool callback_variable_transfer = (handle % 2 == 0);
1275 const unsigned int callback_index = handle / 2;
1281 Assert(sizes_fixed_cumulative.size() > 0,
1283 if (cell_relations.size() > 0)
1285 Assert(dest_data_fixed.size() > 0,
1289 std::vector<char>::const_iterator dest_data_it;
1290 std::vector<char>::const_iterator dest_sizes_cell_it;
1300 if (callback_variable_transfer)
1313 const unsigned int offset_variable_data_sizes =
1314 sizes_fixed_cumulative[sizes_fixed_cumulative.size() - 2];
1321 dest_sizes_cell_it = dest_data_fixed.cbegin() +
1322 offset_variable_data_sizes +
1323 callback_index *
sizeof(
unsigned int);
1326 dest_data_it = dest_data_variable.cbegin();
1333 offset = sizes_fixed_cumulative[callback_index];
1334 size = sizes_fixed_cumulative[callback_index + 1] - offset;
1335 data_increment = sizes_fixed_cumulative.back();
1341 if (cell_relations.begin() != cell_relations.end())
1342 dest_data_it = dest_data_fixed.cbegin() + offset;
1346 auto cell_rel_it = cell_relations.begin();
1347 auto dest_sizes_it = dest_sizes_variable.cbegin();
1348 for (; cell_rel_it != cell_relations.end(); ++cell_rel_it)
1350 const auto &dealii_cell = cell_rel_it->first;
1351 const auto &cell_status = cell_rel_it->second;
1353 if (callback_variable_transfer)
1357 data_increment = *dest_sizes_it;
1361 spacedim>::CELL_INVALID)
1366 if (callback_index == 0)
1369 std::memcpy(&offset,
1370 &(*(dest_sizes_cell_it -
sizeof(
unsigned int))),
1371 sizeof(
unsigned int));
1374 &(*dest_sizes_cell_it),
1375 sizeof(
unsigned int));
1382 dest_data_it += offset;
1383 data_increment -= offset;
1388 if (cell_rel_it != cell_relations.end() - 1)
1389 dest_sizes_cell_it += sizes_fixed_cumulative.back();
1393 switch (cell_status)
1396 spacedim>::CELL_PERSIST:
1398 spacedim>::CELL_COARSEN:
1399 unpack_callback(dealii_cell,
1401 boost::make_iterator_range(dest_data_it,
1402 dest_data_it + size));
1406 spacedim>::CELL_REFINE:
1407 unpack_callback(dealii_cell->parent(),
1409 boost::make_iterator_range(dest_data_it,
1410 dest_data_it + size));
1414 spacedim>::CELL_INVALID:
1423 if (cell_rel_it != cell_relations.end() - 1)
1424 dest_data_it += data_increment;
1430 template <
int dim,
int spacedim>
1433 const
unsigned int global_first_cell,
1434 const
unsigned int global_num_cells,
1435 const
std::
string &filename)
const
1437#ifdef DEAL_II_WITH_MPI
1442 Assert(sizes_fixed_cumulative.size() > 0,
1447 const unsigned int bytes_per_cell = sizes_fixed_cumulative.back();
1453 const std::string fname_fixed = std::string(filename) +
"_fixed.data";
1456 int ierr = MPI_Info_create(&info);
1460 ierr = MPI_File_open(mpi_communicator,
1461 fname_fixed.c_str(),
1462 MPI_MODE_CREATE | MPI_MODE_WRONLY,
1467 ierr = MPI_File_set_size(fh, 0);
1471 ierr = MPI_Barrier(mpi_communicator);
1473 ierr = MPI_Info_free(&info);
1486 sizes_fixed_cumulative.data(),
1487 sizes_fixed_cumulative.size(),
1494 const MPI_Offset size_header =
1495 sizes_fixed_cumulative.size() *
sizeof(
unsigned int);
1499 const MPI_Offset my_global_file_position =
1501 static_cast<MPI_Offset
>(global_first_cell) * bytes_per_cell;
1505 my_global_file_position,
1506 src_data_fixed.data(),
1507 src_data_fixed.size(),
1512 ierr = MPI_File_close(&fh);
1521 if (variable_size_data_stored)
1523 const std::string fname_variable =
1524 std::string(filename) +
"_variable.data";
1527 int ierr = MPI_Info_create(&info);
1531 ierr = MPI_File_open(mpi_communicator,
1532 fname_variable.c_str(),
1533 MPI_MODE_CREATE | MPI_MODE_WRONLY,
1538 ierr = MPI_File_set_size(fh, 0);
1542 ierr = MPI_Barrier(mpi_communicator);
1544 ierr = MPI_Info_free(&info);
1549 const MPI_Offset my_global_file_position =
1550 static_cast<MPI_Offset
>(global_first_cell) *
sizeof(
unsigned int);
1555 static_cast<std::size_t
>(
1556 std::numeric_limits<int>::max()),
1561 my_global_file_position,
1562 src_sizes_variable.data(),
1563 src_sizes_variable.size(),
1572 const std::uint64_t size_on_proc = src_data_variable.size();
1573 std::uint64_t prefix_sum = 0;
1574 ierr = MPI_Exscan(&size_on_proc,
1582 const MPI_Offset my_global_file_position =
1583 static_cast<MPI_Offset
>(global_num_cells) *
sizeof(
unsigned int) +
1589 my_global_file_position,
1590 src_data_variable.data(),
1591 src_data_variable.size(),
1597 ierr = MPI_File_close(&fh);
1601 (void)global_first_cell;
1602 (void)global_num_cells;
1611 template <
int dim,
int spacedim>
1614 const
unsigned int global_first_cell,
1615 const
unsigned int global_num_cells,
1616 const
unsigned int local_num_cells,
1617 const
std::
string &filename,
1618 const
unsigned int n_attached_deserialize_fixed,
1619 const
unsigned int n_attached_deserialize_variable)
1621#ifdef DEAL_II_WITH_MPI
1626 Assert(dest_data_fixed.size() == 0,
1627 ExcMessage(
"Previously loaded data has not been released yet!"));
1629 variable_size_data_stored = (n_attached_deserialize_variable > 0);
1635 const std::string fname_fixed = std::string(filename) +
"_fixed.data";
1638 int ierr = MPI_Info_create(&info);
1642 ierr = MPI_File_open(
1643 mpi_communicator, fname_fixed.c_str(), MPI_MODE_RDONLY, info, &fh);
1646 ierr = MPI_Info_free(&info);
1653 sizes_fixed_cumulative.resize(1 + n_attached_deserialize_fixed +
1654 (variable_size_data_stored ? 1 : 0));
1658 sizes_fixed_cumulative.data(),
1659 sizes_fixed_cumulative.size(),
1665 const unsigned int bytes_per_cell = sizes_fixed_cumulative.back();
1666 dest_data_fixed.resize(
static_cast<size_t>(local_num_cells) *
1670 const MPI_Offset size_header =
1671 sizes_fixed_cumulative.size() *
sizeof(
unsigned int);
1675 const MPI_Offset my_global_file_position =
1677 static_cast<MPI_Offset
>(global_first_cell) * bytes_per_cell;
1680 my_global_file_position,
1681 dest_data_fixed.data(),
1682 dest_data_fixed.size(),
1688 ierr = MPI_File_close(&fh);
1695 if (variable_size_data_stored)
1697 const std::string fname_variable =
1698 std::string(filename) +
"_variable.data";
1701 int ierr = MPI_Info_create(&info);
1705 ierr = MPI_File_open(
1706 mpi_communicator, fname_variable.c_str(), MPI_MODE_RDONLY, info, &fh);
1709 ierr = MPI_Info_free(&info);
1713 dest_sizes_variable.resize(local_num_cells);
1715 const MPI_Offset my_global_file_position_sizes =
1716 static_cast<MPI_Offset
>(global_first_cell) *
sizeof(
unsigned int);
1720 my_global_file_position_sizes,
1721 dest_sizes_variable.data(),
1722 dest_sizes_variable.size(),
1730 const std::uint64_t size_on_proc =
1731 std::accumulate(dest_sizes_variable.begin(),
1732 dest_sizes_variable.end(),
1735 std::uint64_t prefix_sum = 0;
1736 ierr = MPI_Exscan(&size_on_proc,
1744 const MPI_Offset my_global_file_position =
1745 static_cast<MPI_Offset
>(global_num_cells) *
sizeof(
unsigned int) +
1748 dest_data_variable.resize(size_on_proc);
1752 my_global_file_position,
1753 dest_data_variable.data(),
1754 dest_data_variable.size(),
1759 ierr = MPI_File_close(&fh);
1763 (void)global_first_cell;
1764 (void)global_num_cells;
1765 (void)local_num_cells;
1767 (void)n_attached_deserialize_fixed;
1768 (void)n_attached_deserialize_variable;
1776 template <
int dim,
int spacedim>
1780 variable_size_data_stored =
false;
1783 sizes_fixed_cumulative.
clear();
1784 sizes_fixed_cumulative.shrink_to_fit();
1787 src_data_fixed.clear();
1788 src_data_fixed.shrink_to_fit();
1790 dest_data_fixed.clear();
1791 dest_data_fixed.shrink_to_fit();
1794 src_sizes_variable.clear();
1795 src_sizes_variable.shrink_to_fit();
1797 src_data_variable.clear();
1798 src_data_variable.shrink_to_fit();
1800 dest_sizes_variable.clear();
1801 dest_sizes_variable.shrink_to_fit();
1803 dest_data_variable.clear();
1804 dest_data_variable.shrink_to_fit();
1812#include "tria_base.inst"
bool is_locally_owned() const
void add_range(const size_type begin, const size_type end)
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
void save(Archive &ar, const unsigned int version) const
virtual std::size_t memory_consumption() const
virtual void update_reference_cells()
virtual void clear() override
virtual void load(const std::string &filename)=0
typename ::Triangulation< dim, spacedim >::CellStatus CellStatus
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
typename std::pair< cell_iterator, CellStatus > cell_relation_t
const std::set< types::subdomain_id > & level_ghost_owners() const
virtual types::global_cell_index n_global_active_cells() const override
const std::set< types::subdomain_id > & ghost_owners() const
types::subdomain_id locally_owned_subdomain() const override
virtual unsigned int n_global_levels() const override
unsigned int n_locally_owned_active_cells() const
virtual types::coarse_cell_id n_global_coarse_cells() const override
virtual void update_number_cache()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
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)
IndexSet complete_index_set(const IndexSet::size_type N)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
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)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
T max(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 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)
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
void release_all_unused_memory()
constexpr ReferenceCell make_reference_cell_from_int(const std::uint8_t kind)
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
bool is_nan(const double x)
unsigned int global_cell_index
const ::Triangulation< dim, spacedim > & tria