libMesh::GMVIO Class Reference

#include <gmv_io.h>

Inheritance diagram for libMesh::GMVIO:

List of all members.

Public Member Functions

 GMVIO (const MeshBase &)
 GMVIO (MeshBase &)
virtual void write (const std::string &)
virtual void read (const std::string &mesh_file)
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
bool & binary ()
bool & discontinuous ()
bool & partitioning ()
bool & write_subdomain_id_as_material ()
bool & subdivide_second_order ()
bool & p_levels ()
void write_discontinuous_gmv (const std::string &name, const EquationSystems &es, const bool write_partitioning) const
void write_ascii_new_impl (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL)
void add_cell_centered_data (const std::string &cell_centered_data_name, const std::vector< Real > *cell_centered_data_vals)
void copy_nodal_solution (EquationSystems &es)
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL)
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_ascii_old_impl (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL)
void write_binary (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL)
template<typename T >
void to_binary_stream (std::ostream &out, const T i)
void _read_nodes ()
void _read_one_cell ()
ElemType _gmv_elem_to_libmesh_elem (const char *elemname)
void _read_materials ()
void _read_var ()

Private Attributes

bool _binary
bool _discontinuous
bool _partitioning
bool _write_subdomain_id_as_material
bool _subdivide_second_order
bool _p_levels
std::map< std::string, const
std::vector< Real > * > 
_cell_centered_data
unsigned int _next_elem_id
std::map< std::string,
std::vector< Number > > 
_nodal_data

Detailed Description

This class implements writing meshes in the GMV format. For a full description of the GMV format and to obtain the GMV software see the GMV home page

Author:
Benjamin S. Kirk, 2004

Definition at line 53 of file gmv_io.h.


Constructor & Destructor Documentation

libMesh::GMVIO::GMVIO ( const MeshBase mesh  )  [inline, explicit]

Constructor. Takes a reference to a constant mesh object. This constructor will only allow us to write the mesh.

Definition at line 271 of file gmv_io.h.

00271                                   :
00272   MeshOutput<MeshBase>    (mesh),
00273   _binary                 (false),
00274   _discontinuous          (false),
00275   _partitioning           (true),
00276   _write_subdomain_id_as_material (false),
00277   _subdivide_second_order (true),
00278   _p_levels               (true),
00279   _next_elem_id           (0)
00280 {
00281 }

libMesh::GMVIO::GMVIO ( MeshBase mesh  )  [inline, explicit]

Constructor. Takes a writeable reference to a mesh object. This constructor is required to let us read in a mesh.

Definition at line 284 of file gmv_io.h.

00284                             :
00285   MeshInput<MeshBase> (mesh),
00286   MeshOutput<MeshBase>(mesh),
00287   _binary (false),
00288   _discontinuous          (false),
00289   _partitioning           (true),
00290   _write_subdomain_id_as_material (false),
00291   _subdivide_second_order (true),
00292   _p_levels               (true),
00293   _next_elem_id           (0)
00294 {
00295 }


Member Function Documentation

ElemType libMesh::GMVIO::_gmv_elem_to_libmesh_elem ( const char *  elemname  )  [private]

Definition at line 2319 of file gmv_io.C.

References libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMesh::err, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMeshEnums::PRISM15, libMeshEnums::PRISM6, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, and libMeshEnums::TRI6.

Referenced by _read_one_cell().

02320 {
02321   //
02322   // Linear Elements
02323   //
02324   if (!std::strncmp(elemname,"line",4))
02325     return EDGE2;
02326 
02327   if (!std::strncmp(elemname,"tri",3))
02328     return TRI3;
02329 
02330   if (!std::strncmp(elemname,"quad",4))
02331     return QUAD4;
02332 
02333   // FIXME: tet or ptet4?
02334   if ((!std::strncmp(elemname,"tet",3)) ||
02335       (!std::strncmp(elemname,"ptet4",5)))
02336     return TET4;
02337 
02338   // FIXME: hex or phex8?
02339   if ((!std::strncmp(elemname,"hex",3)) ||
02340       (!std::strncmp(elemname,"phex8",5)))
02341     return HEX8;
02342 
02343   // FIXME: prism or pprism6?
02344   if ((!std::strncmp(elemname,"prism",5)) ||
02345       (!std::strncmp(elemname,"pprism6",7)))
02346     return PRISM6;
02347 
02348   //
02349   // Quadratic Elements
02350   //
02351   if (!std::strncmp(elemname,"phex20",6))
02352     return HEX20;
02353 
02354   if (!std::strncmp(elemname,"phex27",6))
02355     return HEX27;
02356 
02357   if (!std::strncmp(elemname,"pprism15",8))
02358     return PRISM15;
02359 
02360   if (!std::strncmp(elemname,"ptet10",6))
02361     return TET10;
02362 
02363   if (!std::strncmp(elemname,"6tri",4))
02364     return TRI6;
02365 
02366   if (!std::strncmp(elemname,"8quad",5))
02367     return QUAD8;
02368 
02369   if (!std::strncmp(elemname,"3line",5))
02370     return EDGE3;
02371 
02372   // Unsupported/Unused types
02373   // if (!std::strncmp(elemname,"vface2d",7))
02374   // if (!std::strncmp(elemname,"vface3d",7))
02375   // if (!std::strncmp(elemname,"pyramid",7))
02376   // if (!std::strncmp(elemname,"ppyrmd5",7))
02377   // if (!std::strncmp(elemname,"ppyrmd13",8))
02378 
02379   // If we didn't return yet, then we didn't find the right cell!
02380   libMesh::err << "Uknown/unsupported element: "
02381                 << elemname
02382                 << " was read."
02383                 << std::endl;
02384   libmesh_error();
02385 }

void libMesh::GMVIO::_read_materials (  )  [private]

Definition at line 2163 of file gmv_io.C.

References libMesh::MeshBase::elem(), libMesh::MeshInput< MeshBase >::mesh(), and libMesh::DofObject::processor_id().

Referenced by read().

02164 {
02165 #ifdef LIBMESH_HAVE_GMV
02166 
02167   // LibMesh assigns materials on a per-cell basis
02168   libmesh_assert_equal_to (GMV::gmv_data.datatype, CELL);
02169 
02170   //   // Material names: LibMesh has no use for these currently...
02171   //   libMesh::out << "Number of material names="
02172   //        << GMV::gmv_data.num
02173   //        << std::endl;
02174 
02175   //   for (int i = 0; i < GMV::gmv_data.num; i++)
02176   //     {
02177   //       // Build a 32-char string from the appropriate entries
02178   //       std::string mat_string(&GMV::gmv_data.chardata1[i*33], 32);
02179 
02180   //       libMesh::out << "Material name " << i << ": " << mat_string << std::endl;
02181   //     }
02182 
02183   //   // Material labels: These correspond to (1-based) CPU IDs, and
02184   //   // there should be 1 of these for each element.
02185   //   libMesh::out << "Number of material labels = "
02186   //        << GMV::gmv_data.nlongdata1
02187   //        << std::endl;
02188 
02189   for (int i = 0; i < GMV::gmv_data.nlongdata1; i++)
02190     {
02191       // Debugging Info
02192       // libMesh::out << "Material ID " << i << ": "
02193       // << GMV::gmv_data.longdata1[i]
02194       // << std::endl;
02195 
02196       MeshInput<MeshBase>::mesh().elem(i)->processor_id() =
02197         GMV::gmv_data.longdata1[i]-1;
02198     }
02199 
02200 #endif
02201 }

void libMesh::GMVIO::_read_nodes (  )  [private]

Helper functions for reading nodes/cells from a GMV file

Definition at line 2206 of file gmv_io.C.

References libMesh::MeshBase::add_point(), and libMesh::MeshInput< MeshBase >::mesh().

Referenced by read().

02207 {
02208 #ifdef LIBMESH_HAVE_GMV
02209   //   // Debug Info
02210   //   libMesh::out << "gmv_data.datatype="
02211   //        <<  GMV::gmv_data.datatype
02212   //        << std::endl;
02213 
02214   // LibMesh writes UNSTRUCT=100 node data
02215   libmesh_assert_equal_to (GMV::gmv_data.datatype, UNSTRUCT);
02216 
02217   // The nodal data is stored in gmv_data.doubledata{1,2,3}
02218   // and is nnodes long
02219   for (int i = 0; i < GMV::gmv_data.num; i++)
02220     {
02221       //       libMesh::out << "(x,y,z)="
02222       //                << "("
02223       //                << GMV::gmv_data.doubledata1[i]
02224       //                << ","
02225       //                << GMV::gmv_data.doubledata2[i]
02226       //                << ","
02227       //                << GMV::gmv_data.doubledata3[i]
02228       //                << ")"
02229       //                << std::endl;
02230 
02231       // Add the point to the Mesh
02232       MeshInput<MeshBase>::mesh().add_point
02233                      ( Point(GMV::gmv_data.doubledata1[i],
02234                              GMV::gmv_data.doubledata2[i],
02235                              GMV::gmv_data.doubledata3[i]), i);
02236     }
02237 #endif
02238 }

void libMesh::GMVIO::_read_one_cell (  )  [private]

Definition at line 2241 of file gmv_io.C.

References _gmv_elem_to_libmesh_elem(), _next_elem_id, libMesh::MeshBase::add_elem(), libMesh::Elem::build(), libMesh::Elem::dim(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::node_ptr(), libMesh::AutoPtr< Tp >::release(), libMesh::DofObject::set_id(), and libMesh::Elem::set_node().

Referenced by read().

02242 {
02243 #ifdef LIBMESH_HAVE_GMV
02244   //   // Debug Info
02245   //   libMesh::out << "gmv_data.datatype="
02246   //        <<  GMV::gmv_data.datatype
02247   //        << std::endl;
02248 
02249   // This is either a REGULAR=111 cell or
02250   // the ENDKEYWORD=207 of the cells
02251 #ifndef NDEBUG
02252     bool recognized =
02253       (GMV::gmv_data.datatype==REGULAR) ||
02254       (GMV::gmv_data.datatype==ENDKEYWORD);
02255 #endif
02256     libmesh_assert (recognized);
02257 
02258   MeshBase& mesh = MeshInput<MeshBase>::mesh();
02259 
02260   if (GMV::gmv_data.datatype == REGULAR)
02261     {
02262       //       libMesh::out << "Name of the cell is: "
02263       //                << GMV::gmv_data.name1
02264       //                << std::endl;
02265 
02266       //       libMesh::out << "Cell has "
02267       //                << GMV::gmv_data.num2
02268       //                << " vertices."
02269       //                << std::endl;
02270 
02271       // We need a mapping from GMV element types to LibMesh
02272       // ElemTypes.  Basically the reverse of the eletypes
02273       // std::map above.
02274       //
02275       // FIXME: Since Quad9's apparently don't exist for GMV, and since
02276       // In general we write linear sub-elements to GMV files, we need
02277       // to be careful to read back in exactly what we wrote out...
02278       ElemType type = this->_gmv_elem_to_libmesh_elem(GMV::gmv_data.name1);
02279 
02280       Elem* elem = Elem::build(type).release();
02281       elem->set_id(_next_elem_id++);
02282 
02283       // Get the ElementDefinition object for this element type
02284       const ElementDefinition& eledef = eletypes[type];
02285 
02286       // Print out the connectivity information for
02287       // this cell.
02288       for (int i=0; i<GMV::gmv_data.num2; i++)
02289         {
02290           //      // Debugging info
02291           //      libMesh::out << "Vertex " << i << " is node "
02292           //                << GMV::gmv_data.longdata1[i]
02293           //                << std::endl;
02294 
02295           // Map index i to GMV's numbering scheme
02296           unsigned mapped_i = eledef.node_map[i];
02297 
02298           // Note: Node numbers (as stored in libmesh) are 1-based
02299           elem->set_node(i) = mesh.node_ptr(GMV::gmv_data.longdata1[mapped_i]-1);
02300         }
02301 
02302       elems_of_dimension[elem->dim()] = true;
02303 
02304       // Add the newly-created element to the mesh
02305       mesh.add_elem(elem);
02306     }
02307 
02308 
02309   if (GMV::gmv_data.datatype == ENDKEYWORD)
02310     {
02311       // There isn't a cell to read, so we just return
02312       return;
02313     }
02314 
02315 #endif
02316 }

void libMesh::GMVIO::_read_var (  )  [private]

Definition at line 2151 of file gmv_io.C.

References _nodal_data.

Referenced by read().

02152 {
02153 #ifdef LIBMESH_HAVE_GMV
02154 
02155   // Copy all the variable's values into a local storage vector.
02156   _nodal_data.insert ( std::make_pair(std::string(GMV::gmv_data.name1),
02157                                       std::vector<Number>(GMV::gmv_data.doubledata1, GMV::gmv_data.doubledata1+GMV::gmv_data.num) ) );
02158 #endif
02159 }

void libMesh::GMVIO::add_cell_centered_data ( const std::string &  cell_centered_data_name,
const std::vector< Real > *  cell_centered_data_vals 
)

Takes a vector of cell-centered data to be plotted. You must ensure that for every active element e, v[e->id()] is a valid number. You can add an arbitrary number of different cell-centered data sets by calling this function multiple times.

.) GMV does not like spaces in the cell_centered_data_name .) No matter what order you add cell-centered data, it will be output alphabetically.

Definition at line 1966 of file gmv_io.C.

References _cell_centered_data.

01968 {
01969   libmesh_assert(cell_centered_data_vals);
01970 
01971   // Make sure there are *at least* enough entries for all the active elements.
01972   // There can also be entries for inactive elements, they will be ignored.
01973   // libmesh_assert_greater_equal (cell_centered_data_vals->size(),
01974   //                               MeshOutput<MeshBase>::mesh().n_active_elem());
01975   this->_cell_centered_data[cell_centered_data_name] = cell_centered_data_vals;
01976 }

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(), write_ascii_new_impl(), and write_ascii_old_impl().

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

Flag indicating whether or not to write a binary file. While binary files may end up being smaller than equivalent ASCII files, they will almost certainly take longer to write. The reason for this is that the ostream::write() function which is used to write "binary" data to streams, only takes a pointer to char as its first argument. This means if you want to write anything other than a buffer of chars, you first have to use a strange memcpy hack to get the data into the desired format. See the templated to_binary_stream() function below.

Definition at line 300 of file gmv_io.h.

References _binary.

Referenced by write(), and write_nodal_data().

00301 {
00302   return _binary;
00303 }

void libMesh::GMVIO::copy_nodal_solution ( EquationSystems es  ) 

If we read in a nodal solution while reading in a mesh, we can attempt to copy that nodal solution into an EquationSystems object.

Definition at line 2390 of file gmv_io.C.

References _nodal_data, libMesh::DofObject::dof_number(), end, libMesh::err, libMesh::FEType::family, libMeshEnums::FIRST, libMesh::EquationSystems::get_system(), libMesh::System::has_variable(), libMeshEnums::LAGRANGE, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_nodes(), libMesh::EquationSystems::n_systems(), libMesh::MeshBase::node_ptr(), libMesh::FEType::order, libMesh::System::solution, libMesh::System::update(), libMesh::System::variable_number(), and libMesh::System::variable_type().

02391 {
02392   // Check for easy return if there isn't any nodal data
02393   if (_nodal_data.empty())
02394     {
02395       libMesh::err << "Unable to copy nodal solution: No nodal "
02396                     << "solution has been read in from file." << std::endl;
02397       return;
02398     }
02399 
02400   // Be sure there is at least one system
02401   libmesh_assert (es.n_systems());
02402 
02403   // Keep track of variable names which have been found and
02404   // copied already.  This could be used to prevent us from
02405   // e.g. copying the same var into 2 different systems ...
02406   // but this seems unlikely.  Also, it is used to tell if
02407   // any variables which were read in were not successfully
02408   // copied to the EquationSystems.
02409   std::set<std::string> vars_copied;
02410 
02411   // For each entry in the nodal data map, try to find a system
02412   // that has the same variable key name.
02413   for (unsigned int sys=0; sys<es.n_systems(); ++sys)
02414     {
02415       // Get a generic refernence to the current System
02416       System& system = es.get_system(sys);
02417 
02418       // And a reference to that system's dof_map
02419       // const DofMap & dof_map = system.get_dof_map();
02420 
02421       // For each var entry in the _nodal_data map, try to find
02422       // that var in the system
02423       std::map<std::string, std::vector<Number> >::iterator it = _nodal_data.begin();
02424       const std::map<std::string, std::vector<Number> >::iterator end = _nodal_data.end();
02425       for (; it != end; ++it)
02426         {
02427           std::string var_name = (*it).first;
02428           // libMesh::out << "Searching for var " << var_name << " in system " << sys << std::endl;
02429 
02430           if (system.has_variable(var_name))
02431             {
02432               // Check if there are as many nodes in the mesh as there are entries
02433               // in the stored nodal data vector
02434               libmesh_assert_equal_to ( (*it).second.size(), MeshInput<MeshBase>::mesh().n_nodes() );
02435 
02436               const unsigned int var_num = system.variable_number(var_name);
02437 
02438               // libMesh::out << "Variable "
02439               //                        << var_name
02440               //                        << " is variable "
02441               //                        << var_num
02442               //                        << " in system " << sys << std::endl;
02443 
02444               // The only type of nodal data we can read in from GMV is for
02445               // linear LAGRANGE type elements.
02446               const FEType& fe_type = system.variable_type(var_num);
02447               if ((fe_type.order != FIRST) || (fe_type.family != LAGRANGE))
02448                 {
02449                   libMesh::err << "Only FIRST-order LAGRANGE variables can be read from GMV files. "
02450                                 << "Skipping variable " << var_name << std::endl;
02451                   //libmesh_error();
02452                   break;
02453                 }
02454 
02455 
02456               // Loop over the stored vector's entries, inserting them into
02457               // the System's solution if appropriate.
02458               for (unsigned int i=0; i<(*it).second.size(); ++i)
02459                 {
02460                   // Since this var came from a GMV file, the index i corresponds to
02461                   // the (single) DOF value of the current variable for node i.
02462                   const unsigned int dof_index =
02463                     MeshInput<MeshBase>::mesh().node_ptr(i)->dof_number(sys,      /*system #*/
02464                                                                         var_num,  /*var # */
02465                                                                         0);       /*component #, always zero for LAGRANGE */
02466 
02467                   // libMesh::out << "Value " << i << ": "
02468                   //                        << (*it).second [i]
02469                   //                        << ", dof index="
02470                   //                        << dof_index << std::endl;
02471 
02472                   // If the dof_index is local to this processor, set the value
02473                   if ((dof_index >= system.solution->first_local_index()) &&
02474                       (dof_index <  system.solution->last_local_index()))
02475                     system.solution->set (dof_index, (*it).second [i]);
02476                 } // end loop over my GMVIO's copy of the solution
02477 
02478               // Add the most recently copied var to the set of copied vars
02479               vars_copied.insert (var_name);
02480             } // end if (system.has_variable)
02481         } // end for loop over _nodal_data
02482 
02483       // Communicate parallel values before going to the next system.
02484       system.solution->close();
02485       system.update();
02486 
02487     } // end loop over all systems
02488 
02489 
02490 
02491   // Warn the user if any GMV variables were not successfully copied over to the EquationSystems object
02492   {
02493     std::map<std::string, std::vector<Number> >::iterator it = _nodal_data.begin();
02494     const std::map<std::string, std::vector<Number> >::iterator end = _nodal_data.end();
02495 
02496     for (; it != end; ++it)
02497       {
02498         if (vars_copied.find( (*it).first ) == vars_copied.end())
02499           {
02500             libMesh::err << "Warning: Variable "
02501                           << (*it).first
02502                           << " was not copied to the EquationSystems object."
02503                           << std::endl;
02504           }
02505       }
02506   }
02507 
02508 }

bool & libMesh::GMVIO::discontinuous (  )  [inline]

Flag indicating whether or not to write the mesh as discontinuous cell patches

Definition at line 308 of file gmv_io.h.

References _discontinuous.

00309 {
00310   return _discontinuous;
00311 }

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

Returns the object as a writeable reference.

Referenced by _read_materials(), _read_nodes(), _read_one_cell(), libMesh::AbaqusIO::assign_boundary_node_ids(), libMesh::AbaqusIO::assign_sideset_ids(), libMesh::AbaqusIO::assign_subdomain_ids(), libMesh::VTKIO::cells_to_vtk(), 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(), libMesh::XdrIO::read(), libMesh::VTKIO::read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), 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(), libMesh::XdrIO::read_serialized_bcs(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::XdrIO::write(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), write_ascii_new_impl(), write_ascii_old_impl(), write_binary(), 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(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::LegacyXdrIO::write_soln().

bool & libMesh::GMVIO::p_levels (  )  [inline]

Flag indicating whether or not to write p level information for p refined meshes

Definition at line 339 of file gmv_io.h.

References _p_levels.

Referenced by write_ascii_new_impl(), write_ascii_old_impl(), and write_binary().

00340 {
00341   return _p_levels;
00342 }

bool & libMesh::GMVIO::partitioning (  )  [inline]

Flag indicating whether or not to write the partitioning information for the mesh.

Definition at line 316 of file gmv_io.h.

References _partitioning.

Referenced by libMesh::UnstructuredMesh::write(), write_ascii_new_impl(), write_ascii_old_impl(), and write_binary().

00317 {
00318   return _partitioning;
00319 }

void libMesh::GMVIO::read ( const std::string &  mesh_file  )  [virtual]

This method implements reading a mesh from a specified file.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 1983 of file gmv_io.C.

References _next_elem_id, _read_materials(), _read_nodes(), _read_one_cell(), _read_var(), libMesh::MeshBase::clear(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Trees::NODES, libMesh::MeshBase::prepare_for_use(), libMesh::processor_id(), and libMesh::MeshBase::set_mesh_dimension().

01984 {
01985   // This is a serial-only process for now;
01986   // the Mesh should be read on processor 0 and
01987   // broadcast later
01988   libmesh_assert_equal_to (libMesh::processor_id(), 0);
01989 
01990   _next_elem_id = 0;
01991 
01992   libmesh_experimental();
01993 
01994 #ifndef LIBMESH_HAVE_GMV
01995 
01996   libMesh::err << "Cannot read GMV file " << name << " without the GMV API." << std::endl;
01997   libmesh_error();
01998 
01999 #else
02000   // We use the file-scope global variable eletypes for mapping nodes
02001   // from GMV to libmesh indices, so initialize that data now.
02002   init_eletypes();
02003 
02004   // Clear the mesh so we are sure to start from a pristeen state.
02005   MeshBase& mesh = MeshInput<MeshBase>::mesh();
02006   mesh.clear();
02007 
02008   // Keep track of what kinds of elements this file contains
02009   elems_of_dimension.clear();
02010   elems_of_dimension.resize(4, false);
02011 
02012   // It is apparently possible for gmv files to contain
02013   // a "fromfile" directive (?) But we currently don't make
02014   // any use of this feature in LibMesh.  Nonzero return val
02015   // from any function usually means an error has occurred.
02016   int ierr = GMV::gmvread_open_fromfileskip(const_cast<char*>(name.c_str()));
02017   if (ierr != 0)
02018     {
02019       libMesh::err << "GMV::gmvread_open_fromfileskip failed!" << std::endl;
02020       libmesh_error();
02021     }
02022 
02023 
02024   // Loop through file until GMVEND.
02025   int iend = 0;
02026   while (iend == 0)
02027     {
02028       GMV::gmvread_data();
02029 
02030       /*  Check for GMVEND.  */
02031       if (GMV::gmv_data.keyword == GMVEND)
02032         {
02033           iend = 1;
02034           GMV::gmvread_close();
02035           break;
02036         }
02037 
02038       /*  Check for GMVERROR.  */
02039       if (GMV::gmv_data.keyword == GMVERROR)
02040         {
02041           libMesh::err << "Encountered GMVERROR while reading!" << std::endl;
02042           libmesh_error();
02043         }
02044 
02045       /*  Process the data.  */
02046       switch (GMV::gmv_data.keyword)
02047         {
02048         case NODES:
02049           {
02050             //libMesh::out << "Reading nodes." << std::endl;
02051 
02052             if (GMV::gmv_data.num2 == NODES)
02053               this->_read_nodes();
02054 
02055             else if (GMV::gmv_data.num2 == NODE_V)
02056               {
02057                 libMesh::err << "Unsupported GMV data type NODE_V!" << std::endl;
02058                 libmesh_error();
02059               }
02060             break;
02061           }
02062 
02063         case CELLS:
02064           {
02065             // Read 1 cell at a time
02066             // libMesh::out << "\nReading one cell." << std::endl;
02067             this->_read_one_cell();
02068             break;
02069           }
02070 
02071         case MATERIAL:
02072           {
02073             // keyword == 6
02074             // These are the materials, which we use to specify the mesh
02075             // partitioning.
02076             this->_read_materials();
02077             break;
02078           }
02079 
02080         case VARIABLE:
02081           {
02082             // keyword == 8
02083             // This is a field variable.
02084 
02085             // Check to see if we're done reading variables and break out.
02086             if (GMV::gmv_data.datatype == ENDKEYWORD)
02087               {
02088                 // libMesh::out << "Done reading GMV variables." << std::endl;
02089                 break;
02090               }
02091 
02092             if (GMV::gmv_data.datatype == NODE)
02093               {
02094                 // libMesh::out << "Reading node field data for variable "
02095                 //        << GMV::gmv_data.name1 << std::endl;
02096                 this->_read_var();
02097                 break;
02098               }
02099 
02100             else
02101               {
02102                 libMesh::err << "Warning: Skipping variable: "
02103                               << GMV::gmv_data.name1
02104                               << " which is of unsupported GMV datatype "
02105                               << GMV::gmv_data.datatype
02106                               << ".  Nodal field data is currently the only type currently supported."
02107                               << std::endl;
02108                 break;
02109               }
02110 
02111           }
02112 
02113         default:
02114           {
02115             libMesh::err << "Encountered unknown GMV keyword "
02116                           << GMV::gmv_data.keyword
02117                           << std::endl;
02118             libmesh_error();
02119           }
02120         } // end switch
02121     } // end while
02122 
02123   // Set the mesh dimension to the largest encountered for an element
02124   for (unsigned int i=0; i!=4; ++i)
02125     if (elems_of_dimension[i])
02126       mesh.set_mesh_dimension(i);
02127 
02128 #if LIBMESH_DIM < 3
02129   if (mesh.mesh_dimension() > LIBMESH_DIM)
02130     {
02131       libMesh::err << "Cannot open dimension " <<
02132                       mesh.mesh_dimension() <<
02133                       " mesh file when configured without " <<
02134                       mesh.mesh_dimension() << "D support." <<
02135                       std::endl;
02136       libmesh_error();
02137     }
02138 #endif
02139 
02140   // Done reading in the mesh, now call find_neighbors, etc.
02141   // mesh.find_neighbors();
02142 
02143   // Pass true flag to skip renumbering nodes and elements
02144   mesh.prepare_for_use(true);
02145 #endif
02146 }

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().

bool & libMesh::GMVIO::subdivide_second_order (  )  [inline]

Flag indicating whether or not to subdivide second order elements

Definition at line 331 of file gmv_io.h.

References _subdivide_second_order.

Referenced by write_ascii_new_impl(), and write_ascii_old_impl().

00332 {
00333   return _subdivide_second_order;
00334 }

template<typename T >
void libMesh::GMVIO::to_binary_stream ( std::ostream &  out,
const T  i 
) [inline, private]

Helper function for writing unsigned ints to an ostream in binary format. Implemented via memcpy as suggested in the standard.

Definition at line 347 of file gmv_io.h.

00349 {
00350   static char buf[sizeof(T)];
00351   memcpy(buf, &i, sizeof(T));
00352   out_str.write(buf, sizeof(T));
00353 }

void libMesh::GMVIO::write ( const std::string &  fname  )  [virtual]

This method implements writing a mesh to a specified file.

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 205 of file gmv_io.C.

References binary(), write_ascii_old_impl(), and write_binary().

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

00206 {
00207   if (this->binary())
00208     this->write_binary (fname);
00209   else
00210     this->write_ascii_old_impl  (fname);
00211 }

void libMesh::GMVIO::write_ascii_new_impl ( const std::string &  fname,
const std::vector< Number > *  v = NULL,
const std::vector< std::string > *  solution_names = NULL 
)

This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. This will write an ASCII file. This is the new implementation (without subcells).

Definition at line 231 of file gmv_io.C.

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshOutput< MeshBase >::ascii_precision(), end, libMesh::err, libMesh::DofObject::id(), std::max(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_active_elem(), libMesh::Elem::n_nodes(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::n_partitions(), libMesh::Elem::n_sub_elem(), n_vars, libMesh::Elem::node(), libMesh::out, libMesh::Elem::p_level(), p_levels(), partitioning(), libMesh::MeshBase::point(), libMesh::processor_id(), libMesh::Real, subdivide_second_order(), libMesh::Elem::type(), write_ascii_old_impl(), and write_subdomain_id_as_material().

00234 {
00235 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00236 
00237   libMesh::err << "WARNING:  GMVIO::write_ascii_new_impl() not infinite-element aware!"
00238                 << std::endl;
00239   libmesh_here();
00240 
00241   // Set it to our current precision
00242   this->write_ascii_old_impl (fname, v, solution_names);
00243 
00244 #else
00245 
00246   // Get a reference to the mesh
00247   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
00248 
00249   // This is a parallel_only function
00250   const unsigned int n_active_elem = mesh.n_active_elem();
00251 
00252   if (libMesh::processor_id() != 0)
00253     return;
00254 
00255   // Open the output file stream
00256   std::ofstream out_stream (fname.c_str());
00257 
00258   out_stream << std::setprecision(this->ascii_precision());
00259 
00260   // Make sure it opened correctly
00261   if (!out_stream.good())
00262     libmesh_file_error(fname.c_str());
00263 
00264   unsigned int mesh_max_p_level = 0;
00265 
00266   // Begin interfacing with the GMV data file
00267   {
00268     out_stream << "gmvinput ascii\n\n";
00269 
00270     // write the nodes
00271     out_stream << "nodes " << mesh.n_nodes() << "\n";
00272     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00273       out_stream << mesh.point(n)(0) << " ";
00274     out_stream << "\n";
00275 
00276     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00277 #if LIBMESH_DIM > 1
00278       out_stream << mesh.point(n)(1) << " ";
00279 #else
00280       out_stream << 0. << " ";
00281 #endif
00282     out_stream << "\n";
00283 
00284     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00285 #if LIBMESH_DIM > 2
00286       out_stream << mesh.point(n)(2) << " ";
00287 #else
00288       out_stream << 0. << " ";
00289 #endif
00290     out_stream << "\n\n";
00291   }
00292 
00293   {
00294     // write the connectivity
00295     out_stream << "cells " << n_active_elem << "\n";
00296 
00297     // initialize the eletypes map (eletypes is a file-global variable)
00298     init_eletypes();
00299 
00300     MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00301     const MeshBase::const_element_iterator end = mesh.active_elements_end();
00302 
00303     for ( ; it != end; ++it)
00304       {
00305         const Elem* elem = *it;
00306 
00307         mesh_max_p_level = std::max(mesh_max_p_level,
00308                                     elem->p_level());
00309 
00310         // Make sure we have a valid entry for
00311         // the current element type.
00312         libmesh_assert (eletypes.count(elem->type()));
00313 
00314         const ElementDefinition& ele = eletypes[elem->type()];
00315 
00316         // The element mapper better not require any more nodes
00317         // than are present in the current element!
00318         libmesh_assert_less_equal (ele.node_map.size(), elem->n_nodes());
00319 
00320         out_stream << ele.label << "\n";
00321         for (unsigned int i=0; i < ele.node_map.size(); i++)
00322           out_stream << elem->node(ele.node_map[i])+1 << " ";
00323         out_stream << "\n";
00324       }
00325     out_stream << "\n";
00326   }
00327 
00328   // optionally write the partition information
00329   if (this->partitioning())
00330     {
00331       if (this->write_subdomain_id_as_material())
00332         {
00333           libMesh::out << "Not yet supported in GMVIO::write_ascii_new_impl" << std::endl;
00334           libmesh_error();
00335         }
00336       else // write processor IDs as materials.  This is the default
00337         {
00338           out_stream << "material "
00339               << mesh.n_partitions()
00340             // Note: GMV may give you errors like
00341             // Error, material for cell 1 is greater than 1
00342             // Error, material for cell 2 is greater than 1
00343             // Error, material for cell 3 is greater than 1
00344             // ... because you put the wrong number of partitions here.
00345             // To ensure you write the correct number of materials, call
00346             // mesh.recalculate_n_partitions() before writing out the
00347             // mesh.
00348             // Note: we can't call it now because the Mesh is const here and
00349             // it is a non-const function.
00350               << " 0\n";
00351 
00352           for (unsigned int proc=0; proc<mesh.n_partitions(); proc++)
00353             out_stream << "proc_" << proc << "\n";
00354 
00355           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00356           const MeshBase::const_element_iterator end = mesh.active_elements_end();
00357 
00358           // FIXME - don't we need to use an ElementDefinition here? - RHS
00359           for ( ; it != end; ++it)
00360             out_stream << (*it)->processor_id()+1 << "\n";
00361           out_stream << "\n";
00362         }
00363     }
00364 
00365   // If there are *any* variables at all in the system (including
00366   // p level, or arbitrary cell-based data)
00367   // to write, the gmv file needs to contain the word "variable"
00368   // on a line by itself.
00369   bool write_variable = false;
00370 
00371   // 1.) p-levels
00372   if (this->p_levels() && mesh_max_p_level)
00373     write_variable = true;
00374 
00375   // 2.) solution data
00376   if ((solution_names != NULL) && (v != NULL))
00377     write_variable = true;
00378 
00379   // 3.) cell-centered data
00380   if ( !(this->_cell_centered_data.empty()) )
00381     write_variable = true;
00382 
00383   if (write_variable)
00384     out_stream << "variable\n";
00385 
00386 //   if ((this->p_levels() && mesh_max_p_level) ||
00387 //     ((solution_names != NULL) && (v != NULL)))
00388 //     out_stream << "variable\n";
00389 
00390   // optionally write the polynomial degree information
00391   if (this->p_levels() && mesh_max_p_level)
00392     {
00393       out_stream << "p_level 0\n";
00394 
00395       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00396       const MeshBase::const_element_iterator end = mesh.active_elements_end();
00397 
00398       for ( ; it != end; ++it)
00399         {
00400           const Elem* elem = *it;
00401 
00402           const ElementDefinition& ele = eletypes[elem->type()];
00403 
00404           // The element mapper better not require any more nodes
00405           // than are present in the current element!
00406           libmesh_assert_less_equal (ele.node_map.size(), elem->n_nodes());
00407 
00408           for (unsigned int i=0; i < ele.node_map.size(); i++)
00409             out_stream << elem->p_level() << " ";
00410         }
00411       out_stream << "\n\n";
00412     }
00413 
00414 
00415   // optionally write cell-centered data
00416   if ( !(this->_cell_centered_data.empty()) )
00417     {
00418       std::map<std::string, const std::vector<Real>* >::iterator       it  = this->_cell_centered_data.begin();
00419       const std::map<std::string, const std::vector<Real>* >::iterator end = this->_cell_centered_data.end();
00420 
00421       for (; it != end; ++it)
00422         {
00423           // write out the variable name, followed by a zero.
00424           out_stream << (*it).first << " 0\n";
00425 
00426           const std::vector<Real>* the_array = (*it).second;
00427 
00428           // Loop over active elements, write out cell data.  If second-order cells
00429           // are split into sub-elements, the sub-elements inherit their parent's
00430           // cell-centered data.
00431           MeshBase::const_element_iterator       elem_it  = mesh.active_elements_begin();
00432           const MeshBase::const_element_iterator elem_end = mesh.active_elements_end();
00433 
00434           for (; elem_it != elem_end; ++elem_it)
00435             {
00436               const Elem* e = *elem_it;
00437 
00438               // Use the element's ID to find the value.
00439               libmesh_assert_less (e->id(), the_array->size());
00440               const Real the_value = the_array->operator[](e->id());
00441 
00442               if (this->subdivide_second_order())
00443                 for (unsigned int se=0; se < e->n_sub_elem(); se++)
00444                   out_stream << the_value << " ";
00445               else
00446                 out_stream << the_value << " ";
00447             }
00448 
00449           out_stream << "\n\n";
00450         }
00451     }
00452 
00453 
00454   // optionally write the data
00455   if ((solution_names != NULL) && (v != NULL))
00456     {
00457       const unsigned int n_vars = solution_names->size();
00458 
00459       if (!(v->size() == mesh.n_nodes()*n_vars))
00460         libMesh::err << "ERROR: v->size()=" << v->size()
00461                       << ", mesh.n_nodes()=" << mesh.n_nodes()
00462                       << ", n_vars=" << n_vars
00463                       << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
00464                       << "\n";
00465 
00466       libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
00467 
00468       for (unsigned int c=0; c<n_vars; c++)
00469         {
00470 
00471 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00472 
00473           // in case of complex data, write _three_ data sets
00474           // for each component
00475 
00476           // this is the real part
00477           out_stream << "r_" << (*solution_names)[c] << " 1\n";
00478 
00479           for (unsigned int n=0; n<mesh.n_nodes(); n++)
00480             out_stream << (*v)[n*n_vars + c].real() << " ";
00481 
00482           out_stream << "\n\n";
00483 
00484           // this is the imaginary part
00485           out_stream << "i_" << (*solution_names)[c] << " 1\n";
00486 
00487           for (unsigned int n=0; n<mesh.n_nodes(); n++)
00488             out_stream << (*v)[n*n_vars + c].imag() << " ";
00489 
00490           out_stream << "\n\n";
00491 
00492           // this is the magnitude
00493           out_stream << "a_" << (*solution_names)[c] << " 1\n";
00494           for (unsigned int n=0; n<mesh.n_nodes(); n++)
00495             out_stream << std::abs((*v)[n*n_vars + c]) << " ";
00496 
00497           out_stream << "\n\n";
00498 
00499 #else
00500 
00501           out_stream << (*solution_names)[c] << " 1\n";
00502 
00503           for (unsigned int n=0; n<mesh.n_nodes(); n++)
00504             out_stream << (*v)[n*n_vars + c] << " ";
00505 
00506           out_stream << "\n\n";
00507 
00508 #endif
00509         }
00510 
00511     }
00512 
00513   // If we wrote any variables, we have to close the variable section now
00514   if (write_variable)
00515     out_stream << "endvars\n";
00516 
00517 
00518   // end of the file
00519   out_stream << "\nendgmv\n";
00520 
00521 #endif
00522 }

void libMesh::GMVIO::write_ascii_old_impl ( const std::string &  fname,
const std::vector< Number > *  v = NULL,
const std::vector< std::string > *  solution_names = NULL 
) [private]

This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. This will write an ASCII file. This is the old implementation (using subcells) which was the default in libMesh-0.4.3-rc2.

Note that the prisms are treated as degenerated phex8's.

Note that the prisms are treated as degenerated phex8's.

Definition at line 529 of file gmv_io.C.

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshOutput< MeshBase >::ascii_precision(), libMesh::Elem::build(), end, libMesh::err, libMeshEnums::FIRST, libMesh::Elem::first_order_equivalent_type(), libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMesh::DofObject::id(), libMeshEnums::INFHEX16, libMeshEnums::INFHEX18, libMeshEnums::INFHEX8, libMeshEnums::INFPRISM12, libMeshEnums::INFPRISM6, libMeshEnums::INFQUAD4, libMeshEnums::INFQUAD6, std::max(), libMesh::MeshBase::max_node_id(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_active_sub_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::n_partitions(), libMesh::Elem::n_sub_elem(), n_vars, libMesh::out, p_levels(), partitioning(), libMesh::MeshBase::point(), libMeshEnums::PRISM15, libMeshEnums::PRISM18, libMeshEnums::PRISM6, libMesh::processor_id(), libMeshEnums::PYRAMID5, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMesh::Real, subdivide_second_order(), libMeshEnums::TECPLOT, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, libMeshEnums::TRI6, and write_subdomain_id_as_material().

Referenced by write(), write_ascii_new_impl(), and write_nodal_data().

00532 {
00533   // Get a reference to the mesh
00534   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
00535 
00536   // Use a MeshSerializer object to gather a parallel mesh before outputting it.
00537   // Note that we cast away constness here (which is bad), but the destructor of
00538   // the MeshSerializer object reparallelizes the Mesh, hopefully keeping it
00539   // "logically const" outside the context of this function...
00540   MeshSerializer serialize(const_cast<MeshBase&>(mesh),
00541                            !MeshOutput<MeshBase>::_is_parallel_format);
00542 
00543   // These are parallel_only functions
00544   const dof_id_type n_active_elem = mesh.n_active_elem(),
00545                     n_active_sub_elem = mesh.n_active_sub_elem();
00546 
00547   if (libMesh::processor_id() != 0)
00548     return;
00549 
00550   // Open the output file stream
00551   std::ofstream out_stream (fname.c_str());
00552 
00553   // Set it to our current precision
00554   out_stream << std::setprecision(this->ascii_precision());
00555 
00556   // Make sure it opened correctly
00557   if (!out_stream.good())
00558     libmesh_file_error(fname.c_str());
00559 
00560   // Make sure our nodes are contiguous and serialized
00561   libmesh_assert_equal_to (mesh.n_nodes(), mesh.max_node_id());
00562 
00563   // libmesh_assert (mesh.is_serial());
00564   // if (!mesh.is_serial())
00565   //   {
00566   //     if (libMesh::processor_id() == 0)
00567   //       libMesh::err << "Error: GMVIO cannot yet write a ParallelMesh solution"
00568   //                     << std::endl;
00569   //     return;
00570   //   }
00571 
00572   unsigned int mesh_max_p_level = 0;
00573 
00574   // Begin interfacing with the GMV data file
00575 
00576   // FIXME - if subdivide_second_order() is off,
00577   // we probably should only be writing the
00578   // vertex nodes - RHS
00579   {
00580     // write the nodes
00581 
00582     out_stream << "gmvinput ascii\n\n";
00583     out_stream << "nodes " << mesh.n_nodes() << '\n';
00584     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00585       out_stream << mesh.point(n)(0) << " ";
00586 
00587     out_stream << '\n';
00588 
00589     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00590 #if LIBMESH_DIM > 1
00591       out_stream << mesh.point(n)(1) << " ";
00592 #else
00593       out_stream << 0. << " ";
00594 #endif
00595 
00596     out_stream << '\n';
00597 
00598     for (unsigned int n=0; n<mesh.n_nodes(); n++)
00599 #if LIBMESH_DIM > 2
00600       out_stream << mesh.point(n)(2) << " ";
00601 #else
00602       out_stream << 0. << " ";
00603 #endif
00604 
00605     out_stream << '\n' << '\n';
00606   }
00607 
00608 
00609 
00610   {
00611     // write the connectivity
00612 
00613     out_stream << "cells ";
00614     if (this->subdivide_second_order())
00615       out_stream << n_active_sub_elem;
00616     else
00617       out_stream << n_active_elem;
00618     out_stream << '\n';
00619 
00620     MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00621     const MeshBase::const_element_iterator end = mesh.active_elements_end();
00622 
00623     switch (mesh.mesh_dimension())
00624       {
00625       case 1:
00626         {
00627           // The same temporary storage will be used for each element
00628           std::vector<dof_id_type> conn;
00629 
00630           for ( ; it != end; ++it)
00631             {
00632               mesh_max_p_level = std::max(mesh_max_p_level,
00633                                           (*it)->p_level());
00634 
00635               if (this->subdivide_second_order())
00636                 for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
00637                   {
00638                     out_stream << "line 2\n";
00639                     (*it)->connectivity(se, TECPLOT, conn);
00640                     for (unsigned int i=0; i<conn.size(); i++)
00641                       out_stream << conn[i] << " ";
00642 
00643                     out_stream << '\n';
00644                   }
00645               else
00646                 {
00647                   out_stream << "line 2\n";
00648                   if ((*it)->default_order() == FIRST)
00649                     (*it)->connectivity(0, TECPLOT, conn);
00650                   else
00651                     {
00652                       AutoPtr<Elem> lo_elem = Elem::build(
00653                         Elem::first_order_equivalent_type((*it)->type()));
00654                       for (unsigned int i = 0; i != lo_elem->n_nodes(); ++i)
00655                         lo_elem->set_node(i) = (*it)->get_node(i);
00656                       lo_elem->connectivity(0, TECPLOT, conn);
00657                     }
00658                   for (unsigned int i=0; i<conn.size(); i++)
00659                     out_stream << conn[i] << " ";
00660 
00661                   out_stream << '\n';
00662                 }
00663             }
00664           break;
00665         }
00666 
00667       case 2:
00668         {
00669           // The same temporary storage will be used for each element
00670           std::vector<dof_id_type> conn;
00671 
00672           for ( ; it != end; ++it)
00673             {
00674               mesh_max_p_level = std::max(mesh_max_p_level,
00675                                           (*it)->p_level());
00676 
00677               if (this->subdivide_second_order())
00678                 for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
00679                   {
00680                     // Quad elements
00681                     if (((*it)->type() == QUAD4) ||
00682                         ((*it)->type() == QUAD8) || // Note: QUAD8 will be output as one central quad and
00683                                                     // four surrounding triangles (though they will be written
00684                                                     // to GMV as QUAD4s).
00685                         ((*it)->type() == QUAD9)
00686 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00687                         || ((*it)->type() == INFQUAD4)
00688                         || ((*it)->type() == INFQUAD6)
00689 #endif
00690                         )
00691                       {
00692                         out_stream << "quad 4\n";
00693                         (*it)->connectivity(se, TECPLOT, conn);
00694                         for (unsigned int i=0; i<conn.size(); i++)
00695                           out_stream << conn[i] << " ";
00696                       }
00697 
00698                     // Triangle elements
00699                     else if (((*it)->type() == TRI3) ||
00700                              ((*it)->type() == TRI6))
00701                       {
00702                         out_stream << "tri 3\n";
00703                         (*it)->connectivity(se, TECPLOT, conn);
00704                         for (unsigned int i=0; i<3; i++)
00705                           out_stream << conn[i] << " ";
00706                       }
00707                     else
00708                       libmesh_error();
00709                   }
00710               else // !this->subdivide_second_order()
00711                 {
00712                   // Quad elements
00713                   if (((*it)->type() == QUAD4)
00714 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00715                       || ((*it)->type() == INFQUAD4)
00716 #endif
00717                       )
00718                     {
00719                       (*it)->connectivity(0, TECPLOT, conn);
00720                       out_stream << "quad 4\n";
00721                       for (unsigned int i=0; i<conn.size(); i++)
00722                         out_stream << conn[i] << " ";
00723                     }
00724                   else if (((*it)->type() == QUAD8) ||
00725                            ((*it)->type() == QUAD9)
00726 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
00727                            || ((*it)->type() == INFQUAD6)
00728 #endif
00729                           )
00730                     {
00731                       AutoPtr<Elem> lo_elem = Elem::build(
00732                         Elem::first_order_equivalent_type((*it)->type()));
00733                       for (unsigned int i = 0; i != lo_elem->n_nodes(); ++i)
00734                         lo_elem->set_node(i) = (*it)->get_node(i);
00735                       lo_elem->connectivity(0, TECPLOT, conn);
00736                       out_stream << "quad 4\n";
00737                       for (unsigned int i=0; i<conn.size(); i++)
00738                         out_stream << conn[i] << " ";
00739                     }
00740                   else if ((*it)->type() == TRI3)
00741                     {
00742                       (*it)->connectivity(0, TECPLOT, conn);
00743                       out_stream << "tri 3\n";
00744                       for (unsigned int i=0; i<3; i++)
00745                         out_stream << conn[i] << " ";
00746                     }
00747                   else if ((*it)->type() == TRI6)
00748                     {
00749                       AutoPtr<Elem> lo_elem = Elem::build(
00750                         Elem::first_order_equivalent_type((*it)->type()));
00751                       for (unsigned int i = 0; i != lo_elem->n_nodes(); ++i)
00752                         lo_elem->set_node(i) = (*it)->get_node(i);
00753                       lo_elem->connectivity(0, TECPLOT, conn);
00754                       out_stream << "tri 3\n";
00755                       for (unsigned int i=0; i<3; i++)
00756                         out_stream << conn[i] << " ";
00757                     }
00758 
00759                   out_stream << '\n';
00760                 }
00761             }
00762 
00763           break;
00764         }
00765 
00766 
00767       case 3:
00768         {
00769           // The same temporary storage will be used for each element
00770           std::vector<dof_id_type> conn;
00771 
00772           for ( ; it != end; ++it)
00773             {
00774               mesh_max_p_level = std::max(mesh_max_p_level,
00775                                           (*it)->p_level());
00776 
00777               if (this->subdivide_second_order())
00778                 for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
00779                   {
00780 
00781 #ifndef  LIBMESH_ENABLE_INFINITE_ELEMENTS
00782                     if (((*it)->type() == HEX8)   ||
00783                         ((*it)->type() == HEX27))
00784                       {
00785                         out_stream << "phex8 8\n";
00786                         (*it)->connectivity(se, TECPLOT, conn);
00787                         for (unsigned int i=0; i<conn.size(); i++)
00788                           out_stream << conn[i] << " ";
00789                       }
00790 
00791                     else if ((*it)->type() == HEX20)
00792                       {
00793                         out_stream << "phex20 20\n";
00794                         out_stream << (*it)->node(0)+1  << " "
00795                             << (*it)->node(1)+1  << " "
00796                             << (*it)->node(2)+1  << " "
00797                             << (*it)->node(3)+1  << " "
00798                             << (*it)->node(4)+1  << " "
00799                             << (*it)->node(5)+1  << " "
00800                             << (*it)->node(6)+1  << " "
00801                             << (*it)->node(7)+1  << " "
00802                             << (*it)->node(8)+1  << " "
00803                             << (*it)->node(9)+1  << " "
00804                             << (*it)->node(10)+1 << " "
00805                             << (*it)->node(11)+1 << " "
00806                             << (*it)->node(16)+1 << " "
00807                             << (*it)->node(17)+1 << " "
00808                             << (*it)->node(18)+1 << " "
00809                             << (*it)->node(19)+1 << " "
00810                             << (*it)->node(12)+1 << " "
00811                             << (*it)->node(13)+1 << " "
00812                             << (*it)->node(14)+1 << " "
00813                             << (*it)->node(15)+1 << " ";
00814                       }
00815 #else
00816                     /*
00817                      * In case of infinite elements, HEX20
00818                      * should be handled just like the
00819                      * INFHEX16, since these connect to each other
00820                      */
00821                     if (((*it)->type() == HEX8)     ||
00822                         ((*it)->type() == HEX27)    ||
00823                         ((*it)->type() == INFHEX8)  ||
00824                         ((*it)->type() == INFHEX16) ||
00825                         ((*it)->type() == INFHEX18) ||
00826                         ((*it)->type() == HEX20))
00827                       {
00828                         out_stream << "phex8 8\n";
00829                         (*it)->connectivity(se, TECPLOT, conn);
00830                         for (unsigned int i=0; i<conn.size(); i++)
00831                           out_stream << conn[i] << " ";
00832                       }
00833 #endif
00834 
00835                     else if (((*it)->type() == TET4)  ||
00836                              ((*it)->type() == TET10))
00837                       {
00838                         out_stream << "tet 4\n";
00839                         // Tecplot connectivity returns 8 entries for
00840                         // the Tet, enough to store it as a degenerate Hex.
00841                         // For GMV we only pick out the four relevant node
00842                         // indices.
00843                         (*it)->connectivity(se, TECPLOT, conn);
00844                         out_stream << conn[0] << " "  // libmesh tet node 0
00845                             << conn[2] << " "  // libmesh tet node 2
00846                             << conn[1] << " "  // libmesh tet node 1
00847                             << conn[4] << " "; // libmesh tet node 3
00848                       }
00849 #ifndef  LIBMESH_ENABLE_INFINITE_ELEMENTS
00850                     else if (((*it)->type() == PRISM6)  ||
00851                              ((*it)->type() == PRISM15) ||
00852                              ((*it)->type() == PRISM18) ||
00853                              ((*it)->type() == PYRAMID5))
00854 #else
00855                     else if (((*it)->type() == PRISM6)     ||
00856                              ((*it)->type() == PRISM15)    ||
00857                              ((*it)->type() == PRISM18)    ||
00858                              ((*it)->type() == PYRAMID5)   ||
00859                              ((*it)->type() == INFPRISM6)  ||
00860                              ((*it)->type() == INFPRISM12))
00861 #endif
00862                       {
00867                         out_stream << "phex8 8\n";
00868                         (*it)->connectivity(se, TECPLOT, conn);
00869                         for (unsigned int i=0; i<conn.size(); i++)
00870                           out_stream << conn[i] << " ";
00871                       }
00872 
00873                     else
00874                       {
00875                         libMesh::out << "Encountered an unrecognized element "
00876                                       << "type: " << (*it)->type()
00877                                       << "\nPossibly a dim-1 dimensional "
00878                                       << "element?  Aborting..."
00879                                       << std::endl;
00880                         libmesh_error();
00881                       }
00882 
00883                     out_stream << '\n';
00884                   }
00885               else // !this->subdivide_second_order()
00886                 {
00887                   AutoPtr<Elem> lo_elem = Elem::build(
00888                     Elem::first_order_equivalent_type((*it)->type()));
00889                   for (unsigned int i = 0; i != lo_elem->n_nodes(); ++i)
00890                     lo_elem->set_node(i) = (*it)->get_node(i);
00891                   if ((lo_elem->type() == HEX8)
00892 #ifdef  LIBMESH_ENABLE_INFINITE_ELEMENTS
00893                       || (lo_elem->type() == HEX27)
00894 #endif
00895                      )
00896                     {
00897                       out_stream << "phex8 8\n";
00898                       lo_elem->connectivity(0, TECPLOT, conn);
00899                       for (unsigned int i=0; i<conn.size(); i++)
00900                         out_stream << conn[i] << " ";
00901                     }
00902 
00903                   else if (lo_elem->type() == TET4)
00904                     {
00905                       out_stream << "tet 4\n";
00906                       lo_elem->connectivity(0, TECPLOT, conn);
00907                       out_stream << conn[0] << " "
00908                           << conn[2] << " "
00909                           << conn[1] << " "
00910                           << conn[4] << " ";
00911                     }
00912                   else if ((lo_elem->type() == PRISM6)
00913 #ifdef  LIBMESH_ENABLE_INFINITE_ELEMENTS
00914                            || (lo_elem->type() == INFPRISM6)
00915 #endif
00916                            )
00917                     {
00922                       out_stream << "phex8 8\n";
00923                       lo_elem->connectivity(0, TECPLOT, conn);
00924                       for (unsigned int i=0; i<conn.size(); i++)
00925                         out_stream << conn[i] << " ";
00926                     }
00927 
00928                   else
00929                     {
00930                       libMesh::out << "Encountered an unrecognized element "
00931                                     << "type.  Possibly a dim-1 dimensional "
00932                                     << "element?  Aborting..."
00933                                     << std::endl;
00934                       libmesh_error();
00935                     }
00936 
00937                   out_stream << '\n';
00938                 }
00939             }
00940 
00941           break;
00942         }
00943 
00944       default:
00945         libmesh_error();
00946       }
00947 
00948     out_stream << '\n';
00949   }
00950 
00951 
00952 
00953   // optionally write the partition information
00954   if (this->partitioning())
00955     {
00956       if (this->write_subdomain_id_as_material())
00957         {
00958           // Subdomain IDs can be non-contiguous and need not
00959           // necessarily start at 0.  Furthermore, since the user is
00960           // free to define subdomain_id_type to be a signed type, we
00961           // can't even assume max(subdomain_id) >= # unique subdomain ids.
00962 
00963           // We build a map<subdomain_id, unsigned> to associate to each
00964           // user-selected subdomain ID a unique, contiguous unsigned value
00965           // which we can write to file.
00966           std::map<subdomain_id_type, unsigned> sbdid_map;
00967           typedef std::map<subdomain_id_type, unsigned>::iterator sbdid_map_iter;
00968           {
00969             MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00970             const MeshBase::const_element_iterator end = mesh.active_elements_end();
00971 
00972             for ( ; it != end; ++it)
00973               {
00974                 // Try to insert with dummy value
00975                 sbdid_map.insert( std::make_pair((*it)->subdomain_id(), 0) );
00976               }
00977           }
00978 
00979           // Map is created, iterate through it to set indices.  They will be
00980           // used repeatedly below.
00981           {
00982             unsigned ctr=0;
00983             for (sbdid_map_iter it=sbdid_map.begin(); it != sbdid_map.end(); ++it)
00984               (*it).second = ctr++;
00985           }
00986 
00987           out_stream << "material "
00988               << sbdid_map.size()
00989               << " 0\n";
00990 
00991           for (unsigned int sbdid=0; sbdid<sbdid_map.size(); sbdid++)
00992             out_stream << "proc_" << sbdid << "\n";
00993 
00994           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00995           const MeshBase::const_element_iterator end = mesh.active_elements_end();
00996 
00997           for ( ; it != end; ++it)
00998             {
00999               // Find the unique index for (*it)->subdomain_id(), print that to file
01000               sbdid_map_iter map_iter = sbdid_map.find( (*it)->subdomain_id() );
01001               unsigned gmv_mat_number = (*map_iter).second;
01002 
01003               if (this->subdivide_second_order())
01004                 for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01005                   out_stream << gmv_mat_number+1 << '\n';
01006               else
01007                 out_stream << gmv_mat_number+1 << "\n";
01008             }
01009           out_stream << '\n';
01010 
01011         }
01012       else // write processor IDs as materials.  This is the default
01013         {
01014           out_stream << "material "
01015               << mesh.n_partitions()
01016               << " 0"<< '\n';
01017 
01018           for (unsigned int proc=0; proc<mesh.n_partitions(); proc++)
01019             out_stream << "proc_" << proc << '\n';
01020 
01021           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01022           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01023 
01024           for ( ; it != end; ++it)
01025             if (this->subdivide_second_order())
01026               for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01027                 out_stream << (*it)->processor_id()+1 << '\n';
01028             else
01029               out_stream << (*it)->processor_id()+1 << '\n';
01030 
01031           out_stream << '\n';
01032         }
01033     }
01034 
01035 
01036   // If there are *any* variables at all in the system (including
01037   // p level, or arbitrary cell-based data)
01038   // to write, the gmv file needs to contain the word "variable"
01039   // on a line by itself.
01040   bool write_variable = false;
01041 
01042   // 1.) p-levels
01043   if (this->p_levels() && mesh_max_p_level)
01044     write_variable = true;
01045 
01046   // 2.) solution data
01047   if ((solution_names != NULL) && (v != NULL))
01048     write_variable = true;
01049 
01050   // 3.) cell-centered data
01051   if ( !(this->_cell_centered_data.empty()) )
01052     write_variable = true;
01053 
01054   if (write_variable)
01055     out_stream << "variable\n";
01056 
01057 
01058   // optionally write the p-level information
01059   if (this->p_levels() && mesh_max_p_level)
01060     {
01061       out_stream << "p_level 0\n";
01062 
01063       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01064       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01065 
01066       for ( ; it != end; ++it)
01067         if (this->subdivide_second_order())
01068           for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01069             out_stream << (*it)->p_level() << " ";
01070         else
01071           out_stream << (*it)->p_level() << " ";
01072       out_stream << "\n\n";
01073     }
01074 
01075 
01076 
01077 
01078   // optionally write cell-centered data
01079   if ( !(this->_cell_centered_data.empty()) )
01080     {
01081       std::map<std::string, const std::vector<Real>* >::iterator       it  = this->_cell_centered_data.begin();
01082       const std::map<std::string, const std::vector<Real>* >::iterator end = this->_cell_centered_data.end();
01083 
01084       for (; it != end; ++it)
01085         {
01086           // write out the variable name, followed by a zero.
01087           out_stream << (*it).first << " 0\n";
01088 
01089           const std::vector<Real>* the_array = (*it).second;
01090 
01091           // Loop over active elements, write out cell data.  If second-order cells
01092           // are split into sub-elements, the sub-elements inherit their parent's
01093           // cell-centered data.
01094           MeshBase::const_element_iterator       elem_it  = mesh.active_elements_begin();
01095           const MeshBase::const_element_iterator elem_end = mesh.active_elements_end();
01096 
01097           for (; elem_it != elem_end; ++elem_it)
01098             {
01099               const Elem* e = *elem_it;
01100 
01101               // Use the element's ID to find the value...
01102               libmesh_assert_less (e->id(), the_array->size());
01103               const Real the_value = the_array->operator[](e->id());
01104 
01105               if (this->subdivide_second_order())
01106                 for (unsigned int se=0; se < e->n_sub_elem(); se++)
01107                   out_stream << the_value << " ";
01108               else
01109                 out_stream << the_value << " ";
01110             }
01111 
01112           out_stream << "\n\n";
01113         }
01114     }
01115 
01116 
01117 
01118 
01119   // optionally write the data
01120   if ((solution_names != NULL) &&
01121       (v != NULL))
01122     {
01123       const unsigned int n_vars =
01124         libmesh_cast_int<unsigned int>(solution_names->size());
01125 
01126       if (!(v->size() == mesh.n_nodes()*n_vars))
01127         libMesh::err << "ERROR: v->size()=" << v->size()
01128                       << ", mesh.n_nodes()=" << mesh.n_nodes()
01129                       << ", n_vars=" << n_vars
01130                       << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
01131                       << std::endl;
01132 
01133       libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
01134 
01135       for (unsigned int c=0; c<n_vars; c++)
01136         {
01137 
01138 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
01139 
01140           // in case of complex data, write _tree_ data sets
01141           // for each component
01142 
01143           // this is the real part
01144           out_stream << "r_" << (*solution_names)[c] << " 1\n";
01145 
01146           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01147             out_stream << (*v)[n*n_vars + c].real() << " ";
01148 
01149           out_stream << '\n' << '\n';
01150 
01151 
01152           // this is the imaginary part
01153           out_stream << "i_" << (*solution_names)[c] << " 1\n";
01154 
01155           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01156             out_stream << (*v)[n*n_vars + c].imag() << " ";
01157 
01158           out_stream << '\n' << '\n';
01159 
01160           // this is the magnitude
01161           out_stream << "a_" << (*solution_names)[c] << " 1\n";
01162           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01163             out_stream << std::abs((*v)[n*n_vars + c]) << " ";
01164 
01165           out_stream << '\n' << '\n';
01166 
01167 #else
01168 
01169           out_stream << (*solution_names)[c] << " 1\n";
01170 
01171           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01172             out_stream << (*v)[n*n_vars + c] << " ";
01173 
01174           out_stream << '\n' << '\n';
01175 
01176 #endif
01177         }
01178 
01179     }
01180 
01181   // If we wrote any variables, we have to close the variable section now
01182   if (write_variable)
01183     out_stream << "endvars\n";
01184 
01185 
01186   // end of the file
01187   out_stream << "\nendgmv\n";
01188 }

void libMesh::GMVIO::write_binary ( const std::string &  fname,
const std::vector< Number > *  vec = NULL,
const std::vector< std::string > *  solution_names = NULL 
) [private]

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

Definition at line 1196 of file gmv_io.C.

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), end, libMesh::err, std::max(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::n_processors(), n_vars, libMesh::out, p_levels(), partitioning(), libMesh::MeshBase::point(), libMesh::processor_id(), libMeshEnums::TECPLOT, and write_subdomain_id_as_material().

Referenced by write(), and write_nodal_data().

01199 {
01200   // Get a reference to the mesh
01201   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
01202 
01203   // This is a parallel_only function
01204   const dof_id_type n_active_elem = mesh.n_active_elem();
01205 
01206   if (libMesh::processor_id() != 0)
01207     return;
01208 
01209   std::ofstream out_stream (fname.c_str());
01210 
01211   libmesh_assert (out_stream.good());
01212 
01213   unsigned int mesh_max_p_level = 0;
01214 
01215   char buf[80];
01216 
01217   // Begin interfacing with the GMV data file
01218   {
01219     // write the nodes
01220     std::strcpy(buf, "gmvinput");
01221     out_stream.write(buf, std::strlen(buf));
01222 
01223     std::strcpy(buf, "ieeei4r4");
01224     out_stream.write(buf, std::strlen(buf));
01225   }
01226 
01227 
01228 
01229   // write the nodes
01230   {
01231     std::strcpy(buf, "nodes   ");
01232     out_stream.write(buf, std::strlen(buf));
01233 
01234     unsigned int tempint = mesh.n_nodes();
01235 
01236     std::memcpy(buf, &tempint, sizeof(unsigned int));
01237 
01238     out_stream.write(buf, sizeof(unsigned int));
01239 
01240     // write the x coordinate
01241     float *temp = new float[mesh.n_nodes()];
01242     for (unsigned int v=0; v<mesh.n_nodes(); v++)
01243       temp[v] = static_cast<float>(mesh.point(v)(0));
01244     out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01245 
01246     // write the y coordinate
01247     for (unsigned int v=0; v<mesh.n_nodes(); v++)
01248 #if LIBMESH_DIM > 1
01249       temp[v] = static_cast<float>(mesh.point(v)(1));
01250 #else
01251       temp[v] = 0.;
01252 #endif
01253     out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01254 
01255     // write the z coordinate
01256     for (unsigned int v=0; v<mesh.n_nodes(); v++)
01257 #if LIBMESH_DIM > 2
01258       temp[v] = static_cast<float>(mesh.point(v)(2));
01259 #else
01260       temp[v] = 0.;
01261 #endif
01262     out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01263 
01264     delete [] temp;
01265   }
01266 
01267 
01268   // write the connectivity
01269   {
01270     std::strcpy(buf, "cells   ");
01271     out_stream.write(buf, std::strlen(buf));
01272 
01273     unsigned int tempint = n_active_elem;
01274 
01275     std::memcpy(buf, &tempint, sizeof(unsigned int));
01276 
01277     out_stream.write(buf, sizeof(unsigned int));
01278 
01279     MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01280     const MeshBase::const_element_iterator end = mesh.active_elements_end();
01281 
01282     switch (mesh.mesh_dimension())
01283       {
01284 
01285       case 1:
01286         for ( ; it != end; ++it)
01287           {
01288             mesh_max_p_level = std::max(mesh_max_p_level,
01289                                         (*it)->p_level());
01290 
01291             for(unsigned se = 0; se < (*it)->n_sub_elem(); ++se)
01292               {
01293                 std::strcpy(buf, "line    ");
01294                 out_stream.write(buf, std::strlen(buf));
01295 
01296                 tempint = 2;
01297                 std::memcpy(buf, &tempint, sizeof(unsigned int));
01298                 out_stream.write(buf, sizeof(unsigned int));
01299 
01300                 std::vector<dof_id_type> conn;
01301                 (*it)->connectivity(se,TECPLOT,conn);
01302 
01303                 out_stream.write(reinterpret_cast<char*>(&conn[0]), sizeof(unsigned int)*tempint);
01304               }
01305           }
01306         break;
01307 
01308       case 2:
01309         for ( ; it != end; ++it)
01310           {
01311             mesh_max_p_level = std::max(mesh_max_p_level,
01312                                         (*it)->p_level());
01313 
01314             for(unsigned se = 0; se < (*it)->n_sub_elem(); ++se)
01315               {
01316                 std::strcpy(buf, "quad    ");
01317                 out_stream.write(buf, std::strlen(buf));
01318                 tempint = 4;
01319                 std::memcpy(buf, &tempint, sizeof(unsigned int));
01320                 out_stream.write(buf, sizeof(unsigned int));
01321                 std::vector<dof_id_type> conn;
01322                 (*it)->connectivity(se,TECPLOT,conn);
01323                 out_stream.write(reinterpret_cast<char*>(&conn[0]), sizeof(unsigned int)*tempint);
01324               }
01325           }
01326         break;
01327       case 3:
01328         for ( ; it != end; ++it)
01329           {
01330             mesh_max_p_level = std::max(mesh_max_p_level,
01331                                         (*it)->p_level());
01332 
01333             for(unsigned se = 0; se < (*it)->n_sub_elem(); ++se)
01334               {
01335                 std::strcpy(buf, "phex8   ");
01336                 out_stream.write(buf, std::strlen(buf));
01337                 tempint = 8;
01338                 std::memcpy(buf, &tempint, sizeof(unsigned int));
01339                 out_stream.write(buf, sizeof(unsigned int));
01340                 std::vector<dof_id_type> conn;
01341                 (*it)->connectivity(se,TECPLOT,conn);
01342                 out_stream.write(reinterpret_cast<char*>(&conn[0]), sizeof(unsigned int)*tempint);
01343               }
01344           }
01345         break;
01346       default:
01347         libmesh_error();
01348 
01349       }
01350   }
01351 
01352 
01353 
01354   // optionally write the partition information
01355   if (this->partitioning())
01356     {
01357       if (this->write_subdomain_id_as_material())
01358         {
01359           libMesh::out << "Not yet supported in GMVIO::write_binary" << std::endl;
01360           libmesh_error();
01361         }
01362       else
01363         {
01364           std::strcpy(buf, "material");
01365           out_stream.write(buf, std::strlen(buf));
01366 
01367           unsigned int tmpint = mesh.n_processors();
01368           std::memcpy(buf, &tmpint, sizeof(unsigned int));
01369           out_stream.write(buf, sizeof(unsigned int));
01370 
01371           tmpint = 0; // IDs are cell based
01372           std::memcpy(buf, &tmpint, sizeof(unsigned int));
01373           out_stream.write(buf, sizeof(unsigned int));
01374 
01375 
01376           for (unsigned int proc=0; proc<mesh.n_processors(); proc++)
01377             {
01378               std::sprintf(buf, "proc_%d", proc);
01379               out_stream.write(buf, 8);
01380             }
01381 
01382           std::vector<unsigned int> proc_id (n_active_elem);
01383 
01384           unsigned int n=0;
01385 
01386           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01387           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01388 
01389           for ( ; it != end; ++it)
01390             for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01391               proc_id[n++] = (*it)->processor_id()+1;
01392 
01393 
01394           out_stream.write(reinterpret_cast<char *>(&proc_id[0]),
01395                     sizeof(unsigned int)*proc_id.size());
01396         }
01397     }
01398 
01399   // If there are *any* variables at all in the system (including
01400   // p level, or arbitrary cell-based data)
01401   // to write, the gmv file needs to contain the word "variable"
01402   // on a line by itself.
01403   bool write_variable = false;
01404 
01405   // 1.) p-levels
01406   if (this->p_levels() && mesh_max_p_level)
01407     write_variable = true;
01408 
01409   // 2.) solution data
01410   if ((solution_names != NULL) && (vec != NULL))
01411     write_variable = true;
01412 
01413   //   // 3.) cell-centered data - unsupported
01414   //   if ( !(this->_cell_centered_data.empty()) )
01415   //     write_variable = true;
01416 
01417   if (write_variable)
01418     {
01419       std::strcpy(buf, "variable");
01420       out_stream.write(buf, std::strlen(buf));
01421     }
01422 
01423   // optionally write the partition information
01424   if (this->p_levels() && mesh_max_p_level)
01425     {
01426       unsigned int n_floats = n_active_elem;
01427       for (unsigned int i=0; i != mesh.mesh_dimension(); ++i)
01428         n_floats *= 2;
01429 
01430       float *temp = new float[n_floats];
01431 
01432       std::strcpy(buf, "p_level");
01433       out_stream.write(buf, std::strlen(buf));
01434 
01435       unsigned int tempint = 0; // p levels are cell data
01436 
01437       std::memcpy(buf, &tempint, sizeof(unsigned int));
01438       out_stream.write(buf, sizeof(unsigned int));
01439 
01440       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01441       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01442       unsigned int n=0;
01443 
01444       for (; it != end; ++it)
01445         for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01446           temp[n++] = static_cast<float>( (*it)->p_level() );
01447 
01448       out_stream.write(reinterpret_cast<char *>(temp),
01449                 sizeof(float)*n_floats);
01450 
01451       delete [] temp;
01452     }
01453 
01454 
01455    // optionally write cell-centered data
01456    if ( !(this->_cell_centered_data.empty()) )
01457      {
01458        libMesh::err << "Cell-centered data not (yet) supported in binary I/O mode!" << std::endl;
01459 
01460 //        std::map<std::string, const std::vector<Real>* >::iterator       it  = this->_cell_centered_data.begin();
01461 //        const std::map<std::string, const std::vector<Real>* >::iterator end = this->_cell_centered_data.end();
01462 
01463 //        for (; it != end; ++it)
01464 //      {
01465 //        // Write out the variable name ...
01466 //        std::strcpy(buf, (*it).first.c_str());
01467 //        out_stream.write(buf, std::strlen(buf));
01468 
01469 //        // ... followed by a zero.
01470 //        unsigned int tempint = 0; // 0 signifies cell data
01471 //        std::memcpy(buf, &tempint, sizeof(unsigned int));
01472 //        out_stream.write(buf, sizeof(unsigned int));
01473 
01474 //        // Get a pointer to the array of cell-centered data values
01475 //        const std::vector<Real>* the_array = (*it).second;
01476 
01477 //        // Since the_array might contain zeros (for inactive elements) we need to
01478 //        // make a copy of it containing just values for active elements.
01479 //        const unsigned int n_floats = n_active_elem * (1<<mesh.mesh_dimension());
01480 //        float *temp = new float[n_floats];
01481 
01482 //        MeshBase::const_element_iterator       elem_it  = mesh.active_elements_begin();
01483 //        const MeshBase::const_element_iterator elem_end = mesh.active_elements_end();
01484 //        unsigned int n=0;
01485 
01486 //        for (; elem_it != elem_end; ++elem_it)
01487 //          {
01488 //            // If there's a seg-fault, it will probably be here!
01489 //            const float the_value = static_cast<float>(the_array->operator[]((*elem_it)->id()));
01490 
01491 //            for (unsigned int se=0; se<(*elem_it)->n_sub_elem(); se++)
01492 //              temp[n++] = the_value;
01493 //          }
01494 
01495 
01496 //        // Write "the_array" directly to the file
01497 //        out_stream.write(reinterpret_cast<char *>(temp),
01498 //                  sizeof(float)*n_floats);
01499 
01500 //        delete [] temp;
01501 //      }
01502      }
01503 
01504 
01505 
01506 
01507   // optionally write the data
01508   if ((solution_names != NULL) &&
01509       (vec != NULL))
01510     {
01511       float *temp = new float[mesh.n_nodes()];
01512 
01513       const unsigned int n_vars = solution_names->size();
01514 
01515       for (unsigned int c=0; c<n_vars; c++)
01516         {
01517 
01518 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
01519           // for complex data, write three datasets
01520 
01521 
01522           // Real part
01523           std::strcpy(buf, "r_");
01524           out_stream.write(buf, 2);
01525           std::strcpy(buf, (*solution_names)[c].c_str());
01526           out_stream.write(buf, 6);
01527 
01528           unsigned int tempint = 1; // always do nodal data
01529           std::memcpy(buf, &tempint, sizeof(unsigned int));
01530           out_stream.write(buf, sizeof(unsigned int));
01531 
01532           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01533             temp[n] = static_cast<float>( (*vec)[n*n_vars + c].real() );
01534 
01535           out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01536 
01537 
01538           // imaginary part
01539           std::strcpy(buf, "i_");
01540           out_stream.write(buf, 2);
01541           std::strcpy(buf, (*solution_names)[c].c_str());
01542           out_stream.write(buf, 6);
01543 
01544           std::memcpy(buf, &tempint, sizeof(unsigned int));
01545           out_stream.write(buf, sizeof(unsigned int));
01546 
01547           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01548             temp[n] = static_cast<float>( (*vec)[n*n_vars + c].imag() );
01549 
01550           out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01551 
01552           // magnitude
01553           std::strcpy(buf, "a_");
01554           out_stream.write(buf, 2);
01555           std::strcpy(buf, (*solution_names)[c].c_str());
01556           out_stream.write(buf, 6);
01557 
01558           std::memcpy(buf, &tempint, sizeof(unsigned int));
01559           out_stream.write(buf, sizeof(unsigned int));
01560 
01561           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01562             temp[n] = static_cast<float>(std::abs((*vec)[n*n_vars + c]));
01563 
01564           out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01565 
01566 #else
01567 
01568 
01569           std::strcpy(buf, (*solution_names)[c].c_str());
01570           out_stream.write(buf, 8);
01571 
01572           unsigned int tempint = 1; // always do nodal data
01573           std::memcpy(buf, &tempint, sizeof(unsigned int));
01574           out_stream.write(buf, sizeof(unsigned int));
01575 
01576           for (unsigned int n=0; n<mesh.n_nodes(); n++)
01577             temp[n] = static_cast<float>((*vec)[n*n_vars + c]);
01578 
01579           out_stream.write(reinterpret_cast<char *>(temp), sizeof(float)*mesh.n_nodes());
01580 
01581 
01582 #endif
01583 
01584 
01585         }
01586 
01587       delete [] temp;
01588 
01589     }
01590 
01591   // If we wrote any variables, we have to close the variable section now
01592   if (write_variable)
01593     {
01594       std::strcpy(buf, "endvars ");
01595       out_stream.write(buf, std::strlen(buf));
01596     }
01597 
01598   // end the file
01599   std::strcpy(buf, "endgmv  ");
01600   out_stream.write(buf, std::strlen(buf));
01601 }

void libMesh::GMVIO::write_discontinuous_gmv ( const std::string &  name,
const EquationSystems es,
const bool  write_partitioning 
) const

Writes a GMV file with discontinuous data

Definition at line 1611 of file gmv_io.C.

References _cell_centered_data, _write_subdomain_id_as_material, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMeshEnums::EDGE4, end, libMesh::err, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMeshEnums::INFEDGE2, libMeshEnums::INFHEX16, libMeshEnums::INFHEX18, libMeshEnums::INFHEX8, libMeshEnums::INFPRISM12, libMeshEnums::INFPRISM6, libMeshEnums::INFQUAD4, libMeshEnums::INFQUAD6, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_processors(), n_vars, libMesh::out, libMeshEnums::PRISM15, libMeshEnums::PRISM18, libMeshEnums::PRISM6, libMesh::MeshBase::processor_id(), libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, and libMeshEnums::TRI6.

01614 {
01615   std::vector<std::string> solution_names;
01616   std::vector<Number>      v;
01617 
01618   // Get a reference to the mesh
01619   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
01620 
01621   es.build_variable_names  (solution_names);
01622   es.build_discontinuous_solution_vector (v);
01623 
01624   // These are parallel_only functions
01625   const unsigned int n_active_elem = mesh.n_active_elem();
01626 
01627   if (mesh.processor_id() != 0)
01628     return;
01629 
01630   std::ofstream out_stream(name.c_str());
01631 
01632   libmesh_assert (out_stream.good());
01633 
01634   // Begin interfacing with the GMV data file
01635   {
01636 
01637     // write the nodes
01638     out_stream << "gmvinput ascii" << std::endl << std::endl;
01639 
01640     // Compute the total weight
01641     {
01642       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01643       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01644 
01645       unsigned int tw=0;
01646 
01647       for ( ; it != end; ++it)
01648         tw += (*it)->n_nodes();
01649 
01650       out_stream << "nodes " << tw << std::endl;
01651     }
01652 
01653 
01654 
01655     // Write all the x values
01656     {
01657       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01658       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01659 
01660       for ( ; it != end; ++it)
01661         for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01662           out_stream << (*it)->point(n)(0) << " ";
01663 
01664       out_stream << std::endl;
01665     }
01666 
01667 
01668     // Write all the y values
01669     {
01670       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01671       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01672 
01673       for ( ; it != end; ++it)
01674         for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01675 #if LIBMESH_DIM > 1
01676           out_stream << (*it)->point(n)(1) << " ";
01677 #else
01678           out_stream << 0. << " ";
01679 #endif
01680 
01681       out_stream << std::endl;
01682     }
01683 
01684 
01685     // Write all the z values
01686     {
01687       MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01688       const MeshBase::const_element_iterator end = mesh.active_elements_end();
01689 
01690       for ( ; it != end; ++it)
01691         for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01692 #if LIBMESH_DIM > 2
01693           out_stream << (*it)->point(n)(2) << " ";
01694 #else
01695           out_stream << 0. << " ";
01696 #endif
01697 
01698       out_stream << std::endl << std::endl;
01699     }
01700   }
01701 
01702 
01703 
01704   {
01705     // write the connectivity
01706 
01707     out_stream << "cells " << n_active_elem << std::endl;
01708 
01709     MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01710     const MeshBase::const_element_iterator end = mesh.active_elements_end();
01711 
01712     unsigned int nn=1;
01713 
01714     switch (mesh.mesh_dimension())
01715       {
01716       case 1:
01717         {
01718           for ( ; it != end; ++it)
01719             for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01720               {
01721                 if (((*it)->type() == EDGE2) ||
01722                     ((*it)->type() == EDGE3) ||
01723                     ((*it)->type() == EDGE4)
01724 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01725                     || ((*it)->type() == INFEDGE2)
01726 #endif
01727                     )
01728                   {
01729                     out_stream << "line 2" << std::endl;
01730                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01731                       out_stream << nn++ << " ";
01732 
01733                   }
01734                 else
01735                   {
01736                     libmesh_error();
01737                   }
01738 
01739                 out_stream << std::endl;
01740               }
01741 
01742           break;
01743         }
01744 
01745       case 2:
01746         {
01747           for ( ; it != end; ++it)
01748             for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01749               {
01750                 if (((*it)->type() == QUAD4) ||
01751                     ((*it)->type() == QUAD8) || // Note: QUAD8 will be output as one central quad and
01752                                                 // four surrounding triangles (though they will be written
01753                                                 // to GMV as QUAD4s).
01754                     ((*it)->type() == QUAD9)
01755 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01756                     || ((*it)->type() == INFQUAD4)
01757                     || ((*it)->type() == INFQUAD6)
01758 #endif
01759                     )
01760                   {
01761                     out_stream << "quad 4" << std::endl;
01762                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01763                       out_stream << nn++ << " ";
01764 
01765                   }
01766                 else if (((*it)->type() == TRI3) ||
01767                          ((*it)->type() == TRI6))
01768                   {
01769                     out_stream << "tri 3" << std::endl;
01770                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01771                       out_stream << nn++ << " ";
01772 
01773                   }
01774                 else
01775                   {
01776                     libmesh_error();
01777                   }
01778 
01779                 out_stream << std::endl;
01780               }
01781 
01782           break;
01783         }
01784 
01785 
01786       case 3:
01787         {
01788           for ( ; it != end; ++it)
01789             for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
01790               {
01791                 if (((*it)->type() == HEX8) ||
01792                     ((*it)->type() == HEX20) ||
01793                     ((*it)->type() == HEX27)
01794 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01795                     || ((*it)->type() == INFHEX8)
01796                     || ((*it)->type() == INFHEX16)
01797                     || ((*it)->type() == INFHEX18)
01798 #endif
01799                     )
01800                   {
01801                     out_stream << "phex8 8" << std::endl;
01802                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01803                       out_stream << nn++ << " ";
01804                   }
01805                 else if (((*it)->type() == PRISM6) ||
01806                          ((*it)->type() == PRISM15) ||
01807                          ((*it)->type() == PRISM18)
01808 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
01809                          || ((*it)->type() == INFPRISM6)
01810                          || ((*it)->type() == INFPRISM12)
01811 #endif
01812                          )
01813                   {
01814                     out_stream << "pprism6 6" << std::endl;
01815                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01816                       out_stream << nn++ << " ";
01817                   }
01818                 else if (((*it)->type() == TET4) ||
01819                          ((*it)->type() == TET10))
01820                   {
01821                     out_stream << "tet 4" << std::endl;
01822                     for (unsigned int i=0; i<(*it)->n_nodes(); i++)
01823                       out_stream << nn++ << " ";
01824                   }
01825                 else
01826                   {
01827                     libmesh_error();
01828                   }
01829 
01830                 out_stream << std::endl;
01831               }
01832 
01833           break;
01834         }
01835 
01836       default:
01837         libmesh_error();
01838       }
01839 
01840     out_stream << std::endl;
01841   }
01842 
01843 
01844 
01845   // optionally write the partition information
01846   if (write_partitioning)
01847     {
01848       if (_write_subdomain_id_as_material)
01849         {
01850           libMesh::out << "Not yet supported in GMVIO::write_discontinuous_gmv" << std::endl;
01851           libmesh_error();
01852         }
01853       else
01854         {
01855           out_stream << "material "
01856               << mesh.n_processors()
01857               << " 0"<< std::endl;
01858 
01859           for (unsigned int proc=0; proc<mesh.n_processors(); proc++)
01860             out_stream << "proc_" << proc << std::endl;
01861 
01862           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01863           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01864 
01865           for ( ; it != end; ++it)
01866             out_stream << (*it)->processor_id()+1 << std::endl;
01867 
01868           out_stream << std::endl;
01869         }
01870     }
01871 
01872 
01873   // Writing cell-centered data is not yet supported in discontinuous GMV files.
01874   if ( !(this->_cell_centered_data.empty()) )
01875     {
01876       libMesh::err << "Cell-centered data not (yet) supported for discontinuous GMV files!" << std::endl;
01877     }
01878 
01879 
01880 
01881   // write the data
01882   {
01883     const unsigned int n_vars = solution_names.size();
01884 
01885     //    libmesh_assert_equal_to (v.size(), tw*n_vars);
01886 
01887     out_stream << "variable" << std::endl;
01888 
01889 
01890     for (unsigned int c=0; c<n_vars; c++)
01891       {
01892 
01893 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
01894 
01895         // in case of complex data, write _tree_ data sets
01896         // for each component
01897 
01898         // this is the real part
01899         out_stream << "r_" << solution_names[c] << " 1" << std::endl;
01900         {
01901           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01902           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01903 
01904           for ( ; it != end; ++it)
01905             for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01906               out_stream << v[(n++)*n_vars + c].real() << " ";
01907         }
01908         out_stream << std::endl << std::endl;
01909 
01910 
01911         // this is the imaginary part
01912         out_stream << "i_" << solution_names[c] << " 1" << std::endl;
01913         {
01914           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01915           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01916 
01917           for ( ; it != end; ++it)
01918             for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01919               out_stream << v[(n++)*n_vars + c].imag() << " ";
01920         }
01921         out_stream << std::endl << std::endl;
01922 
01923         // this is the magnitude
01924         out_stream << "a_" << solution_names[c] << " 1" << std::endl;
01925         {
01926           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01927           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01928 
01929           for ( ; it != end; ++it)
01930             for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01931               out_stream << std::abs(v[(n++)*n_vars + c]) << " ";
01932         }
01933         out_stream << std::endl << std::endl;
01934 
01935 #else
01936 
01937         out_stream << solution_names[c] << " 1" << std::endl;
01938         {
01939           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
01940           const MeshBase::const_element_iterator end = mesh.active_elements_end();
01941 
01942           unsigned int nn=0;
01943 
01944           for ( ; it != end; ++it)
01945             for (unsigned int n=0; n<(*it)->n_nodes(); n++)
01946               out_stream << v[(nn++)*n_vars + c] << " ";
01947         }
01948         out_stream << std::endl << std::endl;
01949 
01950 #endif
01951 
01952       }
01953 
01954     out_stream << "endvars" << std::endl;
01955   }
01956 
01957 
01958   // end of the file
01959   out_stream << std::endl << "endgmv" << std::endl;
01960 }

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().

void libMesh::GMVIO::write_nodal_data ( const std::string &  fname,
const std::vector< Number > &  soln,
const std::vector< std::string > &  names 
) [virtual]

This method implements reading a mesh from a specified file. Extension of the MeshInput::read() routine which also takes an optional EquationSystems pointer and tries to read field variables from the GMV file into the EquationSystems object. This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.

Reimplemented from libMesh::MeshOutput< MeshBase >.

Definition at line 215 of file gmv_io.C.

References binary(), write_ascii_old_impl(), and write_binary().

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

00218 {
00219   START_LOG("write_nodal_data()", "GMVIO");
00220 
00221   if (this->binary())
00222     this->write_binary (fname, &soln, &names);
00223   else
00224     this->write_ascii_old_impl  (fname, &soln, &names);
00225 
00226   STOP_LOG("write_nodal_data()", "GMVIO");
00227 }

bool & libMesh::GMVIO::write_subdomain_id_as_material (  )  [inline]

Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's. Allows you to generate exploded views on user-defined subdomains, potentially creating a pretty picture.

Definition at line 323 of file gmv_io.h.

References _write_subdomain_id_as_material.

Referenced by write_ascii_new_impl(), write_ascii_old_impl(), and write_binary().

00324 {
00325   return _write_subdomain_id_as_material;
00326 }


Member Data Documentation

bool libMesh::GMVIO::_binary [private]

Flag to write binary data.

Definition at line 218 of file gmv_io.h.

Referenced by binary().

std::map<std::string, const std::vector<Real>* > libMesh::GMVIO::_cell_centered_data [private]

Storage for arbitrary cell-centered data. Ex: You can use this to plot the effectivity index for a given cell. The map is between the string representing the variable name and a pointer to a vector containing the data.

Definition at line 252 of file gmv_io.h.

Referenced by add_cell_centered_data(), write_ascii_new_impl(), write_ascii_old_impl(), write_binary(), and write_discontinuous_gmv().

Flag to write the mesh as discontinuous patches.

Definition at line 223 of file gmv_io.h.

Referenced by discontinuous().

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().

unsigned int libMesh::GMVIO::_next_elem_id [private]

Definition at line 258 of file gmv_io.h.

Referenced by _read_one_cell(), and read().

std::map<std::string, std::vector<Number> > libMesh::GMVIO::_nodal_data [private]

Definition at line 263 of file gmv_io.h.

Referenced by _read_var(), and copy_nodal_solution().

bool libMesh::GMVIO::_p_levels [private]

Flag to write the mesh p refinement levels.

Definition at line 244 of file gmv_io.h.

Referenced by p_levels().

Flag to write the mesh partitioning.

Definition at line 228 of file gmv_io.h.

Referenced by partitioning().

Flag to subdivide second order elements

Definition at line 239 of file gmv_io.h.

Referenced by subdivide_second_order().

Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's.

Definition at line 234 of file gmv_io.h.

Referenced by write_discontinuous_gmv(), and write_subdomain_id_as_material().


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

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

Hosted By:
SourceForge.net Logo