TTK
Loading...
Searching...
No Matches
ttkOFFReader.cpp
Go to the documentation of this file.
1#include <ttkOFFReader.h>
2
3#include <vtkCellData.h>
4#include <vtkDoubleArray.h>
5#include <vtkIdList.h>
6#include <vtkInformation.h>
7#include <vtkInformationVector.h>
8#include <vtkNew.h>
9#include <vtkObjectFactory.h>
10#include <vtkPointData.h>
11#include <vtkPoints.h>
12#include <vtkStreamingDemandDrivenPipeline.h>
13#include <vtkUnstructuredGrid.h>
14
15#include <iostream>
16#include <sstream>
17#include <string>
18#include <vector>
19
21
22// Public
23// {{{
24
25void ttkOFFReader::PrintSelf(std::ostream &os, vtkIndent indent) {
26 this->Superclass::PrintSelf(os, indent);
27
28 os << indent << "File Name: " << (this->FileName ? this->FileName : "(none)")
29 << std::endl;
30}
31
32// }}}
33// Protected
34// {{{
35
37 this->setDebugMsgPrefix("OFFReader");
38 this->SetNumberOfInputPorts(0);
39 this->SetNumberOfOutputPorts(1);
40}
41
42static int countVertsData(const std::string &line) {
43 std::istringstream ss(line);
44 double buffer;
45 int nbFields = -1;
46 while(!ss.fail()) {
47 ss >> buffer;
48 ++nbFields;
49 }
50 // remove the 3 coordinates
51 return nbFields - 3;
52}
53
54static int countCellsData(const std::string &line) {
55 std::istringstream ss(line);
56 double buffer;
57 int nbFields = -1;
58 int sizeCell;
59 // This line start by the size of the cell
60 ss >> sizeCell;
61 // we read all the fields
62 while(!ss.fail()) {
63 ss >> buffer;
64 ++nbFields;
65 }
66 // remove the vertices of the cell
67 return nbFields - sizeCell;
68}
69
70static int processLineVert(vtkIdType curLine,
71 std::string &line,
72 vtkPoints *points,
73 std::vector<vtkNew<vtkDoubleArray>> &vertScalars,
74 const vtkIdType nbVertsData) {
75
76 double x, y, z;
77 std::istringstream ss(line);
78
79 // Coords
80 ss >> x >> y >> z;
81
82 points->InsertNextPoint(x, y, z);
83
84 // Scalars data
85 for(vtkIdType i = 0; i < nbVertsData; i++) {
86 double scalar;
87 ss >> scalar;
88 vertScalars[i]->SetTuple1(curLine, scalar);
89 }
90
91 return ++curLine;
92}
93
94static int processLineCell(vtkIdType curLine,
95 std::string &line,
96 vtkUnstructuredGrid *mesh,
97 std::vector<vtkNew<vtkDoubleArray>> &cellScalars,
98 const vtkIdType nbVerts,
99 const vtkIdType nbCellsData,
100 const ttk::Debug &dbg) {
101 int nbCellVerts;
102 vtkNew<vtkIdList> cellVerts{};
103 std::istringstream ss(line);
104 ss >> nbCellVerts;
105 for(int j = 0; j < nbCellVerts; j++) {
106 vtkIdType id;
107 ss >> id;
108 cellVerts->InsertNextId(id);
109 }
110 switch(nbCellVerts) {
111 case 2:
112 mesh->InsertNextCell(VTK_LINE, cellVerts);
113 break;
114 case 3:
115 mesh->InsertNextCell(VTK_TRIANGLE, cellVerts);
116 break;
117 case 4:
118 mesh->InsertNextCell(VTK_TETRA, cellVerts);
119 break;
120 default:
121 dbg.printErr("Unsupported cell type having " + std::to_string(nbCellVerts)
122 + " vertices");
123 return -3;
124 }
125
126 // Scalars data
127 for(vtkIdType i = 0; i < nbCellsData; i++) {
128 double scalar;
129 ss >> scalar;
130 // Currline is after having read all the vertices
131 cellScalars[i]->SetTuple1(curLine - nbVerts, scalar);
132 }
133
134 return ++curLine;
135}
136
137int ttkOFFReader::RequestData(vtkInformation *ttkNotUsed(request),
138 vtkInformationVector **ttkNotUsed(inputVector),
139 vtkInformationVector *outputVector) {
140 std::ifstream offFile(FileName, ios::in);
141
142 if(!offFile) {
143 this->printErr("Can't read file: '" + std::string{FileName} + "'");
144 return -1;
145 }
146
147 std::string FileType;
148
149 offFile >> FileType;
150 if(FileType != "OFF") {
151 this->printErr("Bad format for file: '" + std::string{FileName} + "'");
152 return -2;
153 }
154
155 vtkIdType curLine;
156 std::string line;
157 vtkIdType nbVerts{}, nbCells{}, nbVertsData{}, nbCellsData{};
158
159 // init values
160 offFile >> nbVerts >> nbCells;
161 std::getline(offFile, line);
162
163 if(nbVerts == 0) {
164 // empty file
165 return 0;
166 }
167
168 curLine = 0;
169
170 // count numbers of vertices scalars
171 std::getline(offFile, line);
172 nbVertsData = countVertsData(line);
173
174 // allocation verts
175 std::vector<vtkNew<vtkDoubleArray>> vertScalars{};
176 vertScalars.resize(nbVertsData);
177 for(vtkIdType i = 0; i < nbVertsData; i++) {
178 vertScalars[i]->SetNumberOfComponents(1);
179 vertScalars[i]->SetNumberOfTuples(nbVerts);
180 const std::string name = "VertScalarField_" + std::to_string(i);
181 vertScalars[i]->SetName(name.c_str());
182 }
183
184 vtkNew<vtkPoints> const points{};
185
186 // read vertices
187 while(
188 (curLine = processLineVert(curLine, line, points, vertScalars, nbVertsData))
189 < nbVerts) {
190 std::getline(offFile, line);
191 }
192
193 // add verts data to the mesh
194 vtkNew<vtkUnstructuredGrid> mesh{};
195 mesh->SetPoints(points);
196 for(const auto &scalarArray : vertScalars) {
197 mesh->GetPointData()->AddArray(scalarArray);
198 }
199
200 // count numbers of cells scalars
201 std::getline(offFile, line);
202 nbCellsData = countCellsData(line);
203
204 // allocation cells
205 std::vector<vtkNew<vtkDoubleArray>> cellScalars{};
206 cellScalars.resize(nbCellsData);
207 for(vtkIdType i = 0; i < nbCellsData; i++) {
208 cellScalars[i]->SetNumberOfComponents(1);
209 cellScalars[i]->SetNumberOfTuples(nbCells);
210 const std::string name = "CellScalarField_" + std::to_string(i);
211 cellScalars[i]->SetName(name.c_str());
212 }
213
214 // read cells
215 if(nbCells != 0) {
216 while((curLine = processLineCell(
217 curLine, line, mesh, cellScalars, nbVerts, nbCellsData, *this))
218 < nbVerts + nbCells) {
219 std::getline(offFile, line);
220 }
221 }
222
223 // add cell data to the mesh
224 for(const auto &scalarArray : cellScalars) {
225 mesh->GetCellData()->AddArray(scalarArray);
226 }
227
228#ifndef NDEBUG
229 this->printMsg("Read " + std::to_string(mesh->GetNumberOfPoints())
230 + " vertice(s)");
231 this->printMsg("Read " + std::to_string(mesh->GetNumberOfCells())
232 + " cell(s)");
233#endif
234
235 // get the info object
236 vtkInformation *outInfo = outputVector->GetInformationObject(0);
237 outInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(), 1);
238
239 // Set the output
240 vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
241 outInfo->Get(vtkDataObject::DATA_OBJECT()));
242
243 output->ShallowCopy(mesh);
244
245 return 1;
246}
247
248// }}}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttkOFFReader - Object File Format Reader
void PrintSelf(std::ostream &os, vtkIndent indent) override
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Minimalist debugging class.
Definition Debug.h:88
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
vtkStandardNewMacro(ttkOFFReader)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)