libMesh
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes | List of all members
libMesh::GMVIO Class Reference

This class implements writing meshes in the GMV format. More...

#include <gmv_io.h>

Inheritance diagram for libMesh::GMVIO:
[legend]

Public Member Functions

 GMVIO (const MeshBase &)
 Constructor. More...
 
 GMVIO (MeshBase &)
 Constructor. More...
 
virtual void write (const std::string &) override
 This method implements writing a mesh to a specified file. More...
 
virtual void read (const std::string &mesh_file) override
 This method implements reading a mesh from a specified file. More...
 
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
 This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided. More...
 
bool & binary ()
 Flag indicating whether or not to write a binary file. More...
 
bool & discontinuous ()
 Flag indicating whether or not to write the mesh as discontinuous cell patches. More...
 
bool & partitioning ()
 Flag indicating whether or not to write the partitioning information for the mesh. More...
 
bool & write_subdomain_id_as_material ()
 Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's. More...
 
bool & subdivide_second_order ()
 Flag indicating whether or not to subdivide second order elements. More...
 
bool & p_levels ()
 Flag indicating whether or not to write p level information for p refined meshes. More...
 
void write_discontinuous_gmv (const std::string &name, const EquationSystems &es, const bool write_partitioning, const std::set< std::string > *system_names=nullptr) const
 Writes a GMV file with discontinuous data. More...
 
void write_ascii_new_impl (const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
 This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. More...
 
void 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. More...
 
void 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. More...
 
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object. More...
 
virtual void write_discontinuous_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 This method implements writing a mesh with discontinuous data to a specified file where the data is taken from the EquationSystems object. More...
 
virtual void write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &)
 This method may be overridden by "parallel" output formats for writing nodal data. More...
 
virtual void write_nodal_data (const std::string &, const EquationSystems &, const std::set< std::string > *)
 This method should be overridden by "parallel" output formats for writing nodal data. More...
 
virtual void write_nodal_data_discontinuous (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
 This method implements writing a mesh with discontinuous data to a specified file where the nodal data and variables names are provided. More...
 
unsigned intascii_precision ()
 Return/set the precision to use when writing ASCII files. More...
 

Protected Member Functions

MeshBasemesh ()
 
void set_n_partitions (unsigned int n_parts)
 Sets the number of partitions in the mesh. More...
 
void skip_comment_lines (std::istream &in, const char comment_start)
 Reads input from in, skipping all the lines that start with the character comment_start. More...
 
const MeshBasemesh () const
 

Protected Attributes

std::vector< bool > elems_of_dimension
 A vector of bools describing what dimension elements have been encountered when reading a mesh. More...
 
const bool _is_parallel_format
 Flag specifying whether this format is parallel-capable. More...
 
const bool _serial_only_needed_on_proc_0
 Flag specifying whether this format can be written by only serializing the mesh to processor zero. More...
 

Private Member Functions

void write_ascii_old_impl (const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
 This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. More...
 
void write_binary (const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
 This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. More...
 
void _read_nodes ()
 Helper functions for reading nodes/cells from a GMV file. More...
 
void _read_one_cell ()
 
ElemType gmv_elem_to_libmesh_elem (std::string elemname)
 
void _read_materials ()
 
void _read_var ()
 

Static Private Member Functions

static std::map< std::string, ElemTypebuild_reading_element_map ()
 Static function used to build the _reading_element_map. More...
 

Private Attributes

bool _binary
 Flag to write binary data. More...
 
bool _discontinuous
 Flag to write the mesh as discontinuous patches. More...
 
bool _partitioning
 Flag to write the mesh partitioning. More...
 
bool _write_subdomain_id_as_material
 Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's. More...
 
bool _subdivide_second_order
 Flag to subdivide second order elements. More...
 
bool _p_levels
 Flag to write the mesh p refinement levels. More...
 
std::map< std::string, const std::vector< Real > * > _cell_centered_data
 Storage for arbitrary cell-centered data. More...
 
unsigned int _next_elem_id
 
std::map< std::string, std::vector< Number > > _nodal_data
 
MeshBase_obj
 A pointer to a non-const object object. More...
 
const bool _is_parallel_format
 Flag specifying whether this format is parallel-capable. More...
 
unsigned int _ascii_precision
 Precision to use when writing ASCII files. More...
 

Static Private Attributes

static std::map< std::string, ElemType_reading_element_map = GMVIO::build_reading_element_map()
 Static map from string -> ElementType for use during reading. More...
 

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 http://www.generalmeshviewer.com

Author
Benjamin S. Kirk
Date
2004

Definition at line 54 of file gmv_io.h.

Constructor & Destructor Documentation

◆ GMVIO() [1/2]

libMesh::GMVIO::GMVIO ( const MeshBase mesh)
explicit

Constructor.

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

Definition at line 242 of file gmv_io.C.

242  :
244  _binary (false),
245  _discontinuous (false),
246  _partitioning (true),
249  _p_levels (true),
250  _next_elem_id (0)
251 {
252 }

◆ GMVIO() [2/2]

libMesh::GMVIO::GMVIO ( MeshBase mesh)
explicit

Constructor.

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

Definition at line 256 of file gmv_io.C.

256  :
259  _binary (false),
260  _discontinuous (false),
261  _partitioning (true),
264  _p_levels (true),
265  _next_elem_id (0)
266 {
267 }

Member Function Documentation

◆ _read_materials()

void libMesh::GMVIO::_read_materials ( )
private

Definition at line 2042 of file gmv_io.C.

2043 {
2044 #ifdef LIBMESH_HAVE_GMV
2045 
2046  // LibMesh assigns materials on a per-cell basis
2047  libmesh_assert_equal_to (GMVLib::gmv_data.datatype, CELL);
2048 
2049  // Material names: LibMesh has no use for these currently...
2050  // for (int i = 0; i < GMVLib::gmv_data.num; i++)
2051  // {
2052  // // Build a 32-char string from the appropriate entries
2053  // std::string mat_string(&GMVLib::gmv_data.chardata1[i*33], 32);
2054  // }
2055 
2056  for (int i = 0; i < GMVLib::gmv_data.nlongdata1; i++)
2057  MeshInput<MeshBase>::mesh().elem_ref(i).processor_id() =
2058  cast_int<processor_id_type>(GMVLib::gmv_data.longdata1[i]-1);
2059 
2060 #endif
2061 }

Referenced by read().

◆ _read_nodes()

void libMesh::GMVIO::_read_nodes ( )
private

Helper functions for reading nodes/cells from a GMV file.

Definition at line 2066 of file gmv_io.C.

2067 {
2068 #ifdef LIBMESH_HAVE_GMV
2069  // LibMesh writes UNSTRUCT=100 node data
2070  libmesh_assert_equal_to (GMVLib::gmv_data.datatype, UNSTRUCT);
2071 
2072  // The nodal data is stored in gmv_data.doubledata{1,2,3}
2073  // and is nnodes long
2074  for (int i = 0; i < GMVLib::gmv_data.num; i++)
2075  {
2076  // Add the point to the Mesh
2077  MeshInput<MeshBase>::mesh().add_point(Point(GMVLib::gmv_data.doubledata1[i],
2078  GMVLib::gmv_data.doubledata2[i],
2079  GMVLib::gmv_data.doubledata3[i]), i);
2080  }
2081 #endif
2082 }

References libMesh::MeshInput< MT >::mesh().

Referenced by read().

◆ _read_one_cell()

void libMesh::GMVIO::_read_one_cell ( )
private

Definition at line 2085 of file gmv_io.C.

2086 {
2087 #ifdef LIBMESH_HAVE_GMV
2088  // This is either a REGULAR=111 cell or
2089  // the ENDKEYWORD=207 of the cells
2090 #ifndef NDEBUG
2091  bool recognized =
2092  (GMVLib::gmv_data.datatype==REGULAR) ||
2093  (GMVLib::gmv_data.datatype==ENDKEYWORD);
2094 #endif
2095  libmesh_assert (recognized);
2096 
2098 
2099  if (GMVLib::gmv_data.datatype == REGULAR)
2100  {
2101  // We need a mapping from GMV element types to LibMesh
2102  // ElemTypes. Basically the reverse of the eletypes
2103  // std::map above.
2104  //
2105  // FIXME: Since Quad9's apparently don't exist for GMV, and since
2106  // In general we write linear sub-elements to GMV files, we need
2107  // to be careful to read back in exactly what we wrote out...
2108  //
2109  // The gmv_data.name1 field is padded with blank characters when
2110  // it is read in by GMV, so the function that converts it to the
2111  // libmesh element type needs to take that into account.
2112  ElemType type = this->gmv_elem_to_libmesh_elem(GMVLib::gmv_data.name1);
2113 
2114  Elem * elem = Elem::build(type).release();
2115  elem->set_id(_next_elem_id++);
2116 
2117  // Get the ElementDefinition object for this element type
2118  const ElementDefinition & eledef = eletypes[type];
2119 
2120  // Print out the connectivity information for
2121  // this cell.
2122  for (int i=0; i<GMVLib::gmv_data.num2; i++)
2123  {
2124  // Map index i to GMV's numbering scheme
2125  unsigned mapped_i = eledef.node_map[i];
2126 
2127  // Note: Node numbers (as stored in libmesh) are 1-based
2128  elem->set_node(i) = mesh.node_ptr
2129  (cast_int<dof_id_type>(GMVLib::gmv_data.longdata1[mapped_i]-1));
2130  }
2131 
2132  elems_of_dimension[elem->dim()] = true;
2133 
2134  // Add the newly-created element to the mesh
2135  mesh.add_elem(elem);
2136  }
2137 
2138 
2139  if (GMVLib::gmv_data.datatype == ENDKEYWORD)
2140  {
2141  // There isn't a cell to read, so we just return
2142  return;
2143  }
2144 
2145 #endif
2146 }

References _next_elem_id, libMesh::MeshBase::add_elem(), libMesh::Elem::build(), libMesh::Elem::dim(), libMesh::MeshInput< MeshBase >::elems_of_dimension, gmv_elem_to_libmesh_elem(), libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshInput< MT >::mesh(), libMesh::MeshBase::node_ptr(), libMesh::DofObject::set_id(), and libMesh::Elem::set_node().

Referenced by read().

◆ _read_var()

void libMesh::GMVIO::_read_var ( )
private

Definition at line 2030 of file gmv_io.C.

2031 {
2032 #ifdef LIBMESH_HAVE_GMV
2033 
2034  // Copy all the variable's values into a local storage vector.
2035  _nodal_data.insert ( std::make_pair(std::string(GMVLib::gmv_data.name1),
2036  std::vector<Number>(GMVLib::gmv_data.doubledata1, GMVLib::gmv_data.doubledata1+GMVLib::gmv_data.num) ) );
2037 #endif
2038 }

References _nodal_data.

Referenced by read().

◆ add_cell_centered_data()

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 1862 of file gmv_io.C.

1864 {
1865  libmesh_assert(cell_centered_data_vals);
1866 
1867  // Make sure there are *at least* enough entries for all the active elements.
1868  // There can also be entries for inactive elements, they will be ignored.
1869  // libmesh_assert_greater_equal (cell_centered_data_vals->size(),
1870  // MeshOutput<MeshBase>::mesh().n_active_elem());
1871  this->_cell_centered_data[cell_centered_data_name] = cell_centered_data_vals;
1872 }

References _cell_centered_data, and libMesh::libmesh_assert().

◆ ascii_precision()

unsigned int & libMesh::MeshOutput< MeshBase >::ascii_precision ( )
inlineinherited

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

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

Definition at line 257 of file mesh_output.h.

258 {
259  return _ascii_precision;
260 }

◆ binary()

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 are harder to debug if anything goes wrong, since they are not human-readable.

Definition at line 103 of file gmv_io.h.

103 { return _binary; }

References _binary.

Referenced by write(), and write_nodal_data().

◆ build_reading_element_map()

std::map< std::string, ElemType > libMesh::GMVIO::build_reading_element_map ( )
staticprivate

Static function used to build the _reading_element_map.

Definition at line 209 of file gmv_io.C.

210 {
211  std::map<std::string, ElemType> ret;
212 
213  ret["line"] = EDGE2;
214  ret["tri"] = TRI3;
215  ret["tri3"] = TRI3;
216  ret["quad"] = QUAD4;
217  ret["tet"] = TET4;
218  ret["ptet4"] = TET4;
219  ret["hex"] = HEX8;
220  ret["phex8"] = HEX8;
221  ret["prism"] = PRISM6;
222  ret["pprism6"] = PRISM6;
223  ret["phex20"] = HEX20;
224  ret["phex27"] = HEX27;
225  ret["pprism15"] = PRISM15;
226  ret["ptet10"] = TET10;
227  ret["6tri"] = TRI6;
228  ret["8quad"] = QUAD8;
229  ret["3line"] = EDGE3;
230 
231  // Unsupported/Unused types
232  // ret["vface2d"]
233  // ret["vface3d"]
234  // ret["pyramid"]
235  // ret["ppyrmd5"]
236  // ret["ppyrmd13"]
237 
238  return ret;
239 }

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

◆ copy_nodal_solution()

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 2159 of file gmv_io.C.

2160 {
2161  // Check for easy return if there isn't any nodal data
2162  if (_nodal_data.empty())
2163  {
2164  libMesh::err << "Unable to copy nodal solution: No nodal "
2165  << "solution has been read in from file." << std::endl;
2166  return;
2167  }
2168 
2169  // Be sure there is at least one system
2170  libmesh_assert (es.n_systems());
2171 
2172  // Keep track of variable names which have been found and
2173  // copied already. This could be used to prevent us from
2174  // e.g. copying the same var into 2 different systems ...
2175  // but this seems unlikely. Also, it is used to tell if
2176  // any variables which were read in were not successfully
2177  // copied to the EquationSystems.
2178  std::set<std::string> vars_copied;
2179 
2180  // For each entry in the nodal data map, try to find a system
2181  // that has the same variable key name.
2182  for (auto sys : IntRange<unsigned int>(0, es.n_systems()))
2183  {
2184  // Get a generic reference to the current System
2185  System & system = es.get_system(sys);
2186 
2187  // For each var entry in the _nodal_data map, try to find
2188  // that var in the system
2189  for (const auto & pr : _nodal_data)
2190  {
2191  const std::string & var_name = pr.first;
2192 
2193  if (system.has_variable(var_name))
2194  {
2195  // Check if there are as many nodes in the mesh as there are entries
2196  // in the stored nodal data vector
2197  libmesh_assert_equal_to (pr.second.size(), MeshInput<MeshBase>::mesh().n_nodes());
2198 
2199  const unsigned int var_num = system.variable_number(var_name);
2200 
2201  // The only type of nodal data we can read in from GMV is for
2202  // linear LAGRANGE type elements.
2203  const FEType & fe_type = system.variable_type(var_num);
2204  if ((fe_type.order.get_order() != FIRST) || (fe_type.family != LAGRANGE))
2205  {
2206  libMesh::err << "Only FIRST-order LAGRANGE variables can be read from GMV files. "
2207  << "Skipping variable " << var_name << std::endl;
2208  break;
2209  }
2210 
2211 
2212  // Loop over the stored vector's entries, inserting them into
2213  // the System's solution if appropriate.
2214  for (dof_id_type i=0,
2215  sz = cast_int<dof_id_type>(pr.second.size());
2216  i != sz; ++i)
2217  {
2218  // Since this var came from a GMV file, the index i corresponds to
2219  // the (single) DOF value of the current variable for node i.
2220  const unsigned int dof_index =
2221  MeshInput<MeshBase>::mesh().node_ptr(i)->dof_number(sys, // system #
2222  var_num, // var #
2223  0); // component #, always zero for LAGRANGE
2224 
2225  // If the dof_index is local to this processor, set the value
2226  if ((dof_index >= system.solution->first_local_index()) &&
2227  (dof_index < system.solution->last_local_index()))
2228  system.solution->set (dof_index, pr.second [i]);
2229  } // end loop over my GMVIO's copy of the solution
2230 
2231  // Add the most recently copied var to the set of copied vars
2232  vars_copied.insert (var_name);
2233  } // end if (system.has_variable)
2234  } // end for loop over _nodal_data
2235 
2236  // Communicate parallel values before going to the next system.
2237  system.solution->close();
2238  system.update();
2239  } // end loop over all systems
2240 
2241 
2242 
2243  // Warn the user if any GMV variables were not successfully copied over to the EquationSystems object
2244  for (const auto & pr : _nodal_data)
2245  if (vars_copied.find(pr.first) == vars_copied.end())
2246  libMesh::err << "Warning: Variable "
2247  << pr.first
2248  << " was not copied to the EquationSystems object."
2249  << std::endl;
2250 }

References _nodal_data, libMesh::err, libMesh::FEType::family, libMesh::FIRST, libMesh::OrderWrapper::get_order(), libMesh::EquationSystems::get_system(), libMesh::System::has_variable(), libMesh::LAGRANGE, libMesh::libmesh_assert(), libMesh::MeshInput< MT >::mesh(), libMesh::EquationSystems::n_systems(), libMesh::FEType::order, libMesh::System::solution, libMesh::System::update(), libMesh::System::variable_number(), and libMesh::System::variable_type().

◆ discontinuous()

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

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

Definition at line 109 of file gmv_io.h.

109 { return _discontinuous; }

References _discontinuous.

◆ gmv_elem_to_libmesh_elem()

ElemType libMesh::GMVIO::gmv_elem_to_libmesh_elem ( std::string  elemname)
private

Definition at line 2149 of file gmv_io.C.

2150 {
2151  // Erase any whitespace padding in name coming from gmv before performing comparison.
2152  elemname.erase(std::remove_if(elemname.begin(), elemname.end(), isspace), elemname.end());
2153  return libmesh_map_find(_reading_element_map, elemname);
2154 }

References _reading_element_map.

Referenced by _read_one_cell().

◆ mesh() [1/2]

MeshBase & libMesh::MeshInput< MeshBase >::mesh ( )
inlineprotectedinherited
Returns
The object as a writable reference.

Definition at line 169 of file mesh_input.h.

170 {
171  if (_obj == nullptr)
172  libmesh_error_msg("ERROR: _obj should not be nullptr!");
173  return *_obj;
174 }

◆ mesh() [2/2]

const MeshBase & libMesh::MeshOutput< MeshBase >::mesh ( ) const
inlineprotectedinherited
Returns
The object as a read-only reference.

Definition at line 247 of file mesh_output.h.

248 {
250  return *_obj;
251 }

◆ p_levels()

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

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

Definition at line 134 of file gmv_io.h.

134 { return _p_levels; }

References _p_levels.

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

◆ partitioning()

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

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

Definition at line 115 of file gmv_io.h.

115 { return _partitioning; }

References _partitioning.

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

◆ read()

void libMesh::GMVIO::read ( const std::string &  mesh_file)
overridevirtual

This method implements reading a mesh from a specified file.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 1879 of file gmv_io.C.

1880 {
1881  // This is a serial-only process for now;
1882  // the Mesh should be read on processor 0 and
1883  // broadcast later
1884  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
1885 
1886  _next_elem_id = 0;
1887 
1888  libmesh_experimental();
1889 
1890 #ifndef LIBMESH_HAVE_GMV
1891 
1892  libmesh_error_msg("Cannot read GMV file " << name << " without the GMV API.");
1893 
1894 #else
1895  // We use the file-scope global variable eletypes for mapping nodes
1896  // from GMV to libmesh indices, so initialize that data now.
1897  init_eletypes();
1898 
1899  // Clear the mesh so we are sure to start from a pristine state.
1901  mesh.clear();
1902 
1903  // Keep track of what kinds of elements this file contains
1904  elems_of_dimension.clear();
1905  elems_of_dimension.resize(4, false);
1906 
1907  // It is apparently possible for gmv files to contain
1908  // a "fromfile" directive (?) But we currently don't make
1909  // any use of this feature in LibMesh. Nonzero return val
1910  // from any function usually means an error has occurred.
1911  int ierr = GMVLib::gmvread_open_fromfileskip(const_cast<char *>(name.c_str()));
1912  if (ierr != 0)
1913  libmesh_error_msg("GMVLib::gmvread_open_fromfileskip failed!");
1914 
1915 
1916  // Loop through file until GMVEND.
1917  int iend = 0;
1918  while (iend == 0)
1919  {
1920  GMVLib::gmvread_data();
1921 
1922  // Check for GMVEND.
1923  if (GMVLib::gmv_data.keyword == GMVEND)
1924  {
1925  iend = 1;
1926  GMVLib::gmvread_close();
1927  break;
1928  }
1929 
1930  // Check for GMVERROR.
1931  if (GMVLib::gmv_data.keyword == GMVERROR)
1932  libmesh_error_msg("Encountered GMVERROR while reading!");
1933 
1934  // Process the data.
1935  switch (GMVLib::gmv_data.keyword)
1936  {
1937  case NODES:
1938  {
1939  if (GMVLib::gmv_data.num2 == NODES)
1940  this->_read_nodes();
1941 
1942  else if (GMVLib::gmv_data.num2 == NODE_V)
1943  libmesh_error_msg("Unsupported GMV data type NODE_V!");
1944 
1945  break;
1946  }
1947 
1948  case CELLS:
1949  {
1950  // Read 1 cell at a time
1951  this->_read_one_cell();
1952  break;
1953  }
1954 
1955  case MATERIAL:
1956  {
1957  // keyword == 6
1958  // These are the materials, which we use to specify the mesh
1959  // partitioning.
1960  this->_read_materials();
1961  break;
1962  }
1963 
1964  case VARIABLE:
1965  {
1966  // keyword == 8
1967  // This is a field variable.
1968 
1969  // Check to see if we're done reading variables and break out.
1970  if (GMVLib::gmv_data.datatype == ENDKEYWORD)
1971  {
1972  break;
1973  }
1974 
1975  if (GMVLib::gmv_data.datatype == NODE)
1976  {
1977  this->_read_var();
1978  break;
1979  }
1980 
1981  else
1982  {
1983  libMesh::err << "Warning: Skipping variable: "
1984  << GMVLib::gmv_data.name1
1985  << " which is of unsupported GMV datatype "
1986  << GMVLib::gmv_data.datatype
1987  << ". Nodal field data is currently the only type currently supported."
1988  << std::endl;
1989  break;
1990  }
1991 
1992  }
1993 
1994  default:
1995  libmesh_error_msg("Encountered unknown GMV keyword " << GMVLib::gmv_data.keyword);
1996 
1997  } // end switch
1998  } // end while
1999 
2000  // Set the mesh dimension to the largest encountered for an element
2001  for (unsigned char i=0; i!=4; ++i)
2002  if (elems_of_dimension[i])
2004 
2005 #if LIBMESH_DIM < 3
2006  if (mesh.mesh_dimension() > LIBMESH_DIM)
2007  libmesh_error_msg("Cannot open dimension " \
2008  << mesh.mesh_dimension() \
2009  << " mesh file when configured without " \
2010  << mesh.mesh_dimension() \
2011  << "D support.");
2012 #endif
2013 
2014  // Done reading in the mesh, now call find_neighbors, etc.
2015  // mesh.find_neighbors();
2016 
2017  // Don't allow renumbering of nodes and elements when calling
2018  // prepare_for_use(). It is usually confusing for users when the
2019  // Mesh gets renumbered automatically, since sometimes there are
2020  // other parts of the simulation which rely on the original
2021  // ordering.
2022  mesh.allow_renumbering(false);
2024 #endif
2025 }

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

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

◆ set_n_partitions()

void libMesh::MeshInput< MeshBase >::set_n_partitions ( unsigned int  n_parts)
inlineprotectedinherited

Sets the number of partitions in the mesh.

Typically this gets done by the partitioner, but some parallel file formats begin "pre-partitioned".

Definition at line 91 of file mesh_input.h.

91 { this->mesh().set_n_partitions() = n_parts; }

◆ skip_comment_lines()

void libMesh::MeshInput< MeshBase >::skip_comment_lines ( std::istream &  in,
const char  comment_start 
)
protectedinherited

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

Definition at line 179 of file mesh_input.h.

181 {
182  char c, line[256];
183 
184  while (in.get(c), c==comment_start)
185  in.getline (line, 255);
186 
187  // put back first character of
188  // first non-comment line
189  in.putback (c);
190 }

◆ subdivide_second_order()

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

Flag indicating whether or not to subdivide second order elements.

Definition at line 128 of file gmv_io.h.

128 { return _subdivide_second_order; }

References _subdivide_second_order.

Referenced by write_ascii_new_impl(), and write_ascii_old_impl().

◆ write()

void libMesh::GMVIO::write ( const std::string &  fname)
overridevirtual

This method implements writing a mesh to a specified file.

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 271 of file gmv_io.C.

272 {
273  if (this->binary())
274  this->write_binary (fname);
275  else
276  this->write_ascii_old_impl (fname);
277 }

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

Referenced by main(), and libMesh::NameBasedIO::write().

◆ write_ascii_new_impl()

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

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 295 of file gmv_io.C.

298 {
299 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
300 
301  libMesh::err << "WARNING: GMVIO::write_ascii_new_impl() not infinite-element aware!"
302  << std::endl;
303  libmesh_here();
304 
305  // Set it to our current precision
306  this->write_ascii_old_impl (fname, v, solution_names);
307 
308 #else
309 
310  // Get a reference to the mesh
312 
313  // This is a parallel_only function
314  const unsigned int n_active_elem = mesh.n_active_elem();
315 
316  if (MeshOutput<MeshBase>::mesh().processor_id() != 0)
317  return;
318 
319  // Open the output file stream
320  std::ofstream out_stream (fname.c_str());
321 
322  out_stream << std::setprecision(this->ascii_precision());
323 
324  // Make sure it opened correctly
325  if (!out_stream.good())
326  libmesh_file_error(fname.c_str());
327 
328  unsigned int mesh_max_p_level = 0;
329 
330  // Begin interfacing with the GMV data file
331  {
332  out_stream << "gmvinput ascii\n\n";
333 
334  // write the nodes
335  out_stream << "nodes " << mesh.n_nodes() << "\n";
336  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
337  out_stream << mesh.point(n)(0) << " ";
338  out_stream << "\n";
339 
340  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
341 #if LIBMESH_DIM > 1
342  out_stream << mesh.point(n)(1) << " ";
343 #else
344  out_stream << 0. << " ";
345 #endif
346  out_stream << "\n";
347 
348  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
349 #if LIBMESH_DIM > 2
350  out_stream << mesh.point(n)(2) << " ";
351 #else
352  out_stream << 0. << " ";
353 #endif
354  out_stream << "\n\n";
355  }
356 
357  {
358  // write the connectivity
359  out_stream << "cells " << n_active_elem << "\n";
360 
361  // initialize the eletypes map (eletypes is a file-global variable)
362  init_eletypes();
363 
364  for (const auto & elem : mesh.active_element_ptr_range())
365  {
366  mesh_max_p_level = std::max(mesh_max_p_level,
367  elem->p_level());
368 
369  // Make sure we have a valid entry for
370  // the current element type.
371  libmesh_assert (eletypes.count(elem->type()));
372 
373  const ElementDefinition & ele = eletypes[elem->type()];
374 
375  // The element mapper better not require any more nodes
376  // than are present in the current element!
377  libmesh_assert_less_equal (ele.node_map.size(), elem->n_nodes());
378 
379  out_stream << ele.label << "\n";
380  for (auto i : index_range(ele.node_map))
381  out_stream << elem->node_id(ele.node_map[i])+1 << " ";
382  out_stream << "\n";
383  }
384  out_stream << "\n";
385  }
386 
387  // optionally write the partition information
388  if (this->partitioning())
389  {
390  if (this->write_subdomain_id_as_material())
391  libmesh_error_msg("Not yet supported in GMVIO::write_ascii_new_impl");
392 
393  else // write processor IDs as materials. This is the default
394  {
395  out_stream << "material "
396  << mesh.n_partitions()
397  // Note: GMV may give you errors like
398  // Error, material for cell 1 is greater than 1
399  // Error, material for cell 2 is greater than 1
400  // Error, material for cell 3 is greater than 1
401  // ... because you put the wrong number of partitions here.
402  // To ensure you write the correct number of materials, call
403  // mesh.recalculate_n_partitions() before writing out the
404  // mesh.
405  // Note: we can't call it now because the Mesh is const here and
406  // it is a non-const function.
407  << " 0\n";
408 
409  for (auto proc : IntRange<unsigned int>(0, mesh.n_partitions()))
410  out_stream << "proc_" << proc << "\n";
411 
412  // FIXME - don't we need to use an ElementDefinition here? - RHS
413  for (const auto & elem : mesh.active_element_ptr_range())
414  out_stream << elem->processor_id()+1 << "\n";
415  out_stream << "\n";
416  }
417  }
418 
419  // If there are *any* variables at all in the system (including
420  // p level, or arbitrary cell-based data)
421  // to write, the gmv file needs to contain the word "variable"
422  // on a line by itself.
423  bool write_variable = false;
424 
425  // 1.) p-levels
426  if (this->p_levels() && mesh_max_p_level)
427  write_variable = true;
428 
429  // 2.) solution data
430  if ((solution_names != nullptr) && (v != nullptr))
431  write_variable = true;
432 
433  // 3.) cell-centered data
434  if (!(this->_cell_centered_data.empty()))
435  write_variable = true;
436 
437  if (write_variable)
438  out_stream << "variable\n";
439 
440  // if ((this->p_levels() && mesh_max_p_level) ||
441  // ((solution_names != nullptr) && (v != nullptr)))
442  // out_stream << "variable\n";
443 
444  // optionally write the polynomial degree information
445  if (this->p_levels() && mesh_max_p_level)
446  {
447  out_stream << "p_level 0\n";
448 
449  for (const auto & elem : mesh.active_element_ptr_range())
450  {
451  const ElementDefinition & ele = eletypes[elem->type()];
452 
453  // The element mapper better not require any more nodes
454  // than are present in the current element!
455  libmesh_assert_less_equal (ele.node_map.size(), elem->n_nodes());
456 
457  for (std::size_t i=0, enms=ele.node_map.size(); i < enms; i++)
458  out_stream << elem->p_level() << " ";
459  }
460  out_stream << "\n\n";
461  }
462 
463 
464  // optionally write cell-centered data
465  if (!(this->_cell_centered_data.empty()))
466  {
467  for (auto & pr : this->_cell_centered_data)
468  {
469  // write out the variable name, followed by a zero.
470  out_stream << pr.first << " 0\n";
471 
472  const std::vector<Real> * the_array = pr.second;
473 
474  // Loop over active elements, write out cell data. If second-order cells
475  // are split into sub-elements, the sub-elements inherit their parent's
476  // cell-centered data.
477  for (const auto & elem : mesh.active_element_ptr_range())
478  {
479  // Use the element's ID to find the value.
480  libmesh_assert_less (elem->id(), the_array->size());
481  const Real the_value = the_array->operator[](elem->id());
482 
483  if (this->subdivide_second_order())
484  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
485  out_stream << the_value << " ";
486  else
487  out_stream << the_value << " ";
488  }
489 
490  out_stream << "\n\n";
491  }
492  }
493 
494 
495  // optionally write the data
496  if ((solution_names != nullptr) && (v != nullptr))
497  {
498  const unsigned int n_vars = solution_names->size();
499 
500  if (!(v->size() == mesh.n_nodes()*n_vars))
501  libMesh::err << "ERROR: v->size()=" << v->size()
502  << ", mesh.n_nodes()=" << mesh.n_nodes()
503  << ", n_vars=" << n_vars
504  << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
505  << "\n";
506 
507  libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
508 
509  for (unsigned int c=0; c<n_vars; c++)
510  {
511 
512 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
513 
514  // in case of complex data, write _three_ data sets
515  // for each component
516 
517  // this is the real part
518  out_stream << "r_" << (*solution_names)[c] << " 1\n";
519 
520  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
521  out_stream << (*v)[n*n_vars + c].real() << " ";
522 
523  out_stream << "\n\n";
524 
525  // this is the imaginary part
526  out_stream << "i_" << (*solution_names)[c] << " 1\n";
527 
528  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
529  out_stream << (*v)[n*n_vars + c].imag() << " ";
530 
531  out_stream << "\n\n";
532 
533  // this is the magnitude
534  out_stream << "a_" << (*solution_names)[c] << " 1\n";
535  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
536  out_stream << std::abs((*v)[n*n_vars + c]) << " ";
537 
538  out_stream << "\n\n";
539 
540 #else
541 
542  out_stream << (*solution_names)[c] << " 1\n";
543 
544  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
545  out_stream << (*v)[n*n_vars + c] << " ";
546 
547  out_stream << "\n\n";
548 
549 #endif
550  }
551 
552  }
553 
554  // If we wrote any variables, we have to close the variable section now
555  if (write_variable)
556  out_stream << "endvars\n";
557 
558 
559  // end of the file
560  out_stream << "\nendgmv\n";
561 
562 #endif
563 }

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_element_ptr_range(), libMesh::MeshOutput< MeshBase >::ascii_precision(), libMesh::err, libMesh::index_range(), libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::n_partitions(), n_vars, p_levels(), partitioning(), libMesh::MeshBase::point(), libMesh::Real, subdivide_second_order(), write_ascii_old_impl(), and write_subdomain_id_as_material().

◆ write_ascii_old_impl()

void libMesh::GMVIO::write_ascii_old_impl ( const std::string &  fname,
const std::vector< Number > *  v = nullptr,
const std::vector< std::string > *  solution_names = nullptr 
)
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.

Definition at line 570 of file gmv_io.C.

573 {
574  // Get a reference to the mesh
576 
577  // Use a MeshSerializer object to gather a parallel mesh before outputting it.
578  // Note that we cast away constness here (which is bad), but the destructor of
579  // the MeshSerializer object reparallelizes the Mesh, hopefully keeping it
580  // "logically const" outside the context of this function...
581  MeshSerializer serialize(const_cast<MeshBase &>(mesh),
583 
584  // These are parallel_only functions
585  const dof_id_type n_active_elem = mesh.n_active_elem(),
586  n_active_sub_elem = mesh.n_active_sub_elem();
587 
588  if (MeshOutput<MeshBase>::mesh().processor_id() != 0)
589  return;
590 
591  // Open the output file stream
592  std::ofstream out_stream (fname.c_str());
593 
594  // Set it to our current precision
595  out_stream << std::setprecision(this->ascii_precision());
596 
597  // Make sure it opened correctly
598  if (!out_stream.good())
599  libmesh_file_error(fname.c_str());
600 
601  // Make sure our nodes are contiguous and serialized
602  libmesh_assert_equal_to (mesh.n_nodes(), mesh.max_node_id());
603 
604  // libmesh_assert (mesh.is_serial());
605  // if (!mesh.is_serial())
606  // {
607  // if (MeshOutput<MeshBase>::mesh().processor_id() == 0)
608  // libMesh::err << "Error: GMVIO cannot yet write a DistributedMesh solution"
609  // << std::endl;
610  // return;
611  // }
612 
613  unsigned int mesh_max_p_level = 0;
614 
615  // Begin interfacing with the GMV data file
616 
617  // FIXME - if subdivide_second_order() is off,
618  // we probably should only be writing the
619  // vertex nodes - RHS
620  {
621  // write the nodes
622 
623  out_stream << "gmvinput ascii\n\n";
624  out_stream << "nodes " << mesh.n_nodes() << '\n';
625  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
626  out_stream << mesh.point(n)(0) << " ";
627 
628  out_stream << '\n';
629 
630  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
631 #if LIBMESH_DIM > 1
632  out_stream << mesh.point(n)(1) << " ";
633 #else
634  out_stream << 0. << " ";
635 #endif
636 
637  out_stream << '\n';
638 
639  for (auto n : IntRange<dof_id_type>(0, mesh.n_nodes()))
640 #if LIBMESH_DIM > 2
641  out_stream << mesh.point(n)(2) << " ";
642 #else
643  out_stream << 0. << " ";
644 #endif
645 
646  out_stream << '\n' << '\n';
647  }
648 
649 
650 
651  {
652  // write the connectivity
653 
654  out_stream << "cells ";
655  if (this->subdivide_second_order())
656  out_stream << n_active_sub_elem;
657  else
658  out_stream << n_active_elem;
659  out_stream << '\n';
660 
661  // The same temporary storage will be used for each element
662  std::vector<dof_id_type> conn;
663 
664  for (const auto & elem : mesh.active_element_ptr_range())
665  {
666  mesh_max_p_level = std::max(mesh_max_p_level,
667  elem->p_level());
668 
669  switch (elem->dim())
670  {
671  case 1:
672  {
673  if (this->subdivide_second_order())
674  for (auto se : IntRange<unsigned int>(0, elem->n_sub_elem()))
675  {
676  out_stream << "line 2\n";
677  elem->connectivity(se, TECPLOT, conn);
678  for (const auto & idx : conn)
679  out_stream << idx << " ";
680 
681  out_stream << '\n';
682  }
683  else
684  {
685  out_stream << "line 2\n";
686  if (elem->default_order() == FIRST)
687  elem->connectivity(0, TECPLOT, conn);
688  else
689  {
690  std::unique_ptr<Elem> lo_elem = Elem::build(Elem::first_order_equivalent_type(elem->type()));
691  for (auto i : lo_elem->node_index_range())
692  lo_elem->set_node(i) = elem->node_ptr(i);
693  lo_elem->connectivity(0, TECPLOT, conn);
694  }
695  for (const auto & idx : conn)
696  out_stream << idx << " ";
697 
698  out_stream << '\n';
699  }
700 
701  break;
702  }
703 
704  case 2:
705  {
706  if (this->subdivide_second_order())
707  for (auto se : IntRange<unsigned int>(0, elem->n_sub_elem()))
708  {
709  // Quad elements
710  if ((elem->type() == QUAD4) ||
711  (elem->type() == QUAD8) || // Note: QUAD8 will be output as one central quad and
712  // four surrounding triangles (though they will be written
713  // to GMV as QUAD4s).
714  (elem->type() == QUAD9)
715 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
716  || (elem->type() == INFQUAD4)
717  || (elem->type() == INFQUAD6)
718 #endif
719  )
720  {
721  out_stream << "quad 4\n";
722  elem->connectivity(se, TECPLOT, conn);
723  for (const auto & idx : conn)
724  out_stream << idx << " ";
725  }
726 
727  // Triangle elements
728  else if ((elem->type() == TRI3) ||
729  (elem->type() == TRI6))
730  {
731  out_stream << "tri 3\n";
732  elem->connectivity(se, TECPLOT, conn);
733  for (unsigned int i=0; i<3; i++)
734  out_stream << conn[i] << " ";
735  }
736  else
737  libmesh_error_msg("Unsupported element type: " << Utility::enum_to_string(elem->type()));
738  }
739  else // !this->subdivide_second_order()
740  {
741  // Quad elements
742  if ((elem->type() == QUAD4)
743 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
744  || (elem->type() == INFQUAD4)
745 #endif
746  )
747  {
748  elem->connectivity(0, TECPLOT, conn);
749  out_stream << "quad 4\n";
750  for (const auto & idx : conn)
751  out_stream << idx << " ";
752  }
753  else if ((elem->type() == QUAD8) ||
754  (elem->type() == QUAD9)
755 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
756  || (elem->type() == INFQUAD6)
757 #endif
758  )
759  {
760  std::unique_ptr<Elem> lo_elem = Elem::build(Elem::first_order_equivalent_type(elem->type()));
761  for (auto i : lo_elem->node_index_range())
762  lo_elem->set_node(i) = elem->node_ptr(i);
763  lo_elem->connectivity(0, TECPLOT, conn);
764  out_stream << "quad 4\n";
765  for (const auto & idx : conn)
766  out_stream << idx << " ";
767  }
768  else if (elem->type() == TRI3)
769  {
770  elem->connectivity(0, TECPLOT, conn);
771  out_stream << "tri 3\n";
772  for (unsigned int i=0; i<3; i++)
773  out_stream << conn[i] << " ";
774  }
775  else if (elem->type() == TRI6)
776  {
777  std::unique_ptr<Elem> lo_elem = Elem::build(Elem::first_order_equivalent_type(elem->type()));
778  for (auto i : lo_elem->node_index_range())
779  lo_elem->set_node(i) = elem->node_ptr(i);
780  lo_elem->connectivity(0, TECPLOT, conn);
781  out_stream << "tri 3\n";
782  for (unsigned int i=0; i<3; i++)
783  out_stream << conn[i] << " ";
784  }
785 
786  out_stream << '\n';
787  }
788 
789  break;
790  }
791 
792  case 3:
793  {
794  if (this->subdivide_second_order())
795  for (auto se : IntRange<unsigned int>(0, elem->n_sub_elem()))
796  {
797 
798 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
799  if ((elem->type() == HEX8) ||
800  (elem->type() == HEX27))
801  {
802  out_stream << "phex8 8\n";
803  elem->connectivity(se, TECPLOT, conn);
804  for (const auto & idx : conn)
805  out_stream << idx << " ";
806  }
807 
808  else if (elem->type() == HEX20)
809  {
810  out_stream << "phex20 20\n";
811  out_stream << elem->node_id(0)+1 << " "
812  << elem->node_id(1)+1 << " "
813  << elem->node_id(2)+1 << " "
814  << elem->node_id(3)+1 << " "
815  << elem->node_id(4)+1 << " "
816  << elem->node_id(5)+1 << " "
817  << elem->node_id(6)+1 << " "
818  << elem->node_id(7)+1 << " "
819  << elem->node_id(8)+1 << " "
820  << elem->node_id(9)+1 << " "
821  << elem->node_id(10)+1 << " "
822  << elem->node_id(11)+1 << " "
823  << elem->node_id(16)+1 << " "
824  << elem->node_id(17)+1 << " "
825  << elem->node_id(18)+1 << " "
826  << elem->node_id(19)+1 << " "
827  << elem->node_id(12)+1 << " "
828  << elem->node_id(13)+1 << " "
829  << elem->node_id(14)+1 << " "
830  << elem->node_id(15)+1 << " ";
831  }
832 
833  // According to our copy of gmvread.c, this is the
834  // mapping for the Hex27 element. Unfortunately,
835  // I tried it and Paraview does not seem to be
836  // able to read Hex27 elements. Since this is
837  // unlikely to change any time soon, we'll
838  // continue to write out Hex27 elements as 8 Hex8
839  // sub-elements.
840  //
841  // TODO:
842  // 1.) If we really wanted to use this code for
843  // something, we'd want to avoid repeating the
844  // hard-coded node ordering from the HEX20 case.
845  // These should both be able to use
846  // ElementDefinitions.
847  // 2.) You would also need to change
848  // Hex27::n_sub_elem() to return 1, not sure how
849  // much other code that would affect...
850 
851  // else if (elem->type() == HEX27)
852  // {
853  // out_stream << "phex27 27\n";
854  // out_stream << elem->node_id(0)+1 << " "
855  // << elem->node_id(1)+1 << " "
856  // << elem->node_id(2)+1 << " "
857  // << elem->node_id(3)+1 << " "
858  // << elem->node_id(4)+1 << " "
859  // << elem->node_id(5)+1 << " "
860  // << elem->node_id(6)+1 << " "
861  // << elem->node_id(7)+1 << " "
862  // << elem->node_id(8)+1 << " "
863  // << elem->node_id(9)+1 << " "
864  // << elem->node_id(10)+1 << " "
865  // << elem->node_id(11)+1 << " "
866  // << elem->node_id(16)+1 << " "
867  // << elem->node_id(17)+1 << " "
868  // << elem->node_id(18)+1 << " "
869  // << elem->node_id(19)+1 << " "
870  // << elem->node_id(12)+1 << " "
871  // << elem->node_id(13)+1 << " "
872  // << elem->node_id(14)+1 << " "
873  // << elem->node_id(15)+1 << " "
874  // // mid-face nodes
875  // << elem->node_id(21)+1 << " " // GMV front
876  // << elem->node_id(22)+1 << " " // GMV right
877  // << elem->node_id(23)+1 << " " // GMV back
878  // << elem->node_id(24)+1 << " " // GMV left
879  // << elem->node_id(20)+1 << " " // GMV bottom
880  // << elem->node_id(25)+1 << " " // GMV top
881  // // center node
882  // << elem->node_id(26)+1 << " ";
883  // }
884 
885 #else // LIBMESH_ENABLE_INFINITE_ELEMENTS
886 
887  // In case of infinite elements, HEX20
888  // should be handled just like the
889  // INFHEX16, since these connect to each other
890  if ((elem->type() == HEX8) ||
891  (elem->type() == HEX27) ||
892  (elem->type() == INFHEX8) ||
893  (elem->type() == INFHEX16) ||
894  (elem->type() == INFHEX18) ||
895  (elem->type() == HEX20))
896  {
897  out_stream << "phex8 8\n";
898  elem->connectivity(se, TECPLOT, conn);
899  for (const auto & idx : conn)
900  out_stream << idx << " ";
901  }
902 #endif
903 
904  else if ((elem->type() == TET4) ||
905  (elem->type() == TET10))
906  {
907  out_stream << "tet 4\n";
908  // Tecplot connectivity returns 8 entries for
909  // the Tet, enough to store it as a degenerate Hex.
910  // For GMV we only pick out the four relevant node
911  // indices.
912  elem->connectivity(se, TECPLOT, conn);
913  out_stream << conn[0] << " " // libmesh tet node 0
914  << conn[2] << " " // libmesh tet node 2
915  << conn[1] << " " // libmesh tet node 1
916  << conn[4] << " "; // libmesh tet node 3
917  }
918 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
919  else if ((elem->type() == PRISM6) ||
920  (elem->type() == PRISM15) ||
921  (elem->type() == PRISM18) ||
922  (elem->type() == PYRAMID5))
923 #else
924  else if ((elem->type() == PRISM6) ||
925  (elem->type() == PRISM15) ||
926  (elem->type() == PRISM18) ||
927  (elem->type() == PYRAMID5) ||
928  (elem->type() == INFPRISM6) ||
929  (elem->type() == INFPRISM12))
930 #endif
931  {
932  // Note that the prisms are treated as
933  // degenerated phex8's.
934  out_stream << "phex8 8\n";
935  elem->connectivity(se, TECPLOT, conn);
936  for (const auto & idx : conn)
937  out_stream << idx << " ";
938  }
939 
940  else
941  libmesh_error_msg("Encountered an unrecognized element " \
942  << "type: " << elem->type() \
943  << "\nPossibly a dim-1 dimensional " \
944  << "element? Aborting...");
945 
946  out_stream << '\n';
947  }
948  else // !this->subdivide_second_order()
949  {
950  std::unique_ptr<Elem> lo_elem = Elem::build(Elem::first_order_equivalent_type(elem->type()));
951  for (auto i : lo_elem->node_index_range())
952  lo_elem->set_node(i) = elem->node_ptr(i);
953  if ((lo_elem->type() == HEX8)
954 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
955  || (lo_elem->type() == HEX27)
956 #endif
957  )
958  {
959  out_stream << "phex8 8\n";
960  lo_elem->connectivity(0, TECPLOT, conn);
961  for (const auto & idx : conn)
962  out_stream << idx << " ";
963  }
964 
965  else if (lo_elem->type() == TET4)
966  {
967  out_stream << "tet 4\n";
968  lo_elem->connectivity(0, TECPLOT, conn);
969  out_stream << conn[0] << " "
970  << conn[2] << " "
971  << conn[1] << " "
972  << conn[4] << " ";
973  }
974  else if ((lo_elem->type() == PRISM6)
975 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
976  || (lo_elem->type() == INFPRISM6)
977 #endif
978  )
979  {
980  // Note that the prisms are treated as
981  // degenerated phex8's.
982  out_stream << "phex8 8\n";
983  lo_elem->connectivity(0, TECPLOT, conn);
984  for (const auto & idx : conn)
985  out_stream << idx << " ";
986  }
987 
988  else
989  libmesh_error_msg("Encountered an unrecognized element " \
990  << "type. Possibly a dim-1 dimensional " \
991  << "element? Aborting...");
992 
993  out_stream << '\n';
994  }
995 
996  break;
997  }
998 
999  default:
1000  libmesh_error_msg("Unsupported element dimension: " <<
1001  elem->dim());
1002  }
1003  }
1004 
1005  out_stream << '\n';
1006  }
1007 
1008 
1009 
1010  // optionally write the partition information
1011  if (this->partitioning())
1012  {
1013  if (this->write_subdomain_id_as_material())
1014  {
1015  // Subdomain IDs can be non-contiguous and need not
1016  // necessarily start at 0. Furthermore, since the user is
1017  // free to define subdomain_id_type to be a signed type, we
1018  // can't even assume max(subdomain_id) >= # unique subdomain ids.
1019 
1020  // We build a map<subdomain_id, unsigned> to associate to each
1021  // user-selected subdomain ID a unique, contiguous unsigned value
1022  // which we can write to file.
1023  std::map<subdomain_id_type, unsigned> sbdid_map;
1024 
1025  // Try to insert with dummy value
1026  for (const auto & elem : mesh.active_element_ptr_range())
1027  sbdid_map.insert(std::make_pair(elem->subdomain_id(), 0));
1028 
1029  // Map is created, iterate through it to set indices. They will be
1030  // used repeatedly below.
1031  {
1032  unsigned ctr=0;
1033  for (auto & pr : sbdid_map)
1034  pr.second = ctr++;
1035  }
1036 
1037  out_stream << "material "
1038  << sbdid_map.size()
1039  << " 0\n";
1040 
1041  for (auto sbdid : IntRange<std::size_t>(0, sbdid_map.size()))
1042  out_stream << "proc_" << sbdid << "\n";
1043 
1044  for (const auto & elem : mesh.active_element_ptr_range())
1045  {
1046  // Find the unique index for elem->subdomain_id(), print that to file
1047  unsigned gmv_mat_number = libmesh_map_find(sbdid_map, elem->subdomain_id());
1048 
1049  if (this->subdivide_second_order())
1050  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1051  out_stream << gmv_mat_number+1 << '\n';
1052  else
1053  out_stream << gmv_mat_number+1 << "\n";
1054  }
1055  out_stream << '\n';
1056 
1057  }
1058  else // write processor IDs as materials. This is the default
1059  {
1060  out_stream << "material "
1061  << mesh.n_partitions()
1062  << " 0"<< '\n';
1063 
1064  for (auto proc : IntRange<unsigned int>(0, mesh.n_partitions()))
1065  out_stream << "proc_" << proc << '\n';
1066 
1067  for (const auto & elem : mesh.active_element_ptr_range())
1068  if (this->subdivide_second_order())
1069  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1070  out_stream << elem->processor_id()+1 << '\n';
1071  else
1072  out_stream << elem->processor_id()+1 << '\n';
1073 
1074  out_stream << '\n';
1075  }
1076  }
1077 
1078 
1079  // If there are *any* variables at all in the system (including
1080  // p level, or arbitrary cell-based data)
1081  // to write, the gmv file needs to contain the word "variable"
1082  // on a line by itself.
1083  bool write_variable = false;
1084 
1085  // 1.) p-levels
1086  if (this->p_levels() && mesh_max_p_level)
1087  write_variable = true;
1088 
1089  // 2.) solution data
1090  if ((solution_names != nullptr) && (v != nullptr))
1091  write_variable = true;
1092 
1093  // 3.) cell-centered data
1094  if (!(this->_cell_centered_data.empty()))
1095  write_variable = true;
1096 
1097  if (write_variable)
1098  out_stream << "variable\n";
1099 
1100 
1101  // optionally write the p-level information
1102  if (this->p_levels() && mesh_max_p_level)
1103  {
1104  out_stream << "p_level 0\n";
1105 
1106  for (const auto & elem : mesh.active_element_ptr_range())
1107  if (this->subdivide_second_order())
1108  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1109  out_stream << elem->p_level() << " ";
1110  else
1111  out_stream << elem->p_level() << " ";
1112  out_stream << "\n\n";
1113  }
1114 
1115 
1116 
1117 
1118  // optionally write cell-centered data
1119  if (!(this->_cell_centered_data.empty()))
1120  {
1121  for (const auto & pr : this->_cell_centered_data)
1122  {
1123  // write out the variable name, followed by a zero.
1124  out_stream << pr.first << " 0\n";
1125 
1126  const std::vector<Real> * the_array = pr.second;
1127 
1128  // Loop over active elements, write out cell data. If second-order cells
1129  // are split into sub-elements, the sub-elements inherit their parent's
1130  // cell-centered data.
1131  for (const auto & elem : mesh.active_element_ptr_range())
1132  {
1133  // Use the element's ID to find the value...
1134  libmesh_assert_less (elem->id(), the_array->size());
1135  const Real the_value = (*the_array)[elem->id()];
1136 
1137  if (this->subdivide_second_order())
1138  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1139  out_stream << the_value << " ";
1140  else
1141  out_stream << the_value << " ";
1142  }
1143 
1144  out_stream << "\n\n";
1145  }
1146  }
1147 
1148 
1149 
1150 
1151  // optionally write the data
1152  if ((solution_names != nullptr) &&
1153  (v != nullptr))
1154  {
1155  const unsigned int n_vars =
1156  cast_int<unsigned int>(solution_names->size());
1157 
1158  if (!(v->size() == mesh.n_nodes()*n_vars))
1159  libMesh::err << "ERROR: v->size()=" << v->size()
1160  << ", mesh.n_nodes()=" << mesh.n_nodes()
1161  << ", n_vars=" << n_vars
1162  << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
1163  << std::endl;
1164 
1165  libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
1166 
1167  for (unsigned int c=0; c<n_vars; c++)
1168  {
1169 
1170 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1171 
1172  // in case of complex data, write _tree_ data sets
1173  // for each component
1174 
1175  // this is the real part
1176  out_stream << "r_" << (*solution_names)[c] << " 1\n";
1177 
1178  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1179  out_stream << (*v)[n*n_vars + c].real() << " ";
1180 
1181  out_stream << '\n' << '\n';
1182 
1183 
1184  // this is the imaginary part
1185  out_stream << "i_" << (*solution_names)[c] << " 1\n";
1186 
1187  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1188  out_stream << (*v)[n*n_vars + c].imag() << " ";
1189 
1190  out_stream << '\n' << '\n';
1191 
1192  // this is the magnitude
1193  out_stream << "a_" << (*solution_names)[c] << " 1\n";
1194  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1195  out_stream << std::abs((*v)[n*n_vars + c]) << " ";
1196 
1197  out_stream << '\n' << '\n';
1198 
1199 #else
1200 
1201  out_stream << (*solution_names)[c] << " 1\n";
1202 
1203  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1204  out_stream << (*v)[n*n_vars + c] << " ";
1205 
1206  out_stream << '\n' << '\n';
1207 
1208 #endif
1209  }
1210 
1211  }
1212 
1213  // If we wrote any variables, we have to close the variable section now
1214  if (write_variable)
1215  out_stream << "endvars\n";
1216 
1217 
1218  // end of the file
1219  out_stream << "\nendgmv\n";
1220 }

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_element_ptr_range(), libMesh::MeshOutput< MeshBase >::ascii_precision(), libMesh::Elem::build(), libMesh::Utility::enum_to_string(), libMesh::err, libMesh::FIRST, libMesh::Elem::first_order_equivalent_type(), libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::MeshTools::Generation::Private::idx(), libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::INFQUAD4, libMesh::INFQUAD6, libMesh::MeshBase::max_node_id(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_active_sub_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::n_partitions(), n_vars, p_levels(), partitioning(), libMesh::MeshBase::point(), libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::Real, subdivide_second_order(), libMesh::TECPLOT, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI6, and write_subdomain_id_as_material().

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

◆ write_binary()

void libMesh::GMVIO::write_binary ( const std::string &  fname,
const std::vector< Number > *  vec = nullptr,
const std::vector< std::string > *  solution_names = nullptr 
)
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 1228 of file gmv_io.C.

1231 {
1232  // We use the file-scope global variable eletypes for mapping nodes
1233  // from GMV to libmesh indices, so initialize that data now.
1234  init_eletypes();
1235 
1236  // Get a reference to the mesh
1238 
1239  // This is a parallel_only function
1240  const dof_id_type n_active_elem = mesh.n_active_elem();
1241 
1242  if (MeshOutput<MeshBase>::mesh().processor_id() != 0)
1243  return;
1244 
1245  std::ofstream out_stream (fname.c_str());
1246 
1247  libmesh_assert (out_stream.good());
1248 
1249  unsigned int mesh_max_p_level = 0;
1250 
1251  std::string buffer;
1252 
1253  // Begin interfacing with the GMV data file
1254  {
1255  // write the nodes
1256  buffer = "gmvinput";
1257  out_stream.write(buffer.c_str(), buffer.size());
1258 
1259  buffer = "ieeei4r4";
1260  out_stream.write(buffer.c_str(), buffer.size());
1261  }
1262 
1263 
1264 
1265  // write the nodes
1266  {
1267  buffer = "nodes ";
1268  out_stream.write(buffer.c_str(), buffer.size());
1269 
1270  unsigned int tempint = mesh.n_nodes();
1271  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1272 
1273  // write the x coordinate
1274  std::vector<float> temp(mesh.n_nodes());
1275  for (auto v : IntRange<unsigned int>(0, mesh.n_nodes()))
1276  temp[v] = static_cast<float>(mesh.point(v)(0));
1277  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1278 
1279  // write the y coordinate
1280  for (auto v : IntRange<unsigned int>(0, mesh.n_nodes()))
1281  {
1282 #if LIBMESH_DIM > 1
1283  temp[v] = static_cast<float>(mesh.point(v)(1));
1284 #else
1285  temp[v] = 0.;
1286 #endif
1287  }
1288  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1289 
1290  // write the z coordinate
1291  for (auto v : IntRange<unsigned int>(0, mesh.n_nodes()))
1292  {
1293 #if LIBMESH_DIM > 2
1294  temp[v] = static_cast<float>(mesh.point(v)(2));
1295 #else
1296  temp[v] = 0.;
1297 #endif
1298  }
1299  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1300  }
1301 
1302 
1303  // write the connectivity
1304  {
1305  buffer = "cells ";
1306  out_stream.write(buffer.c_str(), buffer.size());
1307 
1308  unsigned int tempint = n_active_elem;
1309  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1310 
1311  for (const auto & elem : mesh.active_element_ptr_range())
1312  {
1313  mesh_max_p_level = std::max(mesh_max_p_level,
1314  elem->p_level());
1315 
1316  // The ElementDefinition label contains the GMV name
1317  // and the number of nodes for the respective
1318  // element.
1319  const ElementDefinition & ed = eletypes[elem->type()];
1320 
1321  // The ElementDefinition labels all have a name followed by a
1322  // whitespace and then the number of nodes. To write the
1323  // string in binary, we want to drop everything after that
1324  // space, and then pad the string out to 8 characters.
1325  buffer = ed.label;
1326 
1327  // Erase everything from the first whitespace character to the end of the string.
1328  buffer.erase(buffer.find_first_of(' '), std::string::npos);
1329 
1330  // Append whitespaces until the string is exactly 8 characters long.
1331  while (buffer.size() < 8)
1332  buffer.insert(buffer.end(), ' ');
1333 
1334  // Finally, write the 8 character stream to file.
1335  out_stream.write(buffer.c_str(), buffer.size());
1336 
1337  // Write the number of nodes that this type of element has.
1338  // Note: don't want to use the number of nodes of the element,
1339  // because certain elements, like QUAD9 and HEX27 only support
1340  // being written out as lower-order elements (QUAD8 and HEX20,
1341  // respectively).
1342  tempint = cast_int<unsigned int>(ed.node_map.size());
1343  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1344 
1345  // Write the element connectivity
1346  for (const auto & ed_id : ed.node_map)
1347  {
1348  dof_id_type id = elem->node_id(ed_id) + 1;
1349  out_stream.write(reinterpret_cast<char *>(&id), sizeof(dof_id_type));
1350  }
1351  }
1352  }
1353 
1354 
1355 
1356  // optionally write the partition information
1357  if (this->partitioning())
1358  {
1359  if (this->write_subdomain_id_as_material())
1360  libmesh_error_msg("Not yet supported in GMVIO::write_binary");
1361 
1362  else
1363  {
1364  buffer = "material";
1365  out_stream.write(buffer.c_str(), buffer.size());
1366 
1367  unsigned int tmpint = mesh.n_processors();
1368  out_stream.write(reinterpret_cast<char *>(&tmpint), sizeof(unsigned int));
1369 
1370  tmpint = 0; // IDs are cell based
1371  out_stream.write(reinterpret_cast<char *>(&tmpint), sizeof(unsigned int));
1372 
1373  for (auto proc : IntRange<unsigned int>(0, mesh.n_processors()))
1374  {
1375  // We have to write exactly 8 bytes. This means that
1376  // the processor id can be up to 3 characters long, and
1377  // we pad it with blank characters on the end.
1378  std::ostringstream oss;
1379  oss << "proc_" << std::setw(3) << std::left << proc;
1380  out_stream.write(oss.str().c_str(), oss.str().size());
1381  }
1382 
1383  std::vector<unsigned int> proc_id (n_active_elem);
1384 
1385  unsigned int n=0;
1386 
1387  // We no longer write sub-elems in GMV, so just assign a
1388  // processor id value to each element.
1389  for (const auto & elem : mesh.active_element_ptr_range())
1390  proc_id[n++] = elem->processor_id() + 1;
1391 
1392  out_stream.write(reinterpret_cast<char *>(proc_id.data()),
1393  sizeof(unsigned int)*proc_id.size());
1394  }
1395  }
1396 
1397  // If there are *any* variables at all in the system (including
1398  // p level, or arbitrary cell-based data)
1399  // to write, the gmv file needs to contain the word "variable"
1400  // on a line by itself.
1401  bool write_variable = false;
1402 
1403  // 1.) p-levels
1404  if (this->p_levels() && mesh_max_p_level)
1405  write_variable = true;
1406 
1407  // 2.) solution data
1408  if ((solution_names != nullptr) && (vec != nullptr))
1409  write_variable = true;
1410 
1411  // // 3.) cell-centered data - unsupported
1412  // if (!(this->_cell_centered_data.empty()))
1413  // write_variable = true;
1414 
1415  if (write_variable)
1416  {
1417  buffer = "variable";
1418  out_stream.write(buffer.c_str(), buffer.size());
1419  }
1420 
1421  // optionally write the partition information
1422  if (this->p_levels() && mesh_max_p_level)
1423  {
1424  unsigned int n_floats = n_active_elem;
1425  for (unsigned int i=0, d=mesh.mesh_dimension(); i != d; ++i)
1426  n_floats *= 2;
1427 
1428  std::vector<float> temp(n_floats);
1429 
1430  buffer = "p_level";
1431  out_stream.write(buffer.c_str(), buffer.size());
1432 
1433  unsigned int tempint = 0; // p levels are cell data
1434  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1435 
1436  unsigned int n=0;
1437  for (const auto & elem : mesh.active_element_ptr_range())
1438  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1439  temp[n++] = static_cast<float>( elem->p_level() );
1440 
1441  out_stream.write(reinterpret_cast<char *>(temp.data()),
1442  sizeof(float)*n_floats);
1443  }
1444 
1445 
1446  // optionally write cell-centered data
1447  if (!(this->_cell_centered_data.empty()))
1448  libMesh::err << "Cell-centered data not (yet) supported in binary I/O mode!" << std::endl;
1449 
1450 
1451 
1452 
1453  // optionally write the data
1454  if ((solution_names != nullptr) &&
1455  (vec != nullptr))
1456  {
1457  std::vector<float> temp(mesh.n_nodes());
1458 
1459  const unsigned int n_vars =
1460  cast_int<unsigned int>(solution_names->size());
1461 
1462  for (unsigned int c=0; c<n_vars; c++)
1463  {
1464 
1465 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1466  // for complex data, write three datasets
1467 
1468 
1469  // Real part
1470  buffer = "r_";
1471  out_stream.write(buffer.c_str(), buffer.size());
1472 
1473  buffer = (*solution_names)[c];
1474  out_stream.write(buffer.c_str(), buffer.size());
1475 
1476  unsigned int tempint = 1; // always do nodal data
1477  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1478 
1479  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1480  temp[n] = static_cast<float>( (*vec)[n*n_vars + c].real() );
1481 
1482  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1483 
1484 
1485  // imaginary part
1486  buffer = "i_";
1487  out_stream.write(buffer.c_str(), buffer.size());
1488 
1489  buffer = (*solution_names)[c];
1490  out_stream.write(buffer.c_str(), buffer.size());
1491 
1492  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1493 
1494  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1495  temp[n] = static_cast<float>( (*vec)[n*n_vars + c].imag() );
1496 
1497  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1498 
1499  // magnitude
1500  buffer = "a_";
1501  out_stream.write(buffer.c_str(), buffer.size());
1502  buffer = (*solution_names)[c];
1503  out_stream.write(buffer.c_str(), buffer.size());
1504 
1505  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1506 
1507  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1508  temp[n] = static_cast<float>(std::abs((*vec)[n*n_vars + c]));
1509 
1510  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1511 
1512 #else
1513 
1514  buffer = (*solution_names)[c];
1515  out_stream.write(buffer.c_str(), buffer.size());
1516 
1517  unsigned int tempint = 1; // always do nodal data
1518  out_stream.write(reinterpret_cast<char *>(&tempint), sizeof(unsigned int));
1519 
1520  for (auto n : IntRange<unsigned int>(0, mesh.n_nodes()))
1521  temp[n] = static_cast<float>((*vec)[n*n_vars + c]);
1522 
1523  out_stream.write(reinterpret_cast<char *>(temp.data()), sizeof(float)*mesh.n_nodes());
1524 
1525 #endif
1526  }
1527  }
1528 
1529  // If we wrote any variables, we have to close the variable section now
1530  if (write_variable)
1531  {
1532  buffer = "endvars ";
1533  out_stream.write(buffer.c_str(), buffer.size());
1534  }
1535 
1536  // end the file
1537  buffer = "endgmv ";
1538  out_stream.write(buffer.c_str(), buffer.size());
1539 }

References _cell_centered_data, std::abs(), libMesh::MeshBase::active_element_ptr_range(), libMesh::err, std::imag(), libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_nodes(), libMesh::ParallelObject::n_processors(), n_vars, p_levels(), partitioning(), libMesh::MeshBase::point(), std::real(), and write_subdomain_id_as_material().

Referenced by write(), and write_nodal_data().

◆ write_discontinuous_equation_systems()

void libMesh::MeshOutput< MeshBase >::write_discontinuous_equation_systems ( const std::string &  fname,
const EquationSystems es,
const std::set< std::string > *  system_names = nullptr 
)
virtualinherited

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

Definition at line 87 of file mesh_output.C.

90 {
91  LOG_SCOPE("write_discontinuous_equation_systems()", "MeshOutput");
92 
93  // We may need to gather and/or renumber a DistributedMesh to output
94  // it, making that const qualifier in our constructor a dirty lie
95  MT & my_mesh = const_cast<MT &>(*_obj);
96 
97  // If we're asked to write data that's associated with a different
98  // mesh, output files full of garbage are the result.
99  libmesh_assert_equal_to(&es.get_mesh(), _obj);
100 
101  // A non-renumbered mesh may not have a contiguous numbering, and
102  // that needs to be fixed before we can build a solution vector.
103  if (my_mesh.max_elem_id() != my_mesh.n_elem() ||
104  my_mesh.max_node_id() != my_mesh.n_nodes())
105  {
106  // If we were allowed to renumber then we should have already
107  // been properly renumbered...
108  libmesh_assert(!my_mesh.allow_renumbering());
109 
110  libmesh_do_once(libMesh::out <<
111  "Warning: This MeshOutput subclass only supports meshes which are contiguously renumbered!"
112  << std::endl;);
113 
114  my_mesh.allow_renumbering(true);
115 
116  my_mesh.renumber_nodes_and_elements();
117 
118  // Not sure what good going back to false will do here, the
119  // renumbering horses have already left the barn...
120  my_mesh.allow_renumbering(false);
121  }
122 
123  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
124 
125  // Build the list of variable names that will be written.
126  std::vector<std::string> names;
127  es.build_variable_names (names, nullptr, system_names);
128 
129  if (!_is_parallel_format)
130  {
131  // Build the nodal solution values & get the variable
132  // names from the EquationSystems object
133  std::vector<Number> soln;
134  es.build_discontinuous_solution_vector (soln, system_names);
135 
136  this->write_nodal_data_discontinuous (fname, soln, names);
137  }
138  else // _is_parallel_format
139  {
140  libmesh_not_implemented();
141  }
142 }

◆ write_discontinuous_gmv()

void libMesh::GMVIO::write_discontinuous_gmv ( const std::string &  name,
const EquationSystems es,
const bool  write_partitioning,
const std::set< std::string > *  system_names = nullptr 
) const

Writes a GMV file with discontinuous data.

Definition at line 1549 of file gmv_io.C.

1553 {
1554  std::vector<std::string> solution_names;
1555  std::vector<Number> v;
1556 
1557  // Get a reference to the mesh
1559 
1560  es.build_variable_names (solution_names, nullptr, system_names);
1561  es.build_discontinuous_solution_vector (v, system_names);
1562 
1563  // These are parallel_only functions
1564  const unsigned int n_active_elem = mesh.n_active_elem();
1565 
1566  if (mesh.processor_id() != 0)
1567  return;
1568 
1569  std::ofstream out_stream(name.c_str());
1570 
1571  libmesh_assert (out_stream.good());
1572 
1573  // Begin interfacing with the GMV data file
1574  {
1575 
1576  // write the nodes
1577  out_stream << "gmvinput ascii" << std::endl << std::endl;
1578 
1579  // Compute the total weight
1580  {
1581  unsigned int tw=0;
1582 
1583  for (const auto & elem : mesh.active_element_ptr_range())
1584  tw += elem->n_nodes();
1585 
1586  out_stream << "nodes " << tw << std::endl;
1587  }
1588 
1589 
1590 
1591  // Write all the x values
1592  {
1593  for (const auto & elem : mesh.active_element_ptr_range())
1594  for (const Node & node : elem->node_ref_range())
1595  out_stream << node(0) << " ";
1596 
1597  out_stream << std::endl;
1598  }
1599 
1600 
1601  // Write all the y values
1602  {
1603  for (const auto & elem : mesh.active_element_ptr_range())
1604  for (const Node & node : elem->node_ref_range())
1605 #if LIBMESH_DIM > 1
1606  out_stream << node(1) << " ";
1607 #else
1608  out_stream << 0. << " ";
1609 #endif
1610 
1611  out_stream << std::endl;
1612  }
1613 
1614 
1615  // Write all the z values
1616  {
1617  for (const auto & elem : mesh.active_element_ptr_range())
1618  for (const Node & node : elem->node_ref_range())
1619 #if LIBMESH_DIM > 2
1620  out_stream << node(2) << " ";
1621 #else
1622  out_stream << 0. << " ";
1623 #endif
1624 
1625  out_stream << std::endl << std::endl;
1626  }
1627  }
1628 
1629 
1630 
1631  {
1632  // write the connectivity
1633 
1634  out_stream << "cells " << n_active_elem << std::endl;
1635 
1636  unsigned int nn=1;
1637 
1638  switch (mesh.mesh_dimension())
1639  {
1640  case 1:
1641  {
1642  for (const auto & elem : mesh.active_element_ptr_range())
1643  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1644  {
1645  if ((elem->type() == EDGE2) ||
1646  (elem->type() == EDGE3) ||
1647  (elem->type() == EDGE4)
1648 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1649  || (elem->type() == INFEDGE2)
1650 #endif
1651  )
1652  {
1653  out_stream << "line 2" << std::endl;
1654  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1655  out_stream << nn++ << " ";
1656 
1657  }
1658  else
1659  libmesh_error_msg("Unsupported 1D element type: " << Utility::enum_to_string(elem->type()));
1660 
1661  out_stream << std::endl;
1662  }
1663 
1664  break;
1665  }
1666 
1667  case 2:
1668  {
1669  for (const auto & elem : mesh.active_element_ptr_range())
1670  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1671  {
1672  if ((elem->type() == QUAD4) ||
1673  (elem->type() == QUAD8) || // Note: QUAD8 will be output as one central quad and
1674  // four surrounding triangles (though they will be written
1675  // to GMV as QUAD4s).
1676  (elem->type() == QUAD9)
1677 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1678  || (elem->type() == INFQUAD4)
1679  || (elem->type() == INFQUAD6)
1680 #endif
1681  )
1682  {
1683  out_stream << "quad 4" << std::endl;
1684  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1685  out_stream << nn++ << " ";
1686 
1687  }
1688  else if ((elem->type() == TRI3) ||
1689  (elem->type() == TRI6))
1690  {
1691  out_stream << "tri 3" << std::endl;
1692  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1693  out_stream << nn++ << " ";
1694 
1695  }
1696  else
1697  libmesh_error_msg("Unsupported 2D element type: " << Utility::enum_to_string(elem->type()));
1698 
1699  out_stream << std::endl;
1700  }
1701 
1702  break;
1703  }
1704 
1705 
1706  case 3:
1707  {
1708  for (const auto & elem : mesh.active_element_ptr_range())
1709  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
1710  {
1711  if ((elem->type() == HEX8) ||
1712  (elem->type() == HEX20) ||
1713  (elem->type() == HEX27)
1714 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1715  || (elem->type() == INFHEX8)
1716  || (elem->type() == INFHEX16)
1717  || (elem->type() == INFHEX18)
1718 #endif
1719  )
1720  {
1721  out_stream << "phex8 8" << std::endl;
1722  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1723  out_stream << nn++ << " ";
1724  }
1725  else if ((elem->type() == PRISM6) ||
1726  (elem->type() == PRISM15) ||
1727  (elem->type() == PRISM18)
1728 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1729  || (elem->type() == INFPRISM6)
1730  || (elem->type() == INFPRISM12)
1731 #endif
1732  )
1733  {
1734  out_stream << "pprism6 6" << std::endl;
1735  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1736  out_stream << nn++ << " ";
1737  }
1738  else if ((elem->type() == TET4) ||
1739  (elem->type() == TET10))
1740  {
1741  out_stream << "tet 4" << std::endl;
1742  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1743  out_stream << nn++ << " ";
1744  }
1745  else
1746  libmesh_error_msg("Unsupported 3D element type: " << Utility::enum_to_string(elem->type()));
1747 
1748  out_stream << std::endl;
1749  }
1750 
1751  break;
1752  }
1753 
1754  default:
1755  libmesh_error_msg("Unsupported mesh dimension: " << mesh.mesh_dimension());
1756  }
1757 
1758  out_stream << std::endl;
1759  }
1760 
1761 
1762 
1763  // optionally write the partition information
1764  if (write_partitioning)
1765  {
1767  libmesh_error_msg("Not yet supported in GMVIO::write_discontinuous_gmv");
1768 
1769  else
1770  {
1771  out_stream << "material "
1772  << mesh.n_processors()
1773  << " 0"<< std::endl;
1774 
1775  for (auto proc : IntRange<unsigned int>(0, mesh.n_processors()))
1776  out_stream << "proc_" << proc << std::endl;
1777 
1778  for (const auto & elem : mesh.active_element_ptr_range())
1779  out_stream << elem->processor_id()+1 << std::endl;
1780 
1781  out_stream << std::endl;
1782  }
1783  }
1784 
1785 
1786  // Writing cell-centered data is not yet supported in discontinuous GMV files.
1787  if (!(this->_cell_centered_data.empty()))
1788  {
1789  libMesh::err << "Cell-centered data not (yet) supported for discontinuous GMV files!" << std::endl;
1790  }
1791 
1792 
1793 
1794  // write the data
1795  {
1796  const unsigned int n_vars =
1797  cast_int<unsigned int>(solution_names.size());
1798 
1799  // libmesh_assert_equal_to (v.size(), tw*n_vars);
1800 
1801  out_stream << "variable" << std::endl;
1802 
1803 
1804  for (unsigned int c=0; c<n_vars; c++)
1805  {
1806 
1807 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1808 
1809  // in case of complex data, write _tree_ data sets
1810  // for each component
1811 
1812  // this is the real part
1813  out_stream << "r_" << solution_names[c] << " 1" << std::endl;
1814  for (const auto & elem : mesh.active_element_ptr_range())
1815  for (auto n : elem->node_index_range())
1816  out_stream << v[(n++)*n_vars + c].real() << " ";
1817  out_stream << std::endl << std::endl;
1818 
1819 
1820  // this is the imaginary part
1821  out_stream << "i_" << solution_names[c] << " 1" << std::endl;
1822  for (const auto & elem : mesh.active_element_ptr_range())
1823  for (auto n : elem->node_index_range())
1824  out_stream << v[(n++)*n_vars + c].imag() << " ";
1825  out_stream << std::endl << std::endl;
1826 
1827  // this is the magnitude
1828  out_stream << "a_" << solution_names[c] << " 1" << std::endl;
1829  for (const auto & elem : mesh.active_element_ptr_range())
1830  for (auto n : elem->node_index_range())
1831  out_stream << std::abs(v[(n++)*n_vars + c]) << " ";
1832  out_stream << std::endl << std::endl;
1833 
1834 #else
1835 
1836  out_stream << solution_names[c] << " 1" << std::endl;
1837  {
1838  unsigned int nn=0;
1839 
1840  for (const auto & elem : mesh.active_element_ptr_range())
1841  for (unsigned int i=0, enn=elem->n_nodes(); i<enn; i++)
1842  out_stream << v[(nn++)*n_vars + c] << " ";
1843  }
1844  out_stream << std::endl << std::endl;
1845 
1846 #endif
1847 
1848  }
1849 
1850  out_stream << "endvars" << std::endl;
1851  }
1852 
1853 
1854  // end of the file
1855  out_stream << std::endl << "endgmv" << std::endl;
1856 }

References _cell_centered_data, _write_subdomain_id_as_material, std::abs(), libMesh::MeshBase::active_element_ptr_range(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EDGE2, libMesh::EDGE3, libMesh::EDGE4, libMesh::Utility::enum_to_string(), libMesh::err, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::INFEDGE2, libMesh::INFHEX16, libMesh::INFHEX18, libMesh::INFHEX8, libMesh::INFPRISM12, libMesh::INFPRISM6, libMesh::INFQUAD4, libMesh::INFQUAD6, libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::ParallelObject::n_processors(), n_vars, libMesh::Quality::name(), libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::ParallelObject::processor_id(), libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TET10, libMesh::TET4, libMesh::TRI3, and libMesh::TRI6.

Referenced by libMesh::ErrorVector::plot_error().

◆ write_equation_systems()

void libMesh::MeshOutput< MeshBase >::write_equation_systems ( const std::string &  fname,
const EquationSystems es,
const std::set< std::string > *  system_names = nullptr 
)
virtualinherited

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

Reimplemented in libMesh::NameBasedIO.

Definition at line 31 of file mesh_output.C.

34 {
35  LOG_SCOPE("write_equation_systems()", "MeshOutput");
36 
37  // We may need to gather and/or renumber a DistributedMesh to output
38  // it, making that const qualifier in our constructor a dirty lie
39  MT & my_mesh = const_cast<MT &>(*_obj);
40 
41  // If we're asked to write data that's associated with a different
42  // mesh, output files full of garbage are the result.
43  libmesh_assert_equal_to(&es.get_mesh(), _obj);
44 
45  // A non-renumbered mesh may not have a contiguous numbering, and
46  // that needs to be fixed before we can build a solution vector.
47  if (my_mesh.max_elem_id() != my_mesh.n_elem() ||
48  my_mesh.max_node_id() != my_mesh.n_nodes())
49  {
50  // If we were allowed to renumber then we should have already
51  // been properly renumbered...
52  libmesh_assert(!my_mesh.allow_renumbering());
53 
54  libmesh_do_once(libMesh::out <<
55  "Warning: This MeshOutput subclass only supports meshes which are contiguously renumbered!"
56  << std::endl;);
57 
58  my_mesh.allow_renumbering(true);
59 
60  my_mesh.renumber_nodes_and_elements();
61 
62  // Not sure what good going back to false will do here, the
63  // renumbering horses have already left the barn...
64  my_mesh.allow_renumbering(false);
65  }
66 
68  {
69  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
70 
71  // Build the list of variable names that will be written.
72  std::vector<std::string> names;
73  es.build_variable_names (names, nullptr, system_names);
74 
75  // Build the nodal solution values & get the variable
76  // names from the EquationSystems object
77  std::vector<Number> soln;
78  es.build_solution_vector (soln, system_names);
79 
80  this->write_nodal_data (fname, soln, names);
81  }
82  else // _is_parallel_format
83  this->write_nodal_data (fname, es, system_names);
84 }

◆ write_nodal_data() [1/3]

void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  fname,
const EquationSystems es,
const std::set< std::string > *  system_names 
)
virtualinherited

This method should be overridden by "parallel" output formats for writing nodal data.

Instead of getting a localized copy of the nodal solution vector, it directly uses EquationSystems current_local_solution vectors to look up nodal values.

If not implemented, reorders the solutions into a nodal-only NumericVector and calls the above version of this function.

Reimplemented in libMesh::Nemesis_IO.

Definition at line 158 of file mesh_output.C.

161 {
162  std::vector<std::string> names;
163  es.build_variable_names (names, nullptr, system_names);
164 
165  std::unique_ptr<NumericVector<Number>> parallel_soln =
166  es.build_parallel_solution_vector(system_names);
167 
168  this->write_nodal_data (fname, *parallel_soln, names);
169 }

◆ write_nodal_data() [2/3]

void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  fname,
const NumericVector< Number > &  parallel_soln,
const std::vector< std::string > &  names 
)
virtualinherited

This method may be overridden by "parallel" output formats for writing nodal data.

Instead of getting a localized copy of the nodal solution vector, it is passed a NumericVector of type=PARALLEL which is in node-major order i.e. (u0,v0,w0, u1,v1,w1, u2,v2,w2, u3,v3,w3, ...) and contains n_nodes*n_vars total entries. Then, it is up to the individual I/O class to extract the required solution values from this vector and write them in parallel.

If not implemented, localizes the parallel vector into a std::vector and calls the other version of this function.

Reimplemented in libMesh::Nemesis_IO.

Definition at line 145 of file mesh_output.C.

148 {
149  // This is the fallback implementation for parallel I/O formats that
150  // do not yet implement proper writing in parallel, and instead rely
151  // on the full solution vector being available on all processors.
152  std::vector<Number> soln;
153  parallel_soln.localize(soln);
154  this->write_nodal_data(fname, soln, names);
155 }

◆ write_nodal_data() [3/3]

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

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 281 of file gmv_io.C.

284 {
285  LOG_SCOPE("write_nodal_data()", "GMVIO");
286 
287  if (this->binary())
288  this->write_binary (fname, &soln, &names);
289  else
290  this->write_ascii_old_impl (fname, &soln, &names);
291 }

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

Referenced by libMesh::NameBasedIO::write_nodal_data().

◆ write_nodal_data_discontinuous()

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data_discontinuous ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
)
inlinevirtualinherited

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

Reimplemented in libMesh::ExodusII_IO.

Definition at line 114 of file mesh_output.h.

117  { libmesh_not_implemented(); }

◆ write_subdomain_id_as_material()

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 122 of file gmv_io.h.

References _write_subdomain_id_as_material.

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

Member Data Documentation

◆ _ascii_precision

unsigned int libMesh::MeshOutput< MeshBase >::_ascii_precision
privateinherited

Precision to use when writing ASCII files.

Definition at line 195 of file mesh_output.h.

◆ _binary

bool libMesh::GMVIO::_binary
private

Flag to write binary data.

Definition at line 199 of file gmv_io.h.

Referenced by binary().

◆ _cell_centered_data

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

◆ _discontinuous

bool libMesh::GMVIO::_discontinuous
private

Flag to write the mesh as discontinuous patches.

Definition at line 204 of file gmv_io.h.

Referenced by discontinuous().

◆ _is_parallel_format [1/2]

const bool libMesh::MeshInput< MeshBase >::_is_parallel_format
privateinherited

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 121 of file mesh_input.h.

◆ _is_parallel_format [2/2]

const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format
protectedinherited

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 172 of file mesh_output.h.

◆ _next_elem_id

unsigned int libMesh::GMVIO::_next_elem_id
private

Definition at line 239 of file gmv_io.h.

Referenced by _read_one_cell(), and read().

◆ _nodal_data

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

Definition at line 244 of file gmv_io.h.

Referenced by _read_var(), and copy_nodal_solution().

◆ _obj

MeshBase * libMesh::MeshInput< MeshBase >::_obj
privateinherited

A pointer to a non-const object object.

This allows us to read the object from file.

Definition at line 114 of file mesh_input.h.

◆ _p_levels

bool libMesh::GMVIO::_p_levels
private

Flag to write the mesh p refinement levels.

Definition at line 225 of file gmv_io.h.

Referenced by p_levels().

◆ _partitioning

bool libMesh::GMVIO::_partitioning
private

Flag to write the mesh partitioning.

Definition at line 209 of file gmv_io.h.

Referenced by partitioning().

◆ _reading_element_map

std::map< std::string, ElemType > libMesh::GMVIO::_reading_element_map = GMVIO::build_reading_element_map()
staticprivate

Static map from string -> ElementType for use during reading.

Definition at line 250 of file gmv_io.h.

Referenced by gmv_elem_to_libmesh_elem().

◆ _serial_only_needed_on_proc_0

const bool libMesh::MeshOutput< MeshBase >::_serial_only_needed_on_proc_0
protectedinherited

Flag specifying whether this format can be written by only serializing the mesh to processor zero.

If this is false (default) the mesh will be serialized to all processors

Definition at line 181 of file mesh_output.h.

◆ _subdivide_second_order

bool libMesh::GMVIO::_subdivide_second_order
private

Flag to subdivide second order elements.

Definition at line 220 of file gmv_io.h.

Referenced by subdivide_second_order().

◆ _write_subdomain_id_as_material

bool libMesh::GMVIO::_write_subdomain_id_as_material
private

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

Definition at line 215 of file gmv_io.h.

Referenced by write_discontinuous_gmv(), and write_subdomain_id_as_material().

◆ elems_of_dimension

std::vector<bool> libMesh::MeshInput< MeshBase >::elems_of_dimension
protectedinherited

A vector of bools describing what dimension elements have been encountered when reading a mesh.

Definition at line 97 of file mesh_input.h.


The documentation for this class was generated from the following files:
libMesh::GMVIO::_subdivide_second_order
bool _subdivide_second_order
Flag to subdivide second order elements.
Definition: gmv_io.h:220
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::GMVIO::p_levels
bool & p_levels()
Flag indicating whether or not to write p level information for p refined meshes.
Definition: gmv_io.h:134
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::MeshOutput< MeshBase >::_serial_only_needed_on_proc_0
const bool _serial_only_needed_on_proc_0
Flag specifying whether this format can be written by only serializing the mesh to processor zero.
Definition: mesh_output.h:181
libMesh::OrderWrapper::get_order
int get_order() const
Explicitly request the order as an int.
Definition: fe_type.h:77
libMesh::MeshBase::set_n_partitions
unsigned int & set_n_partitions()
Definition: mesh_base.h:1667
libMesh::GMVIO::write_ascii_old_impl
void write_ascii_old_impl(const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: gmv_io.C:570
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::GMVIO::_write_subdomain_id_as_material
bool _write_subdomain_id_as_material
Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's.
Definition: gmv_io.h:215
n_vars
unsigned int n_vars
Definition: adaptivity_ex3.C:116
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::MeshBase::point
virtual const Point & point(const dof_id_type i) const =0
libMesh::EDGE4
Definition: enum_elem_type.h:37
libMesh::INFHEX8
Definition: enum_elem_type.h:60
libMesh::INFQUAD4
Definition: enum_elem_type.h:58
libMesh::MeshBase::active_element_ptr_range
virtual SimpleRange< element_iterator > active_element_ptr_range()=0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
libMesh::GMVIO::_read_materials
void _read_materials()
Definition: gmv_io.C:2042
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::GMVIO::subdivide_second_order
bool & subdivide_second_order()
Flag indicating whether or not to subdivide second order elements.
Definition: gmv_io.h:128
libMesh::GMVIO::_read_one_cell
void _read_one_cell()
Definition: gmv_io.C:2085
libMesh::MeshBase::n_partitions
unsigned int n_partitions() const
Definition: mesh_base.h:1153
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::MeshBase::max_node_id
virtual dof_id_type max_node_id() const =0
libMesh::ierr
ierr
Definition: petsc_dm_wrapper.C:72
libMesh::INFEDGE2
Definition: enum_elem_type.h:57
libMesh::MeshInput< MeshBase >::_obj
MeshBase * _obj
A pointer to a non-const object object.
Definition: mesh_input.h:114
libMesh::MeshBase::n_active_sub_elem
dof_id_type n_active_sub_elem() const
Same as n_sub_elem(), but only counts active elements.
Definition: mesh_base.C:539
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::GMVIO::_partitioning
bool _partitioning
Flag to write the mesh partitioning.
Definition: gmv_io.h:209
libMesh::MeshOutput< MeshBase >::ascii_precision
unsigned int & ascii_precision()
Return/set the precision to use when writing ASCII files.
Definition: mesh_output.h:257
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::INFPRISM6
Definition: enum_elem_type.h:63
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::GMVIO::_p_levels
bool _p_levels
Flag to write the mesh p refinement levels.
Definition: gmv_io.h:225
libMesh::System::has_variable
bool has_variable(const std::string &var) const
Definition: system.C:1225
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
libMesh::MeshSerializer
Temporarily serialize a DistributedMesh for output; a distributed mesh is allgathered by the MeshSeri...
Definition: mesh_serializer.h:42
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::EquationSystems::build_variable_names
void build_variable_names(std::vector< std::string > &var_names, const FEType *type=nullptr, const std::set< std::string > *system_names=nullptr) const
Fill the input vector var_names with the names of the variables for each system.
Definition: equation_systems.C:476
libMesh::Elem::first_order_equivalent_type
static ElemType first_order_equivalent_type(const ElemType et)
Definition: elem.C:2369
libMesh::MeshOutput< MeshBase >::_is_parallel_format
const bool _is_parallel_format
Flag specifying whether this format is parallel-capable.
Definition: mesh_output.h:172
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::INFHEX18
Definition: enum_elem_type.h:62
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::GMVIO::_reading_element_map
static std::map< std::string, ElemType > _reading_element_map
Static map from string -> ElementType for use during reading.
Definition: gmv_io.h:250
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::MeshTools::Generation::Private::idx
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
Definition: mesh_generation.C:72
libMesh::GMVIO::_binary
bool _binary
Flag to write binary data.
Definition: gmv_io.h:199
libMesh::GMVIO::_read_var
void _read_var()
Definition: gmv_io.C:2030
libMesh::EquationSystems::n_systems
unsigned int n_systems() const
Definition: equation_systems.h:652
libMesh::INFHEX16
Definition: enum_elem_type.h:61
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
libMesh::INFPRISM12
Definition: enum_elem_type.h:64
libMesh::GMVIO::_next_elem_id
unsigned int _next_elem_id
Definition: gmv_io.h:239
libMesh::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
libMesh::GMVIO::_read_nodes
void _read_nodes()
Helper functions for reading nodes/cells from a GMV file.
Definition: gmv_io.C:2066
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::MeshOutput< MeshBase >
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::MeshOutput< MeshBase >::write_nodal_data
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: mesh_output.h:105
libMesh::MeshOutput::mesh
const MT & mesh() const
Definition: mesh_output.h:247
libMesh::FEType::order
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:197
libMesh::GMVIO::write_binary
void write_binary(const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: gmv_io.C:1228
libMesh::INFQUAD6
Definition: enum_elem_type.h:59
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::GMVIO::partitioning
bool & partitioning()
Flag indicating whether or not to write the partitioning information for the mesh.
Definition: gmv_io.h:115
libMesh::EquationSystems::build_discontinuous_solution_vector
void build_discontinuous_solution_vector(std::vector< Number > &soln, const std::set< std::string > *system_names=nullptr, const std::vector< std::string > *var_names=nullptr, bool vertices_only=false) const
Fill the input vector soln with solution values.
Definition: equation_systems.C:1049
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::GMVIO::binary
bool & binary()
Flag indicating whether or not to write a binary file.
Definition: gmv_io.h:103
libMesh::GMVIO::_discontinuous
bool _discontinuous
Flag to write the mesh as discontinuous patches.
Definition: gmv_io.h:204
libMesh::MeshOutput< MeshBase >::write_nodal_data_discontinuous
virtual void write_nodal_data_discontinuous(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
This method implements writing a mesh with discontinuous data to a specified file where the nodal dat...
Definition: mesh_output.h:114
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
libMesh::GMVIO::write_subdomain_id_as_material
bool & write_subdomain_id_as_material()
Flag to write element subdomain_id's as GMV "materials" instead of element processor_id's.
Definition: gmv_io.h:122
libMesh::Trees::NODES
Definition: tree_base.h:55
libMesh::MeshInput< MeshBase >::mesh
MeshBase & mesh()
Definition: mesh_input.h:169
libMesh::MeshOutput< MeshBase >::_obj
const MeshBase *const _obj
A pointer to a constant object.
Definition: mesh_output.h:190
libMesh::GMVIO::gmv_elem_to_libmesh_elem
ElemType gmv_elem_to_libmesh_elem(std::string elemname)
Definition: gmv_io.C:2149
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::err
OStreamProxy err
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::MeshBase::allow_renumbering
void allow_renumbering(bool allow)
If false is passed in then this mesh will no longer be renumbered when being prepared for use.
Definition: mesh_base.h:1025
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::MeshBase::prepare_for_use
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:318
libMesh::MeshBase::clear
virtual void clear()
Deletes all the element and node data that is currently stored.
Definition: mesh_base.C:429
libMesh::MeshBase::set_mesh_dimension
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
Definition: mesh_base.h:218
libMesh::PRISM18
Definition: enum_elem_type.h:52
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::out
OStreamProxy out
libMesh::GMVIO::_nodal_data
std::map< std::string, std::vector< Number > > _nodal_data
Definition: gmv_io.h:244
libMesh::System::update
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:408
libMesh::FIRST
Definition: enum_order.h:42
std::imag
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
Definition: float128_shims.h:83
libMesh::TECPLOT
Definition: enum_io_package.h:39
libMesh::MeshBase::n_active_elem
virtual dof_id_type n_active_elem() const =0
libMesh::MeshOutput< MeshBase >::_ascii_precision
unsigned int _ascii_precision
Precision to use when writing ASCII files.
Definition: mesh_output.h:195
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::MeshInput< MeshBase >
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::MeshInput< MeshBase >::elems_of_dimension
std::vector< bool > elems_of_dimension
A vector of bools describing what dimension elements have been encountered when reading a mesh.
Definition: mesh_input.h:97
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::GMVIO::_cell_centered_data
std::map< std::string, const std::vector< Real > * > _cell_centered_data
Storage for arbitrary cell-centered data.
Definition: gmv_io.h:233
std::real
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)
Definition: float128_shims.h:77