libMesh
|
This class implements reading and writing meshes in the VTK format. More...
#include <vtk_io.h>
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 int & | ascii_precision () |
Return/set the precision to use when writing ASCII files. More... | |
Protected Member Functions | |
MeshBase & | mesh () |
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 MeshBase & | mesh () 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... | |
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.
|
explicit |
|
explicit |
Constructor.
Takes a read-only reference to a mesh object. This is the constructor required to write a mesh.
|
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.
|
staticprivate |
Static function used to construct the _element_maps struct.
Definition at line 121 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.
|
private |
write the cells from the mesh into a vtkUnstructuredGrid
Definition at line 425 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().
|
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.
References _local_node_map, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), and libMesh::MeshBase::n_nodes().
Referenced by write_nodal_data().
vtkUnstructuredGrid * libMesh::VTKIO::get_vtk_grid | ( | ) |
|
inlineprotectedinherited |
Definition at line 169 of file mesh_input.h.
|
inlineprotectedinherited |
Definition at line 247 of file mesh_output.h.
|
private |
write the nodal values of soln to a vtkUnstructuredGrid
Definition at line 501 of file vtk_io.C.
References _local_node_map, _vtk_grid, data, libMesh::index_range(), and libMesh::Quality::name().
Referenced by write_nodal_data().
|
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.
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().
|
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.
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().
void libMesh::VTKIO::set_compression | ( | bool | b | ) |
|
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.
|
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.
|
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().
|
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.
|
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.
|
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.
|
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.
|
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.
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().
|
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.
|
privateinherited |
Precision to use when writing ASCII files.
Definition at line 195 of file mesh_output.h.
|
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().
|
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().
|
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.
|
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.
|
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().
|
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.
|
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.
|
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().
|
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.