ViennaLS
Loading...
Searching...
No Matches
lsVTKReader.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 <vtkCellData.h>
14#include <vtkIdList.h>
15#include <vtkPointData.h>
16#include <vtkPolyData.h>
17#include <vtkSmartPointer.h>
18#include <vtkUnstructuredGrid.h>
19#include <vtkXMLPolyDataReader.h>
20#include <vtkXMLUnstructuredGridReader.h>
21#endif // VIENNALS_USE_VTK
22
23namespace viennals {
24
25using namespace viennacore;
26
28template <class T = double> class VTKReader {
29 SmartPointer<Mesh<T>> mesh = nullptr;
31 std::string fileName;
32
33 unsigned vtk_nodes_for_cell_type[15] = {0, 1, 0, 2, 0, 3, 0, 0,
34 4, 4, 4, 8, 8, 6, 5};
35
36public:
38
39 VTKReader(SmartPointer<Mesh<T>> passedMesh) : mesh(passedMesh) {}
40
41 VTKReader(SmartPointer<Mesh<T>> passedMesh, std::string passedFileName)
42 : mesh(passedMesh), fileName(passedFileName) {}
43
44 VTKReader(SmartPointer<Mesh<>> passedMesh, FileFormatEnum passedFormat,
45 std::string passedFileName)
46 : mesh(passedMesh), fileFormat(passedFormat), fileName(passedFileName) {}
47
49 void setMesh(SmartPointer<Mesh<>> passedMesh) { mesh = passedMesh; }
50
52 void setFileFormat(FileFormatEnum passedFormat) { fileFormat = passedFormat; }
53
55 void setFileName(std::string passedFileName) { fileName = passedFileName; }
56
57 void apply() {
58 // check mesh
59 if (mesh == nullptr) {
60 Logger::getInstance()
61 .addWarning("No mesh was passed to VTKReader. Not reading.")
62 .print();
63 return;
64 }
65 // check filename
66 if (fileName.empty()) {
67 Logger::getInstance()
68 .addWarning("No file name specified for VTKReader. Not reading.")
69 .print();
70 return;
71 }
72
73 if (fileFormat == FileFormatEnum::VTK_AUTO) {
74 auto dotPos = fileName.rfind('.');
75 if (dotPos == std::string::npos) {
76 Logger::getInstance()
77 .addWarning("No valid file format found based on the file ending "
78 "passed to VTKReader. Not reading.")
79 .print();
80 return;
81 }
82 auto ending = fileName.substr(dotPos);
83 if (ending == ".vtk") {
84 fileFormat = FileFormatEnum::VTK_LEGACY;
85 } else if (ending == ".vtp") {
86 fileFormat = FileFormatEnum::VTP;
87 } else if (ending == ".vtu") {
88 fileFormat = FileFormatEnum::VTU;
89 } else {
90 Logger::getInstance()
91 .addWarning("No valid file format found based on the file ending "
92 "passed to VTKReader. Not reading.")
93 .print();
94 return;
95 }
96 }
97
98 // check file format
99 switch (fileFormat) {
101 readVTKLegacy(fileName);
102 break;
103#ifdef VIENNALS_USE_VTK
105 readVTP(fileName);
106 break;
108 readVTU(fileName);
109 break;
110#else
113 Logger::getInstance()
114 .addWarning(
115 "VTKReader was built without VTK support. Only VTK_LEGACY "
116 "can be used. File not read.")
117 .print();
118#endif
119 default:
120 Logger::getInstance()
121 .addWarning("No valid file format set for VTKReader. Not reading.")
122 .print();
123 }
124 }
125
126private:
127#ifdef VIENNALS_USE_VTK
128 void readVTP(std::string filename) {
129 if (mesh == nullptr) {
130 Logger::getInstance()
131 .addWarning("No mesh was passed to VTKReader.")
132 .print();
133 return;
134 }
135
136 mesh->clear();
137 vtkSmartPointer<vtkXMLPolyDataReader> pReader =
138 vtkSmartPointer<vtkXMLPolyDataReader>::New();
139 pReader->SetFileName(filename.c_str());
140 pReader->Update();
141
142 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
143 polyData = pReader->GetOutput();
144
145 mesh->nodes.resize(polyData->GetNumberOfPoints());
146 for (unsigned i = 0; i < mesh->nodes.size(); ++i) {
147 std::array<double, 3> coords;
148 polyData->GetPoint(i, coords.data());
149 mesh->nodes[i] = coords;
150 }
151
152 vtkSmartPointer<vtkCellArray> cellArray =
153 vtkSmartPointer<vtkCellArray>::New();
154 // get vertices
155 {
156 mesh->vertices.reserve(polyData->GetNumberOfVerts());
157 cellArray = polyData->GetVerts();
158 cellArray->InitTraversal();
159 vtkIdList *pointList = vtkIdList::New();
160 while (cellArray->GetNextCell(pointList)) {
161 std::array<unsigned, 1> cell;
162 cell[0] = pointList->GetId(0);
163 mesh->vertices.push_back(cell);
164 }
165 }
166
167 // get lines
168 {
169 mesh->lines.reserve(polyData->GetNumberOfLines());
170 cellArray = polyData->GetLines();
171 cellArray->InitTraversal();
172 vtkIdList *pointList = vtkIdList::New();
173 while (cellArray->GetNextCell(pointList)) {
174 std::array<unsigned, 2> cell;
175 for (unsigned i = 0; i < 2; ++i) {
176 cell[i] = pointList->GetId(i);
177 }
178 mesh->lines.push_back(cell);
179 }
180 }
181
182 // get triangles
183 {
184 mesh->triangles.reserve(polyData->GetNumberOfPolys());
185 cellArray = polyData->GetPolys();
186 cellArray->InitTraversal();
187 vtkIdList *pointList = vtkIdList::New();
188 while (cellArray->GetNextCell(pointList)) {
189 std::array<unsigned, 3> cell;
190 for (unsigned i = 0; i < 3; ++i) {
191 cell[i] = pointList->GetId(i);
192 }
193 mesh->triangles.push_back(cell);
194 }
195 }
196
197 // get point data
198 vtkSmartPointer<vtkPointData> pointData =
199 vtkSmartPointer<vtkPointData>::New();
200 pointData = polyData->GetPointData();
201
202 for (unsigned i = 0;
203 i < static_cast<unsigned>(pointData->GetNumberOfArrays()); ++i) {
204 vtkDataArray *dataArray;
205 dataArray = pointData->GetArray(i);
206 if (dataArray->GetNumberOfComponents() == 1) {
207 mesh->pointData.insertNextScalarData(
209 std::string(pointData->GetArrayName(i)));
210 auto &scalars = *(mesh->pointData.getScalarData(i));
211 scalars.resize(pointData->GetNumberOfTuples());
212 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
213 scalars[j] = dataArray->GetTuple1(j);
214 }
215 } else if (dataArray->GetNumberOfComponents() == 3) {
216 mesh->pointData.insertNextVectorData(
218 std::string(pointData->GetArrayName(i)));
219 auto &vectors = *(mesh->pointData.getVectorData(i));
220 vectors.resize(pointData->GetNumberOfTuples());
221 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
222 std::array<double, 3> vector;
223 dataArray->GetTuple(j, &(vector[0]));
224 vectors[j] = vector;
225 }
226 }
227 }
228
229 // get cell data
230 vtkSmartPointer<vtkCellData> cellData = vtkSmartPointer<vtkCellData>::New();
231 cellData = polyData->GetCellData();
232
233 for (unsigned i = 0;
234 i < static_cast<unsigned>(cellData->GetNumberOfArrays()); ++i) {
235 vtkDataArray *dataArray;
236 dataArray = cellData->GetArray(i);
237 if (cellData->GetNumberOfComponents() == 1) {
238 mesh->cellData.insertNextScalarData(
240 std::string(cellData->GetArrayName(i)));
241 auto &scalars = *(mesh->cellData.getScalarData(i));
242 scalars.resize(cellData->GetNumberOfTuples());
243 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
244 scalars[j] = dataArray->GetTuple1(j);
245 }
246 } else if (cellData->GetNumberOfComponents() == 3) {
247 mesh->cellData.insertNextVectorData(
249 std::string(cellData->GetArrayName(i)));
250 auto &vectors = *(mesh->cellData.getVectorData(i));
251 vectors.resize(cellData->GetNumberOfTuples());
252 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
253 std::array<double, 3> vector;
254 dataArray->GetTuple(j, &(vector[0]));
255 vectors[j] = vector;
256 }
257 }
258 }
259 }
260
261 void readVTU(std::string filename) {
262 if (mesh == nullptr) {
263 Logger::getInstance()
264 .addWarning("No mesh was passed to VTKReader.")
265 .print();
266 return;
267 }
268
269 mesh->clear();
270
271 vtkSmartPointer<vtkXMLUnstructuredGridReader> greader =
272 vtkSmartPointer<vtkXMLUnstructuredGridReader>::New();
273 greader->SetFileName(filename.c_str());
274 greader->Update();
275
276 vtkSmartPointer<vtkUnstructuredGrid> ugrid =
277 vtkSmartPointer<vtkUnstructuredGrid>::New();
278 ugrid = greader->GetOutput();
279
280 // get all points
281 mesh->nodes.resize(ugrid->GetNumberOfPoints());
282 for (unsigned i = 0; i < mesh->nodes.size(); ++i) {
283 std::array<double, 3> coords;
284 ugrid->GetPoint(i, &(coords[0]));
285 mesh->nodes[i] = coords;
286 }
287
288 // get cells
289 for (unsigned i = 0; i < ugrid->GetNumberOfCells(); ++i) {
290 vtkIdList *pointList = vtkIdList::New();
291 ugrid->GetCellPoints(i, pointList);
292
293 switch (ugrid->GetCellType(i)) {
294 case 1: // vert
295 {
296 std::array<unsigned, 1> vert;
297 vert[0] = pointList->GetId(0);
298 mesh->vertices.push_back(vert);
299 } break;
300 case 3: // line
301 {
302 std::array<unsigned, 2> elements;
303 for (unsigned j = 0; j < 2; ++j) {
304 elements[j] = pointList->GetId(j);
305 }
306 mesh->lines.push_back(elements);
307 } break;
308 case 5: // triangle
309 {
310 std::array<unsigned, 3> elements;
311 for (unsigned j = 0; j < 3; ++j) {
312 elements[j] = pointList->GetId(j);
313 }
314 mesh->triangles.push_back(elements);
315 } break;
316 case 10: // tetra
317 {
318 std::array<unsigned, 4> elements;
319 for (unsigned j = 0; j < 4; ++j) {
320 elements[j] = pointList->GetId(j);
321 }
322 mesh->tetras.push_back(elements);
323 } break;
324 case 12: // hexa
325 {
326 std::array<unsigned, 8> elements;
327 for (unsigned j = 0; j < 8; ++j) {
328 elements[j] = pointList->GetId(j);
329 }
330 mesh->hexas.push_back(elements);
331 } break;
332 }
333 }
334
335 // get point data
336 vtkSmartPointer<vtkPointData> pointData =
337 vtkSmartPointer<vtkPointData>::New();
338 pointData = ugrid->GetPointData();
339
340 for (unsigned i = 0;
341 i < static_cast<unsigned>(pointData->GetNumberOfArrays()); ++i) {
342 vtkDataArray *dataArray;
343 dataArray = pointData->GetArray(i);
344 if (dataArray->GetNumberOfComponents() == 1) {
345 mesh->pointData.insertNextScalarData(
347 std::string(pointData->GetArrayName(i)));
348 auto &scalars = *(mesh->pointData.getScalarData(i));
349 scalars.resize(pointData->GetNumberOfTuples());
350 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
351 scalars[j] = dataArray->GetTuple1(j);
352 }
353 } else if (dataArray->GetNumberOfComponents() == 3) {
354 mesh->pointData.insertNextVectorData(
356 std::string(pointData->GetArrayName(i)));
357 auto &vectors = *(mesh->pointData.getVectorData(i));
358 vectors.resize(pointData->GetNumberOfTuples());
359 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
360 std::array<double, 3> vector;
361 dataArray->GetTuple(j, &(vector[0]));
362 vectors[j] = vector;
363 }
364 }
365 }
366
367 // get cell data
368 vtkSmartPointer<vtkCellData> cellData = vtkSmartPointer<vtkCellData>::New();
369 cellData = ugrid->GetCellData();
370
371 for (unsigned i = 0;
372 i < static_cast<unsigned>(cellData->GetNumberOfArrays()); ++i) {
373 vtkDataArray *dataArray;
374 dataArray = cellData->GetArray(i);
375 if (cellData->GetNumberOfComponents() == 1) {
376 mesh->cellData.insertNextScalarData(
378 std::string(cellData->GetArrayName(i)));
379 auto &scalars = *(mesh->cellData.getScalarData(i));
380 scalars.resize(cellData->GetNumberOfTuples());
381 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
382 scalars[j] = dataArray->GetTuple1(j);
383 }
384 } else if (cellData->GetNumberOfComponents() == 3) {
385 mesh->cellData.insertNextVectorData(
387 std::string(cellData->GetArrayName(i)));
388 auto &vectors = *(mesh->cellData.getVectorData(i));
389 vectors.resize(cellData->GetNumberOfTuples());
390 for (unsigned j = 0; j < dataArray->GetNumberOfTuples(); ++j) {
391 std::array<double, 3> vector;
392 dataArray->GetTuple(j, &(vector[0]));
393 vectors[j] = vector;
394 }
395 }
396 }
397 }
398
399#endif // VIENNALS_USE_VTK
400
401 void readVTKLegacy(std::string filename) {
402 if (mesh == nullptr) {
403 Logger::getInstance()
404 .addWarning("No mesh was passed to VTKReader.")
405 .print();
406 return;
407 }
408
409 mesh->clear();
410 // open geometry file
411 std::ifstream f(filename.c_str());
412 if (!f)
413 Logger::getInstance().addError("Could not open geometry file!");
414 std::string temp;
415
416 // Check if geometry is an unstructured grid as required
417 while (std::getline(f, temp)) {
418 if (temp.find("DATASET") != std::string::npos)
419 break;
420 }
421 if (temp.find("UNSTRUCTURED_GRID") == std::string::npos) {
422 Logger::getInstance().addError("DATASET is not an UNSTRUCTURED_GRID!");
423 }
424
425 // Find POINTS in file to know number of nodes to read in
426 while (std::getline(f, temp)) {
427 if (temp.find("POINTS") != std::string::npos)
428 break;
429 }
430 int num_nodes = atoi(&temp[temp.find(" ") + 1]);
431
432 mesh->nodes.resize(num_nodes);
433
434 for (int i = 0; i < num_nodes; i++) {
435 double coords[3];
436
437 for (int j = 0; j < 3; j++)
438 f >> coords[j];
439 for (int j = 0; j < 3; j++) {
440 mesh->nodes[i][j] = coords[j];
441 // int shift_size = shift.size();
442 // if (shift_size > j)
443 // mesh->nodes[i][j] += shift[j]; // Assign desired shift
444 // if (InputTransformationSigns[j])
445 // mesh->nodes[i][j] = -mesh->nodes[i][j]; // Assign sign
446 // transformation, if needed
447 // mesh->nodes[i][j] *= scale; // Scale the geometry according to
448 // parameters file
449 }
450 }
451
452 while (std::getline(f, temp)) {
453 if (temp.find("CELLS") == 0)
454 break;
455 }
456
457 int num_elems = atoi(&temp[temp.find(" ") + 1]);
458
459 std::ifstream f_ct(filename.c_str()); // stream to read cell CELL_TYPES
460 std::ifstream f_m(
461 filename.c_str()); // stream for material numbers if they exist
462
463 // advance to cell types and check if there are the right number
464 while (std::getline(f_ct, temp)) {
465 if (temp.find("CELL_TYPES") == 0)
466 break;
467 }
468 int num_cell_types = atoi(&temp[temp.find(" ") + 1]);
469 // need a cell_type for each cell
470 if (num_elems != num_cell_types) {
471 Logger::getInstance().addError(
472 "Corrupt input geometry! Number of CELLS and CELL_TYPES "
473 "is different!");
474 }
475
476 bool is_material = true;
477 // advance to material if it is specified
478 while (std::getline(f_m, temp)) {
479 if (temp.find("CELL_DATA") != std::string::npos) {
480 std::getline(f_m, temp);
481 if ((temp.find("SCALARS material") != std::string::npos) ||
482 (temp.find("SCALARS Material") != std::string::npos)) {
483 std::getline(f_m, temp);
484 break;
485 }
486 }
487 }
488 if (f_m.eof()) {
489 is_material = false;
490 }
491
492 // maximum cell to read in is a tetra with 4 points
493 std::vector<hrleVectorType<unsigned int, 4>> elements;
494 elements.reserve(num_elems);
495
496 std::vector<double> materials;
497
498 unsigned elems_fake;
499 unsigned cell_type;
500 unsigned cell_material;
501 for (int i = 0; i < num_elems; i++) {
502 f >> elems_fake;
503 f_ct >> cell_type;
504 if (is_material)
505 f_m >> cell_material;
506 else
507 cell_material =
508 1; // if there are no materials specified make all the same
509
510 // check if the correct number of nodes for cell_type is given
511 unsigned number_nodes = vtk_nodes_for_cell_type[cell_type];
512 if (number_nodes == elems_fake || number_nodes == 0) {
513 // check for different types to subdivide them into supported types
514 switch (cell_type) {
515 case 1: {
516 std::array<unsigned, 1> elem;
517 f >> elem[0];
518 mesh->template getElements<1>().push_back(elem);
519 materials.push_back(cell_material);
520 break;
521 }
522 case 3: {
523 std::array<unsigned, 2> elem;
524 for (unsigned j = 0; j < number_nodes; ++j) {
525 f >> elem[j];
526 }
527 mesh->template getElements<2>().push_back(elem);
528 materials.push_back(cell_material);
529 break;
530 }
531 case 5: // triangle for 2D
532 {
533 std::array<unsigned, 3> elem;
534 for (unsigned j = 0; j < number_nodes; ++j) {
535 f >> elem[j];
536 }
537 mesh->template getElements<3>().push_back(elem);
538 materials.push_back(cell_material);
539 break;
540 }
541
542 case 10: // tetra for 3D
543 {
544 std::array<unsigned, 4> elem;
545 for (unsigned j = 0; j < number_nodes; ++j) {
546 f >> elem[j];
547 }
548 mesh->template getElements<4>().push_back(elem);
549 materials.push_back(cell_material);
550 break;
551 }
552
553 case 9: // this is a quad, so just plit it into two triangles
554 {
555 std::array<unsigned, 3> elem;
556 for (unsigned j = 0; j < 3; ++j) {
557 f >> elem[j];
558 }
559 mesh->template getElements<3>().push_back(
560 elem); // push the first three nodes as a triangle
561 materials.push_back(cell_material);
562
563 f >> elem[1]; // replace middle element to create other triangle
564 mesh->template getElements<3>().push_back(elem);
565 materials.push_back(cell_material);
566 break;
567 }
568
569 default:
570 std::ostringstream oss;
571 oss << "VTK Cell type " << cell_type
572 << " is not supported. Cell ignored..." << std::endl;
573 Logger::getInstance().addWarning(oss.str()).print();
574 }
575 } else {
576 std::ostringstream oss;
577 oss << "INVALID CELL TYPE! Expected number of nodes: " << number_nodes
578 << ", Found number of nodes: " << elems_fake
579 << "; Ignoring element...";
580 Logger::getInstance().addError(oss.str());
581 // ignore rest of lines
582 f.ignore(std::numeric_limits<std::streamsize>::max(), '\n');
583 }
584 }
585
586 mesh->cellData.insertNextScalarData(materials, "Material");
587
588 // Now read Cell Data
589 int num_cell_data = 0;
590 while (std::getline(f, temp)) {
591 if (temp.find("CELL_DATA") != std::string::npos) {
592 num_cell_data = atoi(&temp[temp.find(" ") + 1]);
593 break;
594 }
595 }
596 std::cout << "Read cell data: " << num_cell_data << std::endl;
597
598 while (!f.eof()) {
599 std::cout << "reading scalar data" << std::endl;
600 while (std::getline(f, temp)) {
601 if (temp.find("SCALARS") != std::string::npos)
602 break;
603 }
604 if (f.eof()) {
605 break;
606 }
607
608 std::string scalarDataName;
609 {
610 int firstS = temp.find(" ") + 1;
611 int secondS = temp.find(" ", firstS + 1);
612 scalarDataName = temp.substr(firstS, secondS - firstS);
613 }
614 std::vector<double> scalarData;
615
616 // consume one line, which defines the lookup table
617 std::getline(f, temp);
618 if (temp.compare("LOOKUP_TABLE default") != 0) {
619 Logger::getInstance()
620 .addWarning("Wrong lookup table for VTKLegacy: " + temp)
621 .print();
622 }
623
624 // now read scalar values
625 for (int i = 0; i < num_cell_data; ++i) {
626 double data;
627 f >> data;
628 scalarData.push_back(data);
629 }
630
631 mesh->cellData.insertNextScalarData(scalarData, scalarDataName);
632 }
633
634 f_ct.close();
635 f_m.close();
636 f.close();
637 }
638};
639
640} // 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
std::vector< std::array< T, 3 > > VectorDataType
Definition lsPointData.hpp:23
std::vector< T > ScalarDataType
Definition lsPointData.hpp:22
Class handling the import of VTK file types.
Definition lsVTKReader.hpp:28
void setFileName(std::string passedFileName)
set file name for file to read
Definition lsVTKReader.hpp:55
void setMesh(SmartPointer< Mesh<> > passedMesh)
set the mesh the file should be read into
Definition lsVTKReader.hpp:49
VTKReader(SmartPointer< Mesh< T > > passedMesh, std::string passedFileName)
Definition lsVTKReader.hpp:41
VTKReader(SmartPointer< Mesh< T > > passedMesh)
Definition lsVTKReader.hpp:39
VTKReader(SmartPointer< Mesh<> > passedMesh, FileFormatEnum passedFormat, std::string passedFileName)
Definition lsVTKReader.hpp:44
void setFileFormat(FileFormatEnum passedFormat)
set file format for file to read. Defaults to VTK_LEGACY.
Definition lsVTKReader.hpp:52
void apply()
Definition lsVTKReader.hpp:57
VTKReader()
Definition lsVTKReader.hpp:37
mesh
Definition AirGapDeposition.py:36
Definition lsAdvect.hpp:46
FileFormatEnum
Definition lsFileFormats.hpp:4