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