19 #include "libmesh/libmesh_common.h"
20 #include "libmesh/parallel.h"
28 #include "libmesh/libmesh_version.h"
29 #include "libmesh/system.h"
30 #include "libmesh/mesh_base.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/xdr_cxx.h"
33 #include "libmesh/numeric_vector.h"
34 #include "libmesh/dof_map.h"
35 #include "libmesh/auto_ptr.h"
63 const std::size_t max_io_blksize = 256000;
71 struct CompareDofObjectsByID
73 bool operator()(
const DofObject * a,
74 const DofObject * b)
const
79 return a->id() < b->id();
86 template <
typename InValType>
91 std::vector<InValType> & _data;
101 if (_data.empty())
return;
102 _io.
data_stream (_data.data(), cast_int<unsigned int>(_data.size()));
115 const std::string & version,
116 const bool read_header_in,
117 const bool read_additional_data,
118 const bool read_legacy_format)
159 const bool read_ifem_info =
160 (version.rfind(
" with infinite elements") < version.size()) ||
170 this->
comm().broadcast(nv);
175 for (
unsigned int var=0; var<nv; var++)
179 std::string var_name;
182 this->
comm().broadcast(var_name);
185 std::set<subdomain_id_type> domains;
186 if (io.
version() >= LIBMESH_VERSION_ID(0,7,2))
188 std::vector<subdomain_id_type> domain_array;
190 io.
data (domain_array);
191 for (
const auto &
id : domain_array)
194 this->
comm().broadcast(domains);
201 this->
comm().broadcast(order);
210 this->
comm().broadcast(rad_order);
217 this->
comm().broadcast(fam);
219 type.
order = static_cast<Order>(order);
220 type.
family = static_cast<FEFamily>(fam);
227 if (read_legacy_format)
233 libMesh::out <<
"*****************************************************************\n"
234 <<
"* WARNING: reading a potentially incompatible restart file!!! *\n"
236 <<
"*****************************************************************"
246 io.
data (radial_fam);
247 this->
comm().broadcast(radial_fam);
250 this->
comm().broadcast(i_map);
253 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
257 type.
inf_map = static_cast<InfMapType>(i_map);
275 unsigned int nvecs=0;
278 this->
comm().broadcast(nvecs);
285 for (
unsigned int vec=0; vec<nvecs; vec++)
289 std::string vec_name;
292 this->
comm().broadcast(vec_name);
294 if (read_additional_data)
308 #ifdef LIBMESH_ENABLE_DEPRECATED
310 const bool read_additional_data)
312 libmesh_deprecated();
328 std::vector<Number> global_vector;
329 std::vector<Number> reordered_vector;
335 io.
data (global_vector);
336 this->
comm().broadcast(global_vector);
341 reordered_vector.resize(global_vector.size());
346 libmesh_assert_equal_to (global_vector.size(), this->
n_dofs());
350 const unsigned int sys = this->
number();
351 const unsigned int nv = cast_int<unsigned int>
353 libmesh_assert_less_equal (nv, this->
n_vars());
355 for (
unsigned int data_var=0; data_var<nv; data_var++)
360 for (
auto & node : this->
get_mesh().node_ptr_range())
363 libmesh_assert_not_equal_to (node->dof_number(sys, var, index),
366 libmesh_assert_less (cnt, global_vector.size());
368 reordered_vector[node->dof_number(sys, var, index)] =
369 global_vector[cnt++];
373 for (
auto & elem : this->
get_mesh().active_element_ptr_range())
376 libmesh_assert_not_equal_to (elem->dof_number(sys, var, index),
379 libmesh_assert_less (cnt, global_vector.size());
381 reordered_vector[elem->dof_number(sys, var, index)] =
382 global_vector[cnt++];
386 *(this->
solution) = reordered_vector;
395 const std::size_t nvecs = this->
_vectors.size();
401 if (read_additional_data && nvecs &&
404 (
"Additional vectors in file do not match system");
406 std::map<std::string, NumericVector<Number> *>::iterator
416 std::fill (global_vector.begin(), global_vector.end(),
libMesh::zero);
419 io.
data (global_vector);
424 if (read_additional_data && nvecs)
426 this->
comm().broadcast(global_vector);
431 std::fill (reordered_vector.begin(),
432 reordered_vector.end(),
435 reordered_vector.resize(global_vector.size());
437 libmesh_assert_equal_to (global_vector.size(), this->
n_dofs());
441 const unsigned int sys = this->
number();
442 const unsigned int nv = cast_int<unsigned int>
444 libmesh_assert_less_equal (nv, this->
n_vars());
446 for (
unsigned int data_var=0; data_var<nv; data_var++)
450 for (
auto & node : this->
get_mesh().node_ptr_range())
453 libmesh_assert_not_equal_to (node->dof_number(sys, var, index),
456 libmesh_assert_less (cnt, global_vector.size());
458 reordered_vector[node->dof_number(sys, var, index)] =
459 global_vector[cnt++];
463 for (
auto & elem : this->
get_mesh().active_element_ptr_range())
466 libmesh_assert_not_equal_to (elem->dof_number(sys, var, index),
469 libmesh_assert_less (cnt, global_vector.size());
471 reordered_vector[elem->dof_number(sys, var, index)] =
472 global_vector[cnt++];
477 *(pos->second) = reordered_vector;
491 template <
typename InValType>
493 const bool read_additional_data)
529 std::vector<const DofObject *> ordered_nodes, ordered_elements;
531 std::set<const DofObject *, CompareDofObjectsByID>
532 ordered_nodes_set (this->
get_mesh().local_nodes_begin(),
533 this->
get_mesh().local_nodes_end());
535 ordered_nodes.insert(ordered_nodes.end(),
536 ordered_nodes_set.begin(),
537 ordered_nodes_set.end());
540 std::set<const DofObject *, CompareDofObjectsByID>
541 ordered_elements_set (this->
get_mesh().local_elements_begin(),
542 this->
get_mesh().local_elements_end());
544 ordered_elements.insert(ordered_elements.end(),
545 ordered_elements_set.begin(),
546 ordered_elements_set.end());
550 std::vector<InValType> io_buffer;
558 total_read_size += cast_int<dof_id_type>(io_buffer.size());
560 const unsigned int sys_num = this->
number();
561 const unsigned int nv = cast_int<unsigned int>
563 libmesh_assert_less_equal (nv, this->
n_vars());
568 for (
unsigned int data_var=0; data_var<nv; data_var++)
574 for (
const auto & node : ordered_nodes)
577 libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
579 libmesh_assert_less (cnt, io_buffer.size());
580 this->
solution->set(node->dof_number(sys_num, var, comp), io_buffer[cnt++]);
584 for (
const auto & elem : ordered_elements)
587 libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
589 libmesh_assert_less (cnt, io_buffer.size());
590 this->
solution->set(elem->dof_number(sys_num, var, comp), io_buffer[cnt++]);
596 for (
unsigned int data_var=0; data_var<nv; data_var++)
604 std::vector<dof_id_type> SCALAR_dofs;
607 for (
auto dof : SCALAR_dofs)
608 this->
solution->set(dof, io_buffer[cnt++]);
622 const std::size_t nvecs = this->
_vectors.size();
628 if (read_additional_data && nvecs &&
631 (
"Additional vectors in file do not match system");
633 std::map<std::string, NumericVector<Number> *>::const_iterator
647 total_read_size += cast_int<dof_id_type>(io_buffer.size());
652 if (read_additional_data && nvecs)
655 for (
unsigned int data_var=0; data_var<nv; data_var++)
661 for (
const auto & node : ordered_nodes)
664 libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
666 libmesh_assert_less (cnt, io_buffer.size());
667 pos->second->set(node->dof_number(sys_num, var, comp), io_buffer[cnt++]);
671 for (
const auto & elem : ordered_elements)
674 libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
676 libmesh_assert_less (cnt, io_buffer.size());
677 pos->second->set(elem->dof_number(sys_num, var, comp), io_buffer[cnt++]);
683 for (
unsigned int data_var=0; data_var<nv; data_var++)
691 std::vector<dof_id_type> SCALAR_dofs;
694 for (
auto dof : SCALAR_dofs)
695 pos->second->set(dof, io_buffer[cnt++]);
701 pos->second->close();
723 template <
typename InValType>
725 const bool read_additional_data)
738 parallel_object_only();
749 this->read_serialized_vector<InValType>(io, this->
solution.get());
761 const std::size_t nvecs = this->
_vectors.size();
767 if (read_additional_data && nvecs &&
770 (
"Additional vectors in file do not match system");
772 std::map<std::string, NumericVector<Number> *>::const_iterator
781 this->read_serialized_vector<InValType>
782 (io, (read_additional_data && nvecs) ? pos->second :
nullptr);
809 template <
typename iterator_type,
typename InValType>
811 const iterator_type begin,
812 const iterator_type
end,
816 const unsigned int var_to_read)
const
838 vars_to_read.clear();
839 vars_to_read.push_back(var_to_read);
844 num_vecs = cast_int<unsigned int>(vecs.size());
846 io_blksize = cast_int<dof_id_type>(std::min(max_io_blksize, static_cast<std::size_t>(n_objs))),
847 num_blks = cast_int<unsigned int>(std::ceil(static_cast<double>(n_objs)/
848 static_cast<double>(io_blksize)));
852 std::size_t n_read_values=0;
854 std::vector<std::vector<dof_id_type>> xfer_ids(num_blks);
855 std::vector<std::vector<Number>> recv_vals(num_blks);
856 std::vector<Parallel::Request>
857 id_requests(num_blks), val_requests(num_blks);
858 std::vector<Parallel::MessageTag>
859 id_tags(num_blks), val_tags(num_blks);
865 std::vector<std::size_t>
866 xfer_ids_size (num_blks,0),
867 recv_vals_size (num_blks,0);
870 for (iterator_type it=begin; it!=
end; ++it)
874 block =
id/io_blksize;
876 libmesh_assert_less (block, num_blks);
878 xfer_ids_size[block] += 2;
881 for (
const auto & var : vars_to_read)
882 n_comp_tot += (*it)->n_comp(sys_num, var);
884 recv_vals_size[block] += n_comp_tot*num_vecs;
889 std::vector<std::size_t> tot_vals_size(recv_vals_size);
890 this->
comm().sum (tot_vals_size);
902 first_object = blk*io_blksize,
903 last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
906 std::vector<dof_id_type> & ids (xfer_ids[blk]);
907 std::vector<Number> & vals (recv_vals[blk]);
911 ids.clear(); ids.reserve (xfer_ids_size[blk]);
912 vals.resize(recv_vals_size[blk]);
914 if (recv_vals_size[blk] != 0)
915 for (iterator_type it=begin; it!=
end; ++it)
916 if (((*it)->id() >= first_object) &&
917 ((*it)->id() < last_object))
919 ids.push_back((*it)->id());
921 unsigned int n_comp_tot=0;
923 for (
const auto & var : vars_to_read)
924 n_comp_tot += (*it)->n_comp(sys_num, var);
926 ids.push_back (n_comp_tot*num_vecs);
929 #ifdef LIBMESH_HAVE_MPI
930 id_tags[blk] = this->
comm().get_unique_tag(100*num_blks + blk);
931 val_tags[blk] = this->
comm().get_unique_tag(200*num_blks + blk);
934 this->
comm().send (0, ids, id_requests[blk], id_tags[blk]);
937 this->
comm().receive (0, vals, val_requests[blk], val_tags[blk]);
947 std::vector<std::vector<dof_id_type>> recv_ids (this->
n_processors());
948 std::vector<std::vector<Number>> send_vals (this->
n_processors());
949 std::vector<Parallel::Request> reply_requests (this->
n_processors());
950 std::vector<unsigned int> obj_val_offsets;
951 std::vector<Number> input_vals;
952 std::vector<InValType> input_vals_tmp;
959 first_object = blk*io_blksize,
960 last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
961 n_objects_blk = last_object - first_object;
969 input_vals.resize(tot_vals_size[blk]);
970 input_vals_tmp.resize(tot_vals_size[blk]);
973 ThreadedIO<InValType> threaded_io(io, input_vals_tmp);
979 obj_val_offsets.resize (n_objects_blk); std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
983 std::size_t n_vals_blk = 0;
989 #ifdef LIBMESH_HAVE_MPI
991 Parallel::Status id_status (this->
comm().probe (Parallel::any_source, id_tags[blk]));
992 std::vector<dof_id_type> & ids (recv_ids[id_status.source()]);
993 std::size_t & n_vals_proc (recv_vals_size[id_status.source()]);
994 this->
comm().receive (id_status.source(), ids, id_tags[blk]);
997 std::vector<dof_id_type> & ids (recv_ids[0]);
998 std::size_t & n_vals_proc (recv_vals_size[0]);
1006 for (std::size_t
idx=0, sz=ids.size();
idx<sz;
idx+=2)
1009 local_idx = ids[
idx+0]-first_object,
1010 n_vals_tot_allvecs = ids[
idx+1];
1012 libmesh_assert_less (local_idx, n_objects_blk);
1014 obj_val_offsets[local_idx] = n_vals_tot_allvecs;
1015 n_vals_proc += n_vals_tot_allvecs;
1019 n_vals_blk += n_vals_proc;
1026 std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
1027 obj_val_offsets.begin());
1029 libmesh_assert_equal_to (n_vals_blk, obj_val_offsets.back());
1030 libmesh_assert_equal_to (n_vals_blk, tot_vals_size[blk]);
1036 input_vals[i_val] = input_vals_tmp[i_val];
1038 n_read_values += input_vals.size();
1043 const std::vector<dof_id_type> & ids (recv_ids[proc]);
1044 std::vector<Number> & vals (send_vals[proc]);
1045 const std::size_t & n_vals_proc (recv_vals_size[proc]);
1047 vals.clear(); vals.reserve(n_vals_proc);
1049 for (std::size_t
idx=0, sz=ids.size();
idx<sz;
idx+=2)
1052 local_idx = ids[
idx+0]-first_object,
1053 n_vals_tot_allvecs = ids[
idx+1];
1055 std::vector<Number>::const_iterator in_vals(input_vals.begin());
1057 std::advance (in_vals, obj_val_offsets[local_idx-1]);
1059 for (
unsigned int val=0; val<n_vals_tot_allvecs; val++, ++in_vals)
1063 vals.push_back(*in_vals);
1067 #ifdef LIBMESH_HAVE_MPI
1069 this->
comm().send (proc, vals, reply_requests[proc], val_tags[blk]);
1071 recv_vals[blk] = vals;
1078 Parallel::wait (val_requests[blk]);
1080 const std::vector<Number> & vals (recv_vals[blk]);
1081 std::vector<Number>::const_iterator val_it(vals.begin());
1083 if (!recv_vals[blk].empty())
1084 for (iterator_type it=begin; it!=
end; ++it)
1085 if (((*it)->id() >= first_object) &&
1086 ((*it)->id() < last_object))
1088 for (
auto & vec : vecs)
1089 for (
const auto & var : vars_to_read)
1091 const unsigned int n_comp = (*it)->n_comp(sys_num, var);
1093 for (
unsigned int comp=0; comp<n_comp; comp++, ++val_it)
1095 const dof_id_type dof_index = (*it)->dof_number (sys_num, var, comp);
1099 libmesh_assert_greater_equal (dof_index, vec->first_local_index());
1100 libmesh_assert_less (dof_index, vec->last_local_index());
1102 vec->set (dof_index, *val_it);
1110 Parallel::wait(reply_requests);
1113 Parallel::wait(id_requests);
1115 return n_read_values;
1124 unsigned int n_assigned_vals = 0;
1128 std::vector<Number> input_buffer(n_SCALAR_dofs);
1130 io.
data_stream(input_buffer.data(), n_SCALAR_dofs);
1132 #ifdef LIBMESH_HAVE_MPI
1135 const Parallel::MessageTag val_tag = this->
comm().get_unique_tag();
1139 this->
comm().receive(0, input_buffer, val_tag);
1151 std::vector<dof_id_type> SCALAR_dofs;
1157 vec->
set (SCALAR_dofs[i], input_buffer[i]);
1162 return n_assigned_vals;
1166 template <
typename InValType>
1170 parallel_object_only();
1185 unsigned int vector_length=0;
1187 std::size_t n_assigned_vals=0;
1192 io.
data(vector_length,
"# vector length");
1193 this->
comm().broadcast(vector_length);
1195 const unsigned int nv = cast_int<unsigned int>
1201 libmesh_assert_less_equal (nv, this->
n_vars());
1204 if (io.
version() >= LIBMESH_VERSION_ID(0,7,4))
1212 this->
get_mesh().local_nodes_begin(),
1213 this->
get_mesh().local_nodes_end(),
1225 this->
get_mesh().local_elements_begin(),
1226 this->
get_mesh().local_elements_end(),
1236 for (
unsigned int data_var=0; data_var<nv; data_var++)
1247 this->
get_mesh().local_nodes_begin(),
1248 this->
get_mesh().local_nodes_end(),
1261 this->
get_mesh().local_elements_begin(),
1262 this->
get_mesh().local_elements_end(),
1273 for (
unsigned int data_var=0; data_var<nv; data_var++)
1289 this->
comm().sum (n_assigned_vals);
1290 libmesh_assert_equal_to (n_assigned_vals, vector_length);
1293 return vector_length;
1299 const std::string & ,
1300 const bool write_additional_data)
const
1343 std::string comment;
1351 comment =
"# No. of Variables in System \"";
1352 comment += this->
name();
1355 unsigned int nv = this->
n_vars();
1356 io.
data (nv, comment.c_str());
1366 comment =
"# Name, Variable No. ";
1367 std::sprintf(buf,
"%u", var);
1369 comment +=
", System \"";
1370 comment += this->
name();
1374 io.
data (var_name, comment.c_str());
1380 comment =
"# Subdomains, Variable \"";
1383 comment +=
"\", System \"";
1384 comment += this->
name();
1388 std::vector<subdomain_id_type> domain_array;
1389 domain_array.assign(domains.begin(), domains.end());
1390 io.
data (domain_array, comment.c_str());
1398 comment =
"# Approximation Order, Variable \"";
1401 comment +=
"\", System \"";
1402 comment += this->
name();
1405 int order = static_cast<int>(this->
variable_type(var).order);
1406 io.
data (order, comment.c_str());
1410 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1414 comment =
"# Radial Approximation Order, Variable \"";
1417 comment +=
"\", System \"";
1418 comment += this->
name();
1421 int rad_order = static_cast<int>(this->
variable_type(var).radial_order);
1422 io.
data (rad_order, comment.c_str());
1431 comment =
"# FE Family, Variable \"";
1434 comment +=
"\", System \"";
1435 comment += this->
name();
1439 int fam = static_cast<int>(type.
family);
1440 io.
data (fam, comment.c_str());
1442 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1444 comment =
"# Radial FE Family, Variable \"";
1447 comment +=
"\", System \"";
1448 comment += this->
name();
1452 io.
data (radial_fam, comment.c_str());
1454 comment =
"# Infinite Mapping Type, Variable \"";
1457 comment +=
"\", System \"";
1458 comment += this->
name();
1461 int i_map = static_cast<int>(type.
inf_map);
1462 io.
data (i_map, comment.c_str());
1474 comment =
"# No. of Additional Vectors, System \"";
1475 comment += this->
name();
1478 unsigned int nvecs = write_additional_data ? this->
n_vectors () : 0;
1479 io.
data (nvecs, comment.c_str());
1482 if (write_additional_data)
1489 comment =
"# Name of ";
1490 std::sprintf(buf,
"%d", cnt++);
1492 comment +=
"th vector";
1493 std::string vec_name = pr.first;
1495 io.
data (vec_name, comment.c_str());
1504 const bool write_additional_data)
const
1529 std::string comment;
1533 std::vector<Number> io_buffer; io_buffer.reserve(this->
solution->local_size());
1543 std::vector<const DofObject *> ordered_nodes, ordered_elements;
1545 std::set<const DofObject *, CompareDofObjectsByID>
1546 ordered_nodes_set (this->
get_mesh().local_nodes_begin(),
1547 this->
get_mesh().local_nodes_end());
1549 ordered_nodes.insert(ordered_nodes.end(),
1550 ordered_nodes_set.begin(),
1551 ordered_nodes_set.end());
1554 std::set<const DofObject *, CompareDofObjectsByID>
1555 ordered_elements_set (this->
get_mesh().local_elements_begin(),
1556 this->
get_mesh().local_elements_end());
1558 ordered_elements.insert(ordered_elements.end(),
1559 ordered_elements_set.begin(),
1560 ordered_elements_set.end());
1563 const unsigned int sys_num = this->
number();
1564 const unsigned int nv = this->
n_vars();
1567 for (
unsigned int var=0; var<nv; var++)
1571 for (
const auto & node : ordered_nodes)
1574 libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
1577 io_buffer.push_back((*this->
solution)(node->dof_number(sys_num, var, comp)));
1581 for (
const auto & elem : ordered_elements)
1584 libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
1587 io_buffer.push_back((*this->
solution)(elem->dof_number(sys_num, var, comp)));
1598 std::vector<dof_id_type> SCALAR_dofs;
1601 for (
auto dof : SCALAR_dofs)
1602 io_buffer.push_back((*this->
solution)(dof));
1613 comment =
"# System \"";
1614 comment += this->
name();
1615 comment +=
"\" Solution Vector";
1618 io.
data (io_buffer, comment.c_str());
1623 if (write_additional_data)
1628 io_buffer.reserve(pr.second->local_size());
1631 for (
unsigned int var=0; var<nv; var++)
1635 for (
const auto & node : ordered_nodes)
1638 libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
1641 io_buffer.push_back((*pr.second)(node->dof_number(sys_num, var, comp)));
1645 for (
const auto & elem : ordered_elements)
1648 libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
1651 io_buffer.push_back((*pr.second)(elem->dof_number(sys_num, var, comp)));
1662 std::vector<dof_id_type> SCALAR_dofs;
1665 for (
auto dof : SCALAR_dofs)
1666 io_buffer.push_back((*pr.second)(dof));
1677 comment =
"# System \"";
1678 comment += this->
name();
1679 comment +=
"\" Additional Vector \"";
1680 comment += pr.first;
1684 io.
data (io_buffer, comment.c_str());
1704 const bool write_additional_data)
const
1719 parallel_object_only();
1720 std::string comment;
1732 comment =
"# System \"";
1733 comment += this->
name();
1734 comment +=
"\" Solution Vector";
1740 if (write_additional_data)
1750 comment =
"# System \"";
1751 comment += this->
name();
1752 comment +=
"\" Additional Vector \"";
1753 comment += pair.first;
1810 template <
typename iterator_type>
1813 const iterator_type begin,
1814 const iterator_type
end,
1816 const unsigned int var_to_write)
const
1818 parallel_object_only();
1836 std::vector<unsigned int> vars_to_write(1, var_to_write);
1840 vars_to_write.clear(); vars_to_write.reserve(this->
n_vars());
1842 vars_to_write.push_back(var);
1845 const dof_id_type io_blksize = cast_int<dof_id_type>
1846 (std::min(max_io_blksize, static_cast<std::size_t>(n_objs)));
1849 sys_num = this->
number(),
1850 num_vecs = cast_int<unsigned int>(vecs.size()),
1851 num_blks = cast_int<unsigned int>(std::ceil(static_cast<double>(n_objs)/
1852 static_cast<double>(io_blksize)));
1859 std::size_t written_length=0;
1860 std::vector<std::vector<dof_id_type>> xfer_ids(num_blks);
1861 std::vector<std::vector<Number>> send_vals(num_blks);
1862 std::vector<Parallel::Request>
1863 id_requests(num_blks), val_requests(num_blks);
1864 std::vector<Parallel::MessageTag>
1865 id_tags(num_blks), val_tags(num_blks);
1871 std::vector<unsigned int>
1872 xfer_ids_size (num_blks,0),
1873 send_vals_size (num_blks,0);
1875 for (iterator_type it=begin; it!=
end; ++it)
1879 block =
id/io_blksize;
1881 libmesh_assert_less (block, num_blks);
1883 xfer_ids_size[block] += 2;
1885 unsigned int n_comp_tot=0;
1887 for (
const auto & var : vars_to_write)
1888 n_comp_tot += (*it)->n_comp(sys_num, var);
1890 send_vals_size[block] += n_comp_tot*num_vecs;
1897 for (
unsigned int blk=0; blk<num_blks; blk++)
1904 first_object = blk*io_blksize,
1905 last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
1908 std::vector<dof_id_type> & ids (xfer_ids[blk]);
1909 std::vector<Number> & vals (send_vals[blk]);
1913 ids.clear(); ids.reserve (xfer_ids_size[blk]);
1914 vals.clear(); vals.reserve (send_vals_size[blk]);
1916 if (send_vals_size[blk] != 0)
1917 for (iterator_type it=begin; it!=
end; ++it)
1918 if (((*it)->id() >= first_object) &&
1919 ((*it)->id() < last_object))
1921 ids.push_back((*it)->id());
1925 unsigned int n_comp_tot=0;
1927 for (
const auto & var : vars_to_write)
1928 n_comp_tot += (*it)->n_comp(sys_num, var);
1930 ids.push_back (n_comp_tot*num_vecs);
1934 for (
const auto & vec : vecs)
1935 for (
const auto & var : vars_to_write)
1937 const unsigned int n_comp = (*it)->n_comp(sys_num, var);
1939 for (
unsigned int comp=0; comp<n_comp; comp++)
1941 libmesh_assert_greater_equal ((*it)->dof_number(sys_num, var, comp), vec->first_local_index());
1942 libmesh_assert_less ((*it)->dof_number(sys_num, var, comp), vec->last_local_index());
1943 vals.push_back((*vec)((*it)->dof_number(sys_num, var, comp)));
1948 #ifdef LIBMESH_HAVE_MPI
1949 id_tags[blk] = this->
comm().get_unique_tag(100*num_blks + blk);
1950 val_tags[blk] = this->
comm().get_unique_tag(200*num_blks + blk);
1953 this->
comm().send (0, ids, id_requests[blk], id_tags[blk]);
1954 this->
comm().send (0, vals, val_requests[blk], val_tags[blk]);
1961 std::vector<std::vector<dof_id_type>> recv_ids (this->
n_processors());
1962 std::vector<std::vector<Number>> recv_vals (this->
n_processors());
1963 std::vector<unsigned int> obj_val_offsets;
1964 std::vector<Number> output_vals;
1967 ThreadedIO<Number> threaded_io(io, output_vals);
1968 std::unique_ptr<Threads::Thread> async_io;
1970 for (
unsigned int blk=0; blk<num_blks; blk++)
1975 first_object = cast_int<dof_id_type>(blk*io_blksize),
1976 last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
1977 n_objects_blk = last_object - first_object;
1982 obj_val_offsets.resize (n_objects_blk); std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
1984 std::size_t n_val_recvd_blk=0;
1989 #ifdef LIBMESH_HAVE_MPI
1991 Parallel::Status id_status (this->
comm().probe (Parallel::any_source, id_tags[blk]));
1992 std::vector<dof_id_type> & ids (recv_ids[id_status.source()]);
1993 this->
comm().receive (id_status.source(), ids, id_tags[blk]);
1995 std::vector<dof_id_type> & ids (recv_ids[0]);
1996 ids = xfer_ids[blk];
2001 for (std::size_t
idx=0, sz=ids.size();
idx<sz;
idx+=2)
2004 local_idx = ids[
idx+0]-first_object,
2005 n_vals_tot_allvecs = ids[
idx+1];
2007 libmesh_assert_less (local_idx, n_objects_blk);
2008 libmesh_assert_less (local_idx, obj_val_offsets.size());
2010 obj_val_offsets[local_idx] = n_vals_tot_allvecs;
2013 #ifdef LIBMESH_HAVE_MPI
2015 Parallel::Status val_status (this->
comm().probe (Parallel::any_source, val_tags[blk]));
2016 std::vector<Number> & vals (recv_vals[val_status.source()]);
2017 this->
comm().receive (val_status.source(), vals, val_tags[blk]);
2020 std::vector<Number> & vals (recv_vals[0]);
2021 vals = send_vals[blk];
2024 n_val_recvd_blk += vals.size();
2030 std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
2031 obj_val_offsets.begin());
2035 if (async_io.get()) async_io->join();
2039 output_vals.resize(n_val_recvd_blk);
2044 const std::vector<dof_id_type> & ids (recv_ids [proc]);
2045 const std::vector<Number> & vals(recv_vals[proc]);
2046 std::vector<Number>::const_iterator proc_vals(vals.begin());
2048 for (std::size_t
idx=0, sz=ids.size();
idx<sz;
idx+=2)
2051 local_idx = ids[
idx+0]-first_object,
2052 n_vals_tot_allvecs = ids[
idx+1];
2056 std::vector<Number>::iterator out_vals(output_vals.begin());
2058 std::advance (out_vals, obj_val_offsets[local_idx-1]);
2060 for (
unsigned int val=0; val<n_vals_tot_allvecs; val++, ++out_vals, ++proc_vals)
2064 *out_vals = *proc_vals;
2071 async_io = libmesh_make_unique<Threads::Thread>(threaded_io);
2072 written_length += output_vals.size();
2080 Parallel::wait(id_requests);
2081 Parallel::wait(val_requests);
2090 this->
comm().broadcast(written_length);
2092 return written_length;
2098 const unsigned int var,
2101 unsigned int written_length=0;
2102 std::vector<Number> vals;
2107 std::vector<dof_id_type> SCALAR_dofs;
2109 const unsigned int n_scalar_dofs = cast_int<unsigned int>
2110 (SCALAR_dofs.size());
2112 for (
unsigned int i=0; i<n_scalar_dofs; i++)
2114 vals.push_back( vec(SCALAR_dofs[i]) );
2118 #ifdef LIBMESH_HAVE_MPI
2121 const Parallel::MessageTag val_tag(1);
2132 this->
comm().send(0, vals, val_tag);
2141 const unsigned int vals_size =
2142 cast_int<unsigned int>(vals.size());
2144 written_length += vals_size;
2147 return written_length;
2155 parallel_object_only();
2166 written_length += cast_int<dof_id_type>
2169 this->
get_mesh().local_nodes_begin(),
2170 this->
get_mesh().local_nodes_end(),
2175 written_length += cast_int<dof_id_type>
2178 this->
get_mesh().local_elements_begin(),
2179 this->
get_mesh().local_elements_end(),
2192 libmesh_assert_equal_to (written_length, vec_length);
2194 return written_length;
2198 template <
typename InValType>
2202 parallel_object_only();
2220 unsigned int num_vecs=0;
2226 io.
data(vector_length);
2228 libmesh_assert_equal_to (num_vecs, vectors.size());
2232 libmesh_assert_not_equal_to (vectors[0], 0);
2233 libmesh_assert_equal_to (vectors[0]->size(), vector_length);
2246 std::size_t read_length = 0;
2252 this->
get_mesh().local_nodes_begin(),
2253 this->
get_mesh().local_nodes_end(),
2262 this->
get_mesh().local_elements_begin(),
2263 this->
get_mesh().local_elements_end(),
2274 libmesh_assert_not_equal_to (vec, 0);
2284 libmesh_assert_not_equal_to (vec, 0);
2296 parallel_object_only();
2305 std::size_t written_length = 0;
2310 n_vec = cast_int<unsigned int>(vectors.size());
2312 vec_size = vectors.empty() ? 0 : vectors[0]->size();
2314 io.
data(n_vec,
"# number of vectors");
2316 io.
data(vec_size,
"# vector length");
2324 this->
get_mesh().local_nodes_begin(),
2325 this->
get_mesh().local_nodes_end(),
2333 this->
get_mesh().local_elements_begin(),
2334 this->
get_mesh().local_elements_end(),
2343 libmesh_assert_not_equal_to (vec, 0);
2349 return written_length;
2355 template void System::read_parallel_data<Number> (
Xdr & io,
const bool read_additional_data);
2356 template void System::read_serialized_data<Number> (
Xdr & io,
const bool read_additional_data);
2358 template std::size_t System::read_serialized_vectors<Number> (
Xdr & io,
const std::vector<
NumericVector<Number> *> & vectors)
const;
2359 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
2360 template void System::read_parallel_data<Real> (
Xdr & io,
const bool read_additional_data);
2361 template void System::read_serialized_data<Real> (
Xdr & io,
const bool read_additional_data);
2363 template std::size_t System::read_serialized_vectors<Real> (
Xdr & io,
const std::vector<
NumericVector<Number> *> & vectors)
const;