system_io.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #include "libmesh/libmesh_common.h"
20 #include "libmesh/parallel.h"
21 
22 // C++ Includes
23 #include <cstdio> // for std::sprintf
24 #include <set>
25 #include <numeric> // for std::partial_sum
26 
27 // Local Include
29 #include "libmesh/system.h"
30 #include "libmesh/mesh_base.h"
31 //#include "libmesh/mesh_tools.h"
32 #include "libmesh/elem.h"
33 #include "libmesh/xdr_cxx.h"
34 #include "libmesh/numeric_vector.h"
35 #include "libmesh/dof_map.h"
36 
37 
38 
39 // Anonymous namespace for implementation details.
40 namespace {
41 
42 using libMesh::DofObject;
43 using libMesh::Number;
44 using libMesh::cast_int;
45 
46 // Comments:
47 // ---------
48 // - The max_io_blksize governs how many nodes or elements will be
49 // treated as a single block when performing parallel IO on large
50 // systems.
51 // - This parameter only loosely affects the size of the actual IO
52 // buffer as this depends on the number of components a given
53 // variable has for the nodes/elements in the block.
54 // - When reading/writing each processor uses an ID map which is
55 // 3*io_blksize*sizeof(dof_id_type) bytes long, so with unsigned int
56 // and // io_blksize=256000 we would expect that buffer alone to be
57 // ~3Mb.
58 // - In general, an increase in max_io_blksize should increase the
59 // efficiency of large parallel read/writes by reducing the number
60 // of MPI messages at the expense of memory.
61 // - If the library exhausts memory during IO you might reduce this
62 // parameter.
63 
64 const std::size_t max_io_blksize = 256000;
65 
66 
72 struct CompareDofObjectsByID
73 {
74  bool operator()(const DofObject * a,
75  const DofObject * b) const
76  {
77  libmesh_assert (a);
78  libmesh_assert (b);
79 
80  return a->id() < b->id();
81  }
82 };
83 
87 template <typename InValType>
88 class ThreadedIO
89 {
90 private:
91  libMesh::Xdr & _io;
92  std::vector<InValType> & _data;
93 
94 public:
95  ThreadedIO (libMesh::Xdr & io, std::vector<InValType> & data) :
96  _io(io),
97  _data(data)
98  {}
99 
100  void operator()()
101  {
102  if (_data.empty()) return;
103  _io.data_stream (_data.data(), cast_int<unsigned int>(_data.size()));
104  }
105 };
106 }
107 
108 
109 namespace libMesh
110 {
111 
112 
113 // ------------------------------------------------------------
114 // System class implementation
116  const std::string & version,
117  const bool read_header_in,
118  const bool read_additional_data,
119  const bool read_legacy_format)
120 {
121  // This method implements the input of a
122  // System object, embedded in the output of
123  // an EquationSystems<T_sys>. This warrants some
124  // documentation. The output file essentially
125  // consists of 5 sections:
126  //
127  // for this system
128  //
129  // 5.) The number of variables in the system (unsigned int)
130  //
131  // for each variable in the system
132  //
133  // 6.) The name of the variable (string)
134  //
135  // 6.1.) Variable subdomains
136  //
137  // 7.) Combined in an FEType:
138  // - The approximation order(s) of the variable
139  // (Order Enum, cast to int/s)
140  // - The finite element family/ies of the variable
141  // (FEFamily Enum, cast to int/s)
142  //
143  // end variable loop
144  //
145  // 8.) The number of additional vectors (unsigned int),
146  //
147  // for each additional vector in the system object
148  //
149  // 9.) the name of the additional vector (string)
150  //
151  // end system
152  libmesh_assert (io.reading());
153 
154  // Possibly clear data structures and start from scratch.
155  if (read_header_in)
156  this->clear ();
157 
158  // Figure out if we need to read infinite element information.
159  // This will be true if the version string contains " with infinite elements"
160  const bool read_ifem_info =
161  (version.rfind(" with infinite elements") < version.size()) ||
162  libMesh::on_command_line ("--read-ifem-systems");
163 
164 
165  {
166  // 5.)
167  // Read the number of variables in the system
168  unsigned int nv=0;
169  if (this->processor_id() == 0)
170  io.data (nv);
171  this->comm().broadcast(nv);
172 
173  _written_var_indices.clear();
174  _written_var_indices.resize(nv, 0);
175 
176  for (unsigned int var=0; var<nv; var++)
177  {
178  // 6.)
179  // Read the name of the var-th variable
180  std::string var_name;
181  if (this->processor_id() == 0)
182  io.data (var_name);
183  this->comm().broadcast(var_name);
184 
185  // 6.1.)
186  std::set<subdomain_id_type> domains;
187  if (io.version() >= LIBMESH_VERSION_ID(0,7,2))
188  {
189  std::vector<subdomain_id_type> domain_array;
190  if (this->processor_id() == 0)
191  io.data (domain_array);
192  for (const auto & id : domain_array)
193  domains.insert(id);
194  }
195  this->comm().broadcast(domains);
196 
197  // 7.)
198  // Read the approximation order(s) of the var-th variable
199  int order=0;
200  if (this->processor_id() == 0)
201  io.data (order);
202  this->comm().broadcast(order);
203 
204 
205  // do the same for infinite element radial_order
206  int rad_order=0;
207  if (read_ifem_info)
208  {
209  if (this->processor_id() == 0)
210  io.data(rad_order);
211  this->comm().broadcast(rad_order);
212  }
213 
214  // Read the finite element type of the var-th variable
215  int fam=0;
216  if (this->processor_id() == 0)
217  io.data (fam);
218  this->comm().broadcast(fam);
219  FEType type;
220  type.order = static_cast<Order>(order);
221  type.family = static_cast<FEFamily>(fam);
222 
223  // Check for incompatibilities. The shape function indexing was
224  // changed for the monomial and xyz finite element families to
225  // simplify extension to arbitrary p. The consequence is that
226  // old restart files will not be read correctly. This is expected
227  // to be an unlikely occurence, but catch it anyway.
228  if (read_legacy_format)
229  if ((type.family == MONOMIAL || type.family == XYZ) &&
230  ((type.order.get_order() > 2 && this->get_mesh().mesh_dimension() == 2) ||
231  (type.order.get_order() > 1 && this->get_mesh().mesh_dimension() == 3)))
232  {
233  libmesh_here();
234  libMesh::out << "*****************************************************************\n"
235  << "* WARNING: reading a potentially incompatible restart file!!! *\n"
236  << "* contact libmesh-users@lists.sourceforge.net for more details *\n"
237  << "*****************************************************************"
238  << std::endl;
239  }
240 
241  // Read additional information for infinite elements
242  int radial_fam=0;
243  int i_map=0;
244  if (read_ifem_info)
245  {
246  if (this->processor_id() == 0)
247  io.data (radial_fam);
248  this->comm().broadcast(radial_fam);
249  if (this->processor_id() == 0)
250  io.data (i_map);
251  this->comm().broadcast(i_map);
252  }
253 
254 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
255 
256  type.radial_order = static_cast<Order>(rad_order);
257  type.radial_family = static_cast<FEFamily>(radial_fam);
258  type.inf_map = static_cast<InfMapType>(i_map);
259 
260 #endif
261 
262  if (read_header_in)
263  {
264  if (domains.empty())
265  _written_var_indices[var] = this->add_variable (var_name, type);
266  else
267  _written_var_indices[var] = this->add_variable (var_name, type, &domains);
268  }
269  else
270  _written_var_indices[var] = this->variable_number(var_name);
271  }
272  }
273 
274  // 8.)
275  // Read the number of additional vectors.
276  unsigned int nvecs=0;
277  if (this->processor_id() == 0)
278  io.data (nvecs);
279  this->comm().broadcast(nvecs);
280 
281  // If nvecs > 0, this means that write_additional_data
282  // was true when this file was written. We will need to
283  // make use of this fact later.
284  this->_additional_data_written = nvecs;
285 
286  for (unsigned int vec=0; vec<nvecs; vec++)
287  {
288  // 9.)
289  // Read the name of the vec-th additional vector
290  std::string vec_name;
291  if (this->processor_id() == 0)
292  io.data (vec_name);
293  this->comm().broadcast(vec_name);
294 
295  if (read_additional_data)
296  {
297  // Systems now can handle adding post-initialization vectors
298  // libmesh_assert(this->_can_add_vectors);
299  // Some systems may have added their own vectors already
300  // libmesh_assert_equal_to (this->_vectors.count(vec_name), 0);
301 
302  this->add_vector(vec_name);
303  }
304  }
305 }
306 
307 
308 
309 #ifdef LIBMESH_ENABLE_DEPRECATED
311  const bool read_additional_data)
312 {
313  libmesh_deprecated();
314 
315  // This method implements the output of the vectors
316  // contained in this System object, embedded in the
317  // output of an EquationSystems<T_sys>.
318  //
319  // 10.) The global solution vector, re-ordered to be node-major
320  // (More on this later.)
321  //
322  // for each additional vector in the object
323  //
324  // 11.) The global additional vector, re-ordered to be
325  // node-major (More on this later.)
326  libmesh_assert (io.reading());
327 
328  // read and reordering buffers
329  std::vector<Number> global_vector;
330  std::vector<Number> reordered_vector;
331 
332  // 10.)
333  // Read and set the solution vector
334  {
335  if (this->processor_id() == 0)
336  io.data (global_vector);
337  this->comm().broadcast(global_vector);
338 
339  // Remember that the stored vector is node-major.
340  // We need to put it into whatever application-specific
341  // ordering we may have using the dof_map.
342  reordered_vector.resize(global_vector.size());
343 
344  //libMesh::out << "global_vector.size()=" << global_vector.size() << std::endl;
345  //libMesh::out << "this->n_dofs()=" << this->n_dofs() << std::endl;
346 
347  libmesh_assert_equal_to (global_vector.size(), this->n_dofs());
348 
349  dof_id_type cnt=0;
350 
351  const unsigned int sys = this->number();
352  const unsigned int nv = cast_int<unsigned int>
353  (this->_written_var_indices.size());
354  libmesh_assert_less_equal (nv, this->n_vars());
355 
356  for (unsigned int data_var=0; data_var<nv; data_var++)
357  {
358  const unsigned int var = _written_var_indices[data_var];
359 
360  // First reorder the nodal DOF values
361  for (auto & node : this->get_mesh().node_ptr_range())
362  for (unsigned int index=0; index<node->n_comp(sys,var); index++)
363  {
364  libmesh_assert_not_equal_to (node->dof_number(sys, var, index),
366 
367  libmesh_assert_less (cnt, global_vector.size());
368 
369  reordered_vector[node->dof_number(sys, var, index)] =
370  global_vector[cnt++];
371  }
372 
373  // Then reorder the element DOF values
374  for (auto & elem : this->get_mesh().active_element_ptr_range())
375  for (unsigned int index=0; index<elem->n_comp(sys,var); index++)
376  {
377  libmesh_assert_not_equal_to (elem->dof_number(sys, var, index),
379 
380  libmesh_assert_less (cnt, global_vector.size());
381 
382  reordered_vector[elem->dof_number(sys, var, index)] =
383  global_vector[cnt++];
384  }
385  }
386 
387  *(this->solution) = reordered_vector;
388  }
389 
390  // For each additional vector, simply go through the list.
391  // ONLY attempt to do this IF additional data was actually
392  // written to the file for this system (controlled by the
393  // _additional_data_written flag).
394  if (this->_additional_data_written)
395  {
396  const std::size_t nvecs = this->_vectors.size();
397 
398  // If the number of additional vectors written is non-zero, and
399  // the number of additional vectors we have is non-zero, and
400  // they don't match, then something is wrong and we can't be
401  // sure we're reading data into the correct places.
402  if (read_additional_data && nvecs &&
403  nvecs != this->_additional_data_written)
404  libmesh_error_msg
405  ("Additional vectors in file do not match system");
406 
407  std::map<std::string, NumericVector<Number> *>::iterator
408  pos = this->_vectors.begin();
409 
410  for (std::size_t i = 0; i != this->_additional_data_written; ++i)
411  {
412  // 11.)
413  // Read the values of the vec-th additional vector.
414  // Prior do _not_ clear, but fill with zero, since the
415  // additional vectors _have_ to have the same size
416  // as the solution vector
417  std::fill (global_vector.begin(), global_vector.end(), libMesh::zero);
418 
419  if (this->processor_id() == 0)
420  io.data (global_vector);
421 
422  // If read_additional_data==true and we have additional vectors,
423  // then we will keep this vector data; otherwise we are going to
424  // throw it away.
425  if (read_additional_data && nvecs)
426  {
427  this->comm().broadcast(global_vector);
428 
429  // Remember that the stored vector is node-major.
430  // We need to put it into whatever application-specific
431  // ordering we may have using the dof_map.
432  std::fill (reordered_vector.begin(),
433  reordered_vector.end(),
434  libMesh::zero);
435 
436  reordered_vector.resize(global_vector.size());
437 
438  libmesh_assert_equal_to (global_vector.size(), this->n_dofs());
439 
440  dof_id_type cnt=0;
441 
442  const unsigned int sys = this->number();
443  const unsigned int nv = cast_int<unsigned int>
444  (this->_written_var_indices.size());
445  libmesh_assert_less_equal (nv, this->n_vars());
446 
447  for (unsigned int data_var=0; data_var<nv; data_var++)
448  {
449  const unsigned int var = _written_var_indices[data_var];
450  // First reorder the nodal DOF values
451  for (auto & node : this->get_mesh().node_ptr_range())
452  for (unsigned int index=0; index<node->n_comp(sys,var); index++)
453  {
454  libmesh_assert_not_equal_to (node->dof_number(sys, var, index),
456 
457  libmesh_assert_less (cnt, global_vector.size());
458 
459  reordered_vector[node->dof_number(sys, var, index)] =
460  global_vector[cnt++];
461  }
462 
463  // Then reorder the element DOF values
464  for (auto & elem : this->get_mesh().active_element_ptr_range())
465  for (unsigned int index=0; index<elem->n_comp(sys,var); index++)
466  {
467  libmesh_assert_not_equal_to (elem->dof_number(sys, var, index),
469 
470  libmesh_assert_less (cnt, global_vector.size());
471 
472  reordered_vector[elem->dof_number(sys, var, index)] =
473  global_vector[cnt++];
474  }
475  }
476 
477  // use the overloaded operator=(std::vector) to assign the values
478  *(pos->second) = reordered_vector;
479  }
480 
481  // If we've got vectors then we need to be iterating through
482  // those too
483  if (pos != this->_vectors.end())
484  ++pos;
485  }
486  } // end if (_additional_data_written)
487 }
488 #endif
489 
490 
491 
492 template <typename InValType>
494  const bool read_additional_data)
495 {
515  // PerfLog pl("IO Performance",false);
516  // pl.push("read_parallel_data");
517  dof_id_type total_read_size = 0;
518 
519  libmesh_assert (io.reading());
520  libmesh_assert (io.is_open());
521 
522  // build the ordered nodes and element maps.
523  // when writing/reading parallel files we need to iterate
524  // over our nodes/elements in order of increasing global id().
525  // however, this is not guaranteed to be ordering we obtain
526  // by using the node_iterators/element_iterators directly.
527  // so build a set, sorted by id(), that provides the ordering.
528  // further, for memory economy build the set but then transfer
529  // its contents to vectors, which will be sorted.
530  std::vector<const DofObject *> ordered_nodes, ordered_elements;
531  {
532  std::set<const DofObject *, CompareDofObjectsByID>
533  ordered_nodes_set (this->get_mesh().local_nodes_begin(),
534  this->get_mesh().local_nodes_end());
535 
536  ordered_nodes.insert(ordered_nodes.end(),
537  ordered_nodes_set.begin(),
538  ordered_nodes_set.end());
539  }
540  {
541  std::set<const DofObject *, CompareDofObjectsByID>
542  ordered_elements_set (this->get_mesh().local_elements_begin(),
543  this->get_mesh().local_elements_end());
544 
545  ordered_elements.insert(ordered_elements.end(),
546  ordered_elements_set.begin(),
547  ordered_elements_set.end());
548  }
549 
550  // std::vector<Number> io_buffer;
551  std::vector<InValType> io_buffer;
552 
553  // 9.)
554  //
555  // Actually read the solution components
556  // for the ith system to disk
557  io.data(io_buffer);
558 
559  total_read_size += cast_int<dof_id_type>(io_buffer.size());
560 
561  const unsigned int sys_num = this->number();
562  const unsigned int nv = cast_int<unsigned int>
563  (this->_written_var_indices.size());
564  libmesh_assert_less_equal (nv, this->n_vars());
565 
566  dof_id_type cnt=0;
567 
568  // Loop over each non-SCALAR variable and each node, and read out the value.
569  for (unsigned int data_var=0; data_var<nv; data_var++)
570  {
571  const unsigned int var = _written_var_indices[data_var];
572  if (this->variable(var).type().family != SCALAR)
573  {
574  // First read the node DOF values
575  for (const auto & node : ordered_nodes)
576  for (unsigned int comp=0; comp<node->n_comp(sys_num, var); comp++)
577  {
578  libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
580  libmesh_assert_less (cnt, io_buffer.size());
581  this->solution->set(node->dof_number(sys_num, var, comp), io_buffer[cnt++]);
582  }
583 
584  // Then read the element DOF values
585  for (const auto & elem : ordered_elements)
586  for (unsigned int comp=0; comp<elem->n_comp(sys_num, var); comp++)
587  {
588  libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
590  libmesh_assert_less (cnt, io_buffer.size());
591  this->solution->set(elem->dof_number(sys_num, var, comp), io_buffer[cnt++]);
592  }
593  }
594  }
595 
596  // Finally, read the SCALAR variables on the last processor
597  for (unsigned int data_var=0; data_var<nv; data_var++)
598  {
599  const unsigned int var = _written_var_indices[data_var];
600  if (this->variable(var).type().family == SCALAR)
601  {
602  if (this->processor_id() == (this->n_processors()-1))
603  {
604  const DofMap & dof_map = this->get_dof_map();
605  std::vector<dof_id_type> SCALAR_dofs;
606  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
607 
608  for (std::size_t i=0; i<SCALAR_dofs.size(); i++)
609  this->solution->set(SCALAR_dofs[i], io_buffer[cnt++]);
610  }
611  }
612  }
613 
614  // And we're done setting solution entries
615  this->solution->close();
616 
617  // For each additional vector, simply go through the list.
618  // ONLY attempt to do this IF additional data was actually
619  // written to the file for this system (controlled by the
620  // _additional_data_written flag).
621  if (this->_additional_data_written)
622  {
623  const std::size_t nvecs = this->_vectors.size();
624 
625  // If the number of additional vectors written is non-zero, and
626  // the number of additional vectors we have is non-zero, and
627  // they don't match, then something is wrong and we can't be
628  // sure we're reading data into the correct places.
629  if (read_additional_data && nvecs &&
630  nvecs != this->_additional_data_written)
631  libmesh_error_msg
632  ("Additional vectors in file do not match system");
633 
634  std::map<std::string, NumericVector<Number> *>::const_iterator
635  pos = _vectors.begin();
636 
637  for (std::size_t i = 0; i != this->_additional_data_written; ++i)
638  {
639  cnt=0;
640  io_buffer.clear();
641 
642  // 10.)
643  //
644  // Actually read the additional vector components
645  // for the ith system from disk
646  io.data(io_buffer);
647 
648  total_read_size += cast_int<dof_id_type>(io_buffer.size());
649 
650  // If read_additional_data==true and we have additional vectors,
651  // then we will keep this vector data; otherwise we are going to
652  // throw it away.
653  if (read_additional_data && nvecs)
654  {
655  // Loop over each non-SCALAR variable and each node, and read out the value.
656  for (unsigned int data_var=0; data_var<nv; data_var++)
657  {
658  const unsigned int var = _written_var_indices[data_var];
659  if (this->variable(var).type().family != SCALAR)
660  {
661  // First read the node DOF values
662  for (const auto & node : ordered_nodes)
663  for (unsigned int comp=0; comp<node->n_comp(sys_num, var); comp++)
664  {
665  libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
667  libmesh_assert_less (cnt, io_buffer.size());
668  pos->second->set(node->dof_number(sys_num, var, comp), io_buffer[cnt++]);
669  }
670 
671  // Then read the element DOF values
672  for (const auto & elem : ordered_elements)
673  for (unsigned int comp=0; comp<elem->n_comp(sys_num, var); comp++)
674  {
675  libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
677  libmesh_assert_less (cnt, io_buffer.size());
678  pos->second->set(elem->dof_number(sys_num, var, comp), io_buffer[cnt++]);
679  }
680  }
681  }
682 
683  // Finally, read the SCALAR variables on the last processor
684  for (unsigned int data_var=0; data_var<nv; data_var++)
685  {
686  const unsigned int var = _written_var_indices[data_var];
687  if (this->variable(var).type().family == SCALAR)
688  {
689  if (this->processor_id() == (this->n_processors()-1))
690  {
691  const DofMap & dof_map = this->get_dof_map();
692  std::vector<dof_id_type> SCALAR_dofs;
693  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
694 
695  for (std::size_t sd=0; sd<SCALAR_dofs.size(); sd++)
696  pos->second->set(SCALAR_dofs[sd], io_buffer[cnt++]);
697  }
698  }
699  }
700 
701  // And we're done setting entries for this variable
702  pos->second->close();
703  }
704 
705  // If we've got vectors then we need to be iterating through
706  // those too
707  if (pos != this->_vectors.end())
708  ++pos;
709  }
710  }
711 
712  // const Real
713  // dt = pl.get_elapsed_time(),
714  // rate = total_read_size*sizeof(Number)/dt;
715 
716  // libMesh::err << "Read " << total_read_size << " \"Number\" values\n"
717  // << " Elapsed time = " << dt << '\n'
718  // << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
719 
720  // pl.pop("read_parallel_data");
721 }
722 
723 
724 template <typename InValType>
726  const bool read_additional_data)
727 {
728  // This method implements the input of the vectors
729  // contained in this System object, embedded in the
730  // output of an EquationSystems<T_sys>.
731  //
732  // 10.) The global solution vector, re-ordered to be node-major
733  // (More on this later.)
734  //
735  // for each additional vector in the object
736  //
737  // 11.) The global additional vector, re-ordered to be
738  // node-major (More on this later.)
739  parallel_object_only();
740  std::string comment;
741 
742  // PerfLog pl("IO Performance",false);
743  // pl.push("read_serialized_data");
744  // std::size_t total_read_size = 0;
745 
746  // 10.)
747  // Read the global solution vector
748  {
749  // total_read_size +=
750  this->read_serialized_vector<InValType>(io, this->solution.get());
751 
752  // get the comment
753  if (this->processor_id() == 0)
754  io.comment (comment);
755  }
756 
757  // 11.)
758  // Only read additional vectors if data is available, and only use
759  // that data to fill our vectors if the user requested it.
760  if (this->_additional_data_written)
761  {
762  const std::size_t nvecs = this->_vectors.size();
763 
764  // If the number of additional vectors written is non-zero, and
765  // the number of additional vectors we have is non-zero, and
766  // they don't match, then we can't read additional vectors
767  // and be sure we're reading data into the correct places.
768  if (read_additional_data && nvecs &&
769  nvecs != this->_additional_data_written)
770  libmesh_error_msg
771  ("Additional vectors in file do not match system");
772 
773  std::map<std::string, NumericVector<Number> *>::const_iterator
774  pos = _vectors.begin();
775 
776  for (std::size_t i = 0; i != this->_additional_data_written; ++i)
777  {
778  // Read data, but only put it into a vector if we've been
779  // asked to and if we have a corresponding vector to read.
780 
781  // total_read_size +=
782  this->read_serialized_vector<InValType>
783  (io, (read_additional_data && nvecs) ? pos->second : nullptr);
784 
785  // get the comment
786  if (this->processor_id() == 0)
787  io.comment (comment);
788 
789 
790  // If we've got vectors then we need to be iterating through
791  // those too
792  if (pos != this->_vectors.end())
793  ++pos;
794  }
795  }
796 
797  // const Real
798  // dt = pl.get_elapsed_time(),
799  // rate = total_read_size*sizeof(Number)/dt;
800 
801  // libMesh::out << "Read " << total_read_size << " \"Number\" values\n"
802  // << " Elapsed time = " << dt << '\n'
803  // << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
804 
805  // pl.pop("read_serialized_data");
806 }
807 
808 
809 
810 template <typename iterator_type, typename InValType>
812  const iterator_type begin,
813  const iterator_type end,
814  const InValType ,
815  Xdr & io,
816  const std::vector<NumericVector<Number> *> & vecs,
817  const unsigned int var_to_read) const
818 {
819  //-------------------------------------------------------
820  // General order: (IO format 0.7.4 & greater)
821  //
822  // for (objects ...)
823  // for (vecs ....)
824  // for (vars ....)
825  // for (comps ...)
826  //
827  // where objects are nodes or elements, sorted to be
828  // partition independent,
829  // vecs are one or more *identically distributed* solution
830  // coefficient vectors, vars are one or more variables
831  // to write, and comps are all the components for said
832  // vars on the object.
833 
834  // variables to read. Unless specified otherwise, defaults to _written_var_indices.
835  std::vector<unsigned int> vars_to_read (_written_var_indices);
836 
837  if (var_to_read != libMesh::invalid_uint)
838  {
839  vars_to_read.clear();
840  vars_to_read.push_back(var_to_read);
841  }
842 
843  const unsigned int
844  sys_num = this->number(),
845  num_vecs = cast_int<unsigned int>(vecs.size());
846  const dof_id_type
847  io_blksize = cast_int<dof_id_type>(std::min(max_io_blksize, static_cast<std::size_t>(n_objs))),
848  num_blks = cast_int<unsigned int>(std::ceil(static_cast<double>(n_objs)/
849  static_cast<double>(io_blksize)));
850 
851  libmesh_assert_less_equal (_written_var_indices.size(), this->n_vars());
852 
853  std::size_t n_read_values=0;
854 
855  std::vector<std::vector<dof_id_type>> xfer_ids(num_blks); // The global IDs and # of components for the local objects in all blocks
856  std::vector<std::vector<Number>> recv_vals(num_blks); // The raw values for the local objects in all blocks
857  std::vector<Parallel::Request>
858  id_requests(num_blks), val_requests(num_blks);
859  std::vector<Parallel::MessageTag>
860  id_tags(num_blks), val_tags(num_blks);
861 
862  // ------------------------------------------------------
863  // First pass - count the number of objects in each block
864  // traverse all the objects and figure out which block they
865  // will ultimately live in.
866  std::vector<std::size_t>
867  xfer_ids_size (num_blks,0),
868  recv_vals_size (num_blks,0);
869 
870 
871  for (iterator_type it=begin; it!=end; ++it)
872  {
873  const dof_id_type
874  id = (*it)->id(),
875  block = id/io_blksize;
876 
877  libmesh_assert_less (block, num_blks);
878 
879  xfer_ids_size[block] += 2; // for each object, we send its id, as well as the total number of components for all variables
880 
881  dof_id_type n_comp_tot=0;
882  for (const auto & var : vars_to_read)
883  n_comp_tot += (*it)->n_comp(sys_num, var); // for each variable, we will receive the nonzero components
884 
885  recv_vals_size[block] += n_comp_tot*num_vecs;
886  }
887 
888  // knowing the recv_vals_size[block] for each processor allows
889  // us to sum them and find the global size for each block.
890  std::vector<std::size_t> tot_vals_size(recv_vals_size);
891  this->comm().sum (tot_vals_size);
892 
893 
894  //------------------------------------------
895  // Collect the ids & number of values needed
896  // for all local objects, binning them into
897  // 'blocks' that will be sent to processor 0
898  for (dof_id_type blk=0; blk<num_blks; blk++)
899  {
900  // Each processor should build up its transfer buffers for its
901  // local objects in [first_object,last_object).
902  const dof_id_type
903  first_object = blk*io_blksize,
904  last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
905 
906  // convenience
907  std::vector<dof_id_type> & ids (xfer_ids[blk]);
908  std::vector<Number> & vals (recv_vals[blk]);
909 
910  // we now know the number of values we will store for each block,
911  // so we can do efficient preallocation
912  ids.clear(); ids.reserve (xfer_ids_size[blk]);
913  vals.resize(recv_vals_size[blk]);
914 
915  if (recv_vals_size[blk] != 0) // only if there are nonzero values to receive
916  for (iterator_type it=begin; it!=end; ++it)
917  if (((*it)->id() >= first_object) && // object in [first_object,last_object)
918  ((*it)->id() < last_object))
919  {
920  ids.push_back((*it)->id());
921 
922  unsigned int n_comp_tot=0;
923 
924  for (const auto & var : vars_to_read)
925  n_comp_tot += (*it)->n_comp(sys_num, var);
926 
927  ids.push_back (n_comp_tot*num_vecs);
928  }
929 
930 #ifdef LIBMESH_HAVE_MPI
931  id_tags[blk] = this->comm().get_unique_tag(100*num_blks + blk);
932  val_tags[blk] = this->comm().get_unique_tag(200*num_blks + blk);
933 
934  // nonblocking send the data for this block
935  this->comm().send (0, ids, id_requests[blk], id_tags[blk]);
936 
937  // Go ahead and post the receive too
938  this->comm().receive (0, vals, val_requests[blk], val_tags[blk]);
939 #endif
940  }
941 
942  //---------------------------------------------------
943  // Here processor 0 will read and distribute the data.
944  // We have to do this block-wise to ensure that we
945  // do not exhaust memory on processor 0.
946 
947  // give these variables scope outside the block to avoid reallocation
948  std::vector<std::vector<dof_id_type>> recv_ids (this->n_processors());
949  std::vector<std::vector<Number>> send_vals (this->n_processors());
950  std::vector<Parallel::Request> reply_requests (this->n_processors());
951  std::vector<unsigned int> obj_val_offsets; // map to traverse entry-wise rather than processor-wise
952  std::vector<Number> input_vals; // The input buffer for the current block
953  std::vector<InValType> input_vals_tmp; // The input buffer for the current block
954 
955  for (dof_id_type blk=0; blk<num_blks; blk++)
956  {
957  // Each processor should build up its transfer buffers for its
958  // local objects in [first_object,last_object).
959  const dof_id_type
960  first_object = blk*io_blksize,
961  last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
962  n_objects_blk = last_object - first_object;
963 
964  // Processor 0 has a special job. It needs to gather the requested indices
965  // in [first_object,last_object) from all processors, read the data from
966  // disk, and reply
967  if (this->processor_id() == 0)
968  {
969  // we know the input buffer size for this block and can begin reading it now
970  input_vals.resize(tot_vals_size[blk]);
971  input_vals_tmp.resize(tot_vals_size[blk]);
972 
973  // a ThreadedIO object to perform asynchronous file IO
974  ThreadedIO<InValType> threaded_io(io, input_vals_tmp);
975  Threads::Thread async_io(threaded_io);
976 
977  // offset array. this will define where each object's values
978  // map into the actual input_vals buffer. this must get
979  // 0-initialized because 0-component objects are not actually sent
980  obj_val_offsets.resize (n_objects_blk); std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
981  recv_vals_size.resize(this->n_processors()); // reuse this to count how many values are going to each processor
982 
983 #ifndef NDEBUG
984  std::size_t n_vals_blk = 0;
985 #endif
986 
987  // loop over all processors and process their index request
988  for (unsigned int comm_step=0; comm_step<this->n_processors(); comm_step++)
989  {
990 #ifdef LIBMESH_HAVE_MPI
991  // blocking receive indices for this block, imposing no particular order on processor
992  Parallel::Status id_status (this->comm().probe (Parallel::any_source, id_tags[blk]));
993  std::vector<dof_id_type> & ids (recv_ids[id_status.source()]);
994  std::size_t & n_vals_proc (recv_vals_size[id_status.source()]);
995  this->comm().receive (id_status.source(), ids, id_tags[blk]);
996 #else
997  // straight copy without MPI
998  std::vector<dof_id_type> & ids (recv_ids[0]);
999  std::size_t & n_vals_proc (recv_vals_size[0]);
1000  ids = xfer_ids[blk];
1001 #endif
1002 
1003  n_vals_proc = 0;
1004 
1005  // note its possible we didn't receive values for objects in
1006  // this block if they have no components allocated.
1007  for (std::size_t idx=0; idx<ids.size(); idx+=2)
1008  {
1009  const dof_id_type
1010  local_idx = ids[idx+0]-first_object,
1011  n_vals_tot_allvecs = ids[idx+1];
1012 
1013  libmesh_assert_less (local_idx, n_objects_blk);
1014 
1015  obj_val_offsets[local_idx] = n_vals_tot_allvecs;
1016  n_vals_proc += n_vals_tot_allvecs;
1017  }
1018 
1019 #ifndef NDEBUG
1020  n_vals_blk += n_vals_proc;
1021 #endif
1022  }
1023 
1024  // We need the offsets into the input_vals vector for each object.
1025  // fortunately, this is simply the partial sum of the total number
1026  // of components for each object
1027  std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
1028  obj_val_offsets.begin());
1029 
1030  libmesh_assert_equal_to (n_vals_blk, obj_val_offsets.back());
1031  libmesh_assert_equal_to (n_vals_blk, tot_vals_size[blk]);
1032 
1033  // Wait for read completion
1034  async_io.join();
1035  // now copy the values back to the main vector for transfer
1036  for (std::size_t i_val=0; i_val<input_vals.size(); i_val++)
1037  input_vals[i_val] = input_vals_tmp[i_val];
1038 
1039  n_read_values += input_vals.size();
1040 
1041  // pack data replies for each processor
1042  for (processor_id_type proc=0; proc<this->n_processors(); proc++)
1043  {
1044  const std::vector<dof_id_type> & ids (recv_ids[proc]);
1045  std::vector<Number> & vals (send_vals[proc]);
1046  const std::size_t & n_vals_proc (recv_vals_size[proc]);
1047 
1048  vals.clear(); vals.reserve(n_vals_proc);
1049 
1050  for (std::size_t idx=0; idx<ids.size(); idx+=2)
1051  {
1052  const dof_id_type
1053  local_idx = ids[idx+0]-first_object,
1054  n_vals_tot_allvecs = ids[idx+1];
1055 
1056  std::vector<Number>::const_iterator in_vals(input_vals.begin());
1057  if (local_idx != 0)
1058  std::advance (in_vals, obj_val_offsets[local_idx-1]);
1059 
1060  for (unsigned int val=0; val<n_vals_tot_allvecs; val++, ++in_vals)
1061  {
1062  libmesh_assert (in_vals != input_vals.end());
1063  //libMesh::out << "*in_vals=" << *in_vals << '\n';
1064  vals.push_back(*in_vals);
1065  }
1066  }
1067 
1068 #ifdef LIBMESH_HAVE_MPI
1069  // send the relevant values to this processor
1070  this->comm().send (proc, vals, reply_requests[proc], val_tags[blk]);
1071 #else
1072  recv_vals[blk] = vals;
1073 #endif
1074  }
1075  } // end processor 0 read/reply
1076 
1077  // all processors complete the (already posted) read for this block
1078  {
1079  Parallel::wait (val_requests[blk]);
1080 
1081  const std::vector<Number> & vals (recv_vals[blk]);
1082  std::vector<Number>::const_iterator val_it(vals.begin());
1083 
1084  if (!recv_vals[blk].empty()) // nonzero values to receive
1085  for (iterator_type it=begin; it!=end; ++it)
1086  if (((*it)->id() >= first_object) && // object in [first_object,last_object)
1087  ((*it)->id() < last_object))
1088  // unpack & set the values
1089  for (auto & vec : vecs)
1090  for (const auto & var : vars_to_read)
1091  {
1092  const unsigned int n_comp = (*it)->n_comp(sys_num, var);
1093 
1094  for (unsigned int comp=0; comp<n_comp; comp++, ++val_it)
1095  {
1096  const dof_id_type dof_index = (*it)->dof_number (sys_num, var, comp);
1097  libmesh_assert (val_it != vals.end());
1098  if (vec)
1099  {
1100  libmesh_assert_greater_equal (dof_index, vec->first_local_index());
1101  libmesh_assert_less (dof_index, vec->last_local_index());
1102  //libMesh::out << "dof_index, *val_it = \t" << dof_index << ", " << *val_it << '\n';
1103  vec->set (dof_index, *val_it);
1104  }
1105  }
1106  }
1107  }
1108 
1109  // processor 0 needs to make sure all replies have been handed off
1110  if (this->processor_id () == 0)
1111  Parallel::wait(reply_requests);
1112  }
1113 
1114  Parallel::wait(id_requests);
1115 
1116  return n_read_values;
1117 }
1118 
1119 
1120 
1121 unsigned int System::read_SCALAR_dofs (const unsigned int var,
1122  Xdr & io,
1123  NumericVector<Number> * vec) const
1124 {
1125  unsigned int n_assigned_vals = 0; // the number of values assigned, this will be returned.
1126 
1127  // Processor 0 will read the block from the buffer stream and send it to the last processor
1128  const unsigned int n_SCALAR_dofs = this->variable(var).type().order.get_order();
1129  std::vector<Number> input_buffer(n_SCALAR_dofs);
1130  if (this->processor_id() == 0)
1131  io.data_stream(input_buffer.data(), n_SCALAR_dofs);
1132 
1133 #ifdef LIBMESH_HAVE_MPI
1134  if (this->n_processors() > 1)
1135  {
1136  const Parallel::MessageTag val_tag = this->comm().get_unique_tag(321);
1137 
1138  // Post the receive on the last processor
1139  if (this->processor_id() == (this->n_processors()-1))
1140  this->comm().receive(0, input_buffer, val_tag);
1141 
1142  // Send the data to processor 0
1143  if (this->processor_id() == 0)
1144  this->comm().send(this->n_processors()-1, input_buffer, val_tag);
1145  }
1146 #endif
1147 
1148  // Finally, set the SCALAR values
1149  if (this->processor_id() == (this->n_processors()-1))
1150  {
1151  const DofMap & dof_map = this->get_dof_map();
1152  std::vector<dof_id_type> SCALAR_dofs;
1153  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
1154 
1155  for (std::size_t i=0; i<SCALAR_dofs.size(); i++)
1156  {
1157  if (vec)
1158  vec->set (SCALAR_dofs[i], input_buffer[i]);
1159  ++n_assigned_vals;
1160  }
1161  }
1162 
1163  return n_assigned_vals;
1164 }
1165 
1166 
1167 template <typename InValType>
1169  NumericVector<Number> * vec)
1170 {
1171  parallel_object_only();
1172 
1173 #ifndef NDEBUG
1174  // In parallel we better be reading a parallel vector -- if not
1175  // we will not set all of its components below!!
1176  if (this->n_processors() > 1 && vec)
1177  {
1178  libmesh_assert (vec->type() == PARALLEL ||
1179  vec->type() == GHOSTED);
1180  }
1181 #endif
1182 
1183  libmesh_assert (io.reading());
1184 
1185  // vector length
1186  unsigned int vector_length=0; // FIXME? size_t would break binary compatibility...
1187 #ifndef NDEBUG
1188  std::size_t n_assigned_vals=0;
1189 #endif
1190 
1191  // Get the buffer size
1192  if (this->processor_id() == 0)
1193  io.data(vector_length, "# vector length");
1194  this->comm().broadcast(vector_length);
1195 
1196  const unsigned int nv = cast_int<unsigned int>
1197  (this->_written_var_indices.size());
1198  const dof_id_type
1199  n_nodes = this->get_mesh().n_nodes(),
1200  n_elem = this->get_mesh().n_elem();
1201 
1202  libmesh_assert_less_equal (nv, this->n_vars());
1203 
1204  // for newer versions, read variables node/elem major
1205  if (io.version() >= LIBMESH_VERSION_ID(0,7,4))
1206  {
1207  //---------------------------------
1208  // Collect the values for all nodes
1209 #ifndef NDEBUG
1210  n_assigned_vals +=
1211 #endif
1212  this->read_serialized_blocked_dof_objects (n_nodes,
1213  this->get_mesh().local_nodes_begin(),
1214  this->get_mesh().local_nodes_end(),
1215  InValType(),
1216  io,
1217  std::vector<NumericVector<Number> *> (1,vec));
1218 
1219 
1220  //------------------------------------
1221  // Collect the values for all elements
1222 #ifndef NDEBUG
1223  n_assigned_vals +=
1224 #endif
1226  this->get_mesh().local_elements_begin(),
1227  this->get_mesh().local_elements_end(),
1228  InValType(),
1229  io,
1230  std::vector<NumericVector<Number> *> (1,vec));
1231  }
1232 
1233  // for older versions, read variables var-major
1234  else
1235  {
1236  // Loop over each variable in the system, and then each node/element in the mesh.
1237  for (unsigned int data_var=0; data_var<nv; data_var++)
1238  {
1239  const unsigned int var = _written_var_indices[data_var];
1240  if (this->variable(var).type().family != SCALAR)
1241  {
1242  //---------------------------------
1243  // Collect the values for all nodes
1244 #ifndef NDEBUG
1245  n_assigned_vals +=
1246 #endif
1247  this->read_serialized_blocked_dof_objects (n_nodes,
1248  this->get_mesh().local_nodes_begin(),
1249  this->get_mesh().local_nodes_end(),
1250  InValType(),
1251  io,
1252  std::vector<NumericVector<Number> *> (1,vec),
1253  var);
1254 
1255 
1256  //------------------------------------
1257  // Collect the values for all elements
1258 #ifndef NDEBUG
1259  n_assigned_vals +=
1260 #endif
1262  this->get_mesh().local_elements_begin(),
1263  this->get_mesh().local_elements_end(),
1264  InValType(),
1265  io,
1266  std::vector<NumericVector<Number> *> (1,vec),
1267  var);
1268  } // end variable loop
1269  }
1270  }
1271 
1272  //-------------------------------------------
1273  // Finally loop over all the SCALAR variables
1274  for (unsigned int data_var=0; data_var<nv; data_var++)
1275  {
1276  const unsigned int var = _written_var_indices[data_var];
1277  if (this->variable(var).type().family == SCALAR)
1278  {
1279 #ifndef NDEBUG
1280  n_assigned_vals +=
1281 #endif
1282  this->read_SCALAR_dofs (var, io, vec);
1283  }
1284  }
1285 
1286  if (vec)
1287  vec->close();
1288 
1289 #ifndef NDEBUG
1290  this->comm().sum (n_assigned_vals);
1291  libmesh_assert_equal_to (n_assigned_vals, vector_length);
1292 #endif
1293 
1294  return vector_length;
1295 }
1296 
1297 
1298 
1300  const std::string & /* version is currently unused */,
1301  const bool write_additional_data) const
1302 {
1336  libmesh_assert (io.writing());
1337 
1338 
1339  // Only write the header information
1340  // if we are processor 0.
1341  if (this->get_mesh().processor_id() != 0)
1342  return;
1343 
1344  std::string comment;
1345  char buf[80];
1346 
1347  // 5.)
1348  // Write the number of variables in the system
1349 
1350  {
1351  // set up the comment
1352  comment = "# No. of Variables in System \"";
1353  comment += this->name();
1354  comment += "\"";
1355 
1356  unsigned int nv = this->n_vars();
1357  io.data (nv, comment.c_str());
1358  }
1359 
1360 
1361  for (unsigned int var=0; var<this->n_vars(); var++)
1362  {
1363  // 6.)
1364  // Write the name of the var-th variable
1365  {
1366  // set up the comment
1367  comment = "# Name, Variable No. ";
1368  std::sprintf(buf, "%u", var);
1369  comment += buf;
1370  comment += ", System \"";
1371  comment += this->name();
1372  comment += "\"";
1373 
1374  std::string var_name = this->variable_name(var);
1375  io.data (var_name, comment.c_str());
1376  }
1377 
1378  // 6.1.) Variable subdomains
1379  {
1380  // set up the comment
1381  comment = "# Subdomains, Variable \"";
1382  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1383  comment += buf;
1384  comment += "\", System \"";
1385  comment += this->name();
1386  comment += "\"";
1387 
1388  const std::set<subdomain_id_type> & domains = this->variable(var).active_subdomains();
1389  std::vector<subdomain_id_type> domain_array;
1390  domain_array.assign(domains.begin(), domains.end());
1391  io.data (domain_array, comment.c_str());
1392  }
1393 
1394  // 7.)
1395  // Write the approximation order of the var-th variable
1396  // in this system
1397  {
1398  // set up the comment
1399  comment = "# Approximation Order, Variable \"";
1400  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1401  comment += buf;
1402  comment += "\", System \"";
1403  comment += this->name();
1404  comment += "\"";
1405 
1406  int order = static_cast<int>(this->variable_type(var).order);
1407  io.data (order, comment.c_str());
1408  }
1409 
1410 
1411 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1412 
1413  // do the same for radial_order
1414  {
1415  comment = "# Radial Approximation Order, Variable \"";
1416  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1417  comment += buf;
1418  comment += "\", System \"";
1419  comment += this->name();
1420  comment += "\"";
1421 
1422  int rad_order = static_cast<int>(this->variable_type(var).radial_order);
1423  io.data (rad_order, comment.c_str());
1424  }
1425 
1426 #endif
1427 
1428  // Write the Finite Element type of the var-th variable
1429  // in this System
1430  {
1431  // set up the comment
1432  comment = "# FE Family, Variable \"";
1433  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1434  comment += buf;
1435  comment += "\", System \"";
1436  comment += this->name();
1437  comment += "\"";
1438 
1439  const FEType & type = this->variable_type(var);
1440  int fam = static_cast<int>(type.family);
1441  io.data (fam, comment.c_str());
1442 
1443 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1444 
1445  comment = "# Radial FE Family, Variable \"";
1446  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1447  comment += buf;
1448  comment += "\", System \"";
1449  comment += this->name();
1450  comment += "\"";
1451 
1452  int radial_fam = static_cast<int>(type.radial_family);
1453  io.data (radial_fam, comment.c_str());
1454 
1455  comment = "# Infinite Mapping Type, Variable \"";
1456  std::sprintf(buf, "%s", this->variable_name(var).c_str());
1457  comment += buf;
1458  comment += "\", System \"";
1459  comment += this->name();
1460  comment += "\"";
1461 
1462  int i_map = static_cast<int>(type.inf_map);
1463  io.data (i_map, comment.c_str());
1464 #endif
1465  }
1466  } // end of the variable loop
1467 
1468  // 8.)
1469  // Write the number of additional vectors in the System.
1470  // If write_additional_data==false, then write zero for
1471  // the number of additional vectors.
1472  {
1473  {
1474  // set up the comment
1475  comment = "# No. of Additional Vectors, System \"";
1476  comment += this->name();
1477  comment += "\"";
1478 
1479  unsigned int nvecs = write_additional_data ? this->n_vectors () : 0;
1480  io.data (nvecs, comment.c_str());
1481  }
1482 
1483  if (write_additional_data)
1484  {
1485  unsigned int cnt=0;
1486  for (const auto & pr : _vectors)
1487  {
1488  // 9.)
1489  // write the name of the cnt-th additional vector
1490  comment = "# Name of ";
1491  std::sprintf(buf, "%d", cnt++);
1492  comment += buf;
1493  comment += "th vector";
1494  std::string vec_name = pr.first;
1495 
1496  io.data (vec_name, comment.c_str());
1497  }
1498  }
1499  }
1500 }
1501 
1502 
1503 
1505  const bool write_additional_data) const
1506 {
1526  // PerfLog pl("IO Performance",false);
1527  // pl.push("write_parallel_data");
1528  // std::size_t total_written_size = 0;
1529 
1530  std::string comment;
1531 
1532  libmesh_assert (io.writing());
1533 
1534  std::vector<Number> io_buffer; io_buffer.reserve(this->solution->local_size());
1535 
1536  // build the ordered nodes and element maps.
1537  // when writing/reading parallel files we need to iterate
1538  // over our nodes/elements in order of increasing global id().
1539  // however, this is not guaranteed to be ordering we obtain
1540  // by using the node_iterators/element_iterators directly.
1541  // so build a set, sorted by id(), that provides the ordering.
1542  // further, for memory economy build the set but then transfer
1543  // its contents to vectors, which will be sorted.
1544  std::vector<const DofObject *> ordered_nodes, ordered_elements;
1545  {
1546  std::set<const DofObject *, CompareDofObjectsByID>
1547  ordered_nodes_set (this->get_mesh().local_nodes_begin(),
1548  this->get_mesh().local_nodes_end());
1549 
1550  ordered_nodes.insert(ordered_nodes.end(),
1551  ordered_nodes_set.begin(),
1552  ordered_nodes_set.end());
1553  }
1554  {
1555  std::set<const DofObject *, CompareDofObjectsByID>
1556  ordered_elements_set (this->get_mesh().local_elements_begin(),
1557  this->get_mesh().local_elements_end());
1558 
1559  ordered_elements.insert(ordered_elements.end(),
1560  ordered_elements_set.begin(),
1561  ordered_elements_set.end());
1562  }
1563 
1564  const unsigned int sys_num = this->number();
1565  const unsigned int nv = this->n_vars();
1566 
1567  // Loop over each non-SCALAR variable and each node, and write out the value.
1568  for (unsigned int var=0; var<nv; var++)
1569  if (this->variable(var).type().family != SCALAR)
1570  {
1571  // First write the node DOF values
1572  for (const auto & node : ordered_nodes)
1573  for (unsigned int comp=0; comp<node->n_comp(sys_num, var); comp++)
1574  {
1575  libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
1577 
1578  io_buffer.push_back((*this->solution)(node->dof_number(sys_num, var, comp)));
1579  }
1580 
1581  // Then write the element DOF values
1582  for (const auto & elem : ordered_elements)
1583  for (unsigned int comp=0; comp<elem->n_comp(sys_num, var); comp++)
1584  {
1585  libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
1587 
1588  io_buffer.push_back((*this->solution)(elem->dof_number(sys_num, var, comp)));
1589  }
1590  }
1591 
1592  // Finally, write the SCALAR data on the last processor
1593  for (unsigned int var=0; var<this->n_vars(); var++)
1594  if (this->variable(var).type().family == SCALAR)
1595  {
1596  if (this->processor_id() == (this->n_processors()-1))
1597  {
1598  const DofMap & dof_map = this->get_dof_map();
1599  std::vector<dof_id_type> SCALAR_dofs;
1600  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
1601 
1602  for (std::size_t i=0; i<SCALAR_dofs.size(); i++)
1603  io_buffer.push_back((*this->solution)(SCALAR_dofs[i]));
1604  }
1605  }
1606 
1607  // 9.)
1608  //
1609  // Actually write the reordered solution vector
1610  // for the ith system to disk
1611 
1612  // set up the comment
1613  {
1614  comment = "# System \"";
1615  comment += this->name();
1616  comment += "\" Solution Vector";
1617  }
1618 
1619  io.data (io_buffer, comment.c_str());
1620 
1621  // total_written_size += io_buffer.size();
1622 
1623  // Only write additional vectors if wanted
1624  if (write_additional_data)
1625  {
1626  for (auto & pr : _vectors)
1627  {
1628  io_buffer.clear();
1629  io_buffer.reserve(pr.second->local_size());
1630 
1631  // Loop over each non-SCALAR variable and each node, and write out the value.
1632  for (unsigned int var=0; var<nv; var++)
1633  if (this->variable(var).type().family != SCALAR)
1634  {
1635  // First write the node DOF values
1636  for (const auto & node : ordered_nodes)
1637  for (unsigned int comp=0; comp<node->n_comp(sys_num, var); comp++)
1638  {
1639  libmesh_assert_not_equal_to (node->dof_number(sys_num, var, comp),
1641 
1642  io_buffer.push_back((*pr.second)(node->dof_number(sys_num, var, comp)));
1643  }
1644 
1645  // Then write the element DOF values
1646  for (const auto & elem : ordered_elements)
1647  for (unsigned int comp=0; comp<elem->n_comp(sys_num, var); comp++)
1648  {
1649  libmesh_assert_not_equal_to (elem->dof_number(sys_num, var, comp),
1651 
1652  io_buffer.push_back((*pr.second)(elem->dof_number(sys_num, var, comp)));
1653  }
1654  }
1655 
1656  // Finally, write the SCALAR data on the last processor
1657  for (unsigned int var=0; var<this->n_vars(); var++)
1658  if (this->variable(var).type().family == SCALAR)
1659  {
1660  if (this->processor_id() == (this->n_processors()-1))
1661  {
1662  const DofMap & dof_map = this->get_dof_map();
1663  std::vector<dof_id_type> SCALAR_dofs;
1664  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
1665 
1666  for (std::size_t i=0; i<SCALAR_dofs.size(); i++)
1667  io_buffer.push_back((*pr.second)(SCALAR_dofs[i]));
1668  }
1669  }
1670 
1671  // 10.)
1672  //
1673  // Actually write the reordered additional vector
1674  // for this system to disk
1675 
1676  // set up the comment
1677  {
1678  comment = "# System \"";
1679  comment += this->name();
1680  comment += "\" Additional Vector \"";
1681  comment += pr.first;
1682  comment += "\"";
1683  }
1684 
1685  io.data (io_buffer, comment.c_str());
1686 
1687  // total_written_size += io_buffer.size();
1688  }
1689  }
1690 
1691  // const Real
1692  // dt = pl.get_elapsed_time(),
1693  // rate = total_written_size*sizeof(Number)/dt;
1694 
1695  // libMesh::err << "Write " << total_written_size << " \"Number\" values\n"
1696  // << " Elapsed time = " << dt << '\n'
1697  // << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
1698 
1699  // pl.pop("write_parallel_data");
1700 }
1701 
1702 
1703 
1705  const bool write_additional_data) const
1706 {
1720  parallel_object_only();
1721  std::string comment;
1722 
1723  // PerfLog pl("IO Performance",false);
1724  // pl.push("write_serialized_data");
1725  // std::size_t total_written_size = 0;
1726 
1727  // total_written_size +=
1728  this->write_serialized_vector(io, *this->solution);
1729 
1730  // set up the comment
1731  if (this->processor_id() == 0)
1732  {
1733  comment = "# System \"";
1734  comment += this->name();
1735  comment += "\" Solution Vector";
1736 
1737  io.comment (comment);
1738  }
1739 
1740  // Only write additional vectors if wanted
1741  if (write_additional_data)
1742  {
1743  std::map<std::string, NumericVector<Number> *>::const_iterator
1744  pos = _vectors.begin();
1745 
1746  for (; pos != this->_vectors.end(); ++pos)
1747  {
1748  // total_written_size +=
1749  this->write_serialized_vector(io, *pos->second);
1750 
1751  // set up the comment
1752  if (this->processor_id() == 0)
1753  {
1754  comment = "# System \"";
1755  comment += this->name();
1756  comment += "\" Additional Vector \"";
1757  comment += pos->first;
1758  comment += "\"";
1759  io.comment (comment);
1760  }
1761  }
1762  }
1763 
1764  // const Real
1765  // dt = pl.get_elapsed_time(),
1766  // rate = total_written_size*sizeof(Number)/dt;
1767 
1768  // libMesh::out << "Write " << total_written_size << " \"Number\" values\n"
1769  // << " Elapsed time = " << dt << '\n'
1770  // << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
1771 
1772  // pl.pop("write_serialized_data");
1773 
1774 
1775 
1776 
1777  // // test the new method
1778  // {
1779  // std::vector<std::string> names;
1780  // std::vector<NumericVector<Number> *> vectors_to_write;
1781 
1782  // names.push_back("Solution Vector");
1783  // vectors_to_write.push_back(this->solution.get());
1784 
1785  // // Only write additional vectors if wanted
1786  // if (write_additional_data)
1787  // {
1788  // std::map<std::string, NumericVector<Number> *>::const_iterator
1789  // pos = _vectors.begin();
1790 
1791  // for (; pos != this->_vectors.end(); ++pos)
1792  // {
1793  // names.push_back("Additional Vector " + pos->first);
1794  // vectors_to_write.push_back(pos->second);
1795  // }
1796  // }
1797 
1798  // total_written_size =
1799  // this->write_serialized_vectors (io, names, vectors_to_write);
1800 
1801  // const Real
1802  // dt2 = pl.get_elapsed_time(),
1803  // rate2 = total_written_size*sizeof(Number)/(dt2-dt);
1804 
1805  // libMesh::out << "Write (new) " << total_written_size << " \"Number\" values\n"
1806  // << " Elapsed time = " << (dt2-dt) << '\n'
1807  // << " Rate = " << rate2/1.e6 << "(MB/sec)\n\n";
1808 
1809  // }
1810 }
1811 
1812 
1813 
1814 template <typename iterator_type>
1815 std::size_t System::write_serialized_blocked_dof_objects (const std::vector<const NumericVector<Number> *> & vecs,
1816  const dof_id_type n_objs,
1817  const iterator_type begin,
1818  const iterator_type end,
1819  Xdr & io,
1820  const unsigned int var_to_write) const
1821 {
1822  parallel_object_only();
1823 
1824  //-------------------------------------------------------
1825  // General order: (IO format 0.7.4 & greater)
1826  //
1827  // for (objects ...)
1828  // for (vecs ....)
1829  // for (vars ....)
1830  // for (comps ...)
1831  //
1832  // where objects are nodes or elements, sorted to be
1833  // partition independent,
1834  // vecs are one or more *identically distributed* solution
1835  // coefficient vectors, vars are one or more variables
1836  // to write, and comps are all the components for said
1837  // vars on the object.
1838 
1839  // We will write all variables unless requested otherwise.
1840  std::vector<unsigned int> vars_to_write(1, var_to_write);
1841 
1842  if (var_to_write == libMesh::invalid_uint)
1843  {
1844  vars_to_write.clear(); vars_to_write.reserve(this->n_vars());
1845  for (unsigned int var=0; var<this->n_vars(); var++)
1846  vars_to_write.push_back(var);
1847  }
1848 
1849  const dof_id_type io_blksize = cast_int<dof_id_type>
1850  (std::min(max_io_blksize, static_cast<std::size_t>(n_objs)));
1851 
1852  const unsigned int
1853  sys_num = this->number(),
1854  num_vecs = cast_int<unsigned int>(vecs.size()),
1855  num_blks = cast_int<unsigned int>(std::ceil(static_cast<double>(n_objs)/
1856  static_cast<double>(io_blksize)));
1857 
1858  // libMesh::out << "io_blksize = " << io_blksize
1859  // << ", num_objects = " << n_objs
1860  // << ", num_blks = " << num_blks
1861  // << std::endl;
1862 
1863  std::size_t written_length=0; // The numer of values written. This will be returned
1864  std::vector<std::vector<dof_id_type>> xfer_ids(num_blks); // The global IDs and # of components for the local objects in all blocks
1865  std::vector<std::vector<Number>> send_vals(num_blks); // The raw values for the local objects in all blocks
1866  std::vector<Parallel::Request>
1867  id_requests(num_blks), val_requests(num_blks); // send request handle for each block
1868  std::vector<Parallel::MessageTag>
1869  id_tags(num_blks), val_tags(num_blks); // tag number for each block
1870 
1871  // ------------------------------------------------------
1872  // First pass - count the number of objects in each block
1873  // traverse all the objects and figure out which block they
1874  // will ultimately live in.
1875  std::vector<unsigned int>
1876  xfer_ids_size (num_blks,0),
1877  send_vals_size (num_blks,0);
1878 
1879  for (iterator_type it=begin; it!=end; ++it)
1880  {
1881  const dof_id_type
1882  id = (*it)->id(),
1883  block = id/io_blksize;
1884 
1885  libmesh_assert_less (block, num_blks);
1886 
1887  xfer_ids_size[block] += 2; // for each object, we store its id, as well as the total number of components for all variables
1888 
1889  unsigned int n_comp_tot=0;
1890 
1891  for (const auto & var : vars_to_write)
1892  n_comp_tot += (*it)->n_comp(sys_num, var); // for each variable, we will store the nonzero components
1893 
1894  send_vals_size[block] += n_comp_tot*num_vecs;
1895  }
1896 
1897  //-----------------------------------------
1898  // Collect the values for all local objects,
1899  // binning them into 'blocks' that will be
1900  // sent to processor 0
1901  for (unsigned int blk=0; blk<num_blks; blk++)
1902  {
1903  // libMesh::out << "Writing object block " << blk << std::endl;
1904 
1905  // Each processor should build up its transfer buffers for its
1906  // local objects in [first_object,last_object).
1907  const dof_id_type
1908  first_object = blk*io_blksize,
1909  last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
1910 
1911  // convenience
1912  std::vector<dof_id_type> & ids (xfer_ids[blk]);
1913  std::vector<Number> & vals (send_vals[blk]);
1914 
1915  // we now know the number of values we will store for each block,
1916  // so we can do efficient preallocation
1917  ids.clear(); ids.reserve (xfer_ids_size[blk]);
1918  vals.clear(); vals.reserve (send_vals_size[blk]);
1919 
1920  if (send_vals_size[blk] != 0) // only send if we have nonzero components to write
1921  for (iterator_type it=begin; it!=end; ++it)
1922  if (((*it)->id() >= first_object) && // object in [first_object,last_object)
1923  ((*it)->id() < last_object))
1924  {
1925  ids.push_back((*it)->id());
1926 
1927  // count the total number of nonzeros transferred for this object
1928  {
1929  unsigned int n_comp_tot=0;
1930 
1931  for (const auto & var : vars_to_write)
1932  n_comp_tot += (*it)->n_comp(sys_num, var);
1933 
1934  ids.push_back (n_comp_tot*num_vecs); // even if 0 - processor 0 has no way of knowing otherwise...
1935  }
1936 
1937  // pack the values to send
1938  for (const auto & vec : vecs)
1939  for (const auto & var : vars_to_write)
1940  {
1941  const unsigned int n_comp = (*it)->n_comp(sys_num, var);
1942 
1943  for (unsigned int comp=0; comp<n_comp; comp++)
1944  {
1945  libmesh_assert_greater_equal ((*it)->dof_number(sys_num, var, comp), vec->first_local_index());
1946  libmesh_assert_less ((*it)->dof_number(sys_num, var, comp), vec->last_local_index());
1947  vals.push_back((*vec)((*it)->dof_number(sys_num, var, comp)));
1948  }
1949  }
1950  }
1951 
1952 #ifdef LIBMESH_HAVE_MPI
1953  id_tags[blk] = this->comm().get_unique_tag(100*num_blks + blk);
1954  val_tags[blk] = this->comm().get_unique_tag(200*num_blks + blk);
1955 
1956  // nonblocking send the data for this block
1957  this->comm().send (0, ids, id_requests[blk], id_tags[blk]);
1958  this->comm().send (0, vals, val_requests[blk], val_tags[blk]);
1959 #endif
1960  }
1961 
1962 
1963  if (this->processor_id() == 0)
1964  {
1965  std::vector<std::vector<dof_id_type>> recv_ids (this->n_processors());
1966  std::vector<std::vector<Number>> recv_vals (this->n_processors());
1967  std::vector<unsigned int> obj_val_offsets; // map to traverse entry-wise rather than processor-wise
1968  std::vector<Number> output_vals; // The output buffer for the current block
1969 
1970  // a ThreadedIO object to perform asynchronous file IO
1971  ThreadedIO<Number> threaded_io(io, output_vals);
1972  std::unique_ptr<Threads::Thread> async_io;
1973 
1974  for (unsigned int blk=0; blk<num_blks; blk++)
1975  {
1976  // Each processor should build up its transfer buffers for its
1977  // local objects in [first_object,last_object).
1978  const dof_id_type
1979  first_object = cast_int<dof_id_type>(blk*io_blksize),
1980  last_object = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
1981  n_objects_blk = last_object - first_object;
1982 
1983  // offset array. this will define where each object's values
1984  // map into the actual output_vals buffer. this must get
1985  // 0-initialized because 0-component objects are not actually sent
1986  obj_val_offsets.resize (n_objects_blk); std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
1987 
1988  std::size_t n_val_recvd_blk=0;
1989 
1990  // receive this block of data from all processors.
1991  for (unsigned int comm_step=0; comm_step<this->n_processors(); comm_step++)
1992  {
1993 #ifdef LIBMESH_HAVE_MPI
1994  // blocking receive indices for this block, imposing no particular order on processor
1995  Parallel::Status id_status (this->comm().probe (Parallel::any_source, id_tags[blk]));
1996  std::vector<dof_id_type> & ids (recv_ids[id_status.source()]);
1997  this->comm().receive (id_status.source(), ids, id_tags[blk]);
1998 #else
1999  std::vector<dof_id_type> & ids (recv_ids[0]);
2000  ids = xfer_ids[blk];
2001 #endif
2002 
2003  // note its possible we didn't receive values for objects in
2004  // this block if they have no components allocated.
2005  for (std::size_t idx=0; idx<ids.size(); idx+=2)
2006  {
2007  const dof_id_type
2008  local_idx = ids[idx+0]-first_object,
2009  n_vals_tot_allvecs = ids[idx+1];
2010 
2011  libmesh_assert_less (local_idx, n_objects_blk);
2012  libmesh_assert_less (local_idx, obj_val_offsets.size());
2013 
2014  obj_val_offsets[local_idx] = n_vals_tot_allvecs;
2015  }
2016 
2017 #ifdef LIBMESH_HAVE_MPI
2018  // blocking receive values for this block, imposing no particular order on processor
2019  Parallel::Status val_status (this->comm().probe (Parallel::any_source, val_tags[blk]));
2020  std::vector<Number> & vals (recv_vals[val_status.source()]);
2021  this->comm().receive (val_status.source(), vals, val_tags[blk]);
2022 #else
2023  // straight copy without MPI
2024  std::vector<Number> & vals (recv_vals[0]);
2025  vals = send_vals[blk];
2026 #endif
2027 
2028  n_val_recvd_blk += vals.size();
2029  }
2030 
2031  // We need the offsets into the output_vals vector for each object.
2032  // fortunately, this is simply the partial sum of the total number
2033  // of components for each object
2034  std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
2035  obj_val_offsets.begin());
2036 
2037  // wait on any previous asynchronous IO - this *must* complete before
2038  // we start messing with the output_vals buffer!
2039  if (async_io.get()) async_io->join();
2040 
2041  // this is the actual output buffer that will be written to disk.
2042  // at ths point we finally know wha size it will be.
2043  output_vals.resize(n_val_recvd_blk);
2044 
2045  // pack data from all processors into output values
2046  for (unsigned int proc=0; proc<this->n_processors(); proc++)
2047  {
2048  const std::vector<dof_id_type> & ids (recv_ids [proc]);
2049  const std::vector<Number> & vals(recv_vals[proc]);
2050  std::vector<Number>::const_iterator proc_vals(vals.begin());
2051 
2052  for (std::size_t idx=0; idx<ids.size(); idx+=2)
2053  {
2054  const dof_id_type
2055  local_idx = ids[idx+0]-first_object,
2056  n_vals_tot_allvecs = ids[idx+1];
2057 
2058  // put this object's data into the proper location
2059  // in the output buffer
2060  std::vector<Number>::iterator out_vals(output_vals.begin());
2061  if (local_idx != 0)
2062  std::advance (out_vals, obj_val_offsets[local_idx-1]);
2063 
2064  for (unsigned int val=0; val<n_vals_tot_allvecs; val++, ++out_vals, ++proc_vals)
2065  {
2066  libmesh_assert (out_vals != output_vals.end());
2067  libmesh_assert (proc_vals != vals.end());
2068  *out_vals = *proc_vals;
2069  }
2070  }
2071  }
2072 
2073  // output_vals buffer is now filled for this block.
2074  // write it to disk
2075  async_io.reset(new Threads::Thread(threaded_io));
2076  written_length += output_vals.size();
2077  }
2078 
2079  // wait on any previous asynchronous IO - this *must* complete before
2080  // our stuff goes out of scope
2081  async_io->join();
2082  }
2083 
2084  Parallel::wait(id_requests);
2085  Parallel::wait(val_requests);
2086 
2087  // we need some synchronization here. Because this method
2088  // can be called for a range of nodes, then a range of elements,
2089  // we need some mechanism to prevent processors from racing past
2090  // to the next range and overtaking ongoing communication. one
2091  // approach would be to figure out unique tags for each range,
2092  // but for now we just impose a barrier here. And might as
2093  // well have it do some useful work.
2094  this->comm().broadcast(written_length);
2095 
2096  return written_length;
2097 }
2098 
2099 
2100 
2102  const unsigned int var,
2103  Xdr & io) const
2104 {
2105  unsigned int written_length=0;
2106  std::vector<Number> vals; // The raw values for the local objects in the current block
2107  // Collect the SCALARs for the current variable
2108  if (this->processor_id() == (this->n_processors()-1))
2109  {
2110  const DofMap & dof_map = this->get_dof_map();
2111  std::vector<dof_id_type> SCALAR_dofs;
2112  dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
2113  const unsigned int n_scalar_dofs = cast_int<unsigned int>
2114  (SCALAR_dofs.size());
2115 
2116  for (unsigned int i=0; i<n_scalar_dofs; i++)
2117  {
2118  vals.push_back( vec(SCALAR_dofs[i]) );
2119  }
2120  }
2121 
2122 #ifdef LIBMESH_HAVE_MPI
2123  if (this->n_processors() > 1)
2124  {
2125  const Parallel::MessageTag val_tag(1);
2126 
2127  // Post the receive on processor 0
2128  if (this->processor_id() == 0)
2129  {
2130  this->comm().receive(this->n_processors()-1, vals, val_tag);
2131  }
2132 
2133  // Send the data to processor 0
2134  if (this->processor_id() == (this->n_processors()-1))
2135  {
2136  this->comm().send(0, vals, val_tag);
2137  }
2138  }
2139 #endif
2140 
2141  // -------------------------------------------------------
2142  // Write the output on processor 0.
2143  if (this->processor_id() == 0)
2144  {
2145  const unsigned int vals_size =
2146  cast_int<unsigned int>(vals.size());
2147  io.data_stream (vals.data(), vals_size);
2148  written_length += vals_size;
2149  }
2150 
2151  return written_length;
2152 }
2153 
2154 
2155 
2157  const NumericVector<Number> & vec) const
2158 {
2159  parallel_object_only();
2160 
2161  libmesh_assert (io.writing());
2162 
2163  dof_id_type vec_length = vec.size();
2164  if (this->processor_id() == 0) io.data (vec_length, "# vector length");
2165 
2166  dof_id_type written_length = 0;
2167 
2168  //---------------------------------
2169  // Collect the values for all nodes
2170  written_length += cast_int<dof_id_type>
2171  (this->write_serialized_blocked_dof_objects (std::vector<const NumericVector<Number> *>(1,&vec),
2172  this->get_mesh().n_nodes(),
2173  this->get_mesh().local_nodes_begin(),
2174  this->get_mesh().local_nodes_end(),
2175  io));
2176 
2177  //------------------------------------
2178  // Collect the values for all elements
2179  written_length += cast_int<dof_id_type>
2180  (this->write_serialized_blocked_dof_objects (std::vector<const NumericVector<Number> *>(1,&vec),
2181  this->get_mesh().n_elem(),
2182  this->get_mesh().local_elements_begin(),
2183  this->get_mesh().local_elements_end(),
2184  io));
2185 
2186  //-------------------------------------------
2187  // Finally loop over all the SCALAR variables
2188  for (unsigned int var=0; var<this->n_vars(); var++)
2189  if (this->variable(var).type().family == SCALAR)
2190  {
2191  written_length +=
2192  this->write_SCALAR_dofs (vec, var, io);
2193  }
2194 
2195  if (this->processor_id() == 0)
2196  libmesh_assert_equal_to (written_length, vec_length);
2197 
2198  return written_length;
2199 }
2200 
2201 
2202 template <typename InValType>
2204  const std::vector<NumericVector<Number> *> & vectors) const
2205 {
2206  parallel_object_only();
2207 
2208  // Error checking
2209  // #ifndef NDEBUG
2210  // // In parallel we better be reading a parallel vector -- if not
2211  // // we will not set all of its components below!!
2212  // if (this->n_processors() > 1)
2213  // {
2214  // libmesh_assert (vec.type() == PARALLEL ||
2215  // vec.type() == GHOSTED);
2216  // }
2217  // #endif
2218 
2219  libmesh_assert (io.reading());
2220 
2221  if (this->processor_id() == 0)
2222  {
2223  // sizes
2224  unsigned int num_vecs=0;
2225  dof_id_type vector_length=0;
2226 
2227  // Get the number of vectors
2228  io.data(num_vecs);
2229  // Get the buffer size
2230  io.data(vector_length);
2231 
2232  libmesh_assert_equal_to (num_vecs, vectors.size());
2233 
2234  if (num_vecs != 0)
2235  {
2236  libmesh_assert_not_equal_to (vectors[0], 0);
2237  libmesh_assert_equal_to (vectors[0]->size(), vector_length);
2238  }
2239  }
2240 
2241  // no need to actually communicate these.
2242  // this->comm().broadcast(num_vecs);
2243  // this->comm().broadcast(vector_length);
2244 
2245  // Cache these - they are not free!
2246  const dof_id_type
2247  n_nodes = this->get_mesh().n_nodes(),
2248  n_elem = this->get_mesh().n_elem();
2249 
2250  std::size_t read_length = 0;
2251 
2252  //---------------------------------
2253  // Collect the values for all nodes
2254  read_length +=
2255  this->read_serialized_blocked_dof_objects (n_nodes,
2256  this->get_mesh().local_nodes_begin(),
2257  this->get_mesh().local_nodes_end(),
2258  InValType(),
2259  io,
2260  vectors);
2261 
2262  //------------------------------------
2263  // Collect the values for all elements
2264  read_length +=
2266  this->get_mesh().local_elements_begin(),
2267  this->get_mesh().local_elements_end(),
2268  InValType(),
2269  io,
2270  vectors);
2271 
2272  //-------------------------------------------
2273  // Finally loop over all the SCALAR variables
2274  for (std::size_t vec=0; vec<vectors.size(); vec++)
2275  for (unsigned int var=0; var<this->n_vars(); var++)
2276  if (this->variable(var).type().family == SCALAR)
2277  {
2278  libmesh_assert_not_equal_to (vectors[vec], 0);
2279 
2280  read_length +=
2281  this->read_SCALAR_dofs (var, io, vectors[vec]);
2282  }
2283 
2284  //---------------------------------------
2285  // last step - must close all the vectors
2286  for (std::size_t vec=0; vec<vectors.size(); vec++)
2287  {
2288  libmesh_assert_not_equal_to (vectors[vec], 0);
2289  vectors[vec]->close();
2290  }
2291 
2292  return read_length;
2293 }
2294 
2295 
2296 
2298  const std::vector<const NumericVector<Number> *> & vectors) const
2299 {
2300  parallel_object_only();
2301 
2302  libmesh_assert (io.writing());
2303 
2304  // Cache these - they are not free!
2305  const dof_id_type
2306  n_nodes = this->get_mesh().n_nodes(),
2307  n_elem = this->get_mesh().n_elem();
2308 
2309  std::size_t written_length = 0;
2310 
2311  if (this->processor_id() == 0)
2312  {
2313  unsigned int
2314  n_vec = cast_int<unsigned int>(vectors.size());
2315  dof_id_type
2316  vec_size = vectors.empty() ? 0 : vectors[0]->size();
2317  // Set the number of vectors
2318  io.data(n_vec, "# number of vectors");
2319  // Set the buffer size
2320  io.data(vec_size, "# vector length");
2321  }
2322 
2323  //---------------------------------
2324  // Collect the values for all nodes
2325  written_length +=
2326  this->write_serialized_blocked_dof_objects (vectors,
2327  n_nodes,
2328  this->get_mesh().local_nodes_begin(),
2329  this->get_mesh().local_nodes_end(),
2330  io);
2331 
2332  //------------------------------------
2333  // Collect the values for all elements
2334  written_length +=
2335  this->write_serialized_blocked_dof_objects (vectors,
2336  n_elem,
2337  this->get_mesh().local_elements_begin(),
2338  this->get_mesh().local_elements_end(),
2339  io);
2340 
2341  //-------------------------------------------
2342  // Finally loop over all the SCALAR variables
2343  for (std::size_t vec=0; vec<vectors.size(); vec++)
2344  for (unsigned int var=0; var<this->n_vars(); var++)
2345  if (this->variable(var).type().family == SCALAR)
2346  {
2347  libmesh_assert_not_equal_to (vectors[vec], 0);
2348 
2349  written_length +=
2350  this->write_SCALAR_dofs (*vectors[vec], var, io);
2351  }
2352 
2353  return written_length;
2354 }
2355 
2356 
2357 
2358 
2359 template void System::read_parallel_data<Number> (Xdr & io, const bool read_additional_data);
2360 template void System::read_serialized_data<Number> (Xdr & io, const bool read_additional_data);
2361 template numeric_index_type System::read_serialized_vector<Number> (Xdr & io, NumericVector<Number> * vec);
2362 template std::size_t System::read_serialized_vectors<Number> (Xdr & io, const std::vector<NumericVector<Number> *> & vectors) const;
2363 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
2364 template void System::read_parallel_data<Real> (Xdr & io, const bool read_additional_data);
2365 template void System::read_serialized_data<Real> (Xdr & io, const bool read_additional_data);
2366 template numeric_index_type System::read_serialized_vector<Real> (Xdr & io, NumericVector<Number> * vec);
2367 template std::size_t System::read_serialized_vectors<Real> (Xdr & io, const std::vector<NumericVector<Number> *> & vectors) const;
2368 #endif
2369 
2370 } // namespace libMesh
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
bool writing() const
Definition: xdr_cxx.h:117
void data(T &a, const char *comment="")
Definition: xdr_cxx.C:761
FEFamily family
Definition: fe_type.h:204
void send(const unsigned int dest_processor_id, const T &buf, const MessageTag &tag=no_tag) const
void wait(std::vector< Request > &r)
Definition: request.C:213
void write_serialized_data(Xdr &io, const bool write_additional_data=true) const
Definition: system_io.C:1704
virtual void clear()
Definition: system.C:205
int source() const
Definition: status.h:133
void write_parallel_data(Xdr &io, const bool write_additional_data) const
Definition: system_io.C:1504
const Variable & variable(unsigned int var) const
Definition: system.h:2133
const unsigned int invalid_uint
Definition: libmesh.h:245
const unsigned int any_source
Definition: communicator.h:70
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Definition: mesh_tools.C:702
void comment(std::string &)
Definition: xdr_cxx.C:1519
std::size_t read_serialized_blocked_dof_objects(const dof_id_type n_objects, const iterator_type begin, const iterator_type end, const InValType dummy, Xdr &io, const std::vector< NumericVector< Number > *> &vecs, const unsigned int var_to_read=libMesh::invalid_uint) const
Definition: system_io.C:811
virtual numeric_index_type size() const =0
OrderWrapper radial_order
Definition: fe_type.h:237
unsigned int write_SCALAR_dofs(const NumericVector< Number > &vec, const unsigned int var, Xdr &io) const
Definition: system_io.C:2101
uint8_t processor_id_type
Definition: id_types.h:99
int version() const
Definition: xdr_cxx.h:164
const Parallel::Communicator & comm() const
OrderWrapper order
Definition: fe_type.h:198
IterBase * end
void read_legacy_data(Xdr &io, const bool read_additional_data=true)
Definition: system_io.C:310
MessageTag get_unique_tag(int tagvalue) const
Definition: communicator.C:201
const Number zero
Definition: libmesh.h:239
const MeshBase & get_mesh() const
Definition: system.h:2033
dof_id_type n_dofs() const
Definition: system.C:150
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
Definition: dof_map.C:2348
Tnew cast_int(Told oldvar)
Manages the degrees of freedom (DOFs) in a simulation.
Definition: dof_map.h:176
processor_id_type n_processors() const
const dof_id_type n_nodes
Definition: tecplot_io.C:68
unsigned int number() const
Definition: system.h:2025
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Definition: system.C:661
const std::set< subdomain_id_type > & active_subdomains() const
Definition: variable.h:150
dof_id_type numeric_index_type
Definition: id_types.h:92
void read_header(Xdr &io, const std::string &version, const bool read_header=true, const bool read_additional_data=true, const bool read_legacy_format=false)
Definition: system_io.C:115
unsigned int n_vectors() const
Definition: system.h:2233
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1243
numeric_index_type read_serialized_vector(Xdr &io, NumericVector< Number > *vec)
Definition: system_io.C:1168
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
static const dof_id_type invalid_id
Definition: dof_object.h:347
unsigned int _additional_data_written
Definition: system.h:1984
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2153
InfMapType inf_map
Definition: fe_type.h:258
std::size_t read_serialized_vectors(Xdr &io, const std::vector< NumericVector< Number > *> &vectors) const
Definition: system_io.C:2203
virtual void close()=0
C++ interface for the XDR (eXternal Data Representation) format.
Definition: xdr_cxx.h:66
FEFamily radial_family
Definition: fe_type.h:250
int get_order() const
Definition: fe_type.h:78
std::size_t write_serialized_vectors(Xdr &io, const std::vector< const NumericVector< Number > *> &vectors) const
Definition: system_io.C:2297
void read_parallel_data(Xdr &io, const bool read_additional_data)
Definition: system_io.C:493
ParallelType type() const
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2183
std::map< std::string, NumericVector< Number > *> _vectors
Definition: system.h:1935
bool is_open() const
Definition: xdr_cxx.C:341
Status receive(const unsigned int dest_processor_id, T &buf, const MessageTag &tag=any_tag) const
unsigned int mesh_dimension() const
Definition: mesh_base.C:126
bool reading() const
Definition: xdr_cxx.h:111
unsigned int read_SCALAR_dofs(const unsigned int var, Xdr &io, NumericVector< Number > *vec) const
Definition: system_io.C:1121
virtual void set(const numeric_index_type i, const T value)=0
bool on_command_line(std::string arg)
Definition: libmesh.C:876
const std::string & name() const
Definition: system.h:2017
IterBase * data
unsigned int n_vars() const
Definition: system.h:2105
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Definition: system.C:1081
virtual dof_id_type n_elem() const =0
processor_id_type processor_id() const
void data_stream(T *val, const unsigned int len, const unsigned int line_break=libMesh::invalid_uint)
Definition: xdr_cxx.C:825
OStreamProxy out(std::cout)
const DofMap & get_dof_map() const
Definition: system.h:2049
long double min(long double a, double b)
void read_serialized_data(Xdr &io, const bool read_additional_data=true)
Definition: system_io.C:725
void broadcast(T &data, const unsigned int root_id=0) const
virtual dof_id_type n_nodes() const =0
std::size_t write_serialized_blocked_dof_objects(const std::vector< const NumericVector< Number > *> &vecs, const dof_id_type n_objects, const iterator_type begin, const iterator_type end, Xdr &io, const unsigned int var_to_write=libMesh::invalid_uint) const
Definition: system_io.C:1815
std::vector< unsigned int > _written_var_indices
Definition: system.h:1996
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
uint8_t dof_id_type
Definition: id_types.h:64
void write_header(Xdr &io, const std::string &version, const bool write_additional_data) const
Definition: system_io.C:1299
const FEType & type() const
Definition: variable.h:119
dof_id_type write_serialized_vector(Xdr &io, const NumericVector< Number > &vec) const
Definition: system_io.C:2156