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 37 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 compressed=false)
 Provide the mesh to write and set if compression should be used. More...
 
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)
 
static MeshLib::MeshreadVTKFile (std::string const &file_name)
 

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  compressed = false 
)
explicit

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

Definition at line 43 of file VtuInterface.cpp.

45  : _mesh(mesh), _data_mode(dataMode), _use_compressor(compress)
46 {
47  if (_data_mode == vtkXMLWriter::Ascii && compress)
48  {
49  WARN("Ascii data cannot be compressed, ignoring compression flag.");
50  }
51 }
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37
const MeshLib::Mesh * _mesh
Definition: VtuInterface.h:68

References _data_mode, and WARN().

Member Function Documentation

◆ readVTKFile()

MeshLib::Mesh * MeshLib::IO::VtuInterface::readVTKFile ( std::string const &  file_name)
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 79 of file VtuInterface.cpp.

80 {
81  if (!BaseLib::IsFileExisting(file_name))
82  {
83  ERR("File '{:s}' does not exist.", file_name);
84  return nullptr;
85  }
86 
87  vtkSmartPointer<vtkGenericDataObjectReader> reader =
88  vtkSmartPointer<vtkGenericDataObjectReader>::New();
89  reader->SetFileName(file_name.c_str());
90  reader->Update();
91 
92  // check for unstructured grid
93  if (reader->ReadOutputType() != 4)
94  {
95  ERR("Only VTK-files with dataset type \"Unstructured Grid\" are "
96  "currently supported.");
97  return nullptr;
98  }
99 
100  reader->ReadAllFieldsOn();
101  reader->ReadAllScalarsOn();
102  vtkUnstructuredGrid* vtkGrid = reader->GetUnstructuredGridOutput();
103  if (vtkGrid->GetNumberOfPoints() == 0)
104  {
105  ERR("Mesh '{:s}' contains zero points.", file_name);
106  return nullptr;
107  }
108 
109  std::string const mesh_name(
112  mesh_name);
113 }
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, 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:43
std::string extractBaseNameWithoutExtension(std::string const &pathname)
Definition: FileTools.cpp:180

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)
static

Read an unstructured grid from a VTU file.

Returns
The converted mesh or a nullptr if reading failed

Definition at line 53 of file VtuInterface.cpp.

54 {
55  if (!BaseLib::IsFileExisting(file_name))
56  {
57  ERR("File '{:s}' does not exist.", file_name);
58  return nullptr;
59  }
60 
61  vtkSmartPointer<vtkXMLUnstructuredGridReader> reader =
62  vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
63  reader->SetFileName(file_name.c_str());
64  reader->Update();
65 
66  vtkUnstructuredGrid* vtkGrid = reader->GetOutput();
67  if (vtkGrid->GetNumberOfPoints() == 0)
68  {
69  ERR("Mesh '{:s}' contains zero points.", file_name);
70  return nullptr;
71  }
72 
73  std::string const mesh_name(
76  mesh_name);
77 }

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

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

◆ 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 138 of file VtuInterface.cpp.

139 {
140 #ifdef USE_PETSC
141  auto const vtu_file_name =
143  int rank;
144  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
145  int mpi_size;
146  MPI_Comm_size(MPI_COMM_WORLD, &mpi_size);
147  return writeVTU<vtkXMLPUnstructuredGridWriter>(vtu_file_name + ".pvtu",
148  mpi_size, rank);
149 #else
150  return writeVTU<vtkXMLUnstructuredGridWriter>(file_path.string());
151 #endif
152 }
std::string getVtuFileNameForPetscOutputWithoutExtension(std::string const &file_name)

References MeshLib::IO::getVtuFileNameForPetscOutputWithoutExtension().

Referenced by SaveMeshDialog::accept(), convert(), OGSFileConverter::convertMSH2VTU(), main(), ProcessLib::makeOutput(), 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 32 of file VtuInterface-impl.h.

34 {
35  if (!_mesh)
36  {
37  ERR("VtuInterface::write(): No mesh specified.");
38  return false;
39  }
40 
41  vtkNew<MeshLib::VtkMappedMeshSource> vtkSource;
42  vtkSource->SetMesh(_mesh);
43 
44  vtkSmartPointer<UnstructuredGridWriter> vtuWriter =
45  vtkSmartPointer<UnstructuredGridWriter>::New();
46 
47  vtkSource->Update();
48  vtuWriter->SetInputData(vtkSource->GetOutput());
49 
50  if (_use_compressor)
51  {
52  vtuWriter->SetCompressorTypeToZLib();
53  }
54  else
55  {
56  vtuWriter->SetCompressorTypeToNone();
57  }
58 
59  vtuWriter->SetDataMode(_data_mode);
60  if (_data_mode == vtkXMLWriter::Appended)
61  {
62  vtuWriter->SetEncodeAppendedData(1);
63  }
64  if (_data_mode == vtkXMLWriter::Ascii)
65  {
66  vtkSource->Update();
67  vtkSmartPointer<vtkUnstructuredGrid> tempGrid =
68  vtkSmartPointer<vtkUnstructuredGrid>::New();
69  tempGrid->DeepCopy(vtkSource->GetOutput());
70  vtuWriter->SetInputDataObject(tempGrid);
71  }
72 
73  vtuWriter->SetFileName(file_name.c_str());
74 #ifdef USE_PETSC
75  vtuWriter->SetGhostLevel(1);
76  vtuWriter->SetNumberOfPieces(num_partitions);
77  vtuWriter->SetStartPiece(rank);
78  vtuWriter->SetEndPiece(rank);
79 #else
80  // avoid unused parameter warnings.
81  (void)num_partitions;
82  (void)rank;
83 #endif
84 
85 #ifdef VTK_USE_64BIT_IDS
86  vtuWriter->SetHeaderTypeToUInt64();
87  // set SetIdTypeToInt64() as well?
88 #endif
89 
90  return (vtuWriter->Write() > 0);
91 }

References _data_mode, _mesh, _use_compressor, and ERR().

Member Data Documentation

◆ _data_mode

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

Definition at line 69 of file VtuInterface.h.

Referenced by VtuInterface(), and writeVTU().

◆ _mesh

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

Definition at line 68 of file VtuInterface.h.

Referenced by writeVTU().

◆ _use_compressor

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

Definition at line 70 of file VtuInterface.h.

Referenced by writeVTU().


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