libMesh::GmshIO Class Reference

#include <gmsh_io.h>

Inheritance diagram for libMesh::GmshIO:

List of all members.

Public Member Functions

 GmshIO (MeshBase &mesh)
 GmshIO (const MeshBase &mesh)
virtual void read (const std::string &name)
virtual void write (const std::string &name)
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
bool & binary ()
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

virtual void read_mesh (std::istream &in)
virtual void write_mesh (std::ostream &out)
void write_post (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL)

Private Attributes

bool _binary

Detailed Description

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

Author:
John W. Peterson, 2004
Martin Lüthi, 2005 (support for reading meshes and writing results)

Definition at line 51 of file gmsh_io.h.


Constructor & Destructor Documentation

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

Constructor. Takes a non-const Mesh reference which it will fill up with elements via the read() command.

Definition at line 153 of file gmsh_io.h.

00153                               :
00154   MeshInput<MeshBase>  (mesh),
00155   MeshOutput<MeshBase> (mesh),
00156   _binary (false)
00157 {}

libMesh::GmshIO::GmshIO ( 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 145 of file gmsh_io.h.

00145                                     :
00146   MeshOutput<MeshBase> (mesh),
00147   _binary        (false)
00148 {
00149 }


Member Function Documentation

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

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

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

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

bool & libMesh::GmshIO::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 160 of file gmsh_io.h.

References _binary.

Referenced by write_post().

00161 {
00162   return _binary;
00163 }

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

Returns the object as a writeable reference.

Referenced by libMesh::GMVIO::_read_materials(), libMesh::GMVIO::_read_nodes(), libMesh::GMVIO::_read_one_cell(), libMesh::AbaqusIO::assign_boundary_node_ids(), libMesh::AbaqusIO::assign_sideset_ids(), libMesh::AbaqusIO::assign_subdomain_ids(), libMesh::VTKIO::cells_to_vtk(), libMesh::GMVIO::copy_nodal_solution(), libMesh::UNVIO::element_in(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::element_out(), libMesh::UNVIO::node_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::node_out(), libMesh::VTKIO::nodes_to_vtk(), libMesh::XdrIO::read(), libMesh::VTKIO::read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), libMesh::AbaqusIO::read(), libMesh::LegacyXdrIO::read_ascii(), libMesh::AbaqusIO::read_elements(), libMesh::UNVIO::read_implementation(), libMesh::UCDIO::read_implementation(), libMesh::LegacyXdrIO::read_mesh(), 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(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::UNVIO::write_implementation(), libMesh::UCDIO::write_implementation(), libMesh::LegacyXdrIO::write_mesh(), write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::XdrIO::write_parallel(), write_post(), libMesh::XdrIO::write_serialized_bcs(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::LegacyXdrIO::write_soln().

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

Reads in a mesh in the Gmsh *.msh format from the ASCII file given by name.

Note that for this method to work (in 2d and 3d) you have to explicitly set the mesh dimension prior to calling GmshIO::read() and that Mesh::prepare_for_use() must be called after reading the mesh and before using it.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 315 of file gmsh_io.C.

References read_mesh().

00316 {
00317   std::ifstream in (name.c_str());
00318   this->read_mesh (in);
00319 }

void libMesh::GmshIO::read_mesh ( std::istream &  in  )  [private, virtual]

Implementation of the read() function. This function is called by the public interface function and implements reading the file.

Read the element block

If the element dimension is smaller than the mesh dimension, this is a boundary element and will be added to mesh.boundary_info.

Because the elements might not yet exist, the sides are put on hold until the elements are created, and inserted once reading elements is finished

add the boundary element nodes to the set of nodes

If the element yet another dimension, just read in the nodes and throw them away

If any lower dimensional elements have been found in the file, try to add them to the mesh.boundary_info as sides and nodes with the respecitve id's (called "physical" in Gmsh).

Definition at line 322 of file gmsh_io.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::MeshBase::boundary_info, libMesh::Elem::build(), libMesh::Elem::build_side(), libMesh::MeshBase::clear(), end, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), libMesh::Elem::n_sides(), libMesh::Elem::neighbor(), libMesh::MeshBase::node_ptr(), libMesh::processor_id(), libMesh::Real, libMesh::AutoPtr< Tp >::release(), libMesh::MeshBase::reserve_elem(), libMesh::MeshBase::reserve_nodes(), libMesh::DofObject::set_id(), libMesh::Elem::set_node(), side, and libMesh::Elem::subdomain_id().

Referenced by read().

00323 {
00324   // This is a serial-only process for now;
00325   // the Mesh should be read on processor 0 and
00326   // broadcast later
00327   libmesh_assert_equal_to (libMesh::processor_id(), 0);
00328 
00329   libmesh_assert(in.good());
00330 
00331   // initialize the map with element types
00332   init_eletypes();
00333 
00334   // clear any data in the mesh
00335   MeshBase& mesh = MeshInput<MeshBase>::mesh();
00336   mesh.clear();
00337 
00338   const unsigned int dim = mesh.mesh_dimension();
00339 
00340   // some variables
00341   const int  bufLen = 256;
00342   char       buf[bufLen+1];
00343   int        format=0, size=0;
00344   Real       version = 1.0;
00345 
00346   // map to hold the node numbers for translation
00347   // note the the nodes can be non-consecutive
00348   std::map<unsigned int, unsigned int> nodetrans;
00349 
00350   {
00351     while (!in.eof()) {
00352       in >> buf;
00353 
00354       if (!std::strncmp(buf,"$MeshFormat",11))
00355         {
00356           in >> version >> format >> size;
00357           if ((version != 2.0) && (version != 2.1)) {
00358             // Some notes on gmsh mesh versions:
00359             //
00360             // Mesh version 2.0 goes back as far as I know.  It's not explicitly
00361             // mentioned here: http://www.geuz.org/gmsh/doc/VERSIONS.txt
00362             //
00363             // As of gmsh-2.4.0:
00364             // bumped mesh version format to 2.1 (small change in the $PhysicalNames
00365             // section, where the group dimension is now required);
00366             // [Since we don't even parse the PhysicalNames section at the time
00367             //  of this writing, I don't think this change affects us.]
00368             libMesh::err << "Error: Wrong msh file version " << version << "\n";
00369             libmesh_error();
00370           }
00371           if(format){
00372             libMesh::err << "Error: Unknown data format for mesh\n";
00373             libmesh_error();
00374           }
00375         }
00376 
00377       // read the node block
00378       else if (!std::strncmp(buf,"$NOD",4) ||
00379                !std::strncmp(buf,"$NOE",4) ||
00380                !std::strncmp(buf,"$Nodes",6)
00381                )
00382         {
00383           unsigned int numNodes = 0;
00384           in >> numNodes;
00385           mesh.reserve_nodes (numNodes);
00386 
00387           // read in the nodal coordinates and form points.
00388           Real x, y, z;
00389           unsigned int id;
00390 
00391           // add the nodal coordinates to the mesh
00392           for (unsigned int i=0; i<numNodes; ++i)
00393             {
00394               in >> id >> x >> y >> z;
00395               mesh.add_point (Point(x, y, z), i);
00396               nodetrans[id] = i;
00397             }
00398           // read the $ENDNOD delimiter
00399           in >> buf;
00400         }
00401 
00412       else if (!std::strncmp(buf,"$ELM",4) ||
00413                !std::strncmp(buf,"$Elements",9)
00414                )
00415         {
00416           unsigned int numElem = 0;
00417           std::vector< boundaryElementInfo > boundary_elem;
00418 
00419           // read how many elements are there, and reserve space in the mesh
00420           in >> numElem;
00421           mesh.reserve_elem (numElem);
00422 
00423           // read the elements
00424           unsigned int elem_id_counter = 0;
00425           for (unsigned int iel=0; iel<numElem; ++iel)
00426             {
00427               unsigned int id, type, physical, elementary,
00428                 /* partition = 1,*/ nnodes, ntags;
00429               // note - partition was assigned but never used - BSK
00430               if(version <= 1.0)
00431                 {
00432                   in >> id >> type >> physical >> elementary >> nnodes;
00433                 }
00434               else
00435                 {
00436                   in >> id >> type >> ntags;
00437                   elementary = physical = /* partition = */ 1;
00438                   for(unsigned int j = 0; j < ntags; j++)
00439                     {
00440                       int tag;
00441                       in >> tag;
00442                       if(j == 0)
00443                         physical = tag;
00444                       else if(j == 1)
00445                         elementary = tag;
00446                       // else if(j == 2)
00447                       //  partition = tag;
00448                       // ignore any other tags for now
00449                     }
00450                 }
00451 
00452               // consult the import element table which element to build
00453               const elementDefinition& eletype = eletypes_imp[type];
00454               nnodes = eletype.nnodes;
00455 
00456               // only elements that match the mesh dimension are added
00457               // if the element dimension is one less than dim, the nodes and
00458               // sides are added to the mesh.boundary_info
00459               if (eletype.dim == dim)
00460                 {
00461                   // add the elements to the mesh
00462                   Elem* elem = Elem::build(eletype.type).release();
00463                   elem->set_id(elem_id_counter);
00464                   mesh.add_elem(elem);
00465 
00466                   // different to iel, lower dimensional elems aren't added
00467                   elem_id_counter++;
00468 
00469                   // check number of nodes. We cannot do that for version 2.0
00470                   if (version <= 1.0)
00471                     {
00472                       if (elem->n_nodes() != nnodes)
00473                         {
00474                           libMesh::err << "Number of nodes for element " << id
00475                                         << " of type " << eletypes_imp[type].type
00476                                         << " (Gmsh type " << type
00477                                         << ") does not match Libmesh definition. "
00478                                         << "I expected " << elem->n_nodes()
00479                                         << " nodes, but got " << nnodes << "\n";
00480                           libmesh_error();
00481                         }
00482                     }
00483 
00484                   // add node pointers to the elements
00485                   int nod = 0;
00486                   // if there is a node translation table, use it
00487                   if (eletype.nodes.size() > 0)
00488                     for (unsigned int i=0; i<nnodes; i++)
00489                       {
00490                         in >> nod;
00491                         elem->set_node(eletype.nodes[i]) = mesh.node_ptr(nodetrans[nod]);
00492                       }
00493                   else
00494                     {
00495                       for (unsigned int i=0; i<nnodes; i++)
00496                         {
00497                           in >> nod;
00498                           elem->set_node(i) = mesh.node_ptr(nodetrans[nod]);
00499                         }
00500                     }
00501 
00502                     // Finally, set the subdomain ID to physical
00503                     elem->subdomain_id() = static_cast<subdomain_id_type>(physical);
00504                 } // if element.dim == dim
00505               // if this is a boundary
00506               else if (eletype.dim == dim-1)
00507                 {
00512                   boundaryElementInfo binfo;
00513                   std::set<dof_id_type>::iterator iter = binfo.nodes.begin();
00514                   int nod = 0;
00515                   for (unsigned int i=0; i<nnodes; i++)
00516                     {
00517                       in >> nod;
00518                       mesh.boundary_info->add_node(nodetrans[nod], physical);
00519                       binfo.nodes.insert(iter, nodetrans[nod]);
00520                     }
00521                   binfo.id = physical;
00522                   boundary_elem.push_back(binfo);
00523                 }
00528               else
00529                 {
00530                   static bool seen_high_dim_element = false;
00531                   if (!seen_high_dim_element)
00532                     {
00533                       std::cerr << "Warning: can't load an element of dimension "
00534                                 << eletype.dim << " into a mesh of dimension "
00535                                 << dim << std::endl;
00536                       seen_high_dim_element = true;
00537                     }
00538                   int nod = 0;
00539                   for (unsigned int i=0; i<nnodes; i++)
00540                     in >> nod;
00541                 }
00542             }//element loop
00543           // read the $ENDELM delimiter
00544           in >> buf;
00545 
00551           if (boundary_elem.size() > 0)
00552             {
00553               // create a index of the boundary nodes to easily locate which
00554               // element might have that boundary
00555               std::map<dof_id_type, std::vector<unsigned int> > node_index;
00556               for (unsigned int i=0; i<boundary_elem.size(); i++)
00557                 {
00558                   boundaryElementInfo binfo = boundary_elem[i];
00559                   std::set<dof_id_type>::iterator iter = binfo.nodes.begin();
00560                   for (;iter!= binfo.nodes.end(); iter++)
00561                     node_index[*iter].push_back(i);
00562                 }
00563 
00564               MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00565               const MeshBase::const_element_iterator end = mesh.active_elements_end();
00566 
00567               // iterate over all elements and see which boundary element has
00568               // the same set of nodes as on of the boundary elements previously read
00569               for ( ; it != end; ++it)
00570                 {
00571                   const Elem* elem = *it;
00572                   for (unsigned int s=0; s<elem->n_sides(); s++)
00573                     if (elem->neighbor(s) == NULL)
00574                       {
00575                         AutoPtr<Elem> side (elem->build_side(s));
00576                         std::set<dof_id_type> side_nodes;
00577                         std::set<dof_id_type>::iterator iter = side_nodes.begin();
00578 
00579                         // make a set with all nodes from this side
00580                         // this allows for easy comparison
00581                         for (unsigned int ns=0; ns<side->n_nodes(); ns++)
00582                           side_nodes.insert(iter, side->node(ns));
00583 
00584                         // See whether one of the side node occurs in the list
00585                         // of tagged nodes. If we would loop over all side
00586                         // nodes, we would just get multiple hits, so taking
00587                         // node 0 is enough to do the job
00588                         dof_id_type sn = side->node(0);
00589                         if (node_index.count(sn) > 0)
00590                           {
00591                             // Loop over all tagged ("physical") "sides" which
00592                             // contain the node sn (typically just 1 to
00593                             // three). For each of these the set of nodes is
00594                             // compared to the current element's side nodes
00595                             for (std::size_t n=0; n<node_index[sn].size(); n++)
00596                               {
00597                                 unsigned int bidx = node_index[sn][n];
00598                                 if (boundary_elem[bidx].nodes == side_nodes)
00599                                   mesh.boundary_info->add_side(elem, s, boundary_elem[bidx].id);
00600                               }
00601                           }
00602                       } // if elem->neighbor(s) == NULL
00603                 } // element loop
00604             } // if boundary_elem.size() > 0
00605         } // if $ELM
00606 
00607     } // while !in.eof()
00608 
00609   }
00610 
00611 }

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

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

This method implements writing a mesh to a specified file in the Gmsh *.msh format.

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 615 of file gmsh_io.C.

References libMesh::processor_id(), and write_mesh().

00616 {
00617   if (libMesh::processor_id() == 0)
00618     {
00619       // Open the output file stream
00620       std::ofstream out_stream (name.c_str());
00621 
00622       // Make sure it opened correctly
00623       if (!out_stream.good())
00624         libmesh_file_error(name.c_str());
00625 
00626       this->write_mesh (out_stream);
00627     }
00628 }

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::GmshIO::write_mesh ( std::ostream &  out  )  [private, virtual]

This method implements writing a mesh to a specified file. This will write an ASCII *.msh file.

Definition at line 645 of file gmsh_io.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), end, libMesh::DofObject::id(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_active_elem(), libMesh::Elem::n_nodes(), libMesh::MeshBase::n_nodes(), libMesh::Elem::node(), libMesh::MeshBase::node(), libMesh::DofObject::processor_id(), libMesh::Real, libMesh::Elem::subdomain_id(), and libMesh::Elem::type().

Referenced by write().

00646 {
00647   // Be sure that the stream is valid.
00648   libmesh_assert (out_stream.good());
00649 
00650   // initialize the map with element types
00651   init_eletypes();
00652 
00653   // Get a const reference to the mesh
00654   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
00655 
00656   // Note: we are using version 2.0 of the gmsh output format.
00657 
00658   {
00659     // Write the file header.
00660     out_stream << "$MeshFormat\n";
00661     out_stream << "2.0 0 " << sizeof(Real) << '\n';
00662     out_stream << "$EndMeshFormat\n";
00663   }
00664 
00665   {
00666     // write the nodes in (n x y z) format
00667     out_stream << "$Nodes\n";
00668     out_stream << mesh.n_nodes() << '\n';
00669 
00670     for (unsigned int v=0; v<mesh.n_nodes(); v++)
00671       out_stream << mesh.node(v).id()+1 << " "
00672           << mesh.node(v)(0) << " "
00673           << mesh.node(v)(1) << " "
00674           << mesh.node(v)(2) << '\n';
00675     out_stream << "$EndNodes\n";
00676   }
00677 
00678   {
00679     // write the connectivity
00680     out_stream << "$Elements\n";
00681     out_stream << mesh.n_active_elem() << '\n';
00682 
00683     MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00684     const MeshBase::const_element_iterator end = mesh.active_elements_end();
00685 
00686     // loop over the elements
00687     for ( ; it != end; ++it)
00688       {
00689         const Elem* elem = *it;
00690 
00691         // Make sure we have a valid entry for
00692         // the current element type.
00693         libmesh_assert (eletypes_exp.count(elem->type()));
00694 
00695         // consult the export element table
00696         const elementDefinition& eletype = eletypes_exp[elem->type()];
00697 
00698         // The element mapper better not require any more nodes
00699         // than are present in the current element!
00700         libmesh_assert_less_equal (eletype.nodes.size(), elem->n_nodes());
00701 
00702         // elements ids are 1 based in Gmsh
00703         out_stream << elem->id()+1 << " ";
00704 
00705         // element type
00706         out_stream << eletype.exptype;
00707 
00708         // write the number of tags and
00709         // tag1 (physical entity), tag2 (geometric entity), and tag3 (partition entity)
00710         out_stream << " 3 "
00711             << static_cast<unsigned int>(elem->subdomain_id())
00712             << " 1 "
00713             << (elem->processor_id()+1) << " ";
00714 
00715         // if there is a node translation table, use it
00716         if (eletype.nodes.size() > 0)
00717           for (unsigned int i=0; i < elem->n_nodes(); i++)
00718             out_stream << elem->node(eletype.nodes[i])+1 << " "; // gmsh is 1-based
00719         // otherwise keep the same node order
00720         else
00721           for (unsigned int i=0; i < elem->n_nodes(); i++)
00722             out_stream << elem->node(i)+1 << " ";                  // gmsh is 1-based
00723         out_stream << "\n";
00724       } // element loop
00725     out_stream << "$EndElements\n";
00726   }
00727 }

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

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 631 of file gmsh_io.C.

References libMesh::processor_id(), and write_post().

00634 {
00635   START_LOG("write_nodal_data()", "GmshIO");
00636 
00637   //this->_binary = true;
00638   if (libMesh::processor_id() == 0)
00639     this->write_post  (fname, &soln, &names);
00640 
00641   STOP_LOG("write_nodal_data()", "GmshIO");
00642 }

void libMesh::GmshIO::write_post ( 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 or binary *.pos file, depending on the binary flag.

Definition at line 730 of file gmsh_io.C.

References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), binary(), libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMeshEnums::EDGE4, end, libMesh::err, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMesh::libmesh_real(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_nodes(), n_vars, libMesh::Elem::n_vertices(), libMesh::Elem::node(), libMesh::out, libMesh::Elem::point(), libMeshEnums::PRISM15, libMeshEnums::PRISM18, libMeshEnums::PRISM6, libMesh::processor_id(), libMeshEnums::PYRAMID5, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, and libMeshEnums::TRI6.

Referenced by write_nodal_data().

00733 {
00734 
00735   // Should only do this on processor 0!
00736   libmesh_assert_equal_to (libMesh::processor_id(), 0);
00737 
00738   // Create an output stream
00739   std::ofstream out_stream(fname.c_str());
00740 
00741   // Make sure it opened correctly
00742   if (!out_stream.good())
00743     libmesh_file_error(fname.c_str());
00744 
00745   // initialize the map with element types
00746   init_eletypes();
00747 
00748   // create a character buffer
00749   char buf[80];
00750 
00751   // Get a constant reference to the mesh.
00752   const MeshBase& mesh = MeshOutput<MeshBase>::mesh();
00753 
00754   //  write the data
00755   if ((solution_names != NULL) && (v != NULL))
00756     {
00757       const unsigned int n_vars =
00758         libmesh_cast_int<unsigned int>(solution_names->size());
00759 
00760       if (!(v->size() == mesh.n_nodes()*n_vars))
00761         libMesh::err << "ERROR: v->size()=" << v->size()
00762                       << ", mesh.n_nodes()=" << mesh.n_nodes()
00763                       << ", n_vars=" << n_vars
00764                       << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
00765                       << "\n";
00766 
00767       libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
00768 
00769       // write the header
00770       out_stream << "$PostFormat\n";
00771       if (this->binary())
00772         out_stream << "1.2 1 " << sizeof(double) << "\n";
00773       else
00774         out_stream << "1.2 0 " << sizeof(double) << "\n";
00775       out_stream << "$EndPostFormat\n";
00776 
00777       // Loop over the elements to see how much of each type there are
00778       unsigned int n_points=0, n_lines=0, n_triangles=0, n_quadrangles=0,
00779         n_tetrahedra=0, n_hexahedra=0, n_prisms=0, n_pyramids=0;
00780       unsigned int n_scalar=0, n_vector=0, n_tensor=0;
00781       unsigned int nb_text2d=0, nb_text2d_chars=0, nb_text3d=0, nb_text3d_chars=0;
00782 
00783       {
00784         MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00785         const MeshBase::const_element_iterator end = mesh.active_elements_end();
00786 
00787 
00788         for ( ; it != end; ++it)
00789           {
00790             const ElemType elemtype = (*it)->type();
00791 
00792             switch (elemtype)
00793               {
00794               case EDGE2:
00795               case EDGE3:
00796               case EDGE4:
00797                 {
00798                   n_lines += 1;
00799                   break;
00800                 }
00801               case TRI3:
00802               case TRI6:
00803                 {
00804                   n_triangles += 1;
00805                   break;
00806                 }
00807               case QUAD4:
00808               case QUAD8:
00809               case QUAD9:
00810                 {
00811                   n_quadrangles += 1;
00812                   break;
00813                 }
00814               case TET4:
00815               case TET10:
00816                 {
00817                   n_tetrahedra += 1;
00818                   break;
00819                 }
00820               case HEX8:
00821               case HEX20:
00822               case HEX27:
00823                 {
00824                   n_hexahedra += 1;
00825                   break;
00826                 }
00827               case PRISM6:
00828               case PRISM15:
00829               case PRISM18:
00830                 {
00831                   n_prisms += 1;
00832                   break;
00833                 }
00834               case PYRAMID5:
00835                 {
00836                   n_pyramids += 1;
00837                   break;
00838                 }
00839               default:
00840                 {
00841                   libMesh::err << "ERROR: Not existant element type "
00842                                 << (*it)->type() << std::endl;
00843                   libmesh_error();
00844                 }
00845               }
00846           }
00847       }
00848 
00849       // create a view for each variable
00850       for (unsigned int ivar=0; ivar < n_vars; ivar++)
00851         {
00852           std::string varname = (*solution_names)[ivar];
00853 
00854           // at the moment, we just write out scalar quantities
00855           // later this should be made configurable through
00856           // options to the writer class
00857           n_scalar = 1;
00858 
00859           // write the variable as a view, and the number of time steps
00860           out_stream << "$View\n" << varname << " " << 1 << "\n";
00861 
00862           // write how many of each geometry type are written
00863           out_stream << n_points * n_scalar << " "
00864               << n_points * n_vector << " "
00865               << n_points * n_tensor << " "
00866               << n_lines * n_scalar << " "
00867               << n_lines * n_vector << " "
00868               << n_lines * n_tensor << " "
00869               << n_triangles * n_scalar << " "
00870               << n_triangles * n_vector << " "
00871               << n_triangles * n_tensor << " "
00872               << n_quadrangles * n_scalar << " "
00873               << n_quadrangles * n_vector << " "
00874               << n_quadrangles * n_tensor << " "
00875               << n_tetrahedra * n_scalar << " "
00876               << n_tetrahedra * n_vector << " "
00877               << n_tetrahedra * n_tensor << " "
00878               << n_hexahedra * n_scalar << " "
00879               << n_hexahedra * n_vector << " "
00880               << n_hexahedra * n_tensor << " "
00881               << n_prisms * n_scalar << " "
00882               << n_prisms * n_vector << " "
00883               << n_prisms * n_tensor << " "
00884               << n_pyramids * n_scalar << " "
00885               << n_pyramids * n_vector << " "
00886               << n_pyramids * n_tensor << " "
00887               << nb_text2d << " "
00888               << nb_text2d_chars << " "
00889               << nb_text3d << " "
00890               << nb_text3d_chars << "\n";
00891 
00892           // if binary, write a marker to identify the endianness of the file
00893           if (this->binary())
00894             {
00895               const int one = 1;
00896               std::memcpy(buf, &one, sizeof(int));
00897               out_stream.write(buf, sizeof(int));
00898             }
00899 
00900           // the time steps (there is just 1 at the moment)
00901           if (this->binary())
00902             {
00903               double one = 1;
00904               std::memcpy(buf, &one, sizeof(double));
00905               out_stream.write(buf, sizeof(double));
00906             }
00907           else
00908             out_stream << "1\n";
00909 
00910           // Loop over the elements and write out the data
00911           MeshBase::const_element_iterator       it  = mesh.active_elements_begin();
00912           const MeshBase::const_element_iterator end = mesh.active_elements_end();
00913 
00914           for ( ; it != end; ++it)
00915             {
00916               const Elem* elem = *it;
00917 
00918               // this is quite crappy, but I did not invent that file format!
00919               for (unsigned int d=0; d<3; d++)  // loop over the dimensions
00920                 {
00921                   for (unsigned int n=0; n < elem->n_vertices(); n++)   // loop over vertices
00922                     {
00923                       const Point vertex = elem->point(n);
00924                       if (this->binary())
00925                         {
00926                           double tmp = vertex(d);
00927                           std::memcpy(buf, &tmp, sizeof(double));
00928                           out_stream.write(reinterpret_cast<char *>(buf), sizeof(double));
00929                         }
00930                       else
00931                         out_stream << vertex(d) << " ";
00932                     }
00933                   if (!this->binary())
00934                     out_stream << "\n";
00935                 }
00936 
00937               // now finally write out the data
00938               for (unsigned int i=0; i < elem->n_vertices(); i++)   // loop over vertices
00939                 if (this->binary())
00940                   {
00941 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00942                     libMesh::out << "WARNING: Gmsh::write_post does not fully support "
00943                                   << "complex numbers. Will only write the real part of "
00944                                   << "variable " << varname << std::endl;
00945 #endif
00946                     double tmp = libmesh_real((*v)[elem->node(i)*n_vars + ivar]);
00947                     std::memcpy(buf, &tmp, sizeof(double));
00948                     out_stream.write(reinterpret_cast<char *>(buf), sizeof(double));
00949                   }
00950                 else
00951                   {
00952 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00953                     libMesh::out << "WARNING: Gmsh::write_post does not fully support "
00954                                   << "complex numbers. Will only write the real part of "
00955                                   << "variable " << varname << std::endl;
00956 #endif
00957                     out_stream << libmesh_real((*v)[elem->node(i)*n_vars + ivar]) << "\n";
00958                   }
00959             }
00960           if (this->binary())
00961             out_stream << "\n";
00962           out_stream << "$EndView\n";
00963 
00964         } // end variable loop (writing the views)
00965     }
00966 
00967 }


Member Data Documentation

bool libMesh::GmshIO::_binary [private]

Flag to write binary data.

Definition at line 134 of file gmsh_io.h.

Referenced by binary().

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

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

Definition at line 126 of file mesh_output.h.

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


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