Go to the documentation of this file.
19 #ifndef LIBMESH_VTK_IO_H
20 #define LIBMESH_VTK_IO_H
23 #include "libmesh/libmesh_common.h"
24 #include "libmesh/mesh_input.h"
25 #include "libmesh/mesh_output.h"
26 #include "libmesh/utility.h"
28 #ifdef LIBMESH_HAVE_VTK
30 #include "libmesh/ignore_warnings.h"
32 #include "vtkSmartPointer.h"
33 #include "libmesh/restore_warnings.h"
41 class vtkUnstructuredGrid;
93 const std::vector<Number> &,
94 const std::vector<std::string> &)
override;
104 virtual void read (
const std::string &)
override;
113 virtual void write (
const std::string &)
override;
115 #ifdef LIBMESH_HAVE_VTK
143 const std::vector<Real> & local_values);
151 std::size_t variable,
152 const std::vector<Number> & soln,
153 const std::vector<std::string> & names);
198 return libmesh_map_find(
writing_map, libmesh_type);
230 #endif // LIBMESH_VTK_IO_H
virtual void read(const std::string &) override
This method implements reading a mesh from a specified file in VTK format.
void cells_to_vtk()
write the cells from the mesh into a vtkUnstructuredGrid
std::map< vtkIdType, ElemType > reading_map
VTKIO(MeshBase &mesh)
Constructor.
void nodes_to_vtk()
write the nodes from the mesh into a vtkUnstructuredGrid and update the local_node_map.
std::map< ElemType, vtkIdType > writing_map
The libMesh namespace provides an interface to certain functionality in the library.
vtkUnstructuredGrid * get_vtk_grid()
Get a pointer to the VTK unstructured grid data structure.
This class implements reading and writing meshes in the VTK format.
This is the MeshBase class.
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Write the system vectors to vtk.
bool _compress
Flag to indicate whether the output should be compressed.
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.
void associate(ElemType libmesh_type, vtkIdType vtk_type)
This class defines an abstract interface for Mesh output.
void set_compression(bool b)
Setter for compression flag.
void node_values_to_vtk(const std::string &name, const std::vector< Real > &local_values)
write the nodal values of soln to a vtkUnstructuredGrid
Helper object that holds a map from VTK to libMesh element types and vice-versa.
std::map< dof_id_type, dof_id_type > _local_node_map
maps global node id to node id of partition
static ElementMaps build_element_maps()
Static function used to construct the _element_maps struct.
vtkIdType find(ElemType libmesh_type)
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 va...
virtual void write(const std::string &) override
Output the mesh without solutions to a .pvtu file.
ElemType find(vtkIdType vtk_type)
static ElementMaps _element_maps
ElementMaps object that is built statically and used by all instances of this class.
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
ElemType
Defines an enum for geometric element types.