OGS
MeshLib::IO::VtuInterface Class Referencefinal

Detailed Description

Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures. This class is currently not inherited from Writer because VTK will implement writing to a string from 6.2 onwards.

Definition at line 40 of file VtuInterface.h.

#include <VtuInterface.h>

Collaboration diagram for MeshLib::IO::VtuInterface:
[legend]

Public Member Functions

 VtuInterface (const MeshLib::Mesh *mesh, int dataMode=vtkXMLWriter::Appended, bool compress=false)
 Provide the mesh to write and set if compression should be used.
 
bool writeToFile (std::filesystem::path const &file_path)
 
template<typename UnstructuredGridWriter >
bool writeVTU (std::string const &file_name, const int num_partitions=1, const int rank=1)
 

Static Public Member Functions

static MeshLib::MeshreadVTUFile (std::string const &file_name, bool const compute_element_neighbors=false)
 
static vtkSmartPointer< vtkUnstructuredGrid > readVtuFileToVtkUnstructuredGrid (std::string const &file_name)
 
static MeshLib::MeshreadVTKFile (std::string const &file_name, bool const compute_element_neighbors=false)
 

Private Attributes

const MeshLib::Mesh_mesh
 
int _data_mode
 
bool _use_compressor
 

Constructor & Destructor Documentation

◆ VtuInterface()

MeshLib::IO::VtuInterface::VtuInterface ( const MeshLib::Mesh * mesh,
int dataMode = vtkXMLWriter::Appended,
bool compress = false )
explicit

Provide the mesh to write and set if compression should be used.

Definition at line 40 of file VtuInterface.cpp.

42 : _mesh(mesh), _data_mode(dataMode), _use_compressor(compress)
43{
44 if (_data_mode == vtkXMLWriter::Ascii && compress)
45 {
46 WARN("Ascii data cannot be compressed, ignoring compression flag.");
47 }
48}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:42
const MeshLib::Mesh * _mesh

References _data_mode, and WARN().

Member Function Documentation

◆ readVTKFile()

MeshLib::Mesh * MeshLib::IO::VtuInterface::readVTKFile ( std::string const & file_name,
bool const compute_element_neighbors = false )
static

Read an unstructured grid from a legacy VTK file. Other data structures such as PolyData are ignored.

Returns
The converted mesh or a nullptr if reading failed

Definition at line 96 of file VtuInterface.cpp.

98{
99 if (!BaseLib::IsFileExisting(file_name))
100 {
101 ERR("File '{:s}' does not exist.", file_name);
102 return nullptr;
103 }
104
105 vtkSmartPointer<vtkGenericDataObjectReader> reader =
106 vtkSmartPointer<vtkGenericDataObjectReader>::New();
107 reader->SetFileName(file_name.c_str());
108 reader->Update();
109
110 // check for unstructured grid
111 if (reader->ReadOutputType() != 4)
112 {
113 ERR("Only VTK-files with dataset type \"Unstructured Grid\" are "
114 "currently supported.");
115 return nullptr;
116 }
117
118 reader->ReadAllFieldsOn();
119 reader->ReadAllScalarsOn();
120 vtkUnstructuredGrid* vtkGrid = reader->GetUnstructuredGridOutput();
121 if (vtkGrid->GetNumberOfPoints() == 0)
122 {
123 ERR("Mesh '{:s}' contains zero points.", file_name);
124 return nullptr;
125 }
126
127 std::string const mesh_name(
130 vtkGrid, compute_element_neighbors, mesh_name);
131}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:48
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, bool const compute_element_neighbors=false, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
Definition FileTools.cpp:48
std::string extractBaseNameWithoutExtension(std::string const &pathname)

References MeshLib::VtkMeshConverter::convertUnstructuredGrid(), ERR(), BaseLib::extractBaseNameWithoutExtension(), and BaseLib::IsFileExisting().

Referenced by anonymous_namespace{readMeshFromFile.cpp}::readMeshFromFileSerial().

◆ readVTUFile()

MeshLib::Mesh * MeshLib::IO::VtuInterface::readVTUFile ( std::string const & file_name,
bool const compute_element_neighbors = false )
static

Read an unstructured grid from a VTU file.

Returns
The converted mesh or a nullptr if reading failed

Definition at line 80 of file VtuInterface.cpp.

82{
83 vtkSmartPointer<vtkUnstructuredGrid> vtkGrid =
85 if (vtkGrid == nullptr)
86 {
87 return nullptr;
88 }
89
90 std::string const mesh_name(
93 vtkGrid, compute_element_neighbors, mesh_name);
94}
static vtkSmartPointer< vtkUnstructuredGrid > readVtuFileToVtkUnstructuredGrid(std::string const &file_name)

References MeshLib::VtkMeshConverter::convertUnstructuredGrid(), BaseLib::extractBaseNameWithoutExtension(), and readVtuFileToVtkUnstructuredGrid().

Referenced by OGSFileConverter::convertVTU2MSH(), main(), and anonymous_namespace{readMeshFromFile.cpp}::readMeshFromFileSerial().

◆ readVtuFileToVtkUnstructuredGrid()

vtkSmartPointer< vtkUnstructuredGrid > MeshLib::IO::VtuInterface::readVtuFileToVtkUnstructuredGrid ( std::string const & file_name)
static

Definition at line 51 of file VtuInterface.cpp.

52{
53 if (!BaseLib::IsFileExisting(file_name))
54 {
55 ERR("File '{:s}' does not exist.", file_name);
56 return nullptr;
57 }
58
59 vtkSmartPointer<vtkXMLUnstructuredGridReader> reader =
60 vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
61 reader->SetFileName(file_name.c_str());
62 {
63 // Reading the VTU files can trigger floating point exceptions. Because
64 // we are not debugging VTK (or other libraries) at this point, the
65 // floating point exceptions are temporarily disabled and are restored
66 // at the end of the scope.
67 [[maybe_unused]] DisableFPE disable_fpe;
68 reader->Update();
69 }
70
71 vtkUnstructuredGrid* vtkGrid = reader->GetOutput();
72 if (vtkGrid->GetNumberOfPoints() == 0)
73 {
74 ERR("Mesh '{:s}' contains zero points.", file_name);
75 return nullptr;
76 }
77 return vtkGrid;
78}

References ERR(), and BaseLib::IsFileExisting().

Referenced by main(), readGrid(), and readVTUFile().

◆ writeToFile()

bool MeshLib::IO::VtuInterface::writeToFile ( std::filesystem::path const & file_path)

Writes the given mesh to file.

Returns
True on success, false on error

Definition at line 156 of file VtuInterface.cpp.

157{
158#ifdef USE_PETSC
160 if (mpi.size == 1)
161 {
162 return writeVTU<vtkXMLUnstructuredGridWriter>(file_path.string());
163 }
164 auto const vtu_file_name =
166 return writeVTU<vtkXMLPUnstructuredGridWriter>(vtu_file_name + ".pvtu",
167 mpi.size, mpi.rank);
168#endif
169 return writeVTU<vtkXMLUnstructuredGridWriter>(file_path.string());
170}
bool writeVTU(std::string const &file_name, const int num_partitions=1, const int rank=1)
std::string getVtuFileNameForPetscOutputWithoutExtension(std::string const &file_name)

References MeshLib::IO::getVtuFileNameForPetscOutputWithoutExtension(), BaseLib::MPI::Mpi::rank, BaseLib::MPI::Mpi::size, and writeVTU().

Referenced by SaveMeshDialog::accept(), convert(), OGSFileConverter::convertMSH2VTU(), main(), ProcessLib::outputMeshVtk(), writeBoundary(), writeDataToMesh(), writeMeshOutput(), MeshLib::IO::writeMeshToFile(), and MeshLib::IO::writeVtu().

◆ writeVTU()

template<typename UnstructuredGridWriter >
bool MeshLib::IO::VtuInterface::writeVTU ( std::string const & file_name,
const int num_partitions = 1,
const int rank = 1 )

Writes the given mesh to vtu file.

Parameters
file_nameFile name.
num_partitionsNumber of partitions to be merged.
rankthe rank of the mpi process.
Returns
True on success, false on error

Definition at line 40 of file VtuInterface-impl.h.

43{
44 if (!_mesh)
45 {
46 ERR("VtuInterface::write(): No mesh specified.");
47 return false;
48 }
49
50#ifdef USE_PETSC
51 if (_mesh->getProperties().existsPropertyVector<unsigned char>(
52 "vtkGhostType", MeshLib::MeshItemType::Cell, 1))
53 {
54 auto* ghost_cell_property =
55 _mesh->getProperties().getPropertyVector<unsigned char>(
56 "vtkGhostType", MeshLib::MeshItemType::Cell, 1);
57 if (ghost_cell_property)
58 {
60 ghost_cell_property)
61 ->is_for_output = true;
62 }
63 }
64 else
65 {
66 DBUG("No vtkGhostType data in mesh '{}'.", _mesh->getName());
67 }
68#endif
69
70 vtkNew<MeshLib::VtkMappedMeshSource> vtkSource;
71 vtkSource->SetMesh(_mesh);
72
73 vtkSmartPointer<UnstructuredGridWriter> vtuWriter =
74 vtkSmartPointer<UnstructuredGridWriter>::New();
75
76 vtkSource->Update();
77 vtuWriter->SetInputData(vtkSource->GetOutput());
78
80 {
81 vtuWriter->SetCompressorTypeToZLib();
82 }
83 else
84 {
85 vtuWriter->SetCompressorTypeToNone();
86 }
87
88 vtuWriter->SetDataMode(_data_mode);
89 if (_data_mode == vtkXMLWriter::Appended)
90 {
91 vtuWriter->SetEncodeAppendedData(1);
92 }
93 if (_data_mode == vtkXMLWriter::Ascii)
94 {
95 vtkSource->Update();
96 vtkSmartPointer<vtkUnstructuredGrid> tempGrid =
97 vtkSmartPointer<vtkUnstructuredGrid>::New();
98 tempGrid->DeepCopy(vtkSource->GetOutput());
99 vtuWriter->SetInputDataObject(tempGrid);
100 }
101
102 vtuWriter->SetFileName(file_name.c_str());
103
104#ifdef USE_PETSC
105 if constexpr (std::is_same_v<UnstructuredGridWriter,
106 vtkXMLPUnstructuredGridWriter>)
107 {
108 // Set the writer controller to same communicator as OGS
109 vtkSmartPointer<vtkMPICommunicator> vtk_comm =
110 vtkSmartPointer<vtkMPICommunicator>::New();
111 MPI_Comm mpi_comm = BaseLib::MPI::OGS_COMM_WORLD;
112 vtkMPICommunicatorOpaqueComm vtk_opaque_comm(&mpi_comm);
113 vtk_comm->InitializeExternal(&vtk_opaque_comm);
114
115 vtkSmartPointer<vtkMPIController> vtk_mpi_ctrl =
116 vtkSmartPointer<vtkMPIController>::New();
117 vtk_mpi_ctrl->SetCommunicator(vtk_comm);
118
119 vtuWriter->SetController(vtk_mpi_ctrl);
120
121 vtuWriter->SetGhostLevel(1);
122 vtuWriter->SetNumberOfPieces(num_partitions);
123 vtuWriter->SetStartPiece(rank);
124 vtuWriter->SetEndPiece(rank);
125 }
126#endif
127
128#ifdef VTK_USE_64BIT_IDS
129 vtuWriter->SetHeaderTypeToUInt64();
130 // set SetIdTypeToInt64() as well?
131#endif
132
133 return (vtuWriter->Write() > 0);
134}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Properties & getProperties()
Definition Mesh.h:136
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:105
bool existsPropertyVector(std::string_view name) const
PropertyVector< T > const * getPropertyVector(std::string_view name) const
MPI_Comm OGS_COMM_WORLD
Definition MPI.cpp:15

References _data_mode, _mesh, _use_compressor, MeshLib::Cell, DBUG(), ERR(), MeshLib::Properties::existsPropertyVector(), MeshLib::Mesh::getName(), MeshLib::Mesh::getProperties(), MeshLib::Properties::getPropertyVector(), and BaseLib::MPI::OGS_COMM_WORLD.

Referenced by writeToFile().

Member Data Documentation

◆ _data_mode

int MeshLib::IO::VtuInterface::_data_mode
private

Definition at line 79 of file VtuInterface.h.

Referenced by VtuInterface(), and writeVTU().

◆ _mesh

const MeshLib::Mesh* MeshLib::IO::VtuInterface::_mesh
private

Definition at line 78 of file VtuInterface.h.

Referenced by writeVTU().

◆ _use_compressor

bool MeshLib::IO::VtuInterface::_use_compressor
private

Definition at line 80 of file VtuInterface.h.

Referenced by writeVTU().


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