ViennaLS
Loading...
Searching...
No Matches
lsVTKWriter.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <fstream>
4#include <string>
5#include <unordered_map>
6#include <utility>
7
8#include <lsFileFormats.hpp>
9#include <lsMesh.hpp>
10
11#include <vcLogger.hpp>
12#include <vcSmartPointer.hpp>
13
14#ifdef VIENNALS_USE_VTK
15#include <vtkCellArray.h>
16#include <vtkCellData.h>
17#include <vtkFloatArray.h>
18#include <vtkPointData.h>
19#include <vtkPoints.h>
20#include <vtkPolyData.h>
21#include <vtkSmartPointer.h>
22#include <vtkXMLPolyDataWriter.h>
23
24#include <vtkUnstructuredGrid.h>
25#include <vtkXMLUnstructuredGridWriter.h>
26#endif // VIENNALS_USE_VTK
27
28namespace viennals {
29
30using namespace viennacore;
31
33template <class T> class VTKWriter {
34 using MetaDataType = std::unordered_map<std::string, std::vector<double>>;
35
36 SmartPointer<Mesh<T>> mesh = nullptr;
38 std::string fileName;
39 MetaDataType metaData;
40
41#ifdef VIENNALS_USE_VTK
42 template <class In, class Out>
43 void addDataFromMesh(const In &inData, Out outData) const {
44 // now add pointData
45 for (unsigned i = 0; i < inData.getScalarDataSize(); ++i) {
46 vtkSmartPointer<vtkFloatArray> pointData =
47 vtkSmartPointer<vtkFloatArray>::New();
48 pointData->SetNumberOfComponents(1);
49 pointData->SetName(inData.getScalarDataLabel(i).c_str());
50 auto scalars = *(inData.getScalarData(i));
51 for (unsigned j = 0; j < inData.getScalarData(i)->size(); ++j) {
52 pointData->InsertNextValue(scalars[j]);
53 }
54 outData->AddArray(pointData);
55 }
56
57 // now add vector data
58 for (unsigned i = 0; i < inData.getVectorDataSize(); ++i) {
59 vtkSmartPointer<vtkFloatArray> vectorData =
60 vtkSmartPointer<vtkFloatArray>::New();
61 vectorData->SetNumberOfComponents(3);
62 vectorData->SetName(inData.getVectorDataLabel(i).c_str());
63 auto vectors = *(inData.getVectorData(i));
64 for (unsigned j = 0; j < inData.getVectorData(i)->size(); ++j) {
65 vectorData->InsertNextTuple3(vectors[j][0], vectors[j][1],
66 vectors[j][2]);
67 }
68 outData->AddArray(vectorData);
69 }
70 }
71
72 void addMetaDataToVTK(vtkDataSet *data) const {
73 if (metaData.empty()) {
74 return;
75 }
76
77 // add metadata to field data
78 vtkSmartPointer<vtkFieldData> fieldData = data->GetFieldData();
79 for (const auto &meta : metaData) {
80 if (meta.second.empty())
81 continue; // skip empty metadata
82
83 vtkSmartPointer<vtkFloatArray> metaDataArray =
84 vtkSmartPointer<vtkFloatArray>::New();
85 metaDataArray->SetName(meta.first.c_str());
86 metaDataArray->SetNumberOfValues(meta.second.size());
87 for (size_t i = 0; i < meta.second.size(); ++i) {
88 metaDataArray->SetValue(i, meta.second[i]);
89 }
90 fieldData->AddArray(metaDataArray);
91 }
92 }
93#endif // VIENNALS_USE_VTK
94
95public:
96 VTKWriter() = default;
97
98 VTKWriter(SmartPointer<Mesh<T>> passedMesh) : mesh(passedMesh) {}
99
100 VTKWriter(SmartPointer<Mesh<T>> passedMesh, std::string passedFileName)
101 : mesh(passedMesh), fileName(std::move(passedFileName)) {}
102
103 VTKWriter(SmartPointer<Mesh<T>> passedMesh, FileFormatEnum passedFormat,
104 std::string passedFileName)
105 : mesh(passedMesh), fileFormat(passedFormat),
106 fileName(std::move(passedFileName)) {}
107
108 void setMesh(SmartPointer<Mesh<T>> passedMesh) { mesh = passedMesh; }
109
111 void setFileFormat(FileFormatEnum passedFormat) { fileFormat = passedFormat; }
112
114 void setFileName(std::string passedFileName) {
115 fileName = std::move(passedFileName);
116 }
117
118 void setMetaData(const MetaDataType &passedMetaData) {
119 metaData = passedMetaData;
120 }
121
122 void addMetaData(const std::string &key, double value) {
123 metaData[key] = std::vector<double>{value};
124 }
125
126 void addMetaData(const std::string &key, const std::vector<double> &values) {
127 metaData[key] = values;
128 }
129
130 void addMetaData(const MetaDataType &newMetaData) {
131 for (const auto &pair : newMetaData) {
132 metaData[pair.first] = pair.second;
133 }
134 }
135
136 void apply() {
137 // check mesh
138 if (mesh == nullptr) {
139 Logger::getInstance()
140 .addError("No mesh was passed to VTKWriter.")
141 .print();
142 return;
143 }
144 // check filename
145 if (fileName.empty()) {
146 VIENNACORE_LOG_ERROR("No file name specified for VTKWriter.");
147 return;
148 }
149 if (mesh->nodes.empty()) {
150 VIENNACORE_LOG_WARNING("Writing empty mesh.");
151 return;
152 }
153
154 if (fileFormat == FileFormatEnum::VTK_AUTO) {
155 auto dotPos = fileName.rfind('.');
156 if (dotPos == std::string::npos) {
157 fileFormat = FileFormatEnum::VTP;
158 } else {
159 auto ending = fileName.substr(dotPos);
160 if (ending == ".vtk") {
161 fileFormat = FileFormatEnum::VTK_LEGACY;
162 } else if (ending == ".vtp") {
163 fileFormat = FileFormatEnum::VTP;
164 } else if (ending == ".vtu") {
165 fileFormat = FileFormatEnum::VTU;
166 } else {
167 Logger::getInstance()
168 .addError("No valid file format found based on the file ending "
169 "passed to VTKWriter.")
170 .print();
171 return;
172 }
173 }
174 }
175
176 // check file format
177 switch (fileFormat) {
179 writeVTKLegacy(fileName);
180 break;
181#ifdef VIENNALS_USE_VTK
183 writeVTP(fileName);
184 break;
186 writeVTU(fileName);
187 break;
188#else
191 VIENNACORE_LOG_WARNING(
192 "VTKWriter was built without VTK support. Falling back "
193 "to VTK_LEGACY.");
194 writeVTKLegacy(fileName);
195 break;
196#endif
197 default:
198 VIENNACORE_LOG_ERROR("No valid file format set for VTKWriter.");
199 }
200 }
201
202private:
203#ifdef VIENNALS_USE_VTK
204 void writeVTP(std::string filename) const {
205 if (mesh == nullptr) {
206 Logger::getInstance()
207 .addError("No mesh was passed to VTKWriter.")
208 .print();
209 return;
210 }
211
212 if (filename.find(".vtp") != filename.size() - 4)
213 filename += ".vtp";
214 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
215
216 // Points
217 vtkSmartPointer<vtkPoints> polyPoints = vtkSmartPointer<vtkPoints>::New();
218 for (auto it = mesh->getNodes().begin(); it != mesh->getNodes().end();
219 ++it) {
220 polyPoints->InsertNextPoint((*it)[0], (*it)[1], (*it)[2]);
221 }
222 polyData->SetPoints(polyPoints);
223
224 // Vertices
225 if (mesh->vertices.size() > 0) {
226 vtkSmartPointer<vtkCellArray> polyCells =
227 vtkSmartPointer<vtkCellArray>::New();
228 for (auto it = mesh->vertices.begin(); it != mesh->vertices.end(); ++it) {
229 polyCells->InsertNextCell(1);
230 polyCells->InsertCellPoint((*it)[0]);
231 }
232 polyData->SetVerts(polyCells);
233 }
234
235 // Lines
236 if (mesh->lines.size() > 0) {
237 vtkSmartPointer<vtkCellArray> polyCells =
238 vtkSmartPointer<vtkCellArray>::New();
239 for (auto it = mesh->lines.begin(); it != mesh->lines.end(); ++it) {
240 polyCells->InsertNextCell(2);
241 for (unsigned i = 0; i < 2; ++i) {
242 polyCells->InsertCellPoint((*it)[i]);
243 }
244 }
245 polyData->SetLines(polyCells);
246 }
247
248 // Triangles
249 if (mesh->triangles.size() > 0) {
250 vtkSmartPointer<vtkCellArray> polyCells =
251 vtkSmartPointer<vtkCellArray>::New();
252 for (auto it = mesh->triangles.begin(); it != mesh->triangles.end();
253 ++it) {
254 polyCells->InsertNextCell(3);
255 for (unsigned i = 0; i < 3; ++i) {
256 polyCells->InsertCellPoint((*it)[i]);
257 }
258 }
259 polyData->SetPolys(polyCells);
260 }
261
262 addDataFromMesh(mesh->pointData, polyData->GetPointData());
263 addDataFromMesh(mesh->cellData, polyData->GetCellData());
264 addMetaDataToVTK(polyData);
265
266 vtkSmartPointer<vtkXMLPolyDataWriter> pwriter =
267 vtkSmartPointer<vtkXMLPolyDataWriter>::New();
268 pwriter->SetFileName(filename.c_str());
269 pwriter->SetInputData(polyData);
270 pwriter->Write();
271 }
272
273 void writeVTU(std::string filename) const {
274 if (mesh == nullptr) {
275 Logger::getInstance()
276 .addError("No mesh was passed to VTKWriter.")
277 .print();
278 return;
279 }
280
281 if (filename.find(".vtu") != filename.size() - 4)
282 filename += ".vtu";
283
284 vtkSmartPointer<vtkUnstructuredGrid> uGrid =
285 vtkSmartPointer<vtkUnstructuredGrid>::New();
286
287 // Points
288 vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
289 for (auto it = mesh->getNodes().begin(); it != mesh->getNodes().end();
290 ++it) {
291 points->InsertNextPoint((*it)[0], (*it)[1], (*it)[2]);
292 }
293 uGrid->SetPoints(points);
294
295 // Now set all cells
296 vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
297 std::vector<int> cellTypes;
298 cellTypes.reserve(mesh->vertices.size() + mesh->lines.size() +
299 mesh->triangles.size() + mesh->tetras.size() +
300 mesh->hexas.size());
301
302 // Vertices
303 if (mesh->vertices.size() > 0) {
304 for (auto it = mesh->vertices.begin(); it != mesh->vertices.end(); ++it) {
305 cells->InsertNextCell(1);
306 cells->InsertCellPoint((*it)[0]);
307 cellTypes.push_back(1); // vtk Vertex
308 }
309 }
310
311 // Lines
312 if (mesh->lines.size() > 0) {
313 for (auto it = mesh->lines.begin(); it != mesh->lines.end(); ++it) {
314 cells->InsertNextCell(2);
315 for (unsigned i = 0; i < 2; ++i) {
316 cells->InsertCellPoint((*it)[i]);
317 }
318 cellTypes.push_back(3); // vtk Line
319 }
320 }
321
322 // Triangles
323 if (mesh->triangles.size() > 0) {
324 for (auto it = mesh->triangles.begin(); it != mesh->triangles.end();
325 ++it) {
326 cells->InsertNextCell(3);
327 for (unsigned i = 0; i < 3; ++i) {
328 cells->InsertCellPoint((*it)[i]);
329 }
330 cellTypes.push_back(5); // vtk Triangle
331 }
332 }
333
334 // Tetras
335 if (mesh->tetras.size() > 0) {
336 for (auto it = mesh->tetras.begin(); it != mesh->tetras.end(); ++it) {
337 cells->InsertNextCell(4);
338 for (unsigned i = 0; i < 4; ++i) {
339 cells->InsertCellPoint((*it)[i]);
340 }
341 cellTypes.push_back(10); // vtk Tetra
342 }
343 }
344
345 // Hexas
346 if (mesh->hexas.size() > 0) {
347 for (auto it = mesh->hexas.begin(); it != mesh->hexas.end(); ++it) {
348 cells->InsertNextCell(8);
349 for (unsigned i = 0; i < 8; ++i) {
350 cells->InsertCellPoint((*it)[i]);
351 }
352 cellTypes.push_back(12); // vtk Hexahedron
353 }
354 }
355
356 // set cells
357 uGrid->SetCells(&(cellTypes[0]), cells);
358
359 addDataFromMesh(mesh->pointData, uGrid->GetPointData());
360 addDataFromMesh(mesh->cellData, uGrid->GetCellData());
361 addMetaDataToVTK(uGrid);
362
363 // // now add pointData
364 // for (unsigned i = 0; i < mesh->cellData.getScalarDataSize(); ++i) {
365 // vtkSmartPointer<vtkFloatArray> pointData =
366 // vtkSmartPointer<vtkFloatArray>::New();
367 // pointData->SetNumberOfComponents(1);
368 // pointData->SetName(mesh->cellData.getScalarDataLabel(i).c_str());
369 // auto scalars = *(mesh->cellData.getScalarData(i));
370 // for (unsigned j = 0; j < scalars.size(); ++j) {
371 // pointData->InsertNextValue(scalars[j]);
372 // }
373 // uGrid->GetCellData()->AddArray(pointData);
374 // }
375
376 // // now add vector data
377 // for (unsigned i = 0; i < mesh->cellData.getVectorDataSize(); ++i) {
378 // vtkSmartPointer<vtkFloatArray> vectorData =
379 // vtkSmartPointer<vtkFloatArray>::New();
380 // vectorData->SetNumberOfComponents(3);
381 // vectorData->SetName(mesh->cellData.getVectorDataLabel(i).c_str());
382 // auto vectors = *(mesh->cellData.getVectorData(i));
383 // for (unsigned j = 0; j < vectors.size(); ++j) {
384 // vectorData->InsertNextTuple3(vectors[j][0], vectors[j][1],
385 // vectors[j][2]);
386 // }
387 // uGrid->GetCellData()->AddArray(vectorData);
388 // }
389
390 vtkSmartPointer<vtkXMLUnstructuredGridWriter> owriter =
391 vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
392 owriter->SetFileName(filename.c_str());
393 owriter->SetInputData(uGrid);
394 owriter->Write();
395 }
396
397#endif // VIENNALS_USE_VTK
398
399 void writeVTKLegacy(const std::string &filename) {
400 if (mesh == nullptr) {
401 Logger::getInstance()
402 .addError("No mesh was passed to VTKWriter.")
403 .print();
404 return;
405 }
406
407 std::ofstream f(filename.c_str());
408
409 f << "# vtk DataFile Version 2.0" << std::endl;
410 f << ((!mesh->lines.empty()) ? 2 : 3) << "D Surface" << std::endl;
411 f << "ASCII" << std::endl;
412 f << "DATASET UNSTRUCTURED_GRID" << std::endl;
413 f << "POINTS " << mesh->nodes.size() << " float" << std::endl;
414
415 // print nodes
416 for (unsigned int i = 0; i < mesh->nodes.size(); i++) {
417 for (int j = 0; j < 3; j++)
418 f << static_cast<float>(mesh->nodes[i][j]) << " ";
419 f << std::endl;
420 }
421
422 const unsigned numberOfCells = mesh->vertices.size() + mesh->lines.size() +
423 mesh->triangles.size() +
424 mesh->tetras.size() + mesh->hexas.size();
425 const unsigned cellDataSize =
426 2 * mesh->vertices.size() + 3 * mesh->lines.size() +
427 4 * mesh->triangles.size() + 5 * mesh->tetras.size() +
428 9 * mesh->hexas.size();
429
430 f << "CELLS " << numberOfCells << " " << cellDataSize << std::endl;
431
432 // print elements
433 for (unsigned int i = 0; i < mesh->vertices.size(); i++) {
434 f << 1 << " " << mesh->vertices[i][0] << std::endl;
435 }
436 for (unsigned int i = 0; i < mesh->lines.size(); i++) {
437 f << 2 << " ";
438 for (int j = 0; j < 2; j++)
439 f << mesh->lines[i][j] << " ";
440 f << std::endl;
441 }
442
443 for (unsigned int i = 0; i < mesh->triangles.size(); i++) {
444 f << 3 << " ";
445 for (int j = 0; j < 3; j++)
446 f << mesh->triangles[i][j] << " ";
447 f << std::endl;
448 }
449
450 for (unsigned int i = 0; i < mesh->tetras.size(); i++) {
451 f << 4 << " ";
452 for (int j = 0; j < 4; j++)
453 f << mesh->tetras[i][j] << " ";
454 f << std::endl;
455 }
456
457 for (unsigned int i = 0; i < mesh->hexas.size(); i++) {
458 f << 8 << " ";
459 for (int j = 0; j < 8; j++)
460 f << mesh->hexas[i][j] << " ";
461 f << std::endl;
462 }
463
464 f << "CELL_TYPES " << numberOfCells << std::endl;
465 for (unsigned i = 0; i < mesh->vertices.size(); ++i)
466 f << 1 << std::endl;
467
468 for (unsigned i = 0; i < mesh->lines.size(); ++i)
469 f << 3 << std::endl;
470
471 for (unsigned i = 0; i < mesh->triangles.size(); ++i)
472 f << 5 << std::endl;
473
474 for (unsigned i = 0; i < mesh->tetras.size(); ++i)
475 f << 10 << std::endl;
476
477 for (unsigned i = 0; i < mesh->hexas.size(); ++i)
478 f << 12 << std::endl;
479
480 // WRITE POINT DATA
481 if (mesh->pointData.getScalarDataSize() ||
482 mesh->pointData.getVectorDataSize()) {
483 VIENNACORE_LOG_WARNING(
484 "Point data output not supported for legacy VTK output. "
485 "Point data is ignored.");
486 }
487
488 // WRITE SCALAR DATA
489 if (mesh->cellData.getScalarDataSize()) {
490 f << "CELL_DATA " << mesh->cellData.getScalarData(0)->size() << std::endl;
491 for (unsigned i = 0; i < mesh->cellData.getScalarDataSize(); ++i) {
492 auto scalars = *(mesh->cellData.getScalarData(i));
493 f << "SCALARS " << mesh->cellData.getScalarDataLabel(i) << " float"
494 << std::endl;
495 f << "LOOKUP_TABLE default" << std::endl;
496 for (unsigned j = 0; j < scalars.size(); ++j) {
497 f << ((std::abs(scalars[j]) < 1e-6) ? 0.0 : scalars[j]) << std::endl;
498 }
499 }
500 }
501
502 // WRITE VECTOR DATA
503 if (mesh->cellData.getVectorDataSize()) {
504 if (!mesh->cellData.getScalarDataSize())
505 f << "CELL_DATA " << mesh->cellData.getVectorData(0)->size()
506 << std::endl;
507 for (unsigned i = 0; i < mesh->cellData.getVectorDataSize(); ++i) {
508 auto vectors = *(mesh->cellData.getVectorData(i));
509 f << "VECTORS " << mesh->cellData.getVectorDataLabel(i) << " float"
510 << std::endl;
511 for (unsigned j = 0; j < vectors.size(); ++j) {
512 for (unsigned k = 0; k < 3; ++k) {
513 f << vectors[j][k] << " ";
514 }
515 f << std::endl;
516 }
517 }
518 }
519
520 f.close();
521 }
522};
523
524} // namespace viennals
This class holds an explicit mesh, which is always given in 3 dimensions. If it describes a 2D mesh,...
Definition lsMesh.hpp:22
void setFileName(std::string passedFileName)
set file name for file to write
Definition lsVTKWriter.hpp:114
VTKWriter(SmartPointer< Mesh< T > > passedMesh, FileFormatEnum passedFormat, std::string passedFileName)
Definition lsVTKWriter.hpp:103
VTKWriter(SmartPointer< Mesh< T > > passedMesh, std::string passedFileName)
Definition lsVTKWriter.hpp:100
void addMetaData(const std::string &key, double value)
Definition lsVTKWriter.hpp:122
void apply()
Definition lsVTKWriter.hpp:136
void setMetaData(const MetaDataType &passedMetaData)
Definition lsVTKWriter.hpp:118
void setMesh(SmartPointer< Mesh< T > > passedMesh)
Definition lsVTKWriter.hpp:108
void setFileFormat(FileFormatEnum passedFormat)
set file format for file to write. Defaults to VTK_LEGACY.
Definition lsVTKWriter.hpp:111
void addMetaData(const std::string &key, const std::vector< double > &values)
Definition lsVTKWriter.hpp:126
void addMetaData(const MetaDataType &newMetaData)
Definition lsVTKWriter.hpp:130
VTKWriter(SmartPointer< Mesh< T > > passedMesh)
Definition lsVTKWriter.hpp:98
Definition lsAdvect.hpp:37
FileFormatEnum
Definition lsFileFormats.hpp:4
@ VTK_AUTO
Definition lsFileFormats.hpp:8
@ VTK_LEGACY
Definition lsFileFormats.hpp:5
@ VTP
Definition lsFileFormats.hpp:6
@ VTU
Definition lsFileFormats.hpp:7