libMesh::XdrIO Class Reference

#include <xdr_io.h>

Inheritance diagram for libMesh::XdrIO:

List of all members.

Public Member Functions

 XdrIO (MeshBase &, const bool=false)
 XdrIO (const MeshBase &, const bool=false)
virtual ~XdrIO ()
virtual void read (const std::string &)
virtual void write (const std::string &)
bool binary () const
bool & binary ()
bool legacy () const
bool & legacy ()
bool write_parallel () const
void set_write_parallel (bool do_parallel=true)
void set_auto_parallel ()
const std::string & version () const
std::string & version ()
const std::string & boundary_condition_file_name () const
std::string & boundary_condition_file_name ()
const std::string & partition_map_file_name () const
std::string & partition_map_file_name ()
const std::string & subdomain_map_file_name () const
std::string & subdomain_map_file_name ()
const std::string & polynomial_level_file_name () const
std::string & polynomial_level_file_name ()
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL)
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
unsigned int & ascii_precision ()

Protected Member Functions

MeshBasemesh ()
void skip_comment_lines (std::istream &in, const char comment_start)
const MeshBasemesh () const

Protected Attributes

std::vector< bool > elems_of_dimension
const bool _is_parallel_format

Private Member Functions

void write_serialized_connectivity (Xdr &io, const dof_id_type n_elem) const
void write_serialized_nodes (Xdr &io, const dof_id_type n_nodes) const
void write_serialized_bcs (Xdr &io, const std::size_t n_bcs) const
void read_serialized_connectivity (Xdr &io, const dof_id_type n_elem)
void read_serialized_nodes (Xdr &io, const dof_id_type n_nodes)
void read_serialized_bcs (Xdr &io)
void pack_element (std::vector< dof_id_type > &conn, const Elem *elem, const dof_id_type parent_id=DofObject::invalid_id, const dof_id_type parent_pid=DofObject::invalid_id) const

Private Attributes

bool _binary
bool _legacy
bool _write_serial
bool _write_parallel
std::string _version
std::string _bc_file_name
std::string _partition_map_file
std::string _subdomain_map_file
std::string _p_level_file

Static Private Attributes

static const std::size_t io_blksize = 128000

Detailed Description

Author:
Benjamin Kirk, John Peterson, 2004.

Definition at line 50 of file xdr_io.h.


Constructor & Destructor Documentation

libMesh::XdrIO::XdrIO ( MeshBase mesh,
const bool  binary_in = false 
) [explicit]

Constructor. Takes a writeable reference to a mesh object. This is the constructor required to read a mesh. The optional parameter binary can be used to switch between ASCII (false, the default) or binary (true) files.

Definition at line 119 of file xdr_io.C.

00119                                                   :
00120   MeshInput<MeshBase> (mesh,/* is_parallel_format = */ true),
00121   MeshOutput<MeshBase>(mesh,/* is_parallel_format = */ true),
00122   _binary             (binary_in),
00123   _legacy             (false),
00124   _write_serial       (false),
00125   _write_parallel     (false),
00126   _version            ("libMesh-0.7.0+"),
00127   _bc_file_name       ("n/a"),
00128   _partition_map_file ("n/a"),
00129   _subdomain_map_file ("n/a"),
00130   _p_level_file       ("n/a")
00131 {
00132 }

libMesh::XdrIO::XdrIO ( const MeshBase mesh,
const bool  binary_in = false 
) [explicit]

Constructor. Takes a reference to a constant mesh object. This constructor will only allow us to write the mesh. The optional parameter binary can be used to switch between ASCII (false, the default) or binary (true) files.

Definition at line 136 of file xdr_io.C.

00136                                                         :
00137   MeshOutput<MeshBase>(mesh,/* is_parallel_format = */ true),
00138   _binary (binary_in)
00139 {
00140 }

libMesh::XdrIO::~XdrIO (  )  [virtual]

Destructor.

Definition at line 144 of file xdr_io.C.

00145 {
00146 }


Member Function Documentation

unsigned int& libMesh::MeshOutput< MeshBase >::ascii_precision (  )  [inherited]

Return/set the precision to use when writing ASCII files.

By default we use numeric_limits<Real>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.

Referenced by libMesh::TecplotIO::write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().

bool& libMesh::XdrIO::binary (  )  [inline]

Definition at line 95 of file xdr_io.h.

References _binary.

00095 { return _binary; }

bool libMesh::XdrIO::binary (  )  const [inline]

Get/Set the flag indicating if we should read/write binary.

Definition at line 94 of file xdr_io.h.

References _binary.

Referenced by read(), libMesh::UnstructuredMesh::read(), and write().

00094 { return _binary; }

std::string& libMesh::XdrIO::boundary_condition_file_name (  )  [inline]

Definition at line 138 of file xdr_io.h.

References _bc_file_name.

00138 { return _bc_file_name; }

const std::string& libMesh::XdrIO::boundary_condition_file_name (  )  const [inline]

Get/Set the boundary condition file name.

Definition at line 137 of file xdr_io.h.

References _bc_file_name.

Referenced by read(), read_serialized_bcs(), and write().

00137 { return _bc_file_name; }

bool& libMesh::XdrIO::legacy (  )  [inline]

Definition at line 101 of file xdr_io.h.

References _legacy.

00101 { return _legacy; }

bool libMesh::XdrIO::legacy (  )  const [inline]

Get/Set the flag indicating if we should read/write legacy.

Definition at line 100 of file xdr_io.h.

References _legacy.

Referenced by read(), libMesh::UnstructuredMesh::read(), and write().

00100 { return _legacy; }

MeshBase & libMesh::MeshInput< MeshBase >::mesh (  )  [protected, inherited]

Returns the object as a writeable reference.

Referenced by libMesh::GMVIO::_read_materials(), libMesh::GMVIO::_read_nodes(), libMesh::GMVIO::_read_one_cell(), libMesh::AbaqusIO::assign_boundary_node_ids(), libMesh::AbaqusIO::assign_sideset_ids(), libMesh::AbaqusIO::assign_subdomain_ids(), libMesh::VTKIO::cells_to_vtk(), libMesh::GMVIO::copy_nodal_solution(), libMesh::UNVIO::element_in(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::element_out(), libMesh::UNVIO::node_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::node_out(), libMesh::VTKIO::nodes_to_vtk(), read(), libMesh::VTKIO::read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), libMesh::AbaqusIO::read(), libMesh::LegacyXdrIO::read_ascii(), libMesh::AbaqusIO::read_elements(), libMesh::UNVIO::read_implementation(), libMesh::UCDIO::read_implementation(), libMesh::LegacyXdrIO::read_mesh(), libMesh::GmshIO::read_mesh(), libMesh::AbaqusIO::read_nodes(), read_serialized_bcs(), read_serialized_connectivity(), read_serialized_nodes(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), write(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::UNVIO::write_implementation(), libMesh::UCDIO::write_implementation(), libMesh::LegacyXdrIO::write_mesh(), libMesh::GmshIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), write_parallel(), libMesh::GmshIO::write_post(), write_serialized_bcs(), write_serialized_nodes(), and libMesh::LegacyXdrIO::write_soln().

void libMesh::XdrIO::pack_element ( std::vector< dof_id_type > &  conn,
const Elem elem,
const dof_id_type  parent_id = DofObject::invalid_id,
const dof_id_type  parent_pid = DofObject::invalid_id 
) const [private]

Pack an element into a transfer buffer for parallel communication.

Definition at line 1166 of file xdr_io.C.

References libMesh::DofObject::id(), libMesh::invalid_uint, libMesh::Elem::n_nodes(), libMesh::Elem::node(), libMesh::Elem::p_level(), libMesh::DofObject::processor_id(), libMesh::Elem::subdomain_id(), libMesh::Elem::type(), and libMesh::Elem::type_to_n_nodes_map.

01168 {
01169   libmesh_assert(elem);
01170   libmesh_assert_equal_to (elem->n_nodes(), Elem::type_to_n_nodes_map[elem->type()]);
01171 
01172   conn.push_back(elem->n_nodes());
01173 
01174   conn.push_back (elem->type());
01175   conn.push_back (elem->id());
01176 
01177   if (parent_id != libMesh::invalid_uint)
01178     {
01179       conn.push_back (parent_id);
01180       libmesh_assert_not_equal_to (parent_pid, libMesh::invalid_uint);
01181       conn.push_back (parent_pid);
01182     }
01183 
01184   conn.push_back (elem->processor_id());
01185   conn.push_back (elem->subdomain_id());
01186 
01187 #ifdef LIBMESH_ENABLE_AMR
01188   conn.push_back (elem->p_level());
01189 #endif
01190 
01191   for (unsigned int n=0; n<elem->n_nodes(); n++)
01192     conn.push_back (elem->node(n));
01193 }

std::string& libMesh::XdrIO::partition_map_file_name (  )  [inline]

Definition at line 144 of file xdr_io.h.

References _partition_map_file.

00144 { return _partition_map_file; }

const std::string& libMesh::XdrIO::partition_map_file_name (  )  const [inline]

Get/Set the partitioning file name.

Definition at line 143 of file xdr_io.h.

References _partition_map_file.

Referenced by read(), read_serialized_connectivity(), and write().

00143 { return _partition_map_file; }

std::string& libMesh::XdrIO::polynomial_level_file_name (  )  [inline]

Definition at line 156 of file xdr_io.h.

References _p_level_file.

00156 { return _p_level_file; }

const std::string& libMesh::XdrIO::polynomial_level_file_name (  )  const [inline]

Get/Set the polynomial degree file name.

Definition at line 155 of file xdr_io.h.

References _p_level_file.

Referenced by read(), read_serialized_connectivity(), and write().

00155 { return _p_level_file; }

void libMesh::XdrIO::read ( const std::string &  name  )  [virtual]

This method implements reading a mesh from a specified file.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 806 of file xdr_io.C.

References binary(), boundary_condition_file_name(), libMesh::Parallel::Communicator::broadcast(), libMesh::Xdr::close(), libMesh::CommWorld, libMesh::Xdr::data(), libMeshEnums::DECODE, legacy(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshTools::n_elem(), partition_map_file_name(), polynomial_level_file_name(), libMesh::processor_id(), libMeshEnums::READ, read_serialized_bcs(), read_serialized_connectivity(), read_serialized_nodes(), libMesh::MeshBase::reserve_elem(), libMesh::MeshBase::reserve_nodes(), libMesh::Partitioner::set_node_processor_ids(), subdomain_map_file_name(), and version().

Referenced by libMesh::UnstructuredMesh::read().

00807 {
00808   // Only open the file on processor 0 -- this is especially important because
00809   // there may be an underlying bzip/bunzip going on, and multiple simultaneous
00810   // calls will produce a race condition.
00811   Xdr io (libMesh::processor_id() == 0 ? name : "", this->binary() ? DECODE : READ);
00812 
00813   // convenient reference to our mesh
00814   MeshBase &mesh = MeshInput<MeshBase>::mesh();
00815 
00816   // get the version string.
00817   if (libMesh::processor_id() == 0)
00818     io.data (this->version());
00819   CommWorld.broadcast (this->version());
00820 
00821   // note that for "legacy" files the first entry is an
00822   // integer -- not a string at all.
00823   this->legacy() = !(this->version().find("libMesh") < this->version().size());
00824 
00825   // Check for a legacy version format.
00826   if (this->legacy())
00827     {
00828       io.close();
00829       LegacyXdrIO(mesh, this->binary()).read(name);
00830       return;
00831     }
00832 
00833   START_LOG("read()","XdrIO");
00834 
00835   unsigned int n_elem, n_nodes;
00836   if (libMesh::processor_id() == 0)
00837     {
00838       io.data (n_elem);
00839       io.data (n_nodes);
00840       io.data (this->boundary_condition_file_name()); // libMesh::out << "bc_file="  << this->boundary_condition_file_name() << std::endl;
00841       io.data (this->subdomain_map_file_name());      // libMesh::out << "sid_file=" << this->subdomain_map_file_name()      << std::endl;
00842       io.data (this->partition_map_file_name());      // libMesh::out << "pid_file=" << this->partition_map_file_name()      << std::endl;
00843       io.data (this->polynomial_level_file_name());   // libMesh::out << "pl_file="  << this->polynomial_level_file_name()   << std::endl;
00844     }
00845 
00846   //TODO:[BSK] a little extra effort here could change this to two broadcasts...
00847   CommWorld.broadcast (n_elem);
00848   CommWorld.broadcast (n_nodes);
00849   CommWorld.broadcast (this->boundary_condition_file_name());
00850   CommWorld.broadcast (this->subdomain_map_file_name());
00851   CommWorld.broadcast (this->partition_map_file_name());
00852   CommWorld.broadcast (this->polynomial_level_file_name());
00853 
00854   // Tell the mesh how many nodes/elements to expect. Depending on the mesh type,
00855   // this may allow for efficient adding of nodes/elements.
00856   mesh.reserve_elem(n_elem);
00857   mesh.reserve_nodes(n_nodes);
00858 
00859   // read connectivity
00860   this->read_serialized_connectivity (io, n_elem);
00861 
00862   // read the nodal locations
00863   this->read_serialized_nodes (io, n_nodes);
00864 
00865   // read the boundary conditions
00866   this->read_serialized_bcs (io);
00867 
00868   STOP_LOG("read()","XdrIO");
00869 
00870   // set the node processor ids
00871   Partitioner::set_node_processor_ids(mesh);
00872 }

void libMesh::XdrIO::read_serialized_bcs ( Xdr io  )  [private]

Read the boundary conditions for a parallel, distributed mesh

Definition at line 1095 of file xdr_io.C.

References libMesh::BoundaryInfo::add_side(), boundary_condition_file_name(), libMesh::Parallel::Communicator::broadcast(), libMesh::CommWorld, libMesh::Xdr::data(), libMesh::Xdr::data_stream(), end, libMesh::MeshTools::Generation::Private::idx(), io_blksize, libMesh::MeshInput< MeshBase >::mesh(), std::min(), libMesh::processor_id(), and libMesh::Xdr::reading().

Referenced by read().

01096 {
01097   if (this->boundary_condition_file_name() == "n/a") return;
01098 
01099   libmesh_assert (io.reading());
01100 
01101   // convenient reference to our mesh
01102   MeshBase &mesh = MeshInput<MeshBase>::mesh();
01103 
01104   // and our boundary info object
01105   BoundaryInfo &boundary_info = *mesh.boundary_info;
01106 
01107   std::vector<ElemBCData> elem_bc_data;
01108   std::vector<int> input_buffer;
01109 
01110   std::size_t n_bcs=0;
01111   if (libMesh::processor_id() == 0)
01112     io.data (n_bcs);
01113   CommWorld.broadcast (n_bcs);
01114 
01115   for (std::size_t blk=0, first_bc=0, last_bc=0; last_bc<n_bcs; blk++)
01116     {
01117       first_bc = blk*io_blksize;
01118       last_bc  = std::min((blk+1)*io_blksize, n_bcs);
01119 
01120       input_buffer.resize (3*(last_bc - first_bc));
01121 
01122       if (libMesh::processor_id() == 0)
01123         io.data_stream (input_buffer.empty() ? NULL : &input_buffer[0], input_buffer.size());
01124 
01125       CommWorld.broadcast (input_buffer);
01126       elem_bc_data.clear();  elem_bc_data.reserve (input_buffer.size()/3);
01127 
01128       // convert the input_buffer to ElemBCData to facilitate searching
01129       for (std::size_t idx=0; idx<input_buffer.size(); idx+=3)
01130         elem_bc_data.push_back (ElemBCData(input_buffer[idx+0],
01131                                            input_buffer[idx+1],
01132                                            input_buffer[idx+2]));
01133       input_buffer.clear();
01134       // note that while the files *we* write should already be sorted by
01135       // element id this is not necessarily guaranteed.
01136       std::sort (elem_bc_data.begin(), elem_bc_data.end());
01137 
01138       MeshBase::const_element_iterator
01139         it  = mesh.level_elements_begin(0),
01140         end = mesh.level_elements_end(0);
01141 
01142       // Look for BCs in this block for all the level-0 elements we have
01143       // (not just local ones).  Do this by finding all the entries
01144       // in elem_bc_data whose elem_id match the ID of the current element.
01145       // We cannot rely on NULL neighbors at this point since the neighbor
01146       // data structure has not been initialized.
01147       for (std::pair<std::vector<ElemBCData>::iterator,
01148                      std::vector<ElemBCData>::iterator> pos; it!=end; ++it)
01149 #if defined(__SUNPRO_CC) || defined(__PGI)
01150         for (pos = std::equal_range (elem_bc_data.begin(), elem_bc_data.end(), (*it)->id(), CompareIntElemBCData());
01151 #else
01152         for (pos = std::equal_range (elem_bc_data.begin(), elem_bc_data.end(), (*it)->id());
01153 #endif
01154              pos.first != pos.second; ++pos.first)
01155           {
01156             libmesh_assert_equal_to (pos.first->elem_id, (*it)->id());
01157             libmesh_assert_less (pos.first->side, (*it)->n_sides());
01158 
01159             boundary_info.add_side (*it, pos.first->side, pos.first->bc_id);
01160           }
01161     }
01162 }

void libMesh::XdrIO::read_serialized_connectivity ( Xdr io,
const dof_id_type  n_elem 
) [private]

Read the connectivity for a parallel, distributed mesh

Definition at line 876 of file xdr_io.C.

References libMesh::Elem::add_child(), libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::Parallel::Communicator::broadcast(), libMesh::Elem::build(), libMesh::CommWorld, libMesh::Xdr::data(), libMesh::Xdr::data_stream(), libMesh::Elem::dim(), libMesh::MeshBase::elem(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::err, libMesh::Elem::hack_p_level(), libMesh::Elem::INACTIVE, libMesh::invalid_uint, io_blksize, libMesh::Elem::JUST_REFINED, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), std::min(), libMesh::Elem::n_nodes(), partition_map_file_name(), polynomial_level_file_name(), libMesh::DofObject::processor_id(), libMesh::processor_id(), libMesh::Xdr::reading(), libMesh::AutoPtr< Tp >::release(), libMesh::DofObject::set_id(), libMesh::MeshBase::set_mesh_dimension(), libMesh::Elem::set_node(), libMesh::Elem::set_refinement_flag(), libMesh::Elem::subdomain_id(), subdomain_map_file_name(), and libMesh::Elem::type_to_n_nodes_map.

Referenced by read().

00877 {
00878   libmesh_assert (io.reading());
00879 
00880   if (!n_elem) return;
00881 
00882   const bool
00883     read_p_level      = ("." == this->polynomial_level_file_name()),
00884     read_partitioning = ("." == this->partition_map_file_name()),
00885     read_subdomain_id = ("." == this->subdomain_map_file_name());
00886 
00887   // convenient reference to our mesh
00888   MeshBase &mesh = MeshInput<MeshBase>::mesh();
00889 
00890   // Keep track of what kinds of elements this file contains
00891   elems_of_dimension.clear();
00892   elems_of_dimension.resize(4, false);
00893 
00894   std::vector<dof_id_type> conn, input_buffer(100 /* oversized ! */);
00895 
00896   int level=-1;
00897 
00898   for (std::size_t blk=0, first_elem=0, last_elem=0, n_elem_at_level=0, n_processed_at_level=0;
00899        last_elem<n_elem; blk++)
00900     {
00901       first_elem = blk*io_blksize;
00902       last_elem  = std::min((blk+1)*io_blksize, std::size_t(n_elem));
00903 
00904       conn.clear();
00905 
00906       if (libMesh::processor_id() == 0)
00907         for (dof_id_type e=first_elem; e<last_elem; e++, n_processed_at_level++)
00908           {
00909             if (n_processed_at_level == n_elem_at_level)
00910               {
00911                 // get the number of elements to read at this level
00912                 io.data (n_elem_at_level);
00913                 n_processed_at_level = 0;
00914                 level++;
00915               }
00916 
00917             // get the element type,
00918             io.data_stream (&input_buffer[0], 1);
00919 
00920             // maybe the parent
00921             if (level)
00922               io.data_stream (&input_buffer[1], 1);
00923             else
00924               input_buffer[1] = libMesh::invalid_uint;
00925 
00926             // maybe the processor id
00927             if (read_partitioning)
00928               io.data_stream (&input_buffer[2], 1);
00929             else
00930               input_buffer[2] = 0;
00931 
00932             // maybe the subdomain id
00933             if (read_subdomain_id)
00934               io.data_stream (&input_buffer[3], 1);
00935             else
00936               input_buffer[3] = 0;
00937 
00938             // maybe the p level
00939             if (read_p_level)
00940               io.data_stream (&input_buffer[4], 1);
00941             else
00942               input_buffer[4] = 0;
00943 
00944             // and all the nodes
00945             libmesh_assert_less (5+Elem::type_to_n_nodes_map[input_buffer[0]], input_buffer.size());
00946             io.data_stream (&input_buffer[5], Elem::type_to_n_nodes_map[input_buffer[0]]);
00947             conn.insert (conn.end(),
00948                          input_buffer.begin(),
00949                          input_buffer.begin() + 5 + Elem::type_to_n_nodes_map[input_buffer[0]]);
00950           }
00951 
00952       std::size_t conn_size = conn.size();
00953       CommWorld.broadcast(conn_size);
00954       conn.resize (conn_size);
00955       CommWorld.broadcast (conn);
00956 
00957       // All processors now have the connectivity for this block.
00958       std::vector<dof_id_type>::const_iterator it = conn.begin();
00959       for (dof_id_type e=first_elem; e<last_elem; e++)
00960         {
00961           const ElemType elem_type        = static_cast<ElemType>(*it); ++it;
00962           const dof_id_type parent_id    = *it; ++it;
00963           const processor_id_type processor_id = *it; ++it;
00964           const subdomain_id_type subdomain_id = *it; ++it;
00965 #ifdef LIBMESH_ENABLE_AMR
00966           const unsigned int p_level      = *it;
00967 #endif
00968           ++it;
00969 
00970           Elem *parent = (parent_id == libMesh::invalid_uint) ? NULL : mesh.elem(parent_id);
00971 
00972           Elem *elem = Elem::build (elem_type, parent).release();
00973 
00974           elem->set_id() = e;
00975           elem->processor_id() = processor_id;
00976           elem->subdomain_id() = subdomain_id;
00977 #ifdef LIBMESH_ENABLE_AMR
00978           elem->hack_p_level(p_level);
00979 
00980           if (parent)
00981             {
00982               parent->add_child(elem);
00983               parent->set_refinement_flag (Elem::INACTIVE);
00984               elem->set_refinement_flag   (Elem::JUST_REFINED);
00985             }
00986 #endif
00987 
00988           for (unsigned int n=0; n<elem->n_nodes(); n++, ++it)
00989             {
00990               const dof_id_type global_node_number = *it;
00991 
00992               elem->set_node(n) =
00993                 mesh.add_point (Point(), global_node_number);
00994             }
00995 
00996           elems_of_dimension[elem->dim()] = true;
00997           mesh.add_elem(elem);
00998         }
00999     }
01000 
01001   // Set the mesh dimension to the largest encountered for an element
01002   for (unsigned int i=0; i!=4; ++i)
01003     if (elems_of_dimension[i])
01004       mesh.set_mesh_dimension(i);
01005 
01006 #if LIBMESH_DIM < 3
01007   if (mesh.mesh_dimension() > LIBMESH_DIM)
01008     {
01009       libMesh::err << "Cannot open dimension " <<
01010                       mesh.mesh_dimension() <<
01011                       " mesh file when configured without " <<
01012                       mesh.mesh_dimension() << "D support." <<
01013                       std::endl;
01014       libmesh_error();
01015     }
01016 #endif
01017 }

void libMesh::XdrIO::read_serialized_nodes ( Xdr io,
const dof_id_type  n_nodes 
) [private]

Read the nodal locations for a parallel, distributed mesh

Definition at line 1021 of file xdr_io.C.

References libMesh::Parallel::Communicator::broadcast(), libMesh::CommWorld, libMesh::Xdr::data_stream(), end, libMesh::MeshTools::Generation::Private::idx(), io_blksize, libMesh::MeshInput< MeshBase >::mesh(), std::min(), libMesh::processor_id(), and libMesh::Xdr::reading().

Referenced by read().

01022 {
01023   libmesh_assert (io.reading());
01024 
01025   // convenient reference to our mesh
01026   MeshBase &mesh = MeshInput<MeshBase>::mesh();
01027 
01028   if (!mesh.n_nodes()) return;
01029 
01030   // At this point the elements have been read from file and placeholder nodes
01031   // have been assigned.  These nodes, however, do not have the proper (x,y,z)
01032   // locations.  This method will read all the nodes from disk, and each processor
01033   // can then grab the individual values it needs.
01034 
01035   // build up a list of the nodes contained in our local mesh.  These are the nodes
01036   // stored on the local processor whose (x,y,z) values need to be corrected.
01037   std::vector<dof_id_type> needed_nodes; needed_nodes.reserve (mesh.n_nodes());
01038   {
01039     MeshBase::node_iterator
01040       it  = mesh.nodes_begin(),
01041       end = mesh.nodes_end();
01042 
01043     for (; it!=end; ++it)
01044       needed_nodes.push_back((*it)->id());
01045 
01046     std::sort (needed_nodes.begin(), needed_nodes.end());
01047 
01048     // We should not have any duplicate node->id()s
01049     libmesh_assert (std::unique(needed_nodes.begin(), needed_nodes.end()) == needed_nodes.end());
01050   }
01051 
01052   // Get the nodes in blocks.
01053   std::vector<Real> coords;
01054   std::pair<std::vector<dof_id_type>::iterator,
01055             std::vector<dof_id_type>::iterator> pos;
01056   pos.first = needed_nodes.begin();
01057 
01058   for (std::size_t blk=0, first_node=0, last_node=0; last_node<n_nodes; blk++)
01059     {
01060       first_node = blk*io_blksize;
01061       last_node  = std::min((blk+1)*io_blksize, std::size_t(n_nodes));
01062 
01063       coords.resize(3*(last_node - first_node));
01064 
01065       if (libMesh::processor_id() == 0)
01066         io.data_stream (coords.empty() ? NULL : &coords[0], coords.size());
01067 
01068       // For large numbers of processors the majority of processors at any given
01069       // block may not actually need these data.  It may be worth profiling this,
01070       // although it is expected that disk IO will be the bottleneck
01071       CommWorld.broadcast (coords);
01072 
01073       for (std::size_t n=first_node, idx=0; n<last_node; n++, idx+=3)
01074         {
01075           // first see if we need this node.  use pos.first as a smart lower
01076           // bound, this will ensure that the size of the searched range
01077           // decreases as we match nodes.
01078           pos = std::equal_range (pos.first, needed_nodes.end(), n);
01079 
01080           if (pos.first != pos.second) // we need this node.
01081             {
01082               libmesh_assert_equal_to (*pos.first, n);
01083               mesh.node(n) =
01084                 Point (coords[idx+0],
01085                        coords[idx+1],
01086                        coords[idx+2]);
01087 
01088             }
01089         }
01090     }
01091 }

void libMesh::XdrIO::set_auto_parallel (  )  [inline]

Insist that we should write parallel files if and only if the mesh is an already distributed ParallelMesh.

Definition at line 259 of file xdr_io.h.

References _write_parallel, and _write_serial.

00260 {
00261   this->_write_serial   = false;
00262   this->_write_parallel = false;
00263 }

void libMesh::XdrIO::set_write_parallel ( bool  do_parallel = true  )  [inline]

Insist that we should/shouldn't write parallel files.

Definition at line 249 of file xdr_io.h.

References _write_parallel, and _write_serial.

00250 {
00251   this->_write_parallel = do_parallel;
00252 
00253   this->_write_serial = !do_parallel;
00254 }

void libMesh::MeshInput< MeshBase >::skip_comment_lines ( std::istream &  in,
const char  comment_start 
) [protected, inherited]

Reads input from in, skipping all the lines that start with the character comment_start.

Referenced by libMesh::TetGenIO::read(), and libMesh::UCDIO::read_implementation().

std::string& libMesh::XdrIO::subdomain_map_file_name (  )  [inline]

Definition at line 150 of file xdr_io.h.

References _subdomain_map_file.

00150 { return _subdomain_map_file; }

const std::string& libMesh::XdrIO::subdomain_map_file_name (  )  const [inline]

Get/Set the subdomain file name.

Definition at line 149 of file xdr_io.h.

References _subdomain_map_file.

Referenced by read(), read_serialized_connectivity(), and write().

00149 { return _subdomain_map_file; }

std::string& libMesh::XdrIO::version (  )  [inline]

Definition at line 132 of file xdr_io.h.

References _version.

00132 { return _version; }

const std::string& libMesh::XdrIO::version (  )  const [inline]

Get/Set the version string. Vailid version strings:


     "libMesh-0.7.0+"
     "libMesh-0.7.0+ parallel"

     

If "libMesh" is not detected in the version string the LegacyXdrIO class will be used to read older (pre version 0.7.0) mesh files.

Definition at line 131 of file xdr_io.h.

References _version.

Referenced by read(), and write().

00131 { return _version; }

void libMesh::XdrIO::write ( const std::string &  name  )  [virtual]

This method implements writing a mesh to a specified file.

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 150 of file xdr_io.C.

References _write_parallel, libMesh::Parallel::Communicator::barrier(), binary(), boundary_condition_file_name(), libMesh::MeshBase::boundary_info, libMesh::Xdr::close(), libMesh::CommWorld, libMesh::Xdr::data(), libMeshEnums::ENCODE, legacy(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_elem(), libMesh::MeshTools::n_elem(), libMesh::MeshBase::n_nodes(), n_nodes, libMesh::MeshTools::n_p_levels(), libMesh::MeshBase::n_subdomains(), libMesh::out, partition_map_file_name(), polynomial_level_file_name(), libMesh::processor_id(), subdomain_map_file_name(), version(), libMeshEnums::WRITE, write_parallel(), write_serialized_bcs(), write_serialized_connectivity(), and write_serialized_nodes().

00151 {
00152   if (this->legacy())
00153     {
00154       libmesh_deprecated();
00155 
00156       // We don't support writing parallel files in the legacy format
00157       libmesh_assert(!this->_write_parallel);
00158 
00159       LegacyXdrIO(MeshOutput<MeshBase>::mesh(), this->binary()).write(name);
00160       return;
00161     }
00162 
00163   Xdr io ((libMesh::processor_id() == 0) ? name : "", this->binary() ? ENCODE : WRITE);
00164 
00165   START_LOG("write()","XdrIO");
00166 
00167   // convenient reference to our mesh
00168   const MeshBase &mesh = MeshOutput<MeshBase>::mesh();
00169 
00170   dof_id_type
00171     n_elem     = mesh.n_elem(),
00172     n_nodes    = mesh.n_nodes();
00173   std::size_t
00174     n_bcs      = mesh.boundary_info->n_boundary_conds();
00175   unsigned int
00176     n_p_levels = MeshTools::n_p_levels (mesh);
00177 
00178   bool write_parallel_files = this->write_parallel();
00179 
00180   //-------------------------------------------------------------
00181   // For all the optional files -- the default file name is "n/a".
00182   // However, the user may specify an optional external file.
00183 
00184   // If there are BCs and the user has not already provided a
00185   // file name then write to "."
00186   if (n_bcs &&
00187       this->boundary_condition_file_name() == "n/a")
00188     this->boundary_condition_file_name() = ".";
00189 
00190   // If there are more than one subdomains and the user has not specified an
00191   // external file then write the subdomain mapping to the default file "."
00192   if ((mesh.n_subdomains() > 0) &&
00193       (this->subdomain_map_file_name() == "n/a"))
00194     this->subdomain_map_file_name() = ".";
00195 
00196   // In general we don't write the partition information.
00197 
00198   // If we have p levels and the user has not already provided
00199   // a file name then write to "."
00200   if ((n_p_levels > 1) &&
00201       (this->polynomial_level_file_name() == "n/a"))
00202     this->polynomial_level_file_name() = ".";
00203 
00204   // write the header
00205   if (libMesh::processor_id() == 0)
00206     {
00207       std::string full_ver = this->version() + (write_parallel_files ?  " parallel" : "");
00208       io.data (full_ver);
00209 
00210       io.data (n_elem,  "# number of elements");
00211       io.data (n_nodes, "# number of nodes");
00212 
00213       io.data (this->boundary_condition_file_name(), "# boundary condition specification file");
00214       io.data (this->subdomain_map_file_name(),      "# subdomain id specification file");
00215       io.data (this->partition_map_file_name(),      "# processor id specification file");
00216       io.data (this->polynomial_level_file_name(),   "# p-level specification file");
00217     }
00218 
00219   if (write_parallel_files)
00220     {
00221       // Parallel xdr mesh files aren't implemented yet; until they
00222       // are we'll just warn the user and write a serial file.
00223       libMesh::out << "Warning!  Parallel xda/xdr is not yet implemented.\n";
00224       libMesh::out << "Writing a serialized file instead." << std::endl;
00225 
00226       // write connectivity
00227       this->write_serialized_connectivity (io, n_elem);
00228 
00229       // write the nodal locations
00230       this->write_serialized_nodes (io, n_nodes);
00231 
00232       // write the boundary condition information
00233       this->write_serialized_bcs (io, n_bcs);
00234     }
00235   else
00236     {
00237       // write connectivity
00238       this->write_serialized_connectivity (io, n_elem);
00239 
00240       // write the nodal locations
00241       this->write_serialized_nodes (io, n_nodes);
00242 
00243       // write the boundary condition information
00244       this->write_serialized_bcs (io, n_bcs);
00245     }
00246 
00247   STOP_LOG("write()","XdrIO");
00248 
00249   // pause all processes until the writing ends -- this will
00250   // protect for the pathological case where a write is
00251   // followed immediately by a read.  The write must be
00252   // guaranteed to complete first.
00253   io.close();
00254   CommWorld.barrier();
00255 }

virtual void libMesh::MeshOutput< MeshBase >::write_equation_systems ( const std::string &  ,
const EquationSystems ,
const std::set< std::string > *  system_names = NULL 
) [virtual, inherited]

This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object.

Referenced by libMesh::Nemesis_IO::write_timestep().

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
) [inline, virtual, inherited]

This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.

Reimplemented in libMesh::ExodusII_IO, libMesh::GmshIO, libMesh::GMVIO, libMesh::GnuPlotIO, libMesh::MEDITIO, libMesh::Nemesis_IO, libMesh::TecplotIO, libMesh::UCDIO, and libMesh::VTKIO.

Definition at line 98 of file mesh_output.h.

00101   { libmesh_error(); }

bool libMesh::XdrIO::write_parallel (  )  const [inline]

Report whether we should write parallel files.

Definition at line 228 of file xdr_io.h.

References _write_parallel, _write_serial, libMesh::MeshBase::is_serial(), and libMesh::MeshInput< MeshBase >::mesh().

Referenced by write().

00229 {
00230   // We can't insist on both serial and parallel
00231   libmesh_assert (!this->_write_serial || !this->_write_parallel);
00232 
00233   // If we insisted on serial, do that
00234   if (this->_write_serial)
00235     return false;
00236 
00237   // If we insisted on parallel, do that
00238   if (this->_write_parallel)
00239     return true;
00240 
00241   // If we're doing things automatically, check the mesh
00242   const MeshBase &mesh = MeshOutput<MeshBase>::mesh();
00243   return !mesh.is_serial();
00244 }

void libMesh::XdrIO::write_serialized_bcs ( Xdr io,
const std::size_t  n_bcs 
) const [private]

Write the boundary conditions for a parallel, distributed mesh

Definition at line 723 of file xdr_io.C.

References bc_id, libMesh::BoundaryInfo::boundary_ids(), libMesh::MeshBase::boundary_info, libMesh::CommWorld, libMesh::Xdr::data(), libMesh::Xdr::data_stream(), end, libMesh::Parallel::Communicator::gather(), libMesh::MeshTools::Generation::Private::idx(), libMesh::BoundaryInfo::invalid_id, libMesh::MeshBase::local_level_elements_begin(), libMesh::MeshBase::local_level_elements_end(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::n_processors(), libMesh::Elem::n_sides(), libMesh::processor_id(), libMesh::Parallel::Communicator::receive(), libMesh::Parallel::Communicator::send(), and libMesh::Xdr::writing().

Referenced by write().

00724 {
00725   libmesh_assert (io.writing());
00726 
00727   if (!n_bcs) return;
00728 
00729   // convenient reference to our mesh
00730   const MeshBase &mesh = MeshOutput<MeshBase>::mesh();
00731 
00732   // and our boundary info object
00733   const BoundaryInfo &boundary_info = *mesh.boundary_info;
00734 
00735   std::size_t n_bcs_out = n_bcs;
00736   if (libMesh::processor_id() == 0)
00737     io.data (n_bcs_out, "# number of boundary conditions");
00738   n_bcs_out = 0;
00739 
00740   std::vector<dof_id_type> xfer_bcs, recv_bcs;
00741   std::vector<std::size_t> bc_sizes(libMesh::n_processors());;
00742 
00743   // Boundary conditions are only specified for level-0 elements
00744   MeshBase::const_element_iterator
00745     it  = mesh.local_level_elements_begin(0),
00746     end = mesh.local_level_elements_end(0);
00747 
00748   dof_id_type n_local_level_0_elem=0;
00749   for (; it!=end; ++it, n_local_level_0_elem++)
00750     {
00751       const Elem *elem = *it;
00752 
00753       for (unsigned int s=0; s<elem->n_sides(); s++)
00754 // We're supporting boundary ids on internal sides now
00755 //      if (elem->neighbor(s) == NULL)
00756           {
00757             const std::vector<boundary_id_type>& bc_ids =
00758               boundary_info.boundary_ids (elem, s);
00759             for (std::vector<boundary_id_type>::const_iterator id_it=bc_ids.begin(); id_it!=bc_ids.end(); ++id_it)
00760               {
00761                 const boundary_id_type bc_id = *id_it;
00762                 if (bc_id != BoundaryInfo::invalid_id)
00763                   {
00764                     xfer_bcs.push_back (n_local_level_0_elem);
00765                     xfer_bcs.push_back (s) ;
00766                     xfer_bcs.push_back (bc_id);
00767                   }
00768               }
00769           }
00770     }
00771 
00772   xfer_bcs.push_back(n_local_level_0_elem);
00773   std::size_t my_size = xfer_bcs.size();
00774   CommWorld.gather (0, my_size, bc_sizes);
00775 
00776   // All processors send their xfer buffers to processor 0
00777   // Processor 0 will receive all buffers and write out the bcs
00778   if (libMesh::processor_id() == 0)
00779     {
00780       dof_id_type elem_offset = 0;
00781       for (unsigned int pid=0; pid<libMesh::n_processors(); pid++)
00782         {
00783           recv_bcs.resize(bc_sizes[pid]);
00784           if (pid == 0)
00785             recv_bcs = xfer_bcs;
00786           else
00787             CommWorld.receive (pid, recv_bcs);
00788 
00789           const dof_id_type my_n_local_level_0_elem
00790             = recv_bcs.back(); recv_bcs.pop_back();
00791 
00792           for (std::size_t idx=0; idx<recv_bcs.size(); idx += 3, n_bcs_out++)
00793             recv_bcs[idx+0] += elem_offset;
00794 
00795           io.data_stream (recv_bcs.empty() ? NULL : &recv_bcs[0], recv_bcs.size(), 3);
00796           elem_offset += my_n_local_level_0_elem;
00797         }
00798       libmesh_assert_equal_to (n_bcs, n_bcs_out);
00799     }
00800   else
00801     CommWorld.send (0, xfer_bcs);
00802 }

void libMesh::XdrIO::write_serialized_connectivity ( Xdr io,
const dof_id_type  n_elem 
) const [private]

Write the connectivity for a parallel, distributed mesh

Referenced by write().

void libMesh::XdrIO::write_serialized_nodes ( Xdr io,
const dof_id_type  n_nodes 
) const [private]

Write the nodal locations for a parallel, distributed mesh

Definition at line 570 of file xdr_io.C.

References libMesh::Parallel::Communicator_World, libMesh::CommWorld, libMesh::Xdr::data_stream(), end, libMesh::Parallel::Communicator::gather(), libMesh::MeshTools::Generation::Private::idx(), io_blksize, libMesh::MeshBase::local_nodes_begin(), libMesh::MeshBase::local_nodes_end(), libMesh::MeshInput< MeshBase >::mesh(), std::min(), libMesh::MeshBase::n_nodes(), libMesh::n_processors(), libMesh::processor_id(), libMesh::Parallel::Communicator::receive(), libMesh::Parallel::Communicator::send(), and libMesh::Parallel::wait().

Referenced by write().

00571 {
00572   // convenient reference to our mesh
00573   const MeshBase &mesh = MeshOutput<MeshBase>::mesh();
00574   libmesh_assert_equal_to (n_nodes, mesh.n_nodes());
00575 
00576   std::vector<dof_id_type> xfer_ids;
00577   std::vector<Real>         xfer_coords, &coords=xfer_coords;
00578 
00579   std::vector<std::vector<dof_id_type> > recv_ids   (libMesh::n_processors());;
00580   std::vector<std::vector<Real> >         recv_coords(libMesh::n_processors());
00581 
00582   std::size_t n_written=0;
00583 
00584   for (std::size_t blk=0, last_node=0; last_node<n_nodes; blk++)
00585     {
00586       const std::size_t first_node = blk*io_blksize;
00587       last_node = std::min((blk+1)*io_blksize, std::size_t(n_nodes));
00588 
00589       // Build up the xfer buffers on each processor
00590       MeshBase::const_node_iterator
00591         it  = mesh.local_nodes_begin(),
00592         end = mesh.local_nodes_end();
00593 
00594       xfer_ids.clear(); xfer_coords.clear();
00595 
00596       for (; it!=end; ++it)
00597         if (((*it)->id() >= first_node) && // node in [first_node, last_node)
00598             ((*it)->id() <  last_node))
00599           {
00600             xfer_ids.push_back((*it)->id());
00601             const Point &p = **it;
00602             xfer_coords.push_back(p(0));
00603 #if LIBMESH_DIM > 1
00604             xfer_coords.push_back(p(1));
00605 #endif
00606 #if LIBMESH_DIM > 2
00607             xfer_coords.push_back(p(2));
00608 #endif
00609           }
00610 
00611       //-------------------------------------
00612       // Send the xfer buffers to processor 0
00613       std::vector<std::size_t> ids_size, coords_size;
00614 
00615       const std::size_t my_ids_size = xfer_ids.size();
00616 
00617       // explicitly gather ids_size
00618       CommWorld.gather (0, my_ids_size, ids_size);
00619 
00620       // infer coords_size on processor 0
00621       if (libMesh::processor_id() == 0)
00622         {
00623           coords_size.reserve(libMesh::n_processors());
00624           for (std::size_t p=0; p<ids_size.size(); p++)
00625             coords_size.push_back(LIBMESH_DIM*ids_size[p]);
00626         }
00627 
00628       // We will have lots of simultaneous receives if we are
00629       // processor 0, so let's use nonblocking receives.
00630       std::vector<Parallel::Request>
00631         id_request_handles(libMesh::n_processors()-1),
00632         coord_request_handles(libMesh::n_processors()-1);
00633 
00634       Parallel::MessageTag
00635         id_tag    = Parallel::Communicator_World.get_unique_tag(1234),
00636         coord_tag = Parallel::Communicator_World.get_unique_tag(1235);
00637 
00638       // Post the receives -- do this on processor 0 only.
00639       if (libMesh::processor_id() == 0)
00640         {
00641           for (unsigned int pid=0; pid<libMesh::n_processors(); pid++)
00642             {
00643               recv_ids[pid].resize(ids_size[pid]);
00644               recv_coords[pid].resize(coords_size[pid]);
00645 
00646               if (pid == 0)
00647                 {
00648                   recv_ids[0] = xfer_ids;
00649                   recv_coords[0] = xfer_coords;
00650                 }
00651               else
00652                 {
00653                   CommWorld.receive (pid, recv_ids[pid],
00654                                      id_request_handles[pid-1],
00655                                      id_tag);
00656                   CommWorld.receive (pid, recv_coords[pid],
00657                                      coord_request_handles[pid-1],
00658                                      coord_tag);
00659                 }
00660             }
00661         }
00662       else
00663         {
00664           // Send -- do this on all other processors.
00665           CommWorld.send(0, xfer_ids,    id_tag);
00666           CommWorld.send(0, xfer_coords, coord_tag);
00667         }
00668 
00669       // -------------------------------------------------------
00670       // Receive the messages and write the output on processor 0.
00671       if (libMesh::processor_id() == 0)
00672         {
00673           // Wait for all the receives to complete. We have no
00674           // need for the statuses since we already know the
00675           // buffer sizes.
00676           Parallel::wait (id_request_handles);
00677           Parallel::wait (coord_request_handles);
00678 
00679           // Write the coordinates in this block.
00680           std::size_t tot_id_size=0, tot_coord_size=0;
00681           for (unsigned int pid=0; pid<libMesh::n_processors(); pid++)
00682             {
00683               tot_id_size    += recv_ids[pid].size();
00684               tot_coord_size += recv_coords[pid].size();
00685             }
00686 
00687           libmesh_assert_less_equal 
00688             (tot_id_size, std::min(io_blksize, std::size_t(n_nodes)));
00689           libmesh_assert_equal_to (tot_coord_size, LIBMESH_DIM*tot_id_size);
00690 
00691           coords.resize (3*tot_id_size);
00692           for (unsigned int pid=0; pid<libMesh::n_processors(); pid++)
00693             for (std::size_t idx=0; idx<recv_ids[pid].size(); idx++)
00694               {
00695                 const std::size_t local_idx = recv_ids[pid][idx] - first_node;
00696                 libmesh_assert_less ((3*local_idx+2), coords.size());
00697                 libmesh_assert_less ((LIBMESH_DIM*idx+LIBMESH_DIM-1), recv_coords[pid].size());
00698 
00699                 coords[3*local_idx+0] = recv_coords[pid][LIBMESH_DIM*idx+0];
00700 #if LIBMESH_DIM > 1
00701                 coords[3*local_idx+1] = recv_coords[pid][LIBMESH_DIM*idx+1];
00702 #else
00703                 coords[3*local_idx+1] = 0.;
00704 #endif
00705 #if LIBMESH_DIM > 2
00706                 coords[3*local_idx+2] = recv_coords[pid][LIBMESH_DIM*idx+2];
00707 #else
00708                 coords[3*local_idx+2] = 0.;
00709 #endif
00710 
00711                 n_written++;
00712               }
00713 
00714           io.data_stream (coords.empty() ? NULL : &coords[0], coords.size(), 3);
00715         }
00716     }
00717   if (libMesh::processor_id() == 0)
00718     libmesh_assert_equal_to (n_written, n_nodes);
00719 }


Member Data Documentation

std::string libMesh::XdrIO::_bc_file_name [private]

Definition at line 212 of file xdr_io.h.

Referenced by boundary_condition_file_name().

bool libMesh::XdrIO::_binary [private]

Definition at line 207 of file xdr_io.h.

Referenced by binary().

const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format [protected, inherited]

Flag specifying whether this format is parallel-capable. If this is false (default) I/O is only permitted when the mesh has been serialized.

Definition at line 126 of file mesh_output.h.

Referenced by libMesh::PostscriptIO::write(), libMesh::FroIO::write(), libMesh::EnsightIO::write(), and libMesh::DivaIO::write().

bool libMesh::XdrIO::_legacy [private]

Definition at line 208 of file xdr_io.h.

Referenced by legacy().

std::string libMesh::XdrIO::_p_level_file [private]

Definition at line 215 of file xdr_io.h.

Referenced by polynomial_level_file_name().

std::string libMesh::XdrIO::_partition_map_file [private]

Definition at line 213 of file xdr_io.h.

Referenced by partition_map_file_name().

std::string libMesh::XdrIO::_subdomain_map_file [private]

Definition at line 214 of file xdr_io.h.

Referenced by subdomain_map_file_name().

std::string libMesh::XdrIO::_version [private]

Definition at line 211 of file xdr_io.h.

Referenced by version().

Definition at line 210 of file xdr_io.h.

Referenced by set_auto_parallel(), set_write_parallel(), write(), and write_parallel().

Definition at line 209 of file xdr_io.h.

Referenced by set_auto_parallel(), set_write_parallel(), and write_parallel().

const std::size_t libMesh::XdrIO::io_blksize = 128000 [static, private]

Define the block size to use for chunked IO.

Definition at line 220 of file xdr_io.h.

Referenced by read_serialized_bcs(), read_serialized_connectivity(), read_serialized_nodes(), and write_serialized_nodes().


The documentation for this class was generated from the following files:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:43 UTC

Hosted By:
SourceForge.net Logo