OGS
VtuInterface-impl.h
Go to the documentation of this file.
1 
19 #include <vtkNew.h>
20 #include <vtkSmartPointer.h>
21 #include <vtkUnstructuredGrid.h>
22 
23 #include "BaseLib/Logging.h"
25 #include "VtuInterface.h"
26 
27 namespace MeshLib
28 {
29 namespace IO
30 {
31 template <typename UnstructuredGridWriter>
32 bool VtuInterface::writeVTU(std::string const& file_name,
33  const int num_partitions, const int rank)
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 }
92 
93 } // end namespace IO
94 } // end namespace MeshLib
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
Implementation of the VtuInterface class.
bool writeVTU(std::string const &file_name, const int num_partitions=1, const int rank=1)
const MeshLib::Mesh * _mesh
Definition: VtuInterface.h:68