libMesh
Classes | 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::VTKIO Class Reference

This class implements reading and writing meshes in the VTK format. More...

#include <vtk_io.h>

Inheritance diagram for libMesh::VTKIO:
[legend]

Classes

struct  ElementMaps
 Helper object that holds a map from VTK to libMesh element types and vice-versa. More...
 

Public Member Functions

 VTKIO (MeshBase &mesh)
 Constructor. More...
 
 VTKIO (const MeshBase &mesh)
 Constructor. 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...
 
virtual void read (const std::string &) override
 This method implements reading a mesh from a specified file in VTK format. More...
 
virtual void write (const std::string &) override
 Output the mesh without solutions to a .pvtu file. More...
 
void set_compression (bool b)
 Setter for compression flag. More...
 
vtkUnstructuredGrid * get_vtk_grid ()
 Get a pointer to the VTK unstructured grid data structure. More...
 
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object. More...
 
virtual void write_discontinuous_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 This method implements writing a mesh with discontinuous data to a specified file where the data is taken from the EquationSystems object. More...
 
virtual void write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &)
 This method may be overridden by "parallel" output formats for writing nodal data. More...
 
virtual void write_nodal_data (const std::string &, const EquationSystems &, const std::set< std::string > *)
 This method should be overridden by "parallel" output formats for writing nodal data. More...
 
virtual void write_nodal_data_discontinuous (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
 This method implements writing a mesh with discontinuous data to a specified file where the nodal data and variables names are provided. More...
 
unsigned intascii_precision ()
 Return/set the precision to use when writing ASCII files. More...
 

Protected Member Functions

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

Protected Attributes

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

Private Member Functions

void nodes_to_vtk ()
 write the nodes from the mesh into a vtkUnstructuredGrid and update the local_node_map. More...
 
void cells_to_vtk ()
 write the cells from the mesh into a vtkUnstructuredGrid More...
 
void node_values_to_vtk (const std::string &name, const std::vector< Real > &local_values)
 write the nodal values of soln to a vtkUnstructuredGrid More...
 
void get_local_node_values (std::vector< Number > &local_values, std::size_t variable, const std::vector< Number > &soln, const std::vector< std::string > &names)
 Extract the values of soln that correspond to the nodes. More...
 

Static Private Member Functions

static ElementMaps build_element_maps ()
 Static function used to construct the _element_maps struct. More...
 

Private Attributes

vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
 Write the system vectors to vtk. More...
 
bool _compress
 Flag to indicate whether the output should be compressed. More...
 
std::map< dof_id_type, dof_id_type_local_node_map
 maps global node id to node id of partition More...
 
MeshBase_obj
 A pointer to a non-const object object. More...
 
const bool _is_parallel_format
 Flag specifying whether this format is parallel-capable. More...
 
unsigned int _ascii_precision
 Precision to use when writing ASCII files. More...
 

Static Private Attributes

static ElementMaps _element_maps = VTKIO::build_element_maps()
 ElementMaps object that is built statically and used by all instances of this class. More...
 

Detailed Description

This class implements reading and writing meshes in the VTK format.

Format description: cf. VTK home page.

This class will not have any functionality unless VTK is detected during configure and hence LIBMESH_HAVE_VTK is defined.

Author
Wout Ruijter
John W. Peterson
Date
2007

Definition at line 60 of file vtk_io.h.

Constructor & Destructor Documentation

◆ VTKIO() [1/2]

libMesh::VTKIO::VTKIO ( MeshBase mesh)
explicit

Constructor.

Takes a writable reference to a mesh object. This is the constructor required to read a mesh.

Definition at line 81 of file vtk_io.C.

81  :
82  MeshInput<MeshBase> (mesh, /*is_parallel_format=*/true),
83  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
84 #ifdef LIBMESH_HAVE_VTK
85  ,_compress(false)
86 #endif
87 {
88 }
89 
90 
91 
92 // Constructor for writing
93 VTKIO::VTKIO (const MeshBase & mesh) :
94  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
95 #ifdef LIBMESH_HAVE_VTK
96  ,_compress(false)
97 #endif
98 {
99 }
100 
101 
102 
103 // Output the mesh without solutions to a .pvtu file
104 void VTKIO::write (const std::string & name)
105 {
106  std::vector<Number> soln;
107  std::vector<std::string> names;
108  this->write_nodal_data(name, soln, names);
109 }

◆ VTKIO() [2/2]

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

Constructor.

Takes a read-only reference to a mesh object. This is the constructor required to write a mesh.

Member Function Documentation

◆ ascii_precision()

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

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

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

Definition at line 257 of file mesh_output.h.

258 {
259  return _ascii_precision;
260 }

◆ build_element_maps()

VTKIO::ElementMaps libMesh::VTKIO::build_element_maps ( )
staticprivate

Static function used to construct the _element_maps struct.

Definition at line 121 of file vtk_io.C.

122 {
123  // Object to be filled up
124  ElementMaps em;
125 
126  em.associate(EDGE2, VTK_LINE);
127  em.associate(EDGE3, VTK_QUADRATIC_EDGE);
128  em.associate(TRI3, VTK_TRIANGLE);
129  em.associate(TRI6, VTK_QUADRATIC_TRIANGLE);
130  em.associate(QUAD4, VTK_QUAD);
131  em.associate(QUAD8, VTK_QUADRATIC_QUAD);
132  em.associate(TET4, VTK_TETRA);
133  em.associate(TET10, VTK_QUADRATIC_TETRA);
134  em.associate(HEX8, VTK_HEXAHEDRON);
135  em.associate(HEX20, VTK_QUADRATIC_HEXAHEDRON);
136  em.associate(HEX27, VTK_TRIQUADRATIC_HEXAHEDRON);
137  em.associate(PRISM6, VTK_WEDGE);
138  em.associate(PRISM15, VTK_QUADRATIC_WEDGE);
139  em.associate(PRISM18, VTK_BIQUADRATIC_QUADRATIC_WEDGE);
140  em.associate(PYRAMID5, VTK_PYRAMID);
141 
142  // VTK_BIQUADRATIC_QUAD has been around since VTK 5.0
143 #if VTK_MAJOR_VERSION > 5 || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION > 0)
144  em.associate(QUAD9, VTK_BIQUADRATIC_QUAD);
145 #endif
146 
147  // TRI3SUBDIVISION is for writing only
148  em.writing_map[TRI3SUBDIVISION] = VTK_TRIANGLE;
149 
150  return em;
151 }

References libMesh::VTKIO::ElementMaps::associate(), libMesh::EDGE2, libMesh::EDGE3, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::PRISM15, libMesh::PRISM18, libMesh::PRISM6, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TET10, libMesh::TET4, libMesh::TRI3, libMesh::TRI3SUBDIVISION, libMesh::TRI6, and libMesh::VTKIO::ElementMaps::writing_map.

◆ cells_to_vtk()

void libMesh::VTKIO::cells_to_vtk ( )
private

write the cells from the mesh into a vtkUnstructuredGrid

Definition at line 425 of file vtk_io.C.

426 {
427  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
428 
429  vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
430  vtkSmartPointer<vtkIdList> pts = vtkSmartPointer<vtkIdList>::New();
431 
432  std::vector<int> types(mesh.n_active_local_elem());
433 
434  vtkSmartPointer<vtkIntArray> elem_id = vtkSmartPointer<vtkIntArray>::New();
435  elem_id->SetName("libmesh_elem_id");
436  elem_id->SetNumberOfComponents(1);
437 
438  vtkSmartPointer<vtkIntArray> subdomain_id = vtkSmartPointer<vtkIntArray>::New();
439  subdomain_id->SetName("subdomain_id");
440  subdomain_id->SetNumberOfComponents(1);
441 
442  vtkSmartPointer<vtkIntArray> elem_proc_id = vtkSmartPointer<vtkIntArray>::New();
443  elem_proc_id->SetName("processor_id");
444  elem_proc_id->SetNumberOfComponents(1);
445 
446  unsigned active_element_counter = 0;
447  for (const auto & elem : mesh.active_local_element_ptr_range())
448  {
449  pts->SetNumberOfIds(elem->n_nodes());
450 
451  // get the connectivity for this element
452  std::vector<dof_id_type> conn;
453  elem->connectivity(0, VTK, conn);
454 
455  for (unsigned int i=0,
456  n_conn = cast_int<unsigned int>(conn.size());
457  i != n_conn; ++i)
458  {
459  // If the node ID is not found in the _local_node_map, we'll
460  // add it to the _vtk_grid. NOTE[JWP]: none of the examples
461  // I have actually enters this section of code...
462  if (_local_node_map.find(conn[i]) == _local_node_map.end())
463  {
464  dof_id_type global_node_id = elem->node_id(i);
465 
466  const Point & the_node = mesh.point(global_node_id);
467 
468  // InsertNextPoint accepts either a double or float array of length 3.
469  double pt[3] = {0., 0., 0.};
470  for (unsigned int d=0; d<LIBMESH_DIM; ++d)
471  pt[d] = the_node(d);
472 
473  // Insert the point into the _vtk_grid
474  vtkIdType local = _vtk_grid->GetPoints()->InsertNextPoint(pt);
475 
476  // Update the _local_node_map with the ID returned by VTK
477  _local_node_map[global_node_id] =
478  cast_int<dof_id_type>(local);
479  }
480 
481  // Otherwise, the node ID was found in the _local_node_map, so
482  // insert it into the vtkIdList.
483  pts->InsertId(i, _local_node_map[conn[i]]);
484  }
485 
486  vtkIdType vtkcellid = cells->InsertNextCell(pts);
487  types[active_element_counter] = cast_int<int>(_element_maps.find(elem->type()));
488 
489  elem_id->InsertTuple1(vtkcellid, elem->id());
490  subdomain_id->InsertTuple1(vtkcellid, elem->subdomain_id());
491  elem_proc_id->InsertTuple1(vtkcellid, elem->processor_id());
492  ++active_element_counter;
493  } // end loop over active elements
494 
495  _vtk_grid->SetCells(types.data(), cells);
496  _vtk_grid->GetCellData()->AddArray(elem_id);
497  _vtk_grid->GetCellData()->AddArray(subdomain_id);
498  _vtk_grid->GetCellData()->AddArray(elem_proc_id);
499 }

References _element_maps, _local_node_map, _vtk_grid, libMesh::MeshBase::active_local_element_ptr_range(), libMesh::VTKIO::ElementMaps::find(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_active_local_elem(), libMesh::MeshBase::point(), and libMesh::VTK.

Referenced by write_nodal_data().

◆ get_local_node_values()

void libMesh::VTKIO::get_local_node_values ( std::vector< Number > &  local_values,
std::size_t  variable,
const std::vector< Number > &  soln,
const std::vector< std::string > &  names 
)
private

Extract the values of soln that correspond to the nodes.

This method overwrites all values in local_values

Definition at line 520 of file vtk_io.C.

524 {
525  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
526  std::size_t num_vars = names.size();
527  dof_id_type num_nodes = mesh.n_nodes();
528 
529  local_values.clear();
530  local_values.resize(_local_node_map.size(), 0.0);
531 
532  // loop over all nodes and get the solution for the current
533  // variable, if the node is in the current partition
534  for (dof_id_type k=0; k<num_nodes; ++k)
535  {
536  std::map<dof_id_type, dof_id_type>::iterator local_node_it = _local_node_map.find(k);
537  if (local_node_it == _local_node_map.end())
538  continue; // not a local node
539 
540  local_values[local_node_it->second] = soln[k*num_vars + variable];
541  }
542 }

References _local_node_map, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), and libMesh::MeshBase::n_nodes().

Referenced by write_nodal_data().

◆ get_vtk_grid()

vtkUnstructuredGrid * libMesh::VTKIO::get_vtk_grid ( )

Get a pointer to the VTK unstructured grid data structure.

Definition at line 368 of file vtk_io.C.

369 {
370  return _vtk_grid;
371 }

References _vtk_grid.

◆ mesh() [1/2]

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

Definition at line 169 of file mesh_input.h.

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

◆ mesh() [2/2]

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

Definition at line 247 of file mesh_output.h.

248 {
250  return *_obj;
251 }

◆ node_values_to_vtk()

void libMesh::VTKIO::node_values_to_vtk ( const std::string &  name,
const std::vector< Real > &  local_values 
)
private

write the nodal values of soln to a vtkUnstructuredGrid

Definition at line 501 of file vtk_io.C.

503 {
504  vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
505  data->SetName(name.c_str());
506 
507  libmesh_assert_equal_to(_local_node_map.size(), local_values.size());
508 
509  // number of local and ghost nodes
510  data->SetNumberOfValues(_local_node_map.size());
511 
512  // copy values into vtk
513  for (auto i : index_range(local_values)) {
514  data->SetValue(i, local_values[i]);
515  }
516 
517  _vtk_grid->GetPointData()->AddArray(data);
518 }

References _local_node_map, _vtk_grid, data, libMesh::index_range(), and libMesh::Quality::name().

Referenced by write_nodal_data().

◆ nodes_to_vtk()

void libMesh::VTKIO::nodes_to_vtk ( )
private

write the nodes from the mesh into a vtkUnstructuredGrid and update the local_node_map.

Definition at line 382 of file vtk_io.C.

383 {
384  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
385 
386  // containers for points and coordinates of points
387  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
388  vtkSmartPointer<vtkDoubleArray> pcoords = vtkSmartPointer<vtkDoubleArray>::New();
389  // if this grid is to be used in VTK then the dimension of the points should be 3
390  pcoords->SetNumberOfComponents(LIBMESH_DIM);
391  pcoords->Allocate(3*mesh.n_local_nodes());
392  points->SetNumberOfPoints(mesh.n_local_nodes()); // it seems that it needs this to prevent a segfault
393 
394  unsigned int local_node_counter = 0;
395 
396  for (const auto & node_ptr : mesh.local_node_ptr_range())
397  {
398  const Node & node = *node_ptr;
399 
400  double pnt[3] = {0, 0, 0};
401  for (unsigned int i=0; i<LIBMESH_DIM; ++i)
402  pnt[i] = node(i);
403 
404  // Fill mapping between global and local node numbers
405  _local_node_map[node.id()] = local_node_counter;
406 
407  // add point
408 #if VTK_VERSION_LESS_THAN(7,1,0)
409  pcoords->InsertNextTupleValue(pnt);
410 #else
411  pcoords->InsertNextTuple(pnt);
412 #endif
413  ++local_node_counter;
414  }
415 
416  // add coordinates to points
417  points->SetData(pcoords);
418 
419  // add points to grid
420  _vtk_grid->SetPoints(points);
421 }

References _local_node_map, _vtk_grid, libMesh::DofObject::id(), libMesh::MeshBase::local_node_ptr_range(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), and libMesh::MeshBase::n_local_nodes().

Referenced by write_nodal_data().

◆ read()

void libMesh::VTKIO::read ( const std::string &  name)
overridevirtual

This method implements reading a mesh from a specified file in VTK format.

FIXME: This is known to write nonsense on AMR meshes and it strips the imaginary parts of complex Numbers.

As with other Mesh IO classes, this interface is still still "available" when !LIBMESH_HAVE_VTK, however, it will throw a runtime error.

This function is not currently used by anything, so it is commented out, and may eventually be removed entirely.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 155 of file vtk_io.C.

156 {
157  // This is a serial-only process for now;
158  // the Mesh should be read on processor 0 and
159  // broadcast later
160  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
161 
162  // Keep track of what kinds of elements this file contains
163  elems_of_dimension.clear();
164  elems_of_dimension.resize(4, false);
165 
166  // Use a typedef, because these names are just crazy
167  typedef vtkSmartPointer<vtkXMLUnstructuredGridReader> MyReader;
168  MyReader reader = MyReader::New();
169 
170  // Pass the filename along to the reader
171  reader->SetFileName(name.c_str());
172 
173  // Force reading
174  reader->Update();
175 
176  // read in the grid
177  _vtk_grid = reader->GetOutput();
178 
179  // Get a reference to the mesh
180  MeshBase & mesh = MeshInput<MeshBase>::mesh();
181 
182  // Clear out any pre-existing data from the Mesh
183  mesh.clear();
184 
185  // Get the number of points from the _vtk_grid object
186  const unsigned int vtk_num_points = static_cast<unsigned int>(_vtk_grid->GetNumberOfPoints());
187 
188  // always numbered nicely so we can loop like this
189  for (unsigned int i=0; i<vtk_num_points; ++i)
190  {
191  // add to the id map
192  // and add the actual point
193  double pnt[3];
194  _vtk_grid->GetPoint(static_cast<vtkIdType>(i), pnt);
195  Point xyz(pnt[0], pnt[1], pnt[2]);
196  mesh.add_point(xyz, i);
197  }
198 
199  // Get the number of cells from the _vtk_grid object
200  const unsigned int vtk_num_cells = static_cast<unsigned int>(_vtk_grid->GetNumberOfCells());
201 
202  vtkSmartPointer<vtkGenericCell> cell = vtkSmartPointer<vtkGenericCell>::New();
203  for (unsigned int i=0; i<vtk_num_cells; ++i)
204  {
205  _vtk_grid->GetCell(i, cell);
206 
207  // Get the libMesh element type corresponding to this VTK element type.
208  ElemType libmesh_elem_type = _element_maps.find(cell->GetCellType());
209  Elem * elem = Elem::build(libmesh_elem_type).release();
210 
211  // get the straightforward numbering from the VTK cells
212  for (auto j : elem->node_index_range())
213  elem->set_node(j) =
214  mesh.node_ptr(cast_int<dof_id_type>(cell->GetPointId(j)));
215 
216  // then get the connectivity
217  std::vector<dof_id_type> conn;
218  elem->connectivity(0, VTK, conn);
219 
220  // then reshuffle the nodes according to the connectivity, this
221  // two-time-assign would evade the definition of the vtk_mapping
222  for (unsigned int j=0,
223  n_conn = cast_int<unsigned int>(conn.size());
224  j != n_conn; ++j)
225  elem->set_node(j) = mesh.node_ptr(conn[j]);
226 
227  elem->set_id(i);
228 
229  elems_of_dimension[elem->dim()] = true;
230 
231  mesh.add_elem(elem);
232  } // end loop over VTK cells
233 
234  // Set the mesh dimension to the largest encountered for an element
235  for (unsigned char i=0; i!=4; ++i)
236  if (elems_of_dimension[i])
238 
239 #if LIBMESH_DIM < 3
240  if (mesh.mesh_dimension() > LIBMESH_DIM)
241  libmesh_error_msg("Cannot open dimension " \
242  << mesh.mesh_dimension() \
243  << " mesh file when configured without " \
244  << mesh.mesh_dimension() \
245  << "D support.");
246 #endif // LIBMESH_DIM < 3
247 }

References _element_maps, _vtk_grid, libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::Elem::build(), libMesh::MeshBase::clear(), libMesh::Elem::connectivity(), libMesh::Elem::dim(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::VTKIO::ElementMaps::find(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshInput< MT >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Quality::name(), libMesh::Elem::node_index_range(), libMesh::MeshBase::node_ptr(), libMesh::DofObject::set_id(), libMesh::MeshBase::set_mesh_dimension(), libMesh::Elem::set_node(), and libMesh::VTK.

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

◆ set_compression()

void libMesh::VTKIO::set_compression ( bool  b)

Setter for compression flag.

Definition at line 375 of file vtk_io.C.

376 {
377  this->_compress = b;
378 }

References _compress.

◆ set_n_partitions()

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

Sets the number of partitions in the mesh.

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

Definition at line 91 of file mesh_input.h.

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

◆ skip_comment_lines()

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

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

Definition at line 179 of file mesh_input.h.

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

◆ write()

virtual void libMesh::VTKIO::write ( const std::string &  )
overridevirtual

Output the mesh without solutions to a .pvtu file.

As with other Mesh IO classes, this interface is still still "available" when !LIBMESH_HAVE_VTK, however, it will throw a runtime error.

Implements libMesh::MeshOutput< MeshBase >.

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

◆ write_discontinuous_equation_systems()

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

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

Definition at line 87 of file mesh_output.C.

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

◆ write_equation_systems()

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

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

Reimplemented in libMesh::NameBasedIO.

Definition at line 31 of file mesh_output.C.

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

◆ write_nodal_data() [1/3]

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

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

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

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

Reimplemented in libMesh::Nemesis_IO.

Definition at line 158 of file mesh_output.C.

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

◆ write_nodal_data() [2/3]

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

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

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

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

Reimplemented in libMesh::Nemesis_IO.

Definition at line 145 of file mesh_output.C.

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

◆ write_nodal_data() [3/3]

void libMesh::VTKIO::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.

As with other Mesh IO classes, this interface is still still "available" when !LIBMESH_HAVE_VTK, however, it will throw a runtime error.

Reimplemented from libMesh::MeshOutput< MeshBase >.

Definition at line 251 of file vtk_io.C.

254 {
255  // Warn that the .pvtu file extension should be used. Paraview
256  // recognizes this, and it works in both serial and parallel. Only
257  // warn about this once.
258  if (fname.substr(fname.rfind("."), fname.size()) != ".pvtu")
259  libmesh_do_once(libMesh::err << "The .pvtu extension should be used when writing VTK files in libMesh.");
260 
261  // If there are variable names being written, the solution vector
262  // should not be empty, it should have been broadcast to all
263  // processors by the MeshOutput base class, since VTK is a parallel
264  // format. Verify this before going further.
265  if (!names.empty() && soln.empty())
266  libmesh_error_msg("Empty soln vector in VTKIO::write_nodal_data().");
267 
268  // Get a reference to the mesh
269  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
270 
271  // we only use Unstructured grids
272  _vtk_grid = vtkSmartPointer<vtkUnstructuredGrid>::New();
273  vtkSmartPointer<vtkXMLPUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
274 #ifdef LIBMESH_HAVE_MPI
275  // Set VTK to the same communicator as libMesh
276  vtkSmartPointer<vtkMPICommunicator> vtk_comm = vtkSmartPointer<vtkMPICommunicator>::New();
277  MPI_Comm mpi_comm = mesh.comm().get();
278  vtkMPICommunicatorOpaqueComm vtk_opaque_comm(&mpi_comm);
279  vtk_comm->InitializeExternal(&vtk_opaque_comm);
280 
281  vtkSmartPointer<vtkMPIController> vtk_mpi_ctrl = vtkSmartPointer<vtkMPIController>::New();
282  vtk_mpi_ctrl->SetCommunicator(vtk_comm);
283 
284  writer->SetController(vtk_mpi_ctrl);
285 #endif
286 
287  // add nodes to the grid and update _local_node_map
288  _local_node_map.clear();
289  this->nodes_to_vtk();
290 
291  // add cells to the grid
292  this->cells_to_vtk();
293 
294  // add nodal solutions to the grid, if solutions are given
295  if (names.size() > 0)
296  {
297  std::size_t num_vars = names.size();
298  std::vector<Number> local_values;
299 
300 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
301  std::vector<Real> local_real_values;
302 #endif
303 
304  for (std::size_t variable=0; variable<num_vars; ++variable)
305  {
306  get_local_node_values(local_values, variable, soln, names);
307 
308 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
309  // write real part
310  local_real_values.resize(local_values.size());
311  std::transform(local_values.begin(), local_values.end(),
312  local_real_values.begin(),
313  [](Number x) { return x.real(); });
314  node_values_to_vtk(names[variable] + "_real", local_real_values);
315 
316  // write imaginary part
317  local_real_values.resize(local_values.size());
318  std::transform(local_values.begin(), local_values.end(),
319  local_real_values.begin(),
320  [](Number x) { return x.imag(); });
321  node_values_to_vtk(names[variable] + "_imag", local_real_values);
322 #else
323  node_values_to_vtk(names[variable], local_values);
324 #endif
325  }
326  }
327 
328  // Tell the writer how many partitions exist and on which processor
329  // we are currently
330  writer->SetNumberOfPieces(mesh.n_processors());
331  writer->SetStartPiece(mesh.processor_id());
332  writer->SetEndPiece(mesh.processor_id());
333 
334  // partitions overlap by one node
335  // FIXME: According to this document
336  // http://paraview.org/Wiki/images/5/51/SC07_tut107_ParaView_Handouts.pdf
337  // the ghosts are cells rather than nodes.
338  writer->SetGhostLevel(1);
339 
340  // VTK 6 replaces SetInput() with SetInputData(). See
341  // http://www.vtk.org/Wiki/VTK/VTK_6_Migration/Replacement_of_SetInput
342  // for the full explanation.
343 #if VTK_VERSION_LESS_THAN(6,0,0)
344  writer->SetInput(_vtk_grid);
345 #else
346  writer->SetInputData(_vtk_grid);
347 #endif
348 
349  writer->SetFileName(fname.c_str());
350  writer->SetDataModeToAscii();
351 
352  // compress the output, if desired (switches also to binary)
353  if (this->_compress)
354  {
355 #if !VTK_VERSION_LESS_THAN(5,6,0)
356  writer->SetCompressorTypeToZLib();
357 #else
358  libmesh_do_once(libMesh::err << "Compression not implemented with old VTK libs!" << std::endl;);
359 #endif
360  }
361 
362  writer->Write();
363 
364 }

References _compress, _local_node_map, _vtk_grid, cells_to_vtk(), libMesh::ParallelObject::comm(), libMesh::err, get_local_node_values(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::ParallelObject::n_processors(), node_values_to_vtk(), nodes_to_vtk(), and libMesh::ParallelObject::processor_id().

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

◆ write_nodal_data_discontinuous()

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

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

Reimplemented in libMesh::ExodusII_IO.

Definition at line 114 of file mesh_output.h.

117  { libmesh_not_implemented(); }

Member Data Documentation

◆ _ascii_precision

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

Precision to use when writing ASCII files.

Definition at line 195 of file mesh_output.h.

◆ _compress

bool libMesh::VTKIO::_compress
private

Flag to indicate whether the output should be compressed.

Definition at line 174 of file vtk_io.h.

Referenced by set_compression(), and write_nodal_data().

◆ _element_maps

VTKIO::ElementMaps libMesh::VTKIO::_element_maps = VTKIO::build_element_maps()
staticprivate

ElementMaps object that is built statically and used by all instances of this class.

Definition at line 215 of file vtk_io.h.

Referenced by cells_to_vtk(), and read().

◆ _is_parallel_format [1/2]

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

Flag specifying whether this format is parallel-capable.

If this is false (default) I/O is only permitted when the mesh has been serialized.

Definition at line 121 of file mesh_input.h.

◆ _is_parallel_format [2/2]

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

Flag specifying whether this format is parallel-capable.

If this is false (default) I/O is only permitted when the mesh has been serialized.

Definition at line 172 of file mesh_output.h.

◆ _local_node_map

std::map<dof_id_type, dof_id_type> libMesh::VTKIO::_local_node_map
private

maps global node id to node id of partition

Definition at line 179 of file vtk_io.h.

Referenced by cells_to_vtk(), get_local_node_values(), node_values_to_vtk(), nodes_to_vtk(), and write_nodal_data().

◆ _obj

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

A pointer to a non-const object object.

This allows us to read the object from file.

Definition at line 114 of file mesh_input.h.

◆ _serial_only_needed_on_proc_0

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

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

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

Definition at line 181 of file mesh_output.h.

◆ _vtk_grid

vtkSmartPointer<vtkUnstructuredGrid> libMesh::VTKIO::_vtk_grid
private

Write the system vectors to vtk.

This function is not currently used by anything, so it is commented out, and may eventually be removed entirely. pointer to the VTK grid. the vtkSmartPointer will automatically initialize the value to null and keep track of reference counting.

Definition at line 169 of file vtk_io.h.

Referenced by cells_to_vtk(), get_vtk_grid(), node_values_to_vtk(), nodes_to_vtk(), read(), and write_nodal_data().

◆ elems_of_dimension

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

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

Definition at line 97 of file mesh_input.h.


The documentation for this class was generated from the following files:
libMesh::MeshOutput< MeshBase >::MeshOutput
MeshOutput(const bool is_parallel_format=false, const bool serial_only_needed_on_proc_0=false)
Default constructor.
Definition: mesh_output.h:207
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::MeshOutput< MeshBase >::_serial_only_needed_on_proc_0
const bool _serial_only_needed_on_proc_0
Flag specifying whether this format can be written by only serializing the mesh to processor zero.
Definition: mesh_output.h:181
libMesh::MeshBase::set_n_partitions
unsigned int & set_n_partitions()
Definition: mesh_base.h:1667
libMesh::VTKIO::cells_to_vtk
void cells_to_vtk()
write the cells from the mesh into a vtkUnstructuredGrid
Definition: vtk_io.C:425
libMesh::VTKIO::VTKIO
VTKIO(MeshBase &mesh)
Constructor.
Definition: vtk_io.C:81
libMesh::VTKIO::nodes_to_vtk
void nodes_to_vtk()
write the nodes from the mesh into a vtkUnstructuredGrid and update the local_node_map.
Definition: vtk_io.C:382
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::MeshBase::point
virtual const Point & point(const dof_id_type i) const =0
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::MeshInput< MeshBase >::_obj
MeshBase * _obj
A pointer to a non-const object object.
Definition: mesh_input.h:114
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::MeshBase::n_local_nodes
dof_id_type n_local_nodes() const
Definition: mesh_base.h:303
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::MeshOutput< MeshBase >::_is_parallel_format
const bool _is_parallel_format
Flag specifying whether this format is parallel-capable.
Definition: mesh_output.h:172
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::VTKIO::_vtk_grid
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Write the system vectors to vtk.
Definition: vtk_io.h:169
libMesh::MeshBase::local_node_ptr_range
virtual SimpleRange< node_iterator > local_node_ptr_range()=0
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::VTKIO::_compress
bool _compress
Flag to indicate whether the output should be compressed.
Definition: vtk_io.h:174
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::VTKIO::get_local_node_values
void get_local_node_values(std::vector< Number > &local_values, std::size_t variable, const std::vector< Number > &soln, const std::vector< std::string > &names)
Extract the values of soln that correspond to the nodes.
Definition: vtk_io.C:520
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::MeshOutput< MeshBase >::write_nodal_data
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: mesh_output.h:105
libMesh::MeshOutput::mesh
const MT & mesh() const
Definition: mesh_output.h:247
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::VTKIO::node_values_to_vtk
void node_values_to_vtk(const std::string &name, const std::vector< Real > &local_values)
write the nodal values of soln to a vtkUnstructuredGrid
Definition: vtk_io.C:501
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::MeshBase::n_active_local_elem
dof_id_type n_active_local_elem() const
Definition: mesh_base.h:420
libMesh::VTKIO::_local_node_map
std::map< dof_id_type, dof_id_type > _local_node_map
maps global node id to node id of partition
Definition: vtk_io.h:179
libMesh::MeshOutput< MeshBase >::write_nodal_data_discontinuous
virtual void write_nodal_data_discontinuous(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
This method implements writing a mesh with discontinuous data to a specified file where the nodal dat...
Definition: mesh_output.h:114
libMesh::VTK
Definition: enum_io_package.h:42
libMesh::VTKIO::ElementMaps::find
vtkIdType find(ElemType libmesh_type)
Definition: vtk_io.h:196
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::MeshInput< MeshBase >::mesh
MeshBase & mesh()
Definition: mesh_input.h:169
libMesh::MeshOutput< MeshBase >::_obj
const MeshBase *const _obj
A pointer to a constant object.
Definition: mesh_output.h:190
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::MeshBase::add_point
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
libMesh::err
OStreamProxy err
libMesh::TRI3SUBDIVISION
Definition: enum_elem_type.h:69
libMesh::MeshBase::clear
virtual void clear()
Deletes all the element and node data that is currently stored.
Definition: mesh_base.C:429
libMesh::MeshBase::set_mesh_dimension
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
Definition: mesh_base.h:218
libMesh::PRISM18
Definition: enum_elem_type.h:52
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::out
OStreamProxy out
libMesh::VTKIO::_element_maps
static ElementMaps _element_maps
ElementMaps object that is built statically and used by all instances of this class.
Definition: vtk_io.h:215
libMesh::MeshOutput< MeshBase >::_ascii_precision
unsigned int _ascii_precision
Precision to use when writing ASCII files.
Definition: mesh_output.h:195
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::MeshInput< MeshBase >::elems_of_dimension
std::vector< bool > elems_of_dimension
A vector of bools describing what dimension elements have been encountered when reading a mesh.
Definition: mesh_input.h:97
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33