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...
 
bool is_parallel_format () const
 Returns true iff this mesh file format and input class are parallelized, so that all processors can read their share of the data at once. 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
 
virtual bool get_add_sides ()
 

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
 

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 46 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 }
bool _binary
Flag to write binary data.
Definition: gmv_io.h:191
unsigned int _next_elem_id
Definition: gmv_io.h:231
bool _discontinuous
Flag to write the mesh as discontinuous patches.
Definition: gmv_io.h:196
bool _write_subdomain_id_as_material
Flag to write element subdomain_id&#39;s as GMV "materials" instead of element processor_id&#39;s.
Definition: gmv_io.h:207
bool _subdivide_second_order
Flag to subdivide second order elements.
Definition: gmv_io.h:212
bool _p_levels
Flag to write the mesh p refinement levels.
Definition: gmv_io.h:217
bool _partitioning
Flag to write the mesh partitioning.
Definition: gmv_io.h:201

◆ 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 }
bool _binary
Flag to write binary data.
Definition: gmv_io.h:191
unsigned int _next_elem_id
Definition: gmv_io.h:231
bool _discontinuous
Flag to write the mesh as discontinuous patches.
Definition: gmv_io.h:196
bool _write_subdomain_id_as_material
Flag to write element subdomain_id&#39;s as GMV "materials" instead of element processor_id&#39;s.
Definition: gmv_io.h:207
bool _subdivide_second_order
Flag to subdivide second order elements.
Definition: gmv_io.h:212
bool _p_levels
Flag to write the mesh p refinement levels.
Definition: gmv_io.h:217
bool _partitioning
Flag to write the mesh partitioning.
Definition: gmv_io.h:201

Member Function Documentation

◆ _read_materials()

void libMesh::GMVIO::_read_materials ( )
private

Definition at line 2028 of file gmv_io.C.

Referenced by read().

2029 {
2030 #ifdef LIBMESH_HAVE_GMV
2031 
2032  // LibMesh assigns materials on a per-cell basis
2033  libmesh_assert_equal_to (GMVLib::gmv_data.datatype, CELL);
2034 
2035  // Material names: LibMesh has no use for these currently...
2036  // for (int i = 0; i < GMVLib::gmv_data.num; i++)
2037  // {
2038  // // Build a 32-char string from the appropriate entries
2039  // std::string mat_string(&GMVLib::gmv_data.chardata1[i*33], 32);
2040  // }
2041 
2042  for (int i = 0; i < GMVLib::gmv_data.nlongdata1; i++)
2043  MeshInput<MeshBase>::mesh().elem_ref(i).processor_id() =
2044  cast_int<processor_id_type>(GMVLib::gmv_data.longdata1[i]-1);
2045 
2046 #endif
2047 }
This class defines an abstract interface for Mesh input.
Definition: mesh_base.h:56

◆ _read_nodes()

void libMesh::GMVIO::_read_nodes ( )
private

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

Definition at line 2052 of file gmv_io.C.

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

Referenced by read().

2053 {
2054 #ifdef LIBMESH_HAVE_GMV
2055  // LibMesh writes UNSTRUCT=100 node data
2056  libmesh_assert_equal_to (GMVLib::gmv_data.datatype, UNSTRUCT);
2057 
2058  // The nodal data is stored in gmv_data.doubledata{1,2,3}
2059  // and is nnodes long
2060  for (int i = 0; i < GMVLib::gmv_data.num; i++)
2061  {
2062  // Add the point to the Mesh
2063  MeshInput<MeshBase>::mesh().add_point(Point(GMVLib::gmv_data.doubledata1[i],
2064  GMVLib::gmv_data.doubledata2[i],
2065  GMVLib::gmv_data.doubledata3[i]), i);
2066  }
2067 #endif
2068 }
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39

◆ _read_one_cell()

void libMesh::GMVIO::_read_one_cell ( )
private

Definition at line 2071 of file gmv_io.C.

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

Referenced by read().

2072 {
2073 #ifdef LIBMESH_HAVE_GMV
2074  // This is either a REGULAR=111 cell or
2075  // the ENDKEYWORD=207 of the cells
2076 #ifndef NDEBUG
2077  bool recognized =
2078  (GMVLib::gmv_data.datatype==REGULAR) ||
2079  (GMVLib::gmv_data.datatype==ENDKEYWORD);
2080 #endif
2081  libmesh_assert (recognized);
2082 
2084 
2085  if (GMVLib::gmv_data.datatype == REGULAR)
2086  {
2087  // We need a mapping from GMV element types to LibMesh
2088  // ElemTypes. Basically the reverse of the eletypes
2089  // std::map above.
2090  //
2091  // FIXME: Since Quad9's apparently don't exist for GMV, and since
2092  // In general we write linear sub-elements to GMV files, we need
2093  // to be careful to read back in exactly what we wrote out...
2094  //
2095  // The gmv_data.name1 field is padded with blank characters when
2096  // it is read in by GMV, so the function that converts it to the
2097  // libmesh element type needs to take that into account.
2098  ElemType type = this->gmv_elem_to_libmesh_elem(GMVLib::gmv_data.name1);
2099 
2100  auto elem = Elem::build(type);
2101  elem->set_id(_next_elem_id++);
2102 
2103  // Get the ElementDefinition object for this element type
2104  const ElementDefinition & eledef = eletypes[type];
2105 
2106  // Print out the connectivity information for
2107  // this cell.
2108  for (int i=0; i<GMVLib::gmv_data.num2; i++)
2109  {
2110  // Map index i to GMV's numbering scheme
2111  unsigned mapped_i = eledef.node_map[i];
2112 
2113  // Note: Node numbers (as stored in libmesh) are 1-based
2114  elem->set_node(i) = mesh.node_ptr
2115  (cast_int<dof_id_type>(GMVLib::gmv_data.longdata1[mapped_i]-1));
2116  }
2117 
2118  elems_of_dimension[elem->dim()] = true;
2119 
2120  // Add the newly-created element to the mesh
2121  mesh.add_elem(std::move(elem));
2122  }
2123 
2124 
2125  if (GMVLib::gmv_data.datatype == ENDKEYWORD)
2126  {
2127  // There isn't a cell to read, so we just return
2128  return;
2129  }
2130 
2131 #endif
2132 }
ElemType
Defines an enum for geometric element types.
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:107
unsigned int _next_elem_id
Definition: gmv_io.h:231
This is the MeshBase class.
Definition: mesh_base.h:74
ElemType gmv_elem_to_libmesh_elem(std::string elemname)
Definition: gmv_io.C:2135
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:273
libmesh_assert(ctx)
virtual const Node * node_ptr(const dof_id_type i) const =0

◆ _read_var()

void libMesh::GMVIO::_read_var ( )
private

Definition at line 2016 of file gmv_io.C.

References _nodal_data.

Referenced by read().

2017 {
2018 #ifdef LIBMESH_HAVE_GMV
2019 
2020  // Copy all the variable's values into a local storage vector.
2021  _nodal_data.emplace(std::string(GMVLib::gmv_data.name1),
2022  std::vector<Number>(GMVLib::gmv_data.doubledata1, GMVLib::gmv_data.doubledata1+GMVLib::gmv_data.num));
2023 #endif
2024 }
std::map< std::string, std::vector< Number > > _nodal_data
Definition: gmv_io.h:236

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

References _cell_centered_data, and libMesh::libmesh_assert().

1851 {
1852  libmesh_assert(cell_centered_data_vals);
1853 
1854  // Make sure there are *at least* enough entries for all the active elements.
1855  // There can also be entries for inactive elements, they will be ignored.
1856  // libmesh_assert_greater_equal (cell_centered_data_vals->size(),
1857  // MeshOutput<MeshBase>::mesh().n_active_elem());
1858  this->_cell_centered_data[cell_centered_data_name] = cell_centered_data_vals;
1859 }
std::map< std::string, const std::vector< Real > *> _cell_centered_data
Storage for arbitrary cell-centered data.
Definition: gmv_io.h:225
libmesh_assert(ctx)

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

Referenced by libMesh::UNVIO::nodes_out(), libMesh::FroIO::write(), libMesh::MEDITIO::write_ascii(), libMesh::TecplotIO::write_ascii(), write_ascii_new_impl(), and write_ascii_old_impl().

270 {
271  return _ascii_precision;
272 }
unsigned int _ascii_precision
Precision to use when writing ASCII files.
Definition: mesh_output.h:207

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

References _binary.

Referenced by write(), and write_nodal_data().

95 { return _binary; }
bool _binary
Flag to write binary data.
Definition: gmv_io.h:191

◆ 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.

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.

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 }

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

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::make_range(), 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().

2148 {
2149  // Check for easy return if there isn't any nodal data
2150  if (_nodal_data.empty())
2151  {
2152  libMesh::err << "Unable to copy nodal solution: No nodal "
2153  << "solution has been read in from file." << std::endl;
2154  return;
2155  }
2156 
2157  // Be sure there is at least one system
2158  libmesh_assert (es.n_systems());
2159 
2160  // Keep track of variable names which have been found and
2161  // copied already. This could be used to prevent us from
2162  // e.g. copying the same var into 2 different systems ...
2163  // but this seems unlikely. Also, it is used to tell if
2164  // any variables which were read in were not successfully
2165  // copied to the EquationSystems.
2166  std::set<std::string> vars_copied;
2167 
2168  // For each entry in the nodal data map, try to find a system
2169  // that has the same variable key name.
2170  for (auto sys : make_range(es.n_systems()))
2171  {
2172  // Get a generic reference to the current System
2173  System & system = es.get_system(sys);
2174 
2175  // For each var entry in the _nodal_data map, try to find
2176  // that var in the system
2177  for (const auto & [var_name, vec] : _nodal_data)
2178  {
2179  if (system.has_variable(var_name))
2180  {
2181  // Check if there are as many nodes in the mesh as there are entries
2182  // in the stored nodal data vector
2183  libmesh_assert_equal_to (vec.size(), MeshInput<MeshBase>::mesh().n_nodes());
2184 
2185  const unsigned int var_num = system.variable_number(var_name);
2186 
2187  // The only type of nodal data we can read in from GMV is for
2188  // linear LAGRANGE type elements.
2189  const FEType & fe_type = system.variable_type(var_num);
2190  if ((fe_type.order.get_order() != FIRST) || (fe_type.family != LAGRANGE))
2191  {
2192  libMesh::err << "Only FIRST-order LAGRANGE variables can be read from GMV files. "
2193  << "Skipping variable " << var_name << std::endl;
2194  break;
2195  }
2196 
2197 
2198  // Loop over the stored vector's entries, inserting them into
2199  // the System's solution if appropriate.
2200  for (dof_id_type i=0,
2201  sz = cast_int<dof_id_type>(vec.size());
2202  i != sz; ++i)
2203  {
2204  // Since this var came from a GMV file, the index i corresponds to
2205  // the (single) DOF value of the current variable for node i.
2206  const unsigned int dof_index =
2207  MeshInput<MeshBase>::mesh().node_ptr(i)->dof_number(sys, // system #
2208  var_num, // var #
2209  0); // component #, always zero for LAGRANGE
2210 
2211  // If the dof_index is local to this processor, set the value
2212  if ((dof_index >= system.solution->first_local_index()) &&
2213  (dof_index < system.solution->last_local_index()))
2214  system.solution->set (dof_index, vec[i]);
2215  } // end loop over my GMVIO's copy of the solution
2216 
2217  // Add the most recently copied var to the set of copied vars
2218  vars_copied.insert (var_name);
2219  } // end if (system.has_variable)
2220  } // end for loop over _nodal_data
2221 
2222  // Communicate parallel values before going to the next system.
2223  system.solution->close();
2224  system.update();
2225  } // end loop over all systems
2226 
2227 
2228 
2229  // Warn the user if any GMV variables were not successfully copied over to the EquationSystems object
2230  for (const auto & pr : _nodal_data)
2231  if (vars_copied.find(pr.first) == vars_copied.end())
2232  libMesh::err << "Warning: Variable "
2233  << pr.first
2234  << " was not copied to the EquationSystems object."
2235  << std::endl;
2236 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
OStreamProxy err
FEFamily family
The type of finite element.
Definition: fe_type.h:207
unsigned int n_systems() const
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:201
const T_sys & get_system(std::string_view name) const
std::map< std::string, std::vector< Number > > _nodal_data
Definition: gmv_io.h:236
unsigned int variable_number(std::string_view var) const
Definition: system.C:1557
bool has_variable(std::string_view var) const
Definition: system.C:1550
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1573
libmesh_assert(ctx)
int get_order() const
Explicitly request the order as an int.
Definition: fe_type.h:80
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:493
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2427
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
uint8_t dof_id_type
Definition: id_types.h:67

◆ discontinuous()

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

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

Definition at line 101 of file gmv_io.h.

References _discontinuous.

101 { return _discontinuous; }
bool _discontinuous
Flag to write the mesh as discontinuous patches.
Definition: gmv_io.h:196

◆ get_add_sides()

virtual bool libMesh::MeshOutput< MeshBase >::get_add_sides ( )
inlineprotectedvirtualinherited
Returns
Whether or not added sides are expected to be output, to plot SIDE_DISCONTINUOUS data. Subclasses should override this if they are capable of plotting such data.

Reimplemented in libMesh::ExodusII_IO.

Definition at line 176 of file mesh_output.h.

176 { return false; }

◆ gmv_elem_to_libmesh_elem()

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

Definition at line 2135 of file gmv_io.C.

References _reading_element_map.

Referenced by _read_one_cell().

2136 {
2137  // Erase any whitespace padding in name coming from gmv before performing comparison.
2138  elemname.erase(std::remove_if(elemname.begin(), elemname.end(),
2139  [](unsigned char const c){return std::isspace(c);}),
2140  elemname.end());
2141  return libmesh_map_find(_reading_element_map, elemname);
2142 }
static std::map< std::string, ElemType > _reading_element_map
Static map from string -> ElementType for use during reading.
Definition: gmv_io.h:242

◆ is_parallel_format()

bool libMesh::MeshInput< MeshBase >::is_parallel_format ( ) const
inlineinherited

Returns true iff this mesh file format and input class are parallelized, so that all processors can read their share of the data at once.

Definition at line 87 of file mesh_input.h.

References libMesh::MeshInput< MT >::_is_parallel_format.

87 { return this->_is_parallel_format; }
const bool _is_parallel_format
Flag specifying whether this format is parallel-capable.
Definition: mesh_input.h:130

◆ mesh() [1/2]

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

Definition at line 178 of file mesh_input.h.

Referenced by _read_one_cell(), libMesh::VTKIO::cells_to_vtk(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::Nemesis_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::elements_out(), libMesh::VTKIO::get_local_node_values(), libMesh::ExodusII_IO::get_sideset_data_indices(), libMesh::UNVIO::groups_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::nodes_in(), libMesh::UNVIO::nodes_out(), libMesh::VTKIO::nodes_to_vtk(), libMesh::Nemesis_IO::prepare_to_write_nodal_data(), read(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::VTKIO::read(), libMesh::CheckpointIO::read_bcs(), libMesh::CheckpointIO::read_connectivity(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), libMesh::GmshIO::read_mesh(), libMesh::DynaIO::read_mesh(), libMesh::CheckpointIO::read_nodes(), libMesh::CheckpointIO::read_nodesets(), libMesh::CheckpointIO::read_remote_elem(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::ExodusII_IO::read_sideset_data(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::CheckpointIO::read_subdomain_names(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::ExodusII_IO::write(), write_ascii_new_impl(), write_ascii_old_impl(), write_binary(), write_discontinuous_gmv(), libMesh::Nemesis_IO::write_element_data(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO::write_elemsets(), libMesh::UCDIO::write_header(), libMesh::UCDIO::write_implementation(), libMesh::UCDIO::write_interior_elems(), libMesh::GmshIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::UCDIO::write_nodes(), libMesh::CheckpointIO::write_nodesets(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::ExodusII_IO::write_sideset_data(), libMesh::UCDIO::write_soln(), and libMesh::CheckpointIO::write_subdomain_names().

179 {
180  libmesh_error_msg_if(_obj == nullptr, "ERROR: _obj should not be nullptr!");
181  return *_obj;
182 }
MeshBase * _obj
A pointer to a non-const object object.
Definition: mesh_input.h:123

◆ mesh() [2/2]

const MeshBase & libMesh::MeshOutput< MeshBase >::mesh ( ) const
inlineprotectedinherited

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

References _p_levels.

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

126 { return _p_levels; }
bool _p_levels
Flag to write the mesh p refinement levels.
Definition: gmv_io.h:217

◆ partitioning()

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

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

Definition at line 107 of file gmv_io.h.

References _partitioning.

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

107 { return _partitioning; }
bool _partitioning
Flag to write the mesh partitioning.
Definition: gmv_io.h:201

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

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

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

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

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

Referenced by libMesh::Nemesis_IO::read(), and libMesh::XdrIO::read_header().

101 { this->mesh().set_n_partitions() = n_parts; }
unsigned int & set_n_partitions()
Definition: mesh_base.h:1789

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

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

189 {
190  char c, line[256];
191 
192  while (in.get(c), c==comment_start)
193  in.getline (line, 255);
194 
195  // put back first character of
196  // first non-comment line
197  in.putback (c);
198 }

◆ subdivide_second_order()

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

Flag indicating whether or not to subdivide second order elements.

Definition at line 120 of file gmv_io.h.

References _subdivide_second_order.

Referenced by write_ascii_new_impl(), and write_ascii_old_impl().

120 { return _subdivide_second_order; }
bool _subdivide_second_order
Flag to subdivide second order elements.
Definition: gmv_io.h:212

◆ 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.

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

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

272 {
273  if (this->binary())
274  this->write_binary (fname);
275  else
276  this->write_ascii_old_impl (fname);
277 }
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:1221
bool & binary()
Flag indicating whether or not to write a binary file.
Definition: gmv_io.h:95
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:565

◆ 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.

References _cell_centered_data, std::abs(), libMesh::MeshOutput< MeshBase >::ascii_precision(), libMesh::err, libMesh::index_range(), libMesh::libmesh_assert(), libMesh::make_range(), 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().

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 : make_range(mesh.n_nodes()))
337  out_stream << mesh.point(n)(0) << " ";
338  out_stream << "\n";
339 
340  for (auto n : make_range(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 : make_range(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  libmesh_error_msg_if(this->write_subdomain_id_as_material(),
391  "Not yet supported in GMVIO::write_ascii_new_impl");
392 
393  out_stream << "material "
394  << mesh.n_partitions()
395  // Note: GMV may give you errors like
396  // Error, material for cell 1 is greater than 1
397  // Error, material for cell 2 is greater than 1
398  // Error, material for cell 3 is greater than 1
399  // ... because you put the wrong number of partitions here.
400  // To ensure you write the correct number of materials, call
401  // mesh.recalculate_n_partitions() before writing out the
402  // mesh.
403  // Note: we can't call it now because the Mesh is const here and
404  // it is a non-const function.
405  << " 0\n";
406 
407  for (auto proc : make_range(mesh.n_partitions()))
408  out_stream << "proc_" << proc << "\n";
409 
410  // FIXME - don't we need to use an ElementDefinition here? - RHS
411  for (const auto & elem : mesh.active_element_ptr_range())
412  out_stream << elem->processor_id()+1 << "\n";
413  out_stream << "\n";
414  }
415 
416  // If there are *any* variables at all in the system (including
417  // p level, or arbitrary cell-based data)
418  // to write, the gmv file needs to contain the word "variable"
419  // on a line by itself.
420  bool write_variable = false;
421 
422  // 1.) p-levels
423  if (this->p_levels() && mesh_max_p_level)
424  write_variable = true;
425 
426  // 2.) solution data
427  if ((solution_names != nullptr) && (v != nullptr))
428  write_variable = true;
429 
430  // 3.) cell-centered data
431  if (!(this->_cell_centered_data.empty()))
432  write_variable = true;
433 
434  if (write_variable)
435  out_stream << "variable\n";
436 
437  // if ((this->p_levels() && mesh_max_p_level) ||
438  // ((solution_names != nullptr) && (v != nullptr)))
439  // out_stream << "variable\n";
440 
441  // optionally write the polynomial degree information
442  if (this->p_levels() && mesh_max_p_level)
443  {
444  out_stream << "p_level 0\n";
445 
446  for (const auto & elem : mesh.active_element_ptr_range())
447  {
448  const ElementDefinition & ele = eletypes[elem->type()];
449 
450  // The element mapper better not require any more nodes
451  // than are present in the current element!
452  libmesh_assert_less_equal (ele.node_map.size(), elem->n_nodes());
453 
454  for (std::size_t i=0, enms=ele.node_map.size(); i < enms; i++)
455  out_stream << elem->p_level() << " ";
456  }
457  out_stream << "\n\n";
458  }
459 
460 
461  // optionally write cell-centered data
462  if (!(this->_cell_centered_data.empty()))
463  {
464  for (const auto & [var_name, the_array] : this->_cell_centered_data)
465  {
466  // write out the variable name, followed by a zero.
467  out_stream << var_name << " 0\n";
468 
469  // Loop over active elements, write out cell data. If second-order cells
470  // are split into sub-elements, the sub-elements inherit their parent's
471  // cell-centered data.
472  for (const auto & elem : mesh.active_element_ptr_range())
473  {
474  // Use the element's ID to find the value.
475  libmesh_assert_less (elem->id(), the_array->size());
476  const Real the_value = (*the_array)[elem->id()];
477 
478  if (this->subdivide_second_order())
479  for (unsigned int se=0, nse=elem->n_sub_elem(); se<nse; se++)
480  out_stream << the_value << " ";
481  else
482  out_stream << the_value << " ";
483  }
484 
485  out_stream << "\n\n";
486  }
487  }
488 
489 
490  // optionally write the data
491  if ((solution_names != nullptr) && (v != nullptr))
492  {
493  const unsigned int n_vars = solution_names->size();
494 
495  if (!(v->size() == mesh.n_nodes()*n_vars))
496  libMesh::err << "ERROR: v->size()=" << v->size()
497  << ", mesh.n_nodes()=" << mesh.n_nodes()
498  << ", n_vars=" << n_vars
499  << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
500  << "\n";
501 
502  libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
503 
504  for (unsigned int c=0; c<n_vars; c++)
505  {
506 
507 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
508 
509  // in case of complex data, write _three_ data sets
510  // for each component
511 
512  // this is the real part
513  out_stream << "r_" << (*solution_names)[c] << " 1\n";
514 
515  for (auto n : make_range(mesh.n_nodes()))
516  out_stream << (*v)[n*n_vars + c].real() << " ";
517 
518  out_stream << "\n\n";
519 
520  // this is the imaginary part
521  out_stream << "i_" << (*solution_names)[c] << " 1\n";
522 
523  for (auto n : make_range(mesh.n_nodes()))
524  out_stream << (*v)[n*n_vars + c].imag() << " ";
525 
526  out_stream << "\n\n";
527 
528  // this is the magnitude
529  out_stream << "a_" << (*solution_names)[c] << " 1\n";
530  for (auto n : make_range(mesh.n_nodes()))
531  out_stream << std::abs((*v)[n*n_vars + c]) << " ";
532 
533  out_stream << "\n\n";
534 
535 #else
536 
537  out_stream << (*solution_names)[c] << " 1\n";
538 
539  for (auto n : make_range(mesh.n_nodes()))
540  out_stream << (*v)[n*n_vars + c] << " ";
541 
542  out_stream << "\n\n";
543 
544 #endif
545  }
546 
547  }
548 
549  // If we wrote any variables, we have to close the variable section now
550  if (write_variable)
551  out_stream << "endvars\n";
552 
553 
554  // end of the file
555  out_stream << "\nendgmv\n";
556 
557 #endif
558 }
OStreamProxy err
const MT & mesh() const
Definition: mesh_output.h:259
virtual dof_id_type n_active_elem() const =0
unsigned int & ascii_precision()
Return/set the precision to use when writing ASCII files.
Definition: mesh_output.h:269
This class defines an abstract interface for Mesh output.
Definition: mesh_output.h:53
bool & partitioning()
Flag indicating whether or not to write the partitioning information for the mesh.
Definition: gmv_io.h:107
This is the MeshBase class.
Definition: mesh_base.h:74
bool & write_subdomain_id_as_material()
Flag to write element subdomain_id&#39;s as GMV "materials" instead of element processor_id&#39;s.
Definition: gmv_io.h:114
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:57
std::map< std::string, const std::vector< Real > *> _cell_centered_data
Storage for arbitrary cell-centered data.
Definition: gmv_io.h:225
unsigned int n_vars
bool & p_levels()
Flag indicating whether or not to write p level information for p refined meshes. ...
Definition: gmv_io.h:126
libmesh_assert(ctx)
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:565
unsigned int n_partitions() const
Definition: mesh_base.h:1322
bool & subdivide_second_order()
Flag indicating whether or not to subdivide second order elements.
Definition: gmv_io.h:120
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
virtual const Point & point(const dof_id_type i) const =0
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111
virtual dof_id_type n_nodes() const =0

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

References _cell_centered_data, std::abs(), 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::make_range(), 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().

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

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

References _cell_centered_data, std::abs(), libMesh::err, libMesh::libmesh_assert(), libMesh::make_range(), 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(), and write_subdomain_id_as_material().

Referenced by write(), and write_nodal_data().

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

◆ 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 89 of file mesh_output.C.

References libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), libMesh::libmesh_assert(), and libMesh::out.

Referenced by libMesh::ExodusII_IO::write_timestep_discontinuous().

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

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

References _cell_centered_data, _write_subdomain_id_as_material, std::abs(), 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::make_range(), 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().

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

◆ 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.

References libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), libMesh::libmesh_assert(), and libMesh::out.

Referenced by libMesh::Nemesis_IO::write_timestep(), and libMesh::ExodusII_IO::write_timestep().

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-parallel format, non-renumbered mesh may not have a contiguous
46  // numbering, and that needs to be fixed before we can build a solution vector.
47  if (!_is_parallel_format &&
48  (my_mesh.max_elem_id() != my_mesh.n_elem() ||
49  my_mesh.max_node_id() != my_mesh.n_nodes()))
50  {
51  // If we were allowed to renumber then we should have already
52  // been properly renumbered...
53  libmesh_assert(!my_mesh.allow_renumbering());
54 
55  libmesh_do_once(libMesh::out <<
56  "Warning: This MeshOutput subclass only supports meshes which are contiguously renumbered!"
57  << std::endl;);
58 
59  my_mesh.allow_renumbering(true);
60 
61  my_mesh.renumber_nodes_and_elements();
62 
63  // Not sure what good going back to false will do here, the
64  // renumbering horses have already left the barn...
65  my_mesh.allow_renumbering(false);
66  }
67 
69  {
70  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
71 
72  // Build the list of variable names that will be written.
73  std::vector<std::string> names;
74  es.build_variable_names (names, nullptr, system_names);
75 
76  // Build the nodal solution values & get the variable
77  // names from the EquationSystems object
78  std::vector<Number> soln;
79  es.build_solution_vector (soln, system_names,
80  this->get_add_sides());
81 
82  this->write_nodal_data (fname, soln, names);
83  }
84  else // _is_parallel_format
85  this->write_nodal_data (fname, es, system_names);
86 }
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:109
const MeshBase *const _obj
A pointer to a constant object.
Definition: mesh_output.h:202
const bool _is_parallel_format
Flag specifying whether this format is parallel-capable.
Definition: mesh_output.h:184
libmesh_assert(ctx)
OStreamProxy out
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:193

◆ write_nodal_data() [1/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.

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

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

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 }
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:1221
bool & binary()
Flag indicating whether or not to write a binary file.
Definition: gmv_io.h:95
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:565

◆ 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 149 of file mesh_output.C.

References libMesh::NumericVector< T >::localize().

152 {
153  // This is the fallback implementation for parallel I/O formats that
154  // do not yet implement proper writing in parallel, and instead rely
155  // on the full solution vector being available on all processors.
156  std::vector<Number> soln;
157  parallel_soln.localize(soln);
158  this->write_nodal_data(fname, soln, names);
159 }
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:109
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.

◆ write_nodal_data() [3/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 162 of file mesh_output.C.

References libMesh::EquationSystems::build_parallel_solution_vector(), and libMesh::EquationSystems::build_variable_names().

165 {
166  std::vector<std::string> names;
167  es.build_variable_names (names, nullptr, system_names);
168 
169  std::unique_ptr<NumericVector<Number>> parallel_soln =
170  es.build_parallel_solution_vector(system_names);
171 
172  this->write_nodal_data (fname, *parallel_soln, names);
173 }
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:109

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

121  { 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 114 of file gmv_io.h.

References _write_subdomain_id_as_material.

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

bool _write_subdomain_id_as_material
Flag to write element subdomain_id&#39;s as GMV "materials" instead of element processor_id&#39;s.
Definition: gmv_io.h:207

Member Data Documentation

◆ _binary

bool libMesh::GMVIO::_binary
private

Flag to write binary data.

Definition at line 191 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 225 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 196 of file gmv_io.h.

Referenced by discontinuous().

◆ _is_parallel_format

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

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

◆ _next_elem_id

unsigned int libMesh::GMVIO::_next_elem_id
private

Definition at line 231 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 236 of file gmv_io.h.

Referenced by _read_var(), and copy_nodal_solution().

◆ _p_levels

bool libMesh::GMVIO::_p_levels
private

Flag to write the mesh p refinement levels.

Definition at line 217 of file gmv_io.h.

Referenced by p_levels().

◆ _partitioning

bool libMesh::GMVIO::_partitioning
private

Flag to write the mesh partitioning.

Definition at line 201 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 242 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 193 of file mesh_output.h.

◆ _subdivide_second_order

bool libMesh::GMVIO::_subdivide_second_order
private

Flag to subdivide second order elements.

Definition at line 212 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 207 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

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