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