libMesh::VTKIO Class Reference

#include <vtk_io.h>

Inheritance diagram for libMesh::VTKIO:

Classes

struct  ElementMaps
 

Public Member Functions

 VTKIO (MeshBase &mesh)
 
 VTKIO (const MeshBase &mesh)
 
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
 
virtual void read (const std::string &) override
 
virtual void write (const std::string &) override
 
void set_compression (bool b)
 
vtkUnstructuredGrid * get_vtk_grid ()
 
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 
virtual void write_discontinuous_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
 
virtual void write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &)
 
virtual void write_nodal_data_discontinuous (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
 
unsigned int & ascii_precision ()
 

Protected Member Functions

MeshBasemesh ()
 
void set_n_partitions (unsigned int n_parts)
 
void skip_comment_lines (std::istream &in, const char comment_start)
 
const MeshBasemesh () const
 

Protected Attributes

std::vector< bool > elems_of_dimension
 
const bool _is_parallel_format
 
const bool _serial_only_needed_on_proc_0
 

Private Member Functions

void nodes_to_vtk ()
 
void cells_to_vtk ()
 

Static Private Member Functions

static ElementMaps build_element_maps ()
 

Private Attributes

vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
 
bool _compress
 
std::map< dof_id_type, dof_id_type_local_node_map
 

Static Private Attributes

static ElementMaps _element_maps = VTKIO::build_element_maps()
 

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 59 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 75 of file vtk_io.C.

75  :
76  MeshInput<MeshBase> (mesh, /*is_parallel_format=*/true),
77  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
78 #ifdef LIBMESH_HAVE_VTK
79  ,_compress(false)
80 #endif
81 {
82 }
83 
84 
85 
86 // Constructor for writing
87 VTKIO::VTKIO (const MeshBase & mesh) :
88  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
89 #ifdef LIBMESH_HAVE_VTK
90  ,_compress(false)
91 #endif
92 {
93 }
94 
95 
96 
97 // Output the mesh without solutions to a .pvtu file
98 void VTKIO::write (const std::string & name)
99 {
100  std::vector<Number> soln;
101  std::vector<std::string> names;
102  this->write_nodal_data(name, soln, names);
103 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
bool _compress
Definition: vtk_io.h:156
virtual void write(const std::string &) override
MeshOutput(const bool is_parallel_format=false, const bool serial_only_needed_on_proc_0=false)
Definition: mesh_output.h:194
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
Definition: vtk_io.C:245
VTKIO(MeshBase &mesh)
Definition: vtk_io.C:75

◆ 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>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.

Definition at line 244 of file mesh_output.h.

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

245 {
246  return _ascii_precision;
247 }

◆ build_element_maps()

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

Static function used to construct the _element_maps struct.

Definition at line 115 of file vtk_io.C.

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.

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

◆ cells_to_vtk()

void libMesh::VTKIO::cells_to_vtk ( )
private

write the cells from the mesh into a vtkUnstructuredGrid

Definition at line 407 of file vtk_io.C.

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

408 {
409  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
410 
411  vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
412  vtkSmartPointer<vtkIdList> pts = vtkSmartPointer<vtkIdList>::New();
413 
414  std::vector<int> types(mesh.n_active_local_elem());
415 
416  vtkSmartPointer<vtkIntArray> elem_id = vtkSmartPointer<vtkIntArray>::New();
417  elem_id->SetName("libmesh_elem_id");
418  elem_id->SetNumberOfComponents(1);
419 
420  vtkSmartPointer<vtkIntArray> subdomain_id = vtkSmartPointer<vtkIntArray>::New();
421  subdomain_id->SetName("subdomain_id");
422  subdomain_id->SetNumberOfComponents(1);
423 
424  vtkSmartPointer<vtkIntArray> elem_proc_id = vtkSmartPointer<vtkIntArray>::New();
425  elem_proc_id->SetName("processor_id");
426  elem_proc_id->SetNumberOfComponents(1);
427 
428  unsigned active_element_counter = 0;
429  for (const auto & elem : mesh.active_local_element_ptr_range())
430  {
431  pts->SetNumberOfIds(elem->n_nodes());
432 
433  // get the connectivity for this element
434  std::vector<dof_id_type> conn;
435  elem->connectivity(0, VTK, conn);
436 
437  for (unsigned int i=0,
438  n_conn = cast_int<unsigned int>(conn.size());
439  i != n_conn; ++i)
440  {
441  // If the node ID is not found in the _local_node_map, we'll
442  // add it to the _vtk_grid. NOTE[JWP]: none of the examples
443  // I have actually enters this section of code...
444  if (_local_node_map.find(conn[i]) == _local_node_map.end())
445  {
446  dof_id_type global_node_id = elem->node_id(i);
447 
448  const Point & the_node = mesh.point(global_node_id);
449 
450  // InsertNextPoint accepts either a double or float array of length 3.
451  double pt[3] = {0., 0., 0.};
452  for (unsigned int d=0; d<LIBMESH_DIM; ++d)
453  pt[d] = the_node(d);
454 
455  // Insert the point into the _vtk_grid
456  vtkIdType local = _vtk_grid->GetPoints()->InsertNextPoint(pt);
457 
458  // Update the _local_node_map with the ID returned by VTK
459  _local_node_map[global_node_id] =
460  cast_int<dof_id_type>(local);
461  }
462 
463  // Otherwise, the node ID was found in the _local_node_map, so
464  // insert it into the vtkIdList.
465  pts->InsertId(i, _local_node_map[conn[i]]);
466  }
467 
468  vtkIdType vtkcellid = cells->InsertNextCell(pts);
469  types[active_element_counter] = cast_int<int>(_element_maps.find(elem->type()));
470 
471  elem_id->InsertTuple1(vtkcellid, elem->id());
472  subdomain_id->InsertTuple1(vtkcellid, elem->subdomain_id());
473  elem_proc_id->InsertTuple1(vtkcellid, elem->processor_id());
474  ++active_element_counter;
475  } // end loop over active elements
476 
477  _vtk_grid->SetCells(types.data(), cells);
478  _vtk_grid->GetCellData()->AddArray(elem_id);
479  _vtk_grid->GetCellData()->AddArray(subdomain_id);
480  _vtk_grid->GetCellData()->AddArray(elem_proc_id);
481 }
const MT & mesh() const
Definition: mesh_output.h:234
std::map< dof_id_type, dof_id_type > _local_node_map
Definition: vtk_io.h:161
vtkIdType find(ElemType libmesh_type)
Definition: vtk_io.h:178
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Definition: vtk_io.h:151
static ElementMaps _element_maps
Definition: vtk_io.h:207
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
dof_id_type n_active_local_elem() const
Definition: mesh_base.h:403
virtual const Point & point(const dof_id_type i) const =0
uint8_t dof_id_type
Definition: id_types.h:64

◆ get_vtk_grid()

vtkUnstructuredGrid * libMesh::VTKIO::get_vtk_grid ( )

Get a pointer to the VTK unstructured grid data structure.

Definition at line 350 of file vtk_io.C.

References _vtk_grid.

351 {
352  return _vtk_grid;
353 }
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Definition: vtk_io.h:151

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

Referenced by libMesh::GMVIO::_read_one_cell(), cells_to_vtk(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::elements_out(), libMesh::UNVIO::groups_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::nodes_in(), libMesh::UNVIO::nodes_out(), nodes_to_vtk(), libMesh::Nemesis_IO::prepare_to_write_nodal_data(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), libMesh::GMVIO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), read(), libMesh::CheckpointIO::read_bcs(), libMesh::CheckpointIO::read_connectivity(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), libMesh::GmshIO::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::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::CheckpointIO::read_subdomain_names(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::ExodusII_IO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::Nemesis_IO::write_element_data(), libMesh::ExodusII_IO::write_element_data(), libMesh::UCDIO::write_header(), libMesh::UCDIO::write_implementation(), libMesh::UCDIO::write_interior_elems(), libMesh::GmshIO::write_mesh(), write_nodal_data(), libMesh::UCDIO::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::UCDIO::write_soln(), and libMesh::CheckpointIO::write_subdomain_names().

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

◆ mesh() [2/2]

◆ nodes_to_vtk()

void libMesh::VTKIO::nodes_to_vtk ( )
private

write the nodes from the mesh into a vtkUnstructuredGrid

Definition at line 364 of file vtk_io.C.

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

365 {
366  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
367 
368  // containers for points and coordinates of points
369  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
370  vtkSmartPointer<vtkDoubleArray> pcoords = vtkSmartPointer<vtkDoubleArray>::New();
371  // if this grid is to be used in VTK then the dimension of the points should be 3
372  pcoords->SetNumberOfComponents(LIBMESH_DIM);
373  pcoords->Allocate(3*mesh.n_local_nodes());
374  points->SetNumberOfPoints(mesh.n_local_nodes()); // it seems that it needs this to prevent a segfault
375 
376  unsigned int local_node_counter = 0;
377 
378  for (const auto & node_ptr : mesh.local_node_ptr_range())
379  {
380  const Node & node = *node_ptr;
381 
382  double pnt[3] = {0, 0, 0};
383  for (unsigned int i=0; i<LIBMESH_DIM; ++i)
384  pnt[i] = node(i);
385 
386  // Fill mapping between global and local node numbers
387  _local_node_map[node.id()] = local_node_counter;
388 
389  // add point
390 #if VTK_VERSION_LESS_THAN(7,1,0)
391  pcoords->InsertNextTupleValue(pnt);
392 #else
393  pcoords->InsertNextTuple(pnt);
394 #endif
395  ++local_node_counter;
396  }
397 
398  // add coordinates to points
399  points->SetData(pcoords);
400 
401  // add points to grid
402  _vtk_grid->SetPoints(points);
403 }
const MT & mesh() const
Definition: mesh_output.h:234
std::map< dof_id_type, dof_id_type > _local_node_map
Definition: vtk_io.h:161
virtual SimpleRange< node_iterator > local_node_ptr_range()=0
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Definition: vtk_io.h:151
dof_id_type n_local_nodes() const
Definition: mesh_base.h:286

◆ read()

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

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

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

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

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

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::Elem::n_nodes(), libMesh::Quality::name(), 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().

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

◆ set_compression()

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

Setter for compression flag

Definition at line 357 of file vtk_io.C.

References _compress.

358 {
359  this->_compress = b;
360 }
bool _compress
Definition: vtk_io.h:156

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

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

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

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

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

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

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

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

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

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

◆ 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_parallel_solution_vector(), libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), 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-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 
67  MeshSerializer serialize(const_cast<MT &>(*_obj), !_is_parallel_format, _serial_only_needed_on_proc_0);
68 
69  // Build the list of variable names that will be written.
70  std::vector<std::string> names;
71  es.build_variable_names (names, nullptr, system_names);
72 
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  {
84  std::unique_ptr<NumericVector<Number>> parallel_soln =
85  es.build_parallel_solution_vector(system_names);
86 
87  this->write_nodal_data (fname, *parallel_soln, names);
88  }
89 }
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
Definition: mesh_output.h:105
const MeshBase *const _obj
Definition: mesh_output.h:177
OStreamProxy out(std::cout)

◆ write_nodal_data() [1/2]

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 245 of file vtk_io.C.

References _compress, _local_node_map, _vtk_grid, cells_to_vtk(), data, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_nodes(), and nodes_to_vtk().

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

248 {
249  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
250 
251  // Warn that the .pvtu file extension should be used. Paraview
252  // recognizes this, and it works in both serial and parallel. Only
253  // warn about this once.
254  if (fname.substr(fname.rfind("."), fname.size()) != ".pvtu")
255  libmesh_do_once(libMesh::err << "The .pvtu extension should be used when writing VTK files in libMesh.");
256 
257  // If there are variable names being written, the solution vector
258  // should not be empty, it should have been broadcast to all
259  // processors by the MeshOutput base class, since VTK is a parallel
260  // format. Verify this before going further.
261  if (!names.empty() && soln.empty())
262  libmesh_error_msg("Empty soln vector in VTKIO::write_nodal_data().");
263 
264  // we only use Unstructured grids
265  _vtk_grid = vtkSmartPointer<vtkUnstructuredGrid>::New();
266  vtkSmartPointer<vtkXMLPUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
267 
268  // add nodes to the grid and update _local_node_map
269  _local_node_map.clear();
270  this->nodes_to_vtk();
271 
272  // add cells to the grid
273  this->cells_to_vtk();
274 
275  // add nodal solutions to the grid, if solutions are given
276  if (names.size() > 0)
277  {
278  std::size_t num_vars = names.size();
279  dof_id_type num_nodes = mesh.n_nodes();
280 
281  for (std::size_t variable=0; variable<num_vars; ++variable)
282  {
283  vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
284  data->SetName(names[variable].c_str());
285 
286  // number of local and ghost nodes
287  data->SetNumberOfValues(_local_node_map.size());
288 
289  // loop over all nodes and get the solution for the current
290  // variable, if the node is in the current partition
291  for (dof_id_type k=0; k<num_nodes; ++k)
292  {
293  std::map<dof_id_type, dof_id_type>::iterator local_node_it = _local_node_map.find(k);
294  if (local_node_it == _local_node_map.end())
295  continue; // not a local node
296 
297 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
298  libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n"
299  << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT
300  << std::endl);
301  data->SetValue(local_node_it->second, soln[k*num_vars + variable].real());
302 #else
303  data->SetValue(local_node_it->second, soln[k*num_vars + variable]);
304 #endif
305  }
306  _vtk_grid->GetPointData()->AddArray(data);
307  }
308  }
309 
310  // Tell the writer how many partitions exist and on which processor
311  // we are currently
312  writer->SetNumberOfPieces(MeshOutput<MeshBase>::mesh().n_processors());
313  writer->SetStartPiece(MeshOutput<MeshBase>::mesh().processor_id());
314  writer->SetEndPiece(MeshOutput<MeshBase>::mesh().processor_id());
315 
316  // partitions overlap by one node
317  // FIXME: According to this document
318  // http://paraview.org/Wiki/images/5/51/SC07_tut107_ParaView_Handouts.pdf
319  // the ghosts are cells rather than nodes.
320  writer->SetGhostLevel(1);
321 
322  // VTK 6 replaces SetInput() with SetInputData(). See
323  // http://www.vtk.org/Wiki/VTK/VTK_6_Migration/Replacement_of_SetInput
324  // for the full explanation.
325 #if VTK_VERSION_LESS_THAN(6,0,0)
326  writer->SetInput(_vtk_grid);
327 #else
328  writer->SetInputData(_vtk_grid);
329 #endif
330 
331  writer->SetFileName(fname.c_str());
332  writer->SetDataModeToAscii();
333 
334  // compress the output, if desired (switches also to binary)
335  if (this->_compress)
336  {
337 #if !VTK_VERSION_LESS_THAN(5,6,0)
338  writer->SetCompressorTypeToZLib();
339 #else
340  libmesh_do_once(libMesh::err << "Compression not implemented with old VTK libs!" << std::endl;);
341 #endif
342  }
343 
344  writer->Write();
345 
346 }
const MT & mesh() const
Definition: mesh_output.h:234
std::map< dof_id_type, dof_id_type > _local_node_map
Definition: vtk_io.h:161
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Definition: vtk_io.h:151
bool _compress
Definition: vtk_io.h:156
void nodes_to_vtk()
Definition: vtk_io.C:364
OStreamProxy err(std::cerr)
void cells_to_vtk()
Definition: vtk_io.C:407
IterBase * data
virtual dof_id_type n_nodes() const =0
uint8_t dof_id_type
Definition: id_types.h:64

◆ write_nodal_data() [2/2]

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 should 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 150 of file mesh_output.C.

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

153 {
154  // This is the fallback implementation for parallel I/O formats that
155  // do not yet implement proper writing in parallel, and instead rely
156  // on the full solution vector being available on all processors.
157  std::vector<Number> soln;
158  parallel_soln.localize(soln);
159  this->write_nodal_data(fname, soln, names);
160 }
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
Definition: mesh_output.h:105
virtual void localize(std::vector< T > &v_local) const =0

◆ 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

◆ _compress

bool libMesh::VTKIO::_compress
private

Flag to indicate whether the output should be compressed

Definition at line 156 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 207 of file vtk_io.h.

Referenced by cells_to_vtk(), and read().

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

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

◆ _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 161 of file vtk_io.h.

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

◆ _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 168 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 151 of file vtk_io.h.

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

◆ elems_of_dimension


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