$extrastylesheet
system_io.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 #include "libmesh/libmesh_common.h"
00020 #include "libmesh/parallel.h"
00021 
00022 // C++ Includes
00023 #include <cstdio> // for std::sprintf
00024 #include <set>
00025 #include <numeric> // for std::partial_sum
00026 
00027 // Local Include
00028 #include "libmesh/libmesh_version.h"
00029 #include "libmesh/system.h"
00030 #include "libmesh/mesh_base.h"
00031 //#include "libmesh/mesh_tools.h"
00032 #include "libmesh/elem.h"
00033 #include "libmesh/xdr_cxx.h"
00034 #include "libmesh/numeric_vector.h"
00035 #include "libmesh/dof_map.h"
00036 
00037 
00038 
00039 // Anonymous namespace for implementation details.
00040 namespace {
00041 
00042 using libMesh::DofObject;
00043 using libMesh::Number;
00044 using libMesh::cast_int;
00045 
00046 // Comments:
00047 // ---------
00048 // - The max_io_blksize governs how many nodes or elements will be
00049 // treated as a single block when performing parallel IO on large
00050 // systems.
00051 // - This parameter only loosely affects the size of the actual IO
00052 // buffer as this depends on the number of components a given
00053 // variable has for the nodes/elements in the block.
00054 // - When reading/writing each processor uses an ID map which is
00055 // 3*io_blksize*sizeof(dof_id_type) bytes long, so with unsigned int
00056 // and // io_blksize=256000 we would expect that buffer alone to be
00057 // ~3Mb.
00058 // - In general, an increase in max_io_blksize should increase the
00059 // efficiency of large parallel read/writes by reducing the number
00060 // of MPI messages at the expense of memory.
00061 // - If the library exhausts memory during IO you might reduce this
00062 // parameter.
00063 
00064 const std::size_t max_io_blksize = 256000;
00065 
00066 
00072 struct CompareDofObjectsByID
00073 {
00074   bool operator()(const DofObject *a,
00075                   const DofObject *b) const
00076   {
00077     libmesh_assert (a);
00078     libmesh_assert (b);
00079 
00080     return a->id() < b->id();
00081   }
00082 };
00083 
00087 template <typename InValType>
00088 class ThreadedIO
00089 {
00090 private:
00091   libMesh::Xdr &_io;
00092   std::vector<InValType> &_data;
00093 
00094 public:
00095   ThreadedIO (libMesh::Xdr &io, std::vector<InValType> &data) :
00096     _io(io),
00097     _data(data)
00098   {}
00099 
00100   void operator()()
00101   {
00102     if (_data.empty()) return;
00103     _io.data_stream (&_data[0], cast_int<unsigned int>(_data.size()));
00104   }
00105 };
00106 }
00107 
00108 
00109 namespace libMesh
00110 {
00111 
00112 
00113 // ------------------------------------------------------------
00114 // System class implementation
00115 void System::read_header (Xdr& io,
00116                           const std::string &version,
00117                           const bool read_header_in,
00118                           const bool read_additional_data,
00119                           const bool read_legacy_format)
00120 {
00121   // This method implements the input of a
00122   // System object, embedded in the output of
00123   // an EquationSystems<T_sys>.  This warrants some
00124   // documentation.  The output file essentially
00125   // consists of 5 sections:
00126   //
00127   // for this system
00128   //
00129   //   5.) The number of variables in the system (unsigned int)
00130   //
00131   //   for each variable in the system
00132   //
00133   //     6.) The name of the variable (string)
00134   //
00135   //     6.1.) Variable subdmains
00136   //
00137   //     7.) Combined in an FEType:
00138   //         - The approximation order(s) of the variable
00139   //           (Order Enum, cast to int/s)
00140   //         - The finite element family/ies of the variable
00141   //           (FEFamily Enum, cast to int/s)
00142   //
00143   //   end variable loop
00144   //
00145   //   8.) The number of additional vectors (unsigned int),
00146   //
00147   //     for each additional vector in the system object
00148   //
00149   //     9.) the name of the additional vector  (string)
00150   //
00151   // end system
00152   libmesh_assert (io.reading());
00153 
00154   // Possibly clear data structures and start from scratch.
00155   if (read_header_in)
00156     this->clear ();
00157 
00158   // Figure out if we need to read infinite element information.
00159   // This will be true if the version string contains " with infinite elements"
00160   const bool read_ifem_info =
00161     (version.rfind(" with infinite elements") < version.size()) ||
00162     libMesh::on_command_line ("--read_ifem_systems");
00163 
00164 
00165   {
00166     // 5.)
00167     // Read the number of variables in the system
00168     unsigned int nv=0;
00169     if (this->processor_id() == 0)
00170       io.data (nv);
00171     this->comm().broadcast(nv);
00172 
00173     _written_var_indices.clear();
00174     _written_var_indices.resize(nv, 0);
00175 
00176     for (unsigned int var=0; var<nv; var++)
00177       {
00178         // 6.)
00179         // Read the name of the var-th variable
00180         std::string var_name;
00181         if (this->processor_id() == 0)
00182           io.data (var_name);
00183         this->comm().broadcast(var_name);
00184 
00185         // 6.1.)
00186         std::set<subdomain_id_type> domains;
00187         if (io.version() >= LIBMESH_VERSION_ID(0,7,2))
00188           {
00189             std::vector<subdomain_id_type> domain_array;
00190             if (this->processor_id() == 0)
00191               io.data (domain_array);
00192             for (std::vector<subdomain_id_type>::iterator it = domain_array.begin(); it != domain_array.end(); ++it)
00193               domains.insert(*it);
00194           }
00195         this->comm().broadcast(domains);
00196 
00197         // 7.)
00198         // Read the approximation order(s) of the var-th variable
00199         int order=0;
00200         if (this->processor_id() == 0)
00201           io.data (order);
00202         this->comm().broadcast(order);
00203 
00204 
00205         // do the same for infinite element radial_order
00206         int rad_order=0;
00207         if (read_ifem_info)
00208           {
00209             if (this->processor_id() == 0)
00210               io.data(rad_order);
00211             this->comm().broadcast(rad_order);
00212           }
00213 
00214         // Read the finite element type of the var-th variable
00215         int fam=0;
00216         if (this->processor_id() == 0)
00217           io.data (fam);
00218         this->comm().broadcast(fam);
00219         FEType type;
00220         type.order  = static_cast<Order>(order);
00221         type.family = static_cast<FEFamily>(fam);
00222 
00223         // Check for incompatibilities.  The shape function indexing was
00224         // changed for the monomial and xyz finite element families to
00225         // simplify extension to arbitrary p.  The consequence is that
00226         // old restart files will not be read correctly.  This is expected
00227         // to be an unlikely occurance, but catch it anyway.
00228         if (read_legacy_format)
00229           if ((type.family == MONOMIAL || type.family == XYZ) &&
00230               ((type.order > 2 && this->get_mesh().mesh_dimension() == 2) ||
00231                (type.order > 1 && this->get_mesh().mesh_dimension() == 3)))
00232             {
00233               libmesh_here();
00234               libMesh::out << "*****************************************************************\n"
00235                            << "* WARNING: reading a potentially incompatible restart file!!!   *\n"
00236                            << "*  contact libmesh-users@lists.sourceforge.net for more details *\n"
00237                            << "*****************************************************************"
00238                            << std::endl;
00239             }
00240 
00241         // Read additional information for infinite elements
00242         int radial_fam=0;
00243         int i_map=0;
00244         if (read_ifem_info)
00245           {
00246             if (this->processor_id() == 0)
00247               io.data (radial_fam);
00248             this->comm().broadcast(radial_fam);
00249             if (this->processor_id() == 0)
00250               io.data (i_map);
00251             this->comm().broadcast(i_map);
00252           }
00253 
00254 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00255 
00256         type.radial_order  = static_cast<Order>(rad_order);
00257         type.radial_family = static_cast<FEFamily>(radial_fam);
00258         type.inf_map       = static_cast<InfMapType>(i_map);
00259 
00260 #endif
00261 
00262         if (read_header_in)
00263           {
00264             if (domains.empty())
00265               _written_var_indices[var] = this->add_variable (var_name, type);
00266             else
00267               _written_var_indices[var] = this->add_variable (var_name, type, &domains);
00268           }
00269         else
00270           _written_var_indices[var] = this->variable_number(var_name);
00271       }
00272   }
00273 
00274   // 8.)
00275   // Read the number of additional vectors.
00276   unsigned int nvecs=0;
00277   if (this->processor_id() == 0)
00278     io.data (nvecs);
00279   this->comm().broadcast(nvecs);
00280 
00281   // If nvecs > 0, this means that write_additional_data
00282   // was true when this file was written.  We will need to
00283   // make use of this fact later.
00284   if (nvecs > 0)
00285     this->_additional_data_written = true;
00286 
00287   for (unsigned int vec=0; vec<nvecs; vec++)
00288     {
00289       // 9.)
00290       // Read the name of the vec-th additional vector
00291       std::string vec_name;
00292       if (this->processor_id() == 0)
00293         io.data (vec_name);
00294       this->comm().broadcast(vec_name);
00295 
00296       if (read_additional_data)
00297         {
00298           // Systems now can handle adding post-initialization vectors
00299           //  libmesh_assert(this->_can_add_vectors);
00300           // Some systems may have added their own vectors already
00301           //  libmesh_assert_equal_to (this->_vectors.count(vec_name), 0);
00302 
00303           this->add_vector(vec_name);
00304         }
00305     }
00306 }
00307 
00308 
00309 
00310 void System::read_legacy_data (Xdr& io,
00311                                const bool read_additional_data)
00312 {
00313   libmesh_deprecated();
00314 
00315   // This method implements the output of the vectors
00316   // contained in this System object, embedded in the
00317   // output of an EquationSystems<T_sys>.
00318   //
00319   //   10.) The global solution vector, re-ordered to be node-major
00320   //       (More on this later.)
00321   //
00322   //      for each additional vector in the object
00323   //
00324   //      11.) The global additional vector, re-ordered to be
00325   //           node-major (More on this later.)
00326   libmesh_assert (io.reading());
00327 
00328   // read and reordering buffers
00329   std::vector<Number> global_vector;
00330   std::vector<Number> reordered_vector;
00331 
00332   // 10.)
00333   // Read and set the solution vector
00334   {
00335     if (this->processor_id() == 0)
00336       io.data (global_vector);
00337     this->comm().broadcast(global_vector);
00338 
00339     // Remember that the stored vector is node-major.
00340     // We need to put it into whatever application-specific
00341     // ordering we may have using the dof_map.
00342     reordered_vector.resize(global_vector.size());
00343 
00344     //libMesh::out << "global_vector.size()=" << global_vector.size() << std::endl;
00345     //libMesh::out << "this->n_dofs()=" << this->n_dofs() << std::endl;
00346 
00347     libmesh_assert_equal_to (global_vector.size(), this->n_dofs());
00348 
00349     dof_id_type cnt=0;
00350 
00351     const unsigned int sys = this->number();
00352     const unsigned int nv  = cast_int<unsigned int>
00353       (this->_written_var_indices.size());
00354     libmesh_assert_less_equal (nv, this->n_vars());
00355 
00356     for (unsigned int data_var=0; data_var<nv; data_var++)
00357       {
00358         const unsigned int var = _written_var_indices[data_var];
00359 
00360         // First reorder the nodal DOF values
00361         {
00362           MeshBase::node_iterator
00363             it  = this->get_mesh().nodes_begin(),
00364             end = this->get_mesh().nodes_end();
00365 
00366           for (; it != end; ++it)
00367             for (unsigned int index=0; index<(*it)->n_comp(sys,var); index++)
00368               {
00369                 libmesh_assert_not_equal_to ((*it)->dof_number(sys, var, index),
00370                                              DofObject::invalid_id);
00371 
00372                 libmesh_assert_less (cnt, global_vector.size());
00373 
00374                 reordered_vector[(*it)->dof_number(sys, var, index)] =
00375                   global_vector[cnt++];
00376               }
00377         }
00378 
00379         // Then reorder the element DOF values
00380         {
00381           MeshBase::element_iterator
00382             it  = this->get_mesh().active_elements_begin(),
00383             end = this->get_mesh().active_elements_end();
00384 
00385           for (; it != end; ++it)
00386             for (unsigned int index=0; index<(*it)->n_comp(sys,var); index++)
00387               {
00388                 libmesh_assert_not_equal_to ((*it)->dof_number(sys, var, index),
00389                                              DofObject::invalid_id);
00390 
00391                 libmesh_assert_less (cnt, global_vector.size());
00392 
00393                 reordered_vector[(*it)->dof_number(sys, var, index)] =
00394                   global_vector[cnt++];
00395               }
00396         }
00397       }
00398 
00399     *(this->solution) = reordered_vector;
00400   }
00401 
00402   // For each additional vector, simply go through the list.
00403   // ONLY attempt to do this IF additional data was actually
00404   // written to the file for this system (controlled by the
00405   // _additional_data_written flag).
00406   if (this->_additional_data_written)
00407     {
00408       std::map<std::string, NumericVector<Number>* >::iterator
00409         pos = this->_vectors.begin();
00410 
00411       for (; pos != this->_vectors.end(); ++pos)
00412         {
00413           // 11.)
00414           // Read the values of the vec-th additional vector.
00415           // Prior do _not_ clear, but fill with zero, since the
00416           // additional vectors _have_ to have the same size
00417           // as the solution vector
00418           std::fill (global_vector.begin(), global_vector.end(), libMesh::zero);
00419 
00420           if (this->processor_id() == 0)
00421             io.data (global_vector);
00422           this->comm().broadcast(global_vector);
00423 
00424           // If read_additional_data==true, then we will keep this vector, otherwise
00425           // we are going to throw it away.
00426           if (read_additional_data)
00427             {
00428               // Remember that the stored vector is node-major.
00429               // We need to put it into whatever application-specific
00430               // ordering we may have using the dof_map.
00431               std::fill (reordered_vector.begin(),
00432                          reordered_vector.end(),
00433                          libMesh::zero);
00434 
00435               reordered_vector.resize(global_vector.size());
00436 
00437               libmesh_assert_equal_to (global_vector.size(), this->n_dofs());
00438 
00439               dof_id_type cnt=0;
00440 
00441               const unsigned int sys = this->number();
00442               const unsigned int nv  = cast_int<unsigned int>
00443                 (this->_written_var_indices.size());
00444               libmesh_assert_less_equal (nv, this->n_vars());
00445 
00446               for (unsigned int data_var=0; data_var<nv; data_var++)
00447                 {
00448                   const unsigned int var = _written_var_indices[data_var];
00449                   // First reorder the nodal DOF values
00450                   {
00451                     MeshBase::node_iterator
00452                       it  = this->get_mesh().nodes_begin(),
00453                       end = this->get_mesh().nodes_end();
00454 
00455                     for (; it!=end; ++it)
00456                       for (unsigned int index=0; index<(*it)->n_comp(sys,var); index++)
00457                         {
00458                           libmesh_assert_not_equal_to ((*it)->dof_number(sys, var, index),
00459                                                        DofObject::invalid_id);
00460 
00461                           libmesh_assert_less (cnt, global_vector.size());
00462 
00463                           reordered_vector[(*it)->dof_number(sys, var, index)] =
00464                             global_vector[cnt++];
00465                         }
00466                   }
00467 
00468                   // Then reorder the element DOF values
00469                   {
00470                     MeshBase::element_iterator
00471                       it  = this->get_mesh().active_elements_begin(),
00472                       end = this->get_mesh().active_elements_end();
00473 
00474                     for (; it!=end; ++it)
00475                       for (unsigned int index=0; index<(*it)->n_comp(sys,var); index++)
00476                         {
00477                           libmesh_assert_not_equal_to ((*it)->dof_number(sys, var, index),
00478                                                        DofObject::invalid_id);
00479 
00480                           libmesh_assert_less (cnt, global_vector.size());
00481 
00482                           reordered_vector[(*it)->dof_number(sys, var, index)] =
00483                             global_vector[cnt++];
00484                         }
00485                   }
00486                 }
00487 
00488               // use the overloaded operator=(std::vector) to assign the values
00489               *(pos->second) = reordered_vector;
00490             }
00491         }
00492     } // end if (_additional_data_written)
00493 }
00494 
00495 
00496 
00497 template <typename InValType>
00498 void System::read_parallel_data (Xdr &io,
00499                                  const bool read_additional_data)
00500 {
00520   // PerfLog pl("IO Performance",false);
00521   // pl.push("read_parallel_data");
00522   dof_id_type total_read_size = 0;
00523 
00524   libmesh_assert (io.reading());
00525   libmesh_assert (io.is_open());
00526 
00527   // build the ordered nodes and element maps.
00528   // when writing/reading parallel files we need to iterate
00529   // over our nodes/elements in order of increasing global id().
00530   // however, this is not guaranteed to be ordering we obtain
00531   // by using the node_iterators/element_iterators directly.
00532   // so build a set, sorted by id(), that provides the ordering.
00533   // further, for memory economy build the set but then transfer
00534   // its contents to vectors, which will be sorted.
00535   std::vector<const DofObject*> ordered_nodes, ordered_elements;
00536   {
00537     std::set<const DofObject*, CompareDofObjectsByID>
00538       ordered_nodes_set (this->get_mesh().local_nodes_begin(),
00539                          this->get_mesh().local_nodes_end());
00540 
00541     ordered_nodes.insert(ordered_nodes.end(),
00542                          ordered_nodes_set.begin(),
00543                          ordered_nodes_set.end());
00544   }
00545   {
00546     std::set<const DofObject*, CompareDofObjectsByID>
00547       ordered_elements_set (this->get_mesh().local_elements_begin(),
00548                             this->get_mesh().local_elements_end());
00549 
00550     ordered_elements.insert(ordered_elements.end(),
00551                             ordered_elements_set.begin(),
00552                             ordered_elements_set.end());
00553   }
00554 
00555   //  std::vector<Number> io_buffer;
00556   std::vector<InValType> io_buffer;
00557 
00558   // 9.)
00559   //
00560   // Actually read the solution components
00561   // for the ith system to disk
00562   io.data(io_buffer);
00563 
00564   total_read_size += cast_int<dof_id_type>(io_buffer.size());
00565 
00566   const unsigned int sys_num = this->number();
00567   const unsigned int nv      = cast_int<unsigned int>
00568     (this->_written_var_indices.size());
00569   libmesh_assert_less_equal (nv, this->n_vars());
00570 
00571   dof_id_type cnt=0;
00572 
00573   // Loop over each non-SCALAR variable and each node, and read out the value.
00574   for (unsigned int data_var=0; data_var<nv; data_var++)
00575     {
00576       const unsigned int var = _written_var_indices[data_var];
00577       if(this->variable(var).type().family != SCALAR)
00578         {
00579           // First read the node DOF values
00580           for (std::vector<const DofObject*>::const_iterator
00581                  it = ordered_nodes.begin(); it != ordered_nodes.end(); ++it)
00582             for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
00583               {
00584                 libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
00585                                              DofObject::invalid_id);
00586                 libmesh_assert_less (cnt, io_buffer.size());
00587                 this->solution->set((*it)->dof_number(sys_num, var, comp), io_buffer[cnt++]);
00588               }
00589 
00590           // Then read the element DOF values
00591           for (std::vector<const DofObject*>::const_iterator
00592                  it = ordered_elements.begin(); it != ordered_elements.end(); ++it)
00593             for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
00594               {
00595                 libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
00596                                              DofObject::invalid_id);
00597                 libmesh_assert_less (cnt, io_buffer.size());
00598                 this->solution->set((*it)->dof_number(sys_num, var, comp), io_buffer[cnt++]);
00599               }
00600         }
00601     }
00602 
00603   // Finally, read the SCALAR variables on the last processor
00604   for (unsigned int data_var=0; data_var<nv; data_var++)
00605     {
00606       const unsigned int var = _written_var_indices[data_var];
00607       if(this->variable(var).type().family == SCALAR)
00608         {
00609           if (this->processor_id() == (this->n_processors()-1))
00610             {
00611               const DofMap& dof_map = this->get_dof_map();
00612               std::vector<dof_id_type> SCALAR_dofs;
00613               dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
00614 
00615               for(unsigned int i=0; i<SCALAR_dofs.size(); i++)
00616                 {
00617                   this->solution->set( SCALAR_dofs[i], io_buffer[cnt++] );
00618                 }
00619             }
00620         }
00621     }
00622 
00623   // And we're done setting solution entries
00624   this->solution->close();
00625 
00626   // Only read additional vectors if wanted
00627   if (read_additional_data)
00628     {
00629       std::map<std::string, NumericVector<Number>* >::const_iterator
00630         pos = _vectors.begin();
00631 
00632       for(; pos != this->_vectors.end(); ++pos)
00633         {
00634           cnt=0;
00635           io_buffer.clear();
00636 
00637           // 10.)
00638           //
00639           // Actually read the additional vector components
00640           // for the ith system to disk
00641           io.data(io_buffer);
00642 
00643           total_read_size += cast_int<dof_id_type>(io_buffer.size());
00644 
00645           // Loop over each non-SCALAR variable and each node, and read out the value.
00646           for (unsigned int data_var=0; data_var<nv; data_var++)
00647             {
00648               const unsigned int var = _written_var_indices[data_var];
00649               if(this->variable(var).type().family != SCALAR)
00650                 {
00651                   // First read the node DOF values
00652                   for (std::vector<const DofObject*>::const_iterator
00653                          it = ordered_nodes.begin(); it != ordered_nodes.end(); ++it)
00654                     for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
00655                       {
00656                         libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
00657                                                      DofObject::invalid_id);
00658                         libmesh_assert_less (cnt, io_buffer.size());
00659                         pos->second->set((*it)->dof_number(sys_num, var, comp), io_buffer[cnt++]);
00660                       }
00661 
00662                   // Then read the element DOF values
00663                   for (std::vector<const DofObject*>::const_iterator
00664                          it = ordered_elements.begin(); it != ordered_elements.end(); ++it)
00665                     for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
00666                       {
00667                         libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
00668                                                      DofObject::invalid_id);
00669                         libmesh_assert_less (cnt, io_buffer.size());
00670                         pos->second->set((*it)->dof_number(sys_num, var, comp), io_buffer[cnt++]);
00671                       }
00672                 }
00673             }
00674 
00675           // Finally, read the SCALAR variables on the last processor
00676           for (unsigned int data_var=0; data_var<nv; data_var++)
00677             {
00678               const unsigned int var = _written_var_indices[data_var];
00679               if(this->variable(var).type().family == SCALAR)
00680                 {
00681                   if (this->processor_id() == (this->n_processors()-1))
00682                     {
00683                       const DofMap& dof_map = this->get_dof_map();
00684                       std::vector<dof_id_type> SCALAR_dofs;
00685                       dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
00686 
00687                       for(unsigned int i=0; i<SCALAR_dofs.size(); i++)
00688                         {
00689                           pos->second->set( SCALAR_dofs[i], io_buffer[cnt++] );
00690                         }
00691                     }
00692                 }
00693             }
00694 
00695           // And we're done setting entries for this variable
00696           pos->second->close();
00697         }
00698     }
00699 
00700   // const Real
00701   //   dt   = pl.get_elapsed_time(),
00702   //   rate = total_read_size*sizeof(Number)/dt;
00703 
00704   // libMesh::err << "Read " << total_read_size << " \"Number\" values\n"
00705   //     << " Elapsed time = " << dt << '\n'
00706   //     << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
00707 
00708   // pl.pop("read_parallel_data");
00709 }
00710 
00711 
00712 template <typename InValType>
00713 void System::read_serialized_data (Xdr& io,
00714                                    const bool read_additional_data)
00715 {
00716   // This method implements the input of the vectors
00717   // contained in this System object, embedded in the
00718   // output of an EquationSystems<T_sys>.
00719   //
00720   //   10.) The global solution vector, re-ordered to be node-major
00721   //       (More on this later.)
00722   //
00723   //      for each additional vector in the object
00724   //
00725   //      11.) The global additional vector, re-ordered to be
00726   //          node-major (More on this later.)
00727   parallel_object_only();
00728   std::string comment;
00729 
00730   // PerfLog pl("IO Performance",false);
00731   // pl.push("read_serialized_data");
00732   // std::size_t total_read_size = 0;
00733 
00734   // 10.)
00735   // Read the global solution vector
00736   {
00737     // total_read_size +=
00738     this->read_serialized_vector<InValType>(io, *this->solution);
00739 
00740     // get the comment
00741     if (this->processor_id() == 0)
00742       io.comment (comment);
00743   }
00744 
00745   // 11.)
00746   // Only read additional vectors if wanted
00747   if (read_additional_data)
00748     {
00749       std::map<std::string, NumericVector<Number>* >::const_iterator
00750         pos = _vectors.begin();
00751 
00752       for(; pos != this->_vectors.end(); ++pos)
00753         {
00754           // total_read_size +=
00755           this->read_serialized_vector<InValType>(io, *pos->second);
00756 
00757           // get the comment
00758           if (this->processor_id() == 0)
00759             io.comment (comment);
00760 
00761         }
00762     }
00763 
00764   // const Real
00765   //   dt   = pl.get_elapsed_time(),
00766   //   rate = total_read_size*sizeof(Number)/dt;
00767 
00768   // libMesh::out << "Read " << total_read_size << " \"Number\" values\n"
00769   //     << " Elapsed time = " << dt << '\n'
00770   //     << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
00771 
00772   // pl.pop("read_serialized_data");
00773 }
00774 
00775 
00776 
00777 template <typename iterator_type, typename InValType>
00778 std::size_t System::read_serialized_blocked_dof_objects (const dof_id_type n_objs,
00779                                                          const iterator_type begin,
00780                                                          const iterator_type end,
00781                                                          const InValType ,
00782                                                          Xdr &io,
00783                                                          const std::vector<NumericVector<Number>*> &vecs,
00784                                                          const unsigned int var_to_read) const
00785 {
00786   //-------------------------------------------------------
00787   // General order: (IO format 0.7.4 & greater)
00788   //
00789   // for (objects ...)
00790   //   for (vecs ....)
00791   //     for (vars ....)
00792   //       for (comps ...)
00793   //
00794   // where objects are nodes or elements, sorted to be
00795   // partition independent,
00796   // vecs are one or more *identically distributed* solution
00797   // coefficient vectors, vars are one or more variables
00798   // to write, and comps are all the components for said
00799   // vars on the object.
00800 
00801   typedef std::vector<NumericVector<Number>*>::const_iterator vec_iterator_type;
00802 
00803   // variables to read.  Unless specified otherwise, defaults to _written_var_indices.
00804   std::vector<unsigned int> vars_to_read (_written_var_indices);
00805 
00806   if (var_to_read != libMesh::invalid_uint)
00807     {
00808       vars_to_read.clear();
00809       vars_to_read.push_back(var_to_read);
00810     }
00811 
00812   const unsigned int
00813     sys_num    = this->number(),
00814     num_vecs   = cast_int<unsigned int>(vecs.size());
00815   const dof_id_type
00816     io_blksize = cast_int<dof_id_type>(std::min(max_io_blksize, static_cast<std::size_t>(n_objs))),
00817     num_blks   = std::ceil(static_cast<double>(n_objs)/static_cast<double>(io_blksize));
00818 
00819   libmesh_assert_less_equal (_written_var_indices.size(), this->n_vars());
00820 
00821   std::size_t n_read_values=0;
00822 
00823   std::vector<std::vector<dof_id_type> > xfer_ids(num_blks);  // The global IDs and # of components for the local objects in all blocks
00824   std::vector<std::vector<Number> >      recv_vals(num_blks); // The raw values for the local objects in all blocks
00825   std::vector<Parallel::Request>
00826     id_requests(num_blks), val_requests(num_blks);
00827 
00828   // ------------------------------------------------------
00829   // First pass - count the number of objects in each block
00830   // traverse all the objects and figure out which block they
00831   // will utlimately live in.
00832   std::vector<std::size_t>
00833     xfer_ids_size  (num_blks,0),
00834     recv_vals_size (num_blks,0);
00835 
00836 
00837   for (iterator_type it=begin; it!=end; ++it)
00838     {
00839       const dof_id_type
00840         id    = (*it)->id(),
00841         block = id/io_blksize;
00842 
00843       libmesh_assert_less (block, num_blks);
00844 
00845       xfer_ids_size[block] += 2; // for each object, we send its id, as well as the total number of components for all variables
00846 
00847       dof_id_type n_comp_tot=0;
00848       for (std::vector<unsigned int>::const_iterator var_it=vars_to_read.begin();
00849            var_it!=vars_to_read.end(); ++var_it)
00850         n_comp_tot += (*it)->n_comp(sys_num, *var_it); // for each variable, we will receive the nonzero components
00851 
00852       recv_vals_size[block] += n_comp_tot*num_vecs;
00853     }
00854 
00855   // knowing the recv_vals_size[block] for each processor allows
00856   // us to sum them and find the global size for each block.
00857   std::vector<std::size_t> tot_vals_size(recv_vals_size);
00858   this->comm().sum (tot_vals_size);
00859 
00860 
00861   //------------------------------------------
00862   // Collect the ids & number of values needed
00863   // for all local objects, binning them into
00864   // 'blocks' that will be sent to processor 0
00865   for (dof_id_type blk=0; blk<num_blks; blk++)
00866     {
00867       // Each processor should build up its transfer buffers for its
00868       // local objects in [first_object,last_object).
00869       const dof_id_type
00870         first_object = blk*io_blksize,
00871         last_object  = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
00872 
00873       // convenience
00874       std::vector<dof_id_type> &ids  (xfer_ids[blk]);
00875       std::vector<Number>       &vals (recv_vals[blk]);
00876 
00877       // we now know the number of values we will store for each block,
00878       // so we can do efficient preallocation
00879       ids.clear();  ids.reserve (xfer_ids_size[blk]);
00880       vals.resize(recv_vals_size[blk]);
00881 
00882       if (recv_vals_size[blk] != 0) // only if there are nonzero values to receive
00883         for (iterator_type it=begin; it!=end; ++it)
00884           if (((*it)->id() >= first_object) && // object in [first_object,last_object)
00885               ((*it)->id() <   last_object))
00886             {
00887               ids.push_back((*it)->id());
00888 
00889               unsigned int n_comp_tot=0;
00890 
00891               for (std::vector<unsigned int>::const_iterator var_it=vars_to_read.begin();
00892                    var_it!=vars_to_read.end(); ++var_it)
00893                 n_comp_tot += (*it)->n_comp(sys_num,*var_it);
00894 
00895               ids.push_back (n_comp_tot*num_vecs);
00896             }
00897 
00898 #ifdef LIBMESH_HAVE_MPI
00899       Parallel::MessageTag id_tag  = this->comm().get_unique_tag(100*num_blks + blk);
00900       Parallel::MessageTag val_tag = this->comm().get_unique_tag(200*num_blks + blk);
00901 
00902       // nonblocking send the data for this block
00903       this->comm().send (0, ids,  id_requests[blk],  id_tag);
00904 
00905       // Go ahead and post the receive too
00906       this->comm().receive (0, vals, val_requests[blk], val_tag);
00907 #endif
00908     }
00909 
00910   //---------------------------------------------------
00911   // Here processor 0 will read and distribute the data.
00912   // We have to do this block-wise to ensure that we
00913   // do not exhaust memory on processor 0.
00914 
00915   // give these variables scope outside the block to avoid reallocation
00916   std::vector<std::vector<dof_id_type> > recv_ids       (this->n_processors());
00917   std::vector<std::vector<Number> >      send_vals      (this->n_processors());
00918   std::vector<Parallel::Request>         reply_requests (this->n_processors());
00919   std::vector<unsigned int>              obj_val_offsets;          // map to traverse entry-wise rather than processor-wise
00920   std::vector<Number>                    input_vals;               // The input buffer for the current block
00921   std::vector<InValType>                 input_vals_tmp;               // The input buffer for the current block
00922 
00923   for (dof_id_type blk=0; blk<num_blks; blk++)
00924     {
00925       // Each processor should build up its transfer buffers for its
00926       // local objects in [first_object,last_object).
00927       const dof_id_type
00928         first_object  = blk*io_blksize,
00929         last_object   = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
00930         n_objects_blk = last_object - first_object;
00931 
00932       // Processor 0 has a special job.  It needs to gather the requested indices
00933       // in [first_object,last_object) from all processors, read the data from
00934       // disk, and reply
00935       if (this->processor_id() == 0)
00936         {
00937           // we know the input buffer size for this block and can begin reading it now
00938           input_vals.resize(tot_vals_size[blk]);
00939           input_vals_tmp.resize(tot_vals_size[blk]);
00940 
00941           // a ThreadedIO object to perform asychronous file IO
00942           ThreadedIO<InValType> threaded_io(io, input_vals_tmp);
00943           Threads::Thread async_io(threaded_io);
00944 
00945           Parallel::MessageTag id_tag  = this->comm().get_unique_tag(100*num_blks + blk);
00946           Parallel::MessageTag val_tag = this->comm().get_unique_tag(200*num_blks + blk);
00947 
00948           // offset array. this will define where each object's values
00949           // map into the actual input_vals buffer.  this must get
00950           // 0-initialized because 0-component objects are not actually sent
00951           obj_val_offsets.resize (n_objects_blk);  std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
00952           recv_vals_size.resize(this->n_processors()); // reuse this to count how many values are going to each processor
00953 
00954 #ifndef NDEBUG
00955           std::size_t n_vals_blk = 0;
00956 #endif
00957 
00958           // loop over all processors and process their index request
00959           for (unsigned int comm_step=0; comm_step<this->n_processors(); comm_step++)
00960             {
00961 #ifdef LIBMESH_HAVE_MPI
00962               // blocking receive indices for this block, imposing no particular order on processor
00963               Parallel::Status id_status (this->comm().probe (Parallel::any_source, id_tag));
00964               std::vector<dof_id_type> &ids (recv_ids[id_status.source()]);
00965               std::size_t &n_vals_proc (recv_vals_size[id_status.source()]);
00966               this->comm().receive (id_status.source(), ids,  id_tag);
00967 #else
00968               // straight copy without MPI
00969               std::vector<dof_id_type> &ids (recv_ids[0]);
00970               std::size_t &n_vals_proc (recv_vals_size[0]);
00971               ids = xfer_ids[blk];
00972 #endif
00973 
00974               n_vals_proc = 0;
00975 
00976               // note its possible we didn't receive values for objects in
00977               // this block if they have no components allocated.
00978               for (dof_id_type idx=0; idx<ids.size(); idx+=2)
00979                 {
00980                   const dof_id_type
00981                     local_idx          = ids[idx+0]-first_object,
00982                     n_vals_tot_allvecs = ids[idx+1];
00983 
00984                   libmesh_assert_less (local_idx, n_objects_blk);
00985 
00986                   obj_val_offsets[local_idx] = n_vals_tot_allvecs;
00987                   n_vals_proc += n_vals_tot_allvecs;
00988                 }
00989 
00990 #ifndef NDEBUG
00991               n_vals_blk += n_vals_proc;
00992 #endif
00993             }
00994 
00995           // We need the offests into the input_vals vector for each object.
00996           // fortunately, this is simply the partial sum of the total number
00997           // of components for each object
00998           std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
00999                            obj_val_offsets.begin());
01000 
01001           libmesh_assert_equal_to (n_vals_blk, obj_val_offsets.back());
01002           libmesh_assert_equal_to (n_vals_blk, tot_vals_size[blk]);
01003 
01004           // Wait for read completion
01005           async_io.join();
01006           // now copy the values back to the main vector for transfer
01007           for (unsigned int i_val=0; i_val<input_vals.size(); i_val++)
01008             input_vals[i_val] = input_vals_tmp[i_val];
01009 
01010           n_read_values += input_vals.size();
01011 
01012           // pack data replies for each processor
01013           for (processor_id_type proc=0; proc<this->n_processors(); proc++)
01014             {
01015               const std::vector<dof_id_type> &ids (recv_ids[proc]);
01016               std::vector<Number> &vals (send_vals[proc]);
01017               const std::size_t &n_vals_proc (recv_vals_size[proc]);
01018 
01019               vals.clear();  vals.reserve(n_vals_proc);
01020 
01021               for (std::size_t idx=0; idx<ids.size(); idx+=2)
01022                 {
01023                   const dof_id_type
01024                     local_idx          = ids[idx+0]-first_object,
01025                     n_vals_tot_allvecs = ids[idx+1];
01026 
01027                   std::vector<Number>::const_iterator in_vals(input_vals.begin());
01028                   if (local_idx != 0)
01029                     std::advance (in_vals, obj_val_offsets[local_idx-1]);
01030 
01031                   for (unsigned int val=0; val<n_vals_tot_allvecs; val++, ++in_vals)
01032                     {
01033                       libmesh_assert (in_vals != input_vals.end());
01034                       //libMesh::out << "*in_vals=" << *in_vals << '\n';
01035                       vals.push_back(*in_vals);
01036                     }
01037                 }
01038 
01039 #ifdef LIBMESH_HAVE_MPI
01040               // send the relevant values to this processor
01041               this->comm().send (proc, vals, reply_requests[proc], val_tag);
01042 #else
01043               recv_vals[blk] = vals;
01044 #endif
01045             }
01046         } // end processor 0 read/reply
01047 
01048       // all processors complete the (already posted) read for this block
01049       {
01050         Parallel::wait (val_requests[blk]);
01051 
01052         const std::vector<Number> &vals (recv_vals[blk]);
01053         std::vector<Number>::const_iterator val_it(vals.begin());
01054 
01055         if (!recv_vals[blk].empty()) // nonzero values to receive
01056           for (iterator_type it=begin; it!=end; ++it)
01057             if (((*it)->id() >= first_object) && // object in [first_object,last_object)
01058                 ((*it)->id() <   last_object))
01059               // unpack & set the values
01060               for (vec_iterator_type vec_it=vecs.begin(); vec_it!=vecs.end(); ++vec_it)
01061                 {
01062                   NumericVector<Number> &vec(**vec_it);
01063 
01064                   for (std::vector<unsigned int>::const_iterator var_it=vars_to_read.begin();
01065                        var_it!=vars_to_read.end(); ++var_it)
01066                     {
01067                       const unsigned int n_comp = (*it)->n_comp(sys_num,*var_it);
01068 
01069                       for (unsigned int comp=0; comp<n_comp; comp++, ++val_it)
01070                         {
01071                           const dof_id_type dof_index = (*it)->dof_number (sys_num, *var_it, comp);
01072                           libmesh_assert (val_it != vals.end());
01073                           libmesh_assert_greater_equal (dof_index, vec.first_local_index());
01074                           libmesh_assert_less (dof_index, vec.last_local_index());
01075                           //libMesh::out << "dof_index, *val_it = \t" << dof_index << ", " << *val_it << '\n';
01076                           vec.set (dof_index, *val_it);
01077                         }
01078                     }
01079                 }
01080       }
01081 
01082       // processor 0 needs to make sure all replies have been handed off
01083       if (this->processor_id () == 0)
01084         Parallel::wait(reply_requests);
01085     }
01086 
01087   return n_read_values;
01088 }
01089 
01090 
01091 
01092 unsigned int System::read_SCALAR_dofs (const unsigned int var,
01093                                        Xdr &io,
01094                                        NumericVector<Number> &vec) const
01095 {
01096   unsigned int n_assigned_vals = 0; // the number of values assigned, this will be returned.
01097 
01098   // Processor 0 will read the block from the buffer stream and send it to the last processor
01099   const unsigned int n_SCALAR_dofs = this->variable(var).type().order;
01100   std::vector<Number> input_buffer(n_SCALAR_dofs);
01101   if (this->processor_id() == 0)
01102     {
01103       io.data_stream(&input_buffer[0], n_SCALAR_dofs);
01104     }
01105 
01106 #ifdef LIBMESH_HAVE_MPI
01107   if ( this->n_processors() > 1 )
01108     {
01109       const Parallel::MessageTag val_tag = this->comm().get_unique_tag(321);
01110 
01111       // Post the receive on the last processor
01112       if (this->processor_id() == (this->n_processors()-1))
01113         this->comm().receive(0, input_buffer, val_tag);
01114 
01115       // Send the data to processor 0
01116       if (this->processor_id() == 0)
01117         this->comm().send(this->n_processors()-1, input_buffer, val_tag);
01118     }
01119 #endif
01120 
01121   // Finally, set the SCALAR values
01122   if (this->processor_id() == (this->n_processors()-1))
01123     {
01124       const DofMap& dof_map = this->get_dof_map();
01125       std::vector<dof_id_type> SCALAR_dofs;
01126       dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
01127 
01128       for(unsigned int i=0; i<SCALAR_dofs.size(); i++)
01129         {
01130           vec.set (SCALAR_dofs[i], input_buffer[i]);
01131           ++n_assigned_vals;
01132         }
01133     }
01134 
01135   return n_assigned_vals;
01136 }
01137 
01138 
01139 template <typename InValType>
01140 numeric_index_type System::read_serialized_vector (Xdr& io, NumericVector<Number>& vec)
01141 {
01142   parallel_object_only();
01143 
01144 #ifndef NDEBUG
01145   // In parallel we better be reading a parallel vector -- if not
01146   // we will not set all of its components below!!
01147   if (this->n_processors() > 1)
01148     {
01149       libmesh_assert (vec.type() == PARALLEL ||
01150                       vec.type() == GHOSTED);
01151     }
01152 #endif
01153 
01154   libmesh_assert (io.reading());
01155 
01156   // vector length
01157   unsigned int vector_length=0; // FIXME?  size_t would break binary compatibility...
01158 #ifndef NDEBUG
01159   std::size_t n_assigned_vals=0;
01160 #endif
01161 
01162   // Get the buffer size
01163   if (this->processor_id() == 0)
01164     io.data(vector_length, "# vector length");
01165   this->comm().broadcast(vector_length);
01166 
01167   const unsigned int nv = cast_int<unsigned int>
01168     (this->_written_var_indices.size());
01169   const dof_id_type
01170     n_nodes = this->get_mesh().n_nodes(),
01171     n_elem  = this->get_mesh().n_elem();
01172 
01173   libmesh_assert_less_equal (nv, this->n_vars());
01174 
01175   // for newer versions, read variables node/elem major
01176   if (io.version() >= LIBMESH_VERSION_ID(0,7,4))
01177     {
01178       //---------------------------------
01179       // Collect the values for all nodes
01180 #ifndef NDEBUG
01181       n_assigned_vals +=
01182 #endif
01183         this->read_serialized_blocked_dof_objects (n_nodes,
01184                                                    this->get_mesh().local_nodes_begin(),
01185                                                    this->get_mesh().local_nodes_end(),
01186                                                    InValType(),
01187                                                    io,
01188                                                    std::vector<NumericVector<Number>*> (1,&vec));
01189 
01190 
01191       //------------------------------------
01192       // Collect the values for all elements
01193 #ifndef NDEBUG
01194       n_assigned_vals +=
01195 #endif
01196         this->read_serialized_blocked_dof_objects (n_elem,
01197                                                    this->get_mesh().local_elements_begin(),
01198                                                    this->get_mesh().local_elements_end(),
01199                                                    InValType(),
01200                                                    io,
01201                                                    std::vector<NumericVector<Number>*> (1,&vec));
01202     }
01203 
01204   // for older versions, read variables var-major
01205   else
01206     {
01207       // Loop over each variable in the system, and then each node/element in the mesh.
01208       for (unsigned int data_var=0; data_var<nv; data_var++)
01209         {
01210           const unsigned int var = _written_var_indices[data_var];
01211           if(this->variable(var).type().family != SCALAR)
01212             {
01213               //---------------------------------
01214               // Collect the values for all nodes
01215 #ifndef NDEBUG
01216               n_assigned_vals +=
01217 #endif
01218                 this->read_serialized_blocked_dof_objects (n_nodes,
01219                                                            this->get_mesh().local_nodes_begin(),
01220                                                            this->get_mesh().local_nodes_end(),
01221                                                            InValType(),
01222                                                            io,
01223                                                            std::vector<NumericVector<Number>*> (1,&vec),
01224                                                            var);
01225 
01226 
01227               //------------------------------------
01228               // Collect the values for all elements
01229 #ifndef NDEBUG
01230               n_assigned_vals +=
01231 #endif
01232                 this->read_serialized_blocked_dof_objects (n_elem,
01233                                                            this->get_mesh().local_elements_begin(),
01234                                                            this->get_mesh().local_elements_end(),
01235                                                            InValType(),
01236                                                            io,
01237                                                            std::vector<NumericVector<Number>*> (1,&vec),
01238                                                            var);
01239             } // end variable loop
01240         }
01241     }
01242 
01243   //-------------------------------------------
01244   // Finally loop over all the SCALAR variables
01245   for (unsigned int data_var=0; data_var<nv; data_var++)
01246     {
01247       const unsigned int var = _written_var_indices[data_var];
01248       if(this->variable(var).type().family == SCALAR)
01249         {
01250 #ifndef NDEBUG
01251           n_assigned_vals +=
01252 #endif
01253             this->read_SCALAR_dofs (var, io, vec);
01254         }
01255     }
01256 
01257   vec.close();
01258 
01259 #ifndef NDEBUG
01260   this->comm().sum (n_assigned_vals);
01261   libmesh_assert_equal_to (n_assigned_vals, vector_length);
01262 #endif
01263 
01264   return vector_length;
01265 }
01266 
01267 
01268 
01269 void System::write_header (Xdr& io,
01270                            const std::string & /* version is currently unused */,
01271                            const bool write_additional_data) const
01272 {
01306   libmesh_assert (io.writing());
01307 
01308 
01309   // Only write the header information
01310   // if we are processor 0.
01311   if (this->get_mesh().processor_id() != 0)
01312     return;
01313 
01314   std::string comment;
01315   char buf[80];
01316 
01317   // 5.)
01318   // Write the number of variables in the system
01319 
01320   {
01321     // set up the comment
01322     comment = "# No. of Variables in System \"";
01323     comment += this->name();
01324     comment += "\"";
01325 
01326     unsigned int nv = this->n_vars();
01327     io.data (nv, comment.c_str());
01328   }
01329 
01330 
01331   for (unsigned int var=0; var<this->n_vars(); var++)
01332     {
01333       // 6.)
01334       // Write the name of the var-th variable
01335       {
01336         // set up the comment
01337         comment  = "#   Name, Variable No. ";
01338         std::sprintf(buf, "%u", var);
01339         comment += buf;
01340         comment += ", System \"";
01341         comment += this->name();
01342         comment += "\"";
01343 
01344         std::string var_name = this->variable_name(var);
01345         io.data (var_name, comment.c_str());
01346       }
01347 
01348       // 6.1.) Variable subdomains
01349       {
01350         // set up the comment
01351         comment  = "#     Subdomains, Variable \"";
01352         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01353         comment += buf;
01354         comment += "\", System \"";
01355         comment += this->name();
01356         comment += "\"";
01357 
01358         const std::set<subdomain_id_type> & domains = this->variable(var).active_subdomains();
01359         std::vector<subdomain_id_type> domain_array;
01360         domain_array.assign(domains.begin(), domains.end());
01361         io.data (domain_array, comment.c_str());
01362       }
01363 
01364       // 7.)
01365       // Write the approximation order of the var-th variable
01366       // in this system
01367       {
01368         // set up the comment
01369         comment = "#     Approximation Order, Variable \"";
01370         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01371         comment += buf;
01372         comment += "\", System \"";
01373         comment += this->name();
01374         comment += "\"";
01375 
01376         int order = static_cast<int>(this->variable_type(var).order);
01377         io.data (order, comment.c_str());
01378       }
01379 
01380 
01381 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01382 
01383       // do the same for radial_order
01384       {
01385         comment = "#     Radial Approximation Order, Variable \"";
01386         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01387         comment += buf;
01388         comment += "\", System \"";
01389         comment += this->name();
01390         comment += "\"";
01391 
01392         int rad_order = static_cast<int>(this->variable_type(var).radial_order);
01393         io.data (rad_order, comment.c_str());
01394       }
01395 
01396 #endif
01397 
01398       // Write the Finite Element type of the var-th variable
01399       // in this System
01400       {
01401         // set up the comment
01402         comment = "#     FE Family, Variable \"";
01403         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01404         comment += buf;
01405         comment += "\", System \"";
01406         comment += this->name();
01407         comment += "\"";
01408 
01409         const FEType& type = this->variable_type(var);
01410         int fam = static_cast<int>(type.family);
01411         io.data (fam, comment.c_str());
01412 
01413 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01414 
01415         comment = "#     Radial FE Family, Variable \"";
01416         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01417         comment += buf;
01418         comment += "\", System \"";
01419         comment += this->name();
01420         comment += "\"";
01421 
01422         int radial_fam = static_cast<int>(type.radial_family);
01423         io.data (radial_fam, comment.c_str());
01424 
01425         comment = "#     Infinite Mapping Type, Variable \"";
01426         std::sprintf(buf, "%s", this->variable_name(var).c_str());
01427         comment += buf;
01428         comment += "\", System \"";
01429         comment += this->name();
01430         comment += "\"";
01431 
01432         int i_map = static_cast<int>(type.inf_map);
01433         io.data (i_map, comment.c_str());
01434 #endif
01435       }
01436     } // end of the variable loop
01437 
01438   // 8.)
01439   // Write the number of additional vectors in the System.
01440   // If write_additional_data==false, then write zero for
01441   // the number of additional vectors.
01442   {
01443     {
01444       // set up the comment
01445       comment = "# No. of Additional Vectors, System \"";
01446       comment += this->name();
01447       comment += "\"";
01448 
01449       unsigned int nvecs = write_additional_data ? this->n_vectors () : 0;
01450       io.data (nvecs, comment.c_str());
01451     }
01452 
01453     if (write_additional_data)
01454       {
01455         std::map<std::string, NumericVector<Number>* >::const_iterator
01456           vec_pos = this->_vectors.begin();
01457         unsigned int cnt=0;
01458 
01459         for (; vec_pos != this->_vectors.end(); ++vec_pos)
01460           {
01461             // 9.)
01462             // write the name of the cnt-th additional vector
01463             comment =  "# Name of ";
01464             std::sprintf(buf, "%d", cnt++);
01465             comment += buf;
01466             comment += "th vector";
01467             std::string vec_name = vec_pos->first;
01468 
01469             io.data (vec_name, comment.c_str());
01470           }
01471       }
01472   }
01473 }
01474 
01475 
01476 
01477 void System::write_parallel_data (Xdr &io,
01478                                   const bool write_additional_data) const
01479 {
01499   // PerfLog pl("IO Performance",false);
01500   // pl.push("write_parallel_data");
01501   // std::size_t total_written_size = 0;
01502 
01503   std::string comment;
01504 
01505   libmesh_assert (io.writing());
01506 
01507   std::vector<Number> io_buffer; io_buffer.reserve(this->solution->local_size());
01508 
01509   // build the ordered nodes and element maps.
01510   // when writing/reading parallel files we need to iterate
01511   // over our nodes/elements in order of increasing global id().
01512   // however, this is not guaranteed to be ordering we obtain
01513   // by using the node_iterators/element_iterators directly.
01514   // so build a set, sorted by id(), that provides the ordering.
01515   // further, for memory economy build the set but then transfer
01516   // its contents to vectors, which will be sorted.
01517   std::vector<const DofObject*> ordered_nodes, ordered_elements;
01518   {
01519     std::set<const DofObject*, CompareDofObjectsByID>
01520       ordered_nodes_set (this->get_mesh().local_nodes_begin(),
01521                          this->get_mesh().local_nodes_end());
01522 
01523     ordered_nodes.insert(ordered_nodes.end(),
01524                          ordered_nodes_set.begin(),
01525                          ordered_nodes_set.end());
01526   }
01527   {
01528     std::set<const DofObject*, CompareDofObjectsByID>
01529       ordered_elements_set (this->get_mesh().local_elements_begin(),
01530                             this->get_mesh().local_elements_end());
01531 
01532     ordered_elements.insert(ordered_elements.end(),
01533                             ordered_elements_set.begin(),
01534                             ordered_elements_set.end());
01535   }
01536 
01537   const unsigned int sys_num = this->number();
01538   const unsigned int nv      = this->n_vars();
01539 
01540   // Loop over each non-SCALAR variable and each node, and write out the value.
01541   for (unsigned int var=0; var<nv; var++)
01542     if (this->variable(var).type().family != SCALAR)
01543       {
01544         // First write the node DOF values
01545         for (std::vector<const DofObject*>::const_iterator
01546                it = ordered_nodes.begin(); it != ordered_nodes.end(); ++it)
01547           for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
01548             {
01549               //libMesh::out << "(*it)->id()=" << (*it)->id() << std::endl;
01550               libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
01551                                            DofObject::invalid_id);
01552 
01553               io_buffer.push_back((*this->solution)((*it)->dof_number(sys_num, var, comp)));
01554             }
01555 
01556         // Then write the element DOF values
01557         for (std::vector<const DofObject*>::const_iterator
01558                it = ordered_elements.begin(); it != ordered_elements.end(); ++it)
01559           for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
01560             {
01561               libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
01562                                            DofObject::invalid_id);
01563 
01564               io_buffer.push_back((*this->solution)((*it)->dof_number(sys_num, var, comp)));
01565             }
01566       }
01567 
01568   // Finally, write the SCALAR data on the last processor
01569   for (unsigned int var=0; var<this->n_vars(); var++)
01570     if(this->variable(var).type().family == SCALAR)
01571       {
01572         if (this->processor_id() == (this->n_processors()-1))
01573           {
01574             const DofMap& dof_map = this->get_dof_map();
01575             std::vector<dof_id_type> SCALAR_dofs;
01576             dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
01577 
01578             for(unsigned int i=0; i<SCALAR_dofs.size(); i++)
01579               {
01580                 io_buffer.push_back( (*this->solution)(SCALAR_dofs[i]) );
01581               }
01582           }
01583       }
01584 
01585   // 9.)
01586   //
01587   // Actually write the reordered solution vector
01588   // for the ith system to disk
01589 
01590   // set up the comment
01591   {
01592     comment = "# System \"";
01593     comment += this->name();
01594     comment += "\" Solution Vector";
01595   }
01596 
01597   io.data (io_buffer, comment.c_str());
01598 
01599   // total_written_size += io_buffer.size();
01600 
01601   // Only write additional vectors if wanted
01602   if (write_additional_data)
01603     {
01604       std::map<std::string, NumericVector<Number>* >::const_iterator
01605         pos = _vectors.begin();
01606 
01607       for(; pos != this->_vectors.end(); ++pos)
01608         {
01609           io_buffer.clear(); io_buffer.reserve( pos->second->local_size());
01610 
01611           // Loop over each non-SCALAR variable and each node, and write out the value.
01612           for (unsigned int var=0; var<nv; var++)
01613             if(this->variable(var).type().family != SCALAR)
01614               {
01615                 // First write the node DOF values
01616                 for (std::vector<const DofObject*>::const_iterator
01617                        it = ordered_nodes.begin(); it != ordered_nodes.end(); ++it)
01618                   for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
01619                     {
01620                       libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
01621                                                    DofObject::invalid_id);
01622 
01623                       io_buffer.push_back((*pos->second)((*it)->dof_number(sys_num, var, comp)));
01624                     }
01625 
01626                 // Then write the element DOF values
01627                 for (std::vector<const DofObject*>::const_iterator
01628                        it = ordered_elements.begin(); it != ordered_elements.end(); ++it)
01629                   for (unsigned int comp=0; comp<(*it)->n_comp(sys_num, var); comp++)
01630                     {
01631                       libmesh_assert_not_equal_to ((*it)->dof_number(sys_num, var, comp),
01632                                                    DofObject::invalid_id);
01633 
01634                       io_buffer.push_back((*pos->second)((*it)->dof_number(sys_num, var, comp)));
01635                     }
01636               }
01637 
01638           // Finally, write the SCALAR data on the last processor
01639           for (unsigned int var=0; var<this->n_vars(); var++)
01640             if(this->variable(var).type().family == SCALAR)
01641               {
01642                 if (this->processor_id() == (this->n_processors()-1))
01643                   {
01644                     const DofMap& dof_map = this->get_dof_map();
01645                     std::vector<dof_id_type> SCALAR_dofs;
01646                     dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
01647 
01648                     for(unsigned int i=0; i<SCALAR_dofs.size(); i++)
01649                       {
01650                         io_buffer.push_back( (*pos->second)(SCALAR_dofs[i]) );
01651                       }
01652                   }
01653               }
01654 
01655           // 10.)
01656           //
01657           // Actually write the reordered additional vector
01658           // for this system to disk
01659 
01660           // set up the comment
01661           {
01662             comment = "# System \"";
01663             comment += this->name();
01664             comment += "\" Additional Vector \"";
01665             comment += pos->first;
01666             comment += "\"";
01667           }
01668 
01669           io.data (io_buffer, comment.c_str());
01670 
01671           // total_written_size += io_buffer.size();
01672         }
01673     }
01674 
01675   // const Real
01676   //   dt   = pl.get_elapsed_time(),
01677   //   rate = total_written_size*sizeof(Number)/dt;
01678 
01679   // libMesh::err << "Write " << total_written_size << " \"Number\" values\n"
01680   //     << " Elapsed time = " << dt << '\n'
01681   //     << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
01682 
01683   // pl.pop("write_parallel_data");
01684 }
01685 
01686 
01687 
01688 void System::write_serialized_data (Xdr& io,
01689                                     const bool write_additional_data) const
01690 {
01704   parallel_object_only();
01705   std::string comment;
01706 
01707   // PerfLog pl("IO Performance",false);
01708   // pl.push("write_serialized_data");
01709   // std::size_t total_written_size = 0;
01710 
01711   // total_written_size +=
01712   this->write_serialized_vector(io, *this->solution);
01713 
01714   // set up the comment
01715   if (this->processor_id() == 0)
01716     {
01717       comment = "# System \"";
01718       comment += this->name();
01719       comment += "\" Solution Vector";
01720 
01721       io.comment (comment);
01722     }
01723 
01724   // Only write additional vectors if wanted
01725   if (write_additional_data)
01726     {
01727       std::map<std::string, NumericVector<Number>* >::const_iterator
01728         pos = _vectors.begin();
01729 
01730       for(; pos != this->_vectors.end(); ++pos)
01731         {
01732           // total_written_size +=
01733           this->write_serialized_vector(io, *pos->second);
01734 
01735           // set up the comment
01736           if (this->processor_id() == 0)
01737             {
01738               comment = "# System \"";
01739               comment += this->name();
01740               comment += "\" Additional Vector \"";
01741               comment += pos->first;
01742               comment += "\"";
01743               io.comment (comment);
01744             }
01745         }
01746     }
01747 
01748   // const Real
01749   //   dt   = pl.get_elapsed_time(),
01750   //   rate = total_written_size*sizeof(Number)/dt;
01751 
01752   // libMesh::out << "Write " << total_written_size << " \"Number\" values\n"
01753   //     << " Elapsed time = " << dt << '\n'
01754   //     << " Rate = " << rate/1.e6 << "(MB/sec)\n\n";
01755 
01756   // pl.pop("write_serialized_data");
01757 
01758 
01759 
01760 
01761   // // test the new method
01762   // {
01763   //   std::vector<std::string> names;
01764   //   std::vector<NumericVector<Number>*> vectors_to_write;
01765 
01766   //   names.push_back("Solution Vector");
01767   //   vectors_to_write.push_back(this->solution.get());
01768 
01769   //   // Only write additional vectors if wanted
01770   //   if (write_additional_data)
01771   //     {
01772   // std::map<std::string, NumericVector<Number>* >::const_iterator
01773   //   pos = _vectors.begin();
01774 
01775   // for(; pos != this->_vectors.end(); ++pos)
01776   //   {
01777   //     names.push_back("Additional Vector " + pos->first);
01778   //     vectors_to_write.push_back(pos->second);
01779   //   }
01780   //     }
01781 
01782   //   total_written_size =
01783   //     this->write_serialized_vectors (io, names, vectors_to_write);
01784 
01785   //   const Real
01786   //     dt2   = pl.get_elapsed_time(),
01787   //     rate2 = total_written_size*sizeof(Number)/(dt2-dt);
01788 
01789   //   libMesh::out << "Write (new) " << total_written_size << " \"Number\" values\n"
01790   //       << " Elapsed time = " << (dt2-dt) << '\n'
01791   //       << " Rate = " << rate2/1.e6 << "(MB/sec)\n\n";
01792 
01793   // }
01794 }
01795 
01796 
01797 
01798 template <typename iterator_type>
01799 std::size_t System::write_serialized_blocked_dof_objects (const std::vector<const NumericVector<Number>*> &vecs,
01800                                                           const dof_id_type n_objs,
01801                                                           const iterator_type begin,
01802                                                           const iterator_type end,
01803                                                           Xdr &io,
01804                                                           const unsigned int var_to_write) const
01805 {
01806   //-------------------------------------------------------
01807   // General order: (IO format 0.7.4 & greater)
01808   //
01809   // for (objects ...)
01810   //   for (vecs ....)
01811   //     for (vars ....)
01812   //       for (comps ...)
01813   //
01814   // where objects are nodes or elements, sorted to be
01815   // partition independent,
01816   // vecs are one or more *identically distributed* solution
01817   // coefficient vectors, vars are one or more variables
01818   // to write, and comps are all the components for said
01819   // vars on the object.
01820 
01821   typedef std::vector<const NumericVector<Number>*>::const_iterator vec_iterator_type;
01822 
01823   // We will write all variables unless requested otherwise.
01824   std::vector<unsigned int> vars_to_write(1, var_to_write);
01825 
01826   if (var_to_write == libMesh::invalid_uint)
01827     {
01828       vars_to_write.clear();  vars_to_write.reserve(this->n_vars());
01829       for (unsigned int var=0; var<this->n_vars(); var++)
01830         vars_to_write.push_back(var);
01831     }
01832 
01833   const dof_id_type io_blksize = cast_int<dof_id_type>
01834     (std::min(max_io_blksize, static_cast<std::size_t>(n_objs)));
01835 
01836   const unsigned int
01837     sys_num    = this->number(),
01838     num_vecs   = cast_int<unsigned int>(vecs.size()),
01839     num_blks   = std::ceil(static_cast<double>(n_objs)/static_cast<double>(io_blksize));
01840 
01841   // libMesh::out << "io_blksize = "    << io_blksize
01842   //     << ", num_objects = " << n_objs
01843   //     << ", num_blks = "    << num_blks
01844   //     << std::endl;
01845 
01846   dof_id_type written_length=0;                                   // The numer of values written.  This will be returned
01847   std::vector<std::vector<dof_id_type> > xfer_ids(num_blks);      // The global IDs and # of components for the local objects in all blocks
01848   std::vector<std::vector<Number> >       send_vals(num_blks);     // The raw values for the local objects in all blocks
01849   std::vector<Parallel::Request>
01850     id_requests(num_blks), val_requests(num_blks);                 // send request handle for each block
01851 
01852   // ------------------------------------------------------
01853   // First pass - count the number of objects in each block
01854   // traverse all the objects and figure out which block they
01855   // will utlimately live in.
01856   std::vector<unsigned int>
01857     xfer_ids_size  (num_blks,0),
01858     send_vals_size (num_blks,0);
01859 
01860   for (iterator_type it=begin; it!=end; ++it)
01861     {
01862       const dof_id_type
01863         id    = (*it)->id(),
01864         block = id/io_blksize;
01865 
01866       libmesh_assert_less (block, num_blks);
01867 
01868       xfer_ids_size[block] += 2; // for each object, we store its id, as well as the total number of components for all variables
01869 
01870       unsigned int n_comp_tot=0;
01871 
01872       for (std::vector<unsigned int>::const_iterator var_it=vars_to_write.begin();
01873            var_it!=vars_to_write.end(); ++var_it)
01874         n_comp_tot += (*it)->n_comp(sys_num, *var_it); // for each variable, we will store the nonzero components
01875 
01876       send_vals_size[block] += n_comp_tot*num_vecs;
01877     }
01878 
01879   //-----------------------------------------
01880   // Collect the values for all local objects,
01881   // binning them into 'blocks' that will be
01882   // sent to processor 0
01883   for (unsigned int blk=0; blk<num_blks; blk++)
01884     {
01885       // libMesh::out << "Writing object block " << blk << std::endl;
01886 
01887       // Each processor should build up its transfer buffers for its
01888       // local objects in [first_object,last_object).
01889       const dof_id_type
01890         first_object = blk*io_blksize,
01891         last_object  = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs);
01892 
01893       // convenience
01894       std::vector<dof_id_type> &ids  (xfer_ids[blk]);
01895       std::vector<Number>      &vals (send_vals[blk]);
01896 
01897       // we now know the number of values we will store for each block,
01898       // so we can do efficient preallocation
01899       ids.clear();   ids.reserve  (xfer_ids_size[blk]);
01900       vals.clear();  vals.reserve (send_vals_size[blk]);
01901 
01902       if (send_vals_size[blk] != 0) // only send if we have nonzero components to write
01903         for (iterator_type it=begin; it!=end; ++it)
01904           if (((*it)->id() >= first_object) && // object in [first_object,last_object)
01905               ((*it)->id() <   last_object))
01906             {
01907               ids.push_back((*it)->id());
01908 
01909               // count the total number of nonzeros transferred for this object
01910               {
01911                 unsigned int n_comp_tot=0;
01912 
01913                 for (std::vector<unsigned int>::const_iterator var_it=vars_to_write.begin();
01914                      var_it!=vars_to_write.end(); ++var_it)
01915                   n_comp_tot += (*it)->n_comp(sys_num,*var_it);
01916 
01917                 ids.push_back (n_comp_tot*num_vecs); // even if 0 - processor 0 has no way of knowing otherwise...
01918               }
01919 
01920               // pack the values to send
01921               for (vec_iterator_type vec_it=vecs.begin(); vec_it!=vecs.end(); ++vec_it)
01922                 {
01923                   const NumericVector<Number> &vec(**vec_it);
01924 
01925                   for (std::vector<unsigned int>::const_iterator var_it=vars_to_write.begin();
01926                        var_it!=vars_to_write.end(); ++var_it)
01927                     {
01928                       const unsigned int n_comp = (*it)->n_comp(sys_num,*var_it);
01929 
01930                       for (unsigned int comp=0; comp<n_comp; comp++)
01931                         {
01932                           libmesh_assert_greater_equal ((*it)->dof_number(sys_num, *var_it, comp), vec.first_local_index());
01933                           libmesh_assert_less ((*it)->dof_number(sys_num, *var_it, comp), vec.last_local_index());
01934                           vals.push_back(vec((*it)->dof_number(sys_num, *var_it, comp)));
01935                         }
01936                     }
01937                 }
01938             }
01939 
01940 #ifdef LIBMESH_HAVE_MPI
01941       Parallel::MessageTag id_tag  = this->comm().get_unique_tag(100*num_blks + blk);
01942       Parallel::MessageTag val_tag = this->comm().get_unique_tag(200*num_blks + blk);
01943 
01944       // nonblocking send the data for this block
01945       this->comm().send (0, ids,  id_requests[blk],  id_tag);
01946       this->comm().send (0, vals, val_requests[blk], val_tag);
01947 #endif
01948     }
01949 
01950 
01951   if (this->processor_id() == 0)
01952     {
01953       std::vector<std::vector<dof_id_type> > recv_ids  (this->n_processors());
01954       std::vector<std::vector<Number> >      recv_vals (this->n_processors());
01955       std::vector<unsigned int> obj_val_offsets;          // map to traverse entry-wise rather than processor-wise
01956       std::vector<Number>       output_vals;              // The output buffer for the current block
01957 
01958       // a ThreadedIO object to perform asychronous file IO
01959       ThreadedIO<Number> threaded_io(io, output_vals);
01960       UniquePtr<Threads::Thread> async_io;
01961 
01962       for (unsigned int blk=0; blk<num_blks; blk++)
01963         {
01964           // Each processor should build up its transfer buffers for its
01965           // local objects in [first_object,last_object).
01966           const dof_id_type
01967             first_object  = cast_int<dof_id_type>(blk*io_blksize),
01968             last_object   = std::min(cast_int<dof_id_type>((blk+1)*io_blksize), n_objs),
01969             n_objects_blk = last_object - first_object;
01970 
01971           // offset array. this will define where each object's values
01972           // map into the actual output_vals buffer.  this must get
01973           // 0-initialized because 0-component objects are not actually sent
01974           obj_val_offsets.resize (n_objects_blk);  std::fill (obj_val_offsets.begin(), obj_val_offsets.end(), 0);
01975 
01976           std::size_t n_val_recvd_blk=0;
01977 
01978           // tags to select data received
01979           Parallel::MessageTag id_tag  (this->comm().get_unique_tag(100*num_blks + blk));
01980           Parallel::MessageTag val_tag (this->comm().get_unique_tag(200*num_blks + blk));
01981 
01982           // receive this block of data from all processors.
01983           for (unsigned int comm_step=0; comm_step<this->n_processors(); comm_step++)
01984             {
01985 #ifdef LIBMESH_HAVE_MPI
01986               // blocking receive indices for this block, imposing no particular order on processor
01987               Parallel::Status id_status (this->comm().probe (Parallel::any_source, id_tag));
01988               std::vector<dof_id_type> &ids (recv_ids[id_status.source()]);
01989               this->comm().receive (id_status.source(), ids,  id_tag);
01990 #else
01991               std::vector<dof_id_type> &ids (recv_ids[0]);
01992 #endif
01993 
01994               // note its possible we didn't receive values for objects in
01995               // this block if they have no components allocated.
01996               for (dof_id_type idx=0; idx<ids.size(); idx+=2)
01997                 {
01998                   const dof_id_type
01999                     local_idx          = ids[idx+0]-first_object,
02000                     n_vals_tot_allvecs = ids[idx+1];
02001 
02002                   libmesh_assert_less (local_idx, n_objects_blk);
02003                   libmesh_assert_less (local_idx, obj_val_offsets.size());
02004 
02005                   obj_val_offsets[local_idx] = n_vals_tot_allvecs;
02006                 }
02007 
02008 #ifdef LIBMESH_HAVE_MPI
02009               // blocking receive values for this block, imposing no particular order on processor
02010               Parallel::Status val_status  (this->comm().probe (Parallel::any_source, val_tag));
02011               std::vector<Number> &vals    (recv_vals[val_status.source()]);
02012               this->comm().receive (val_status.source(), vals, val_tag);
02013 #else
02014               // straight copy without MPI
02015               std::vector<Number> &vals (recv_vals[0]);
02016               vals = send_vals[blk];
02017 #endif
02018 
02019               n_val_recvd_blk += vals.size();
02020             }
02021 
02022           // We need the offests into the output_vals vector for each object.
02023           // fortunately, this is simply the partial sum of the total number
02024           // of components for each object
02025           std::partial_sum(obj_val_offsets.begin(), obj_val_offsets.end(),
02026                            obj_val_offsets.begin());
02027 
02028           // wait on any previous asynchronous IO - this *must* complete before
02029           // we start messing with the output_vals buffer!
02030           if (async_io.get()) async_io->join();
02031 
02032           // this is the actual output buffer that will be written to disk.
02033           // at ths point we finally know wha size it will be.
02034           output_vals.resize(n_val_recvd_blk);
02035 
02036           // pack data from all processors into output values
02037           for (unsigned int proc=0; proc<this->n_processors(); proc++)
02038             {
02039               const std::vector<dof_id_type> &ids (recv_ids [proc]);
02040               const std::vector<Number>      &vals(recv_vals[proc]);
02041               std::vector<Number>::const_iterator proc_vals(vals.begin());
02042 
02043               for (dof_id_type idx=0; idx<ids.size(); idx+=2)
02044                 {
02045                   const dof_id_type
02046                     local_idx          = ids[idx+0]-first_object,
02047                     n_vals_tot_allvecs = ids[idx+1];
02048 
02049                   // put this object's data into the proper location
02050                   // in  the output buffer
02051                   std::vector<Number>::iterator out_vals(output_vals.begin());
02052                   if (local_idx != 0)
02053                     std::advance (out_vals, obj_val_offsets[local_idx-1]);
02054 
02055                   for (unsigned int val=0; val<n_vals_tot_allvecs; val++, ++out_vals, ++proc_vals)
02056                     {
02057                       libmesh_assert (out_vals  != output_vals.end());
02058                       libmesh_assert (proc_vals != vals.end());
02059                       *out_vals = *proc_vals;
02060                     }
02061                 }
02062             }
02063 
02064           // output_vals buffer is now filled for this block.
02065           // write it to disk
02066           async_io.reset(new Threads::Thread(threaded_io));
02067           written_length += output_vals.size();
02068         }
02069 
02070       // wait on any previous asynchronous IO - this *must* complete before
02071       // our stuff goes out of scope
02072       async_io->join();
02073     }
02074 
02075   Parallel::wait(id_requests);
02076   Parallel::wait(val_requests);
02077 
02078   // we need some synchronization here.  Because this method
02079   // can be called for a range of nodes, then a range of elements,
02080   // we need some mechanism to prevent processors from racing past
02081   // to the next range and overtaking ongoing communication.  one
02082   // approach would be to figure out unique tags for each range,
02083   // but for now we just impose a barrier here.  And might as
02084   // well have it do some useful work.
02085   this->comm().broadcast(written_length);
02086 
02087   return written_length;
02088 }
02089 
02090 
02091 
02092 unsigned int System::write_SCALAR_dofs (const NumericVector<Number> &vec,
02093                                         const unsigned int var,
02094                                         Xdr &io) const
02095 {
02096   unsigned int written_length=0;
02097   std::vector<Number> vals; // The raw values for the local objects in the current block
02098   // Collect the SCALARs for the current variable
02099   if (this->processor_id() == (this->n_processors()-1))
02100     {
02101       const DofMap& dof_map = this->get_dof_map();
02102       std::vector<dof_id_type> SCALAR_dofs;
02103       dof_map.SCALAR_dof_indices(SCALAR_dofs, var);
02104       const unsigned int n_scalar_dofs = cast_int<unsigned int>
02105         (SCALAR_dofs.size());
02106 
02107       for(unsigned int i=0; i<n_scalar_dofs; i++)
02108         {
02109           vals.push_back( vec(SCALAR_dofs[i]) );
02110         }
02111     }
02112 
02113 #ifdef LIBMESH_HAVE_MPI
02114   if ( this->n_processors() > 1 )
02115     {
02116       const Parallel::MessageTag val_tag(1);
02117 
02118       // Post the receive on processor 0
02119       if ( this->processor_id() == 0 )
02120         {
02121           this->comm().receive(this->n_processors()-1, vals, val_tag);
02122         }
02123 
02124       // Send the data to processor 0
02125       if (this->processor_id() == (this->n_processors()-1))
02126         {
02127           this->comm().send(0, vals, val_tag);
02128         }
02129     }
02130 #endif
02131 
02132   // -------------------------------------------------------
02133   // Write the output on processor 0.
02134   if (this->processor_id() == 0)
02135     {
02136       const unsigned int vals_size =
02137         cast_int<unsigned int>(vals.size());
02138       io.data_stream (&vals[0], vals_size);
02139       written_length += vals_size;
02140     }
02141 
02142   return written_length;
02143 }
02144 
02145 
02146 
02147 dof_id_type System::write_serialized_vector (Xdr& io, const NumericVector<Number>& vec) const
02148 {
02149   parallel_object_only();
02150 
02151   libmesh_assert (io.writing());
02152 
02153   dof_id_type vec_length = vec.size();
02154   if (this->processor_id() == 0) io.data (vec_length, "# vector length");
02155 
02156   dof_id_type written_length = 0;
02157 
02158   //---------------------------------
02159   // Collect the values for all nodes
02160   written_length += cast_int<dof_id_type>
02161     (this->write_serialized_blocked_dof_objects (std::vector<const NumericVector<Number>*>(1,&vec),
02162                                                  this->get_mesh().n_nodes(),
02163                                                  this->get_mesh().local_nodes_begin(),
02164                                                  this->get_mesh().local_nodes_end(),
02165                                                  io));
02166 
02167   //------------------------------------
02168   // Collect the values for all elements
02169   written_length += cast_int<dof_id_type>
02170     (this->write_serialized_blocked_dof_objects (std::vector<const NumericVector<Number>*>(1,&vec),
02171                                                  this->get_mesh().n_elem(),
02172                                                  this->get_mesh().local_elements_begin(),
02173                                                  this->get_mesh().local_elements_end(),
02174                                                  io));
02175 
02176   //-------------------------------------------
02177   // Finally loop over all the SCALAR variables
02178   for (unsigned int var=0; var<this->n_vars(); var++)
02179     if(this->variable(var).type().family == SCALAR)
02180       {
02181         written_length +=
02182           this->write_SCALAR_dofs (vec, var, io);
02183       }
02184 
02185   if (this->processor_id() == 0)
02186     libmesh_assert_equal_to (written_length, vec_length);
02187 
02188   return written_length;
02189 }
02190 
02191 
02192 template <typename InValType>
02193 std::size_t System::read_serialized_vectors (Xdr &io,
02194                                              const std::vector<NumericVector<Number>*> &vectors) const
02195 {
02196   parallel_object_only();
02197 
02198   // Error checking
02199   // #ifndef NDEBUG
02200   //   // In parallel we better be reading a parallel vector -- if not
02201   //   // we will not set all of its components below!!
02202   //   if (this->n_processors() > 1)
02203   //     {
02204   //       libmesh_assert (vec.type() == PARALLEL ||
02205   //       vec.type() == GHOSTED);
02206   //     }
02207   // #endif
02208 
02209   libmesh_assert (io.reading());
02210 
02211   if (this->processor_id() == 0)
02212     {
02213       // sizes
02214       unsigned int num_vecs=0;
02215       dof_id_type vector_length=0;
02216 
02217       // Get the number of vectors
02218       io.data(num_vecs);
02219       // Get the buffer size
02220       io.data(vector_length);
02221 
02222       libmesh_assert_equal_to (num_vecs, vectors.size());
02223 
02224       if (num_vecs != 0)
02225         {
02226           libmesh_assert_not_equal_to (vectors[0], 0);
02227           libmesh_assert_equal_to     (vectors[0]->size(), vector_length);
02228         }
02229     }
02230 
02231   // no need to actually communicate these.
02232   // this->comm().broadcast(num_vecs);
02233   // this->comm().broadcast(vector_length);
02234 
02235   // Cache these - they are not free!
02236   const dof_id_type
02237     n_nodes = this->get_mesh().n_nodes(),
02238     n_elem  = this->get_mesh().n_elem();
02239 
02240   std::size_t read_length = 0.;
02241 
02242   //---------------------------------
02243   // Collect the values for all nodes
02244   read_length +=
02245     this->read_serialized_blocked_dof_objects (n_nodes,
02246                                                this->get_mesh().local_nodes_begin(),
02247                                                this->get_mesh().local_nodes_end(),
02248                                                InValType(),
02249                                                io,
02250                                                vectors);
02251 
02252   //------------------------------------
02253   // Collect the values for all elements
02254   read_length +=
02255     this->read_serialized_blocked_dof_objects (n_elem,
02256                                                this->get_mesh().local_elements_begin(),
02257                                                this->get_mesh().local_elements_end(),
02258                                                InValType(),
02259                                                io,
02260                                                vectors);
02261 
02262   //-------------------------------------------
02263   // Finally loop over all the SCALAR variables
02264   for (unsigned int vec=0; vec<vectors.size(); vec++)
02265     for (unsigned int var=0; var<this->n_vars(); var++)
02266       if(this->variable(var).type().family == SCALAR)
02267         {
02268           libmesh_assert_not_equal_to (vectors[vec], 0);
02269 
02270           read_length +=
02271             this->read_SCALAR_dofs (var, io, *vectors[vec]);
02272         }
02273 
02274   //---------------------------------------
02275   // last step - must close all the vectors
02276   for (unsigned int vec=0; vec<vectors.size(); vec++)
02277     {
02278       libmesh_assert_not_equal_to (vectors[vec], 0);
02279       vectors[vec]->close();
02280     }
02281 
02282   return read_length;
02283 }
02284 
02285 
02286 
02287 std::size_t System::write_serialized_vectors (Xdr &io,
02288                                               const std::vector<const NumericVector<Number>*> &vectors) const
02289 {
02290   parallel_object_only();
02291 
02292   libmesh_assert (io.writing());
02293 
02294   // Cache these - they are not free!
02295   const dof_id_type
02296     n_nodes       = this->get_mesh().n_nodes(),
02297     n_elem        = this->get_mesh().n_elem();
02298 
02299   std::size_t written_length = 0.;
02300 
02301   if (this->processor_id() == 0)
02302     {
02303       unsigned int
02304         n_vec    = cast_int<unsigned int>(vectors.size());
02305       dof_id_type
02306         vec_size = vectors.empty() ? 0 : vectors[0]->size();
02307       // Set the number of vectors
02308       io.data(n_vec, "# number of vectors");
02309       // Set the buffer size
02310       io.data(vec_size, "# vector length");
02311     }
02312 
02313   //---------------------------------
02314   // Collect the values for all nodes
02315   written_length +=
02316     this->write_serialized_blocked_dof_objects (vectors,
02317                                                 n_nodes,
02318                                                 this->get_mesh().local_nodes_begin(),
02319                                                 this->get_mesh().local_nodes_end(),
02320                                                 io);
02321 
02322   //------------------------------------
02323   // Collect the values for all elements
02324   written_length +=
02325     this->write_serialized_blocked_dof_objects (vectors,
02326                                                 n_elem,
02327                                                 this->get_mesh().local_elements_begin(),
02328                                                 this->get_mesh().local_elements_end(),
02329                                                 io);
02330 
02331   //-------------------------------------------
02332   // Finally loop over all the SCALAR variables
02333   for (unsigned int vec=0; vec<vectors.size(); vec++)
02334     for (unsigned int var=0; var<this->n_vars(); var++)
02335       if(this->variable(var).type().family == SCALAR)
02336         {
02337           libmesh_assert_not_equal_to (vectors[vec], 0);
02338 
02339           written_length +=
02340             this->write_SCALAR_dofs (*vectors[vec], var, io);
02341         }
02342 
02343   return written_length;
02344 }
02345 
02346 
02347 
02348 
02349 template void System::read_parallel_data<Number> (Xdr &io, const bool read_additional_data);
02350 template void System::read_serialized_data<Number> (Xdr& io, const bool read_additional_data);
02351 template numeric_index_type System::read_serialized_vector<Number> (Xdr& io, NumericVector<Number>& vec);
02352 template std::size_t System::read_serialized_vectors<Number> (Xdr &io, const std::vector<NumericVector<Number>*> &vectors) const;
02353 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
02354 template void System::read_parallel_data<Real> (Xdr &io, const bool read_additional_data);
02355 template void System::read_serialized_data<Real> (Xdr& io, const bool read_additional_data);
02356 template numeric_index_type System::read_serialized_vector<Real> (Xdr& io, NumericVector<Number>& vec);
02357 template std::size_t System::read_serialized_vectors<Real> (Xdr &io, const std::vector<NumericVector<Number>*> &vectors) const;
02358 #endif
02359 
02360 } // namespace libMesh