TTK
Loading...
Searching...
No Matches
ttkUtils.cpp
Go to the documentation of this file.
1#include <ttkUtils.h>
2
3#include <vtkAbstractArray.h>
4#include <vtkCellArray.h>
5#include <vtkDoubleArray.h>
6#include <vtkFieldData.h>
7#include <vtkIdTypeArray.h>
8#include <vtkNew.h>
9#include <vtkPoints.h>
10#include <vtkPolyData.h>
11#include <vtkSmartPointer.h>
12#include <vtkStringArray.h>
13#include <vtkUnstructuredGrid.h>
14
15#ifdef TTK_ENABLE_MPI_TIME
16#include <mpi.h>
17#endif
18int ttkUtils::replaceVariable(const std::string &iString,
19 vtkFieldData *fieldData,
20 std::string &oString,
21 std::string &errorMsg) {
22 std::string varName = iString;
23 int varIndex = -1;
24 bool varIndexDefined = false;
25
26 // Check if varIndex is specified
27 size_t const indexDelimiter0 = iString.find('[');
28 size_t const indexDelimiter1 = iString.find(']');
29 if(indexDelimiter0 != std::string::npos
30 && indexDelimiter1 != std::string::npos) {
31 if(indexDelimiter0 > indexDelimiter1
32 || iString.find('[', indexDelimiter0 + 1) != std::string::npos
33 || iString.find('}', indexDelimiter1 + 1) != std::string::npos) {
34 errorMsg = "Invalid Syntax:\n" + iString;
35 return 0;
36 }
37
38 varName = iString.substr(0, indexDelimiter0);
39 varIndex = stoi(iString.substr(
40 indexDelimiter0 + 1, indexDelimiter1 - indexDelimiter0 - 1));
41 varIndexDefined = true;
42 }
43
44 // Search in candidates
45 auto column = fieldData->GetAbstractArray(varName.data());
46 if(column == nullptr) {
47 errorMsg = "FieldData does not contain array '" + varName + "'";
48 return 0;
49 }
50
51 size_t const n = column->GetNumberOfTuples();
52 size_t const m = column->GetNumberOfComponents();
53 int const s = n * m;
54
55 if(!varIndexDefined) {
56 if(s > 0) {
57 oString = column->GetVariantValue(0).ToString();
58 for(int i = 1; i < s; i++)
59 oString += "," + column->GetVariantValue(i).ToString();
60 }
61 } else {
62 if(varIndex < 0 || varIndex >= s) {
63 errorMsg = "Index " + std::to_string(varIndex) + "/" + std::to_string(s)
64 + " for FieldData Array '" + varName + "' out of range";
65 return 0;
66 }
67 oString = column->GetVariantValue(varIndex).ToString();
68 }
69
70 return 1;
71}
72
73int ttkUtils::replaceVariables(const std::string &iString,
74 vtkFieldData *fieldData,
75 std::string &oString,
76 std::string &errorMsg) {
77 oString = iString;
78
79 while(oString.find('{') != std::string::npos
80 && oString.find('}') != std::string::npos) {
81 size_t o = oString.find('{');
82 size_t const c = oString.find('}');
83 // {...{....{...}...}..}
84 // | |
85 // o c
86
87 size_t oNext = oString.find('{', o + 1);
88 while(oNext != std::string::npos && oNext < c) {
89 o = oNext;
90 oNext = oString.find('{', o + 1);
91 }
92
93 // {...{....{var}...}..}
94 // | |
95 // o c
96 std::string const var = oString.substr(o + 1, c - 1 - o);
97
98 std::string rVar;
99 if(!replaceVariable(var, fieldData, rVar, errorMsg))
100 return 0;
101
102 oString = oString.substr(0, o).append(rVar).append(
103 oString.substr(c + 1, oString.length() - c - 1));
104 }
105
106 if(oString.find('{') != std::string::npos
107 || oString.find('}') != std::string::npos) {
108 errorMsg = "Invalid Syntax:\n" + iString;
109 return 0;
110 }
111
112 return 1;
113}
114
115int ttkUtils::stringListToVector(const std::string &iString,
116 std::vector<std::string> &v) {
117 // ...,...,....,
118 // | |
119 // i j
120 size_t i = 0;
121 size_t j = iString.find(',');
122 while(j != std::string::npos) {
123 v.push_back(iString.substr(i, j - i));
124 i = j + 1;
125 j = iString.find(',', i);
126 }
127 if(iString.length() > i)
128 v.push_back(iString.substr(i, iString.length() - i));
129
130 return 1;
131}
132
133int ttkUtils::stringListToDoubleVector(const std::string &iString,
134 std::vector<double> &v) {
135 std::vector<std::string> stringVector;
136 if(!ttkUtils::stringListToVector(iString, stringVector))
137 return 0;
138
139 size_t const n = stringVector.size();
140 v.resize(n);
141 // try {
142 for(size_t i = 0; i < n; i++)
143 v[i] = stod(stringVector[i]);
144 // } catch(std::invalid_argument &e) {
145 // return 0;
146 // }
147
148 return 1;
149}
150
152 ttkUtils::csvToVtkArray(const std::string &line) {
153 size_t const firstComma = line.find(',', 0);
154
155 if(firstComma == std::string::npos)
156 return nullptr;
157
158 std::string const arrayName = line.substr(0, firstComma);
159
160 std::vector<std::string> valuesAsString;
162 line.substr(firstComma + 1, std::string::npos), valuesAsString);
163 size_t const nValues = valuesAsString.size();
164 if(nValues < 1)
165 return nullptr;
166
167 // Check if all elements are numbers
168 bool isNumeric = true;
169
170 std::vector<double> valuesAsDouble(nValues);
171 try {
172 for(size_t i = 0; i < nValues; i++)
173 valuesAsDouble[i] = std::stod(valuesAsString[i]);
174 } catch(const std::invalid_argument &) {
175 isNumeric = false;
176 } catch(const std::out_of_range &) {
177 isNumeric = false;
178 }
179
180 if(isNumeric) {
182 array->SetName(arrayName.data());
183 array->SetNumberOfValues(nValues);
184 for(size_t i = 0; i < nValues; i++)
185 array->SetValue(i, valuesAsDouble[i]);
186 return array;
187 } else {
189 array->SetName(arrayName.data());
190 array->SetNumberOfValues(nValues);
191 for(size_t i = 0; i < nValues; i++)
192 array->SetValue(i, valuesAsString[i]);
193 return array;
194 }
195}
196
198 ttkUtils::csvToDoubleArray(const std::string &line) {
199 size_t const firstComma = line.find(',', 0);
200
201 if(firstComma == std::string::npos)
202 return nullptr;
203
204 std::string const arrayName = line.substr(0, firstComma);
205 std::string const valuesAsString
206 = line.substr(firstComma + 1, std::string::npos);
207
208 std::vector<double> values;
209 ttkUtils::stringListToDoubleVector(valuesAsString, values);
210 size_t const n = values.size();
211
213 array->SetName(arrayName.data());
214 array->SetNumberOfComponents(1);
215 array->SetNumberOfTuples(n);
216 auto arrayData = reinterpret_cast<double *>(GetVoidPointer(array));
217 for(size_t i = 0; i < n; i++)
218 arrayData[i] = values[i];
219
220 return array;
221}
222
226void *ttkUtils::GetVoidPointer(vtkDataArray *array, vtkIdType start) {
227 void *outPtr = nullptr;
228 if(array == nullptr)
229 return outPtr;
230
231 switch(array->GetDataType()) {
232 vtkTemplateMacro(
233 auto *aosArray = vtkAOSDataArrayTemplate<VTK_TT>::FastDownCast(array);
234 if(aosArray) { outPtr = aosArray->GetVoidPointer(start); });
235 }
236 return outPtr;
237}
238
239void *ttkUtils::GetVoidPointer(vtkPoints *points, vtkIdType start) {
240 return GetVoidPointer(points->GetData(), start);
241}
242
244 vtkIdType idx) {
245 auto slicedArray
246 = vtkSmartPointer<vtkAbstractArray>::Take(array->NewInstance());
247 slicedArray->SetName(array->GetName());
248 slicedArray->SetNumberOfComponents(array->GetNumberOfComponents());
249 slicedArray->SetNumberOfTuples(1);
250 slicedArray->SetTuple(0, idx, array);
251 return slicedArray;
252}
253
254void *ttkUtils::WriteVoidPointer(vtkDataArray *array,
255 vtkIdType valueIdx,
256 vtkIdType numValues) {
257 void *outPtr = nullptr;
258 switch(array->GetDataType()) {
259 vtkTemplateMacro(auto *aosArray
260 = vtkAOSDataArrayTemplate<VTK_TT>::FastDownCast(array);
261 if(aosArray) {
262 outPtr = aosArray->WriteVoidPointer(valueIdx, numValues);
263 });
264 }
265 return outPtr;
266}
267
268void *ttkUtils::WritePointer(vtkDataArray *array,
269 vtkIdType valueIdx,
270 vtkIdType numValues) {
271 void *outPtr = nullptr;
272 switch(array->GetDataType()) {
273 vtkTemplateMacro(
274 auto *aosArray = vtkAOSDataArrayTemplate<VTK_TT>::FastDownCast(array);
275 if(aosArray) { outPtr = aosArray->WritePointer(valueIdx, numValues); });
276 }
277 return outPtr;
278}
279
280void ttkUtils::SetVoidArray(vtkDataArray *array,
281 void *data,
282 vtkIdType size,
283 int save) {
284 switch(array->GetDataType()) {
285 vtkTemplateMacro(
286 auto *aosArray = vtkAOSDataArrayTemplate<VTK_TT>::FastDownCast(array);
287 if(aosArray) { aosArray->SetVoidArray(data, size, save); } else {
288 std::cerr << "SetVoidArray on incompatible vtkDataArray:" << endl;
289 array->Print(std::cerr);
290 });
291 }
292}
293
294[[deprecated]] void ttkUtils::FillCellArrayFromSingle(vtkIdType const *cells,
295 vtkIdType ncells,
296 vtkCellArray *cellArray) {
297 size_t curPos = 0;
298 vtkNew<vtkIdList> verts;
299 for(vtkIdType cid = 0; cid < ncells; cid++) {
300 const vtkIdType nbVerts = cells[curPos];
301 verts->SetNumberOfIds(nbVerts);
302 curPos++;
303 for(vtkIdType v = 0; v < nbVerts; v++) {
304 verts->SetId(v, cells[curPos]);
305 curPos++;
306 }
307 cellArray->InsertNextCell(verts);
308 }
309}
310
311void ttkUtils::FillCellArrayFromDual(vtkIdType const *cells_co,
312 vtkIdType const *cells_off,
313 vtkIdType ncells,
314 vtkCellArray *cellArray) {
315 size_t curPos = 0;
316 vtkNew<vtkIdList> verts;
317 for(vtkIdType cid = 0; cid < ncells; cid++) {
318 const vtkIdType nbVerts = cells_off[cid + 1] - cells_off[cid];
319 verts->SetNumberOfIds(nbVerts);
320 for(vtkIdType v = 0; v < nbVerts; v++) {
321 verts->SetId(v, cells_co[curPos]);
322 curPos++;
323 }
324 cellArray->InsertNextCell(verts);
325 }
326}
327
328int ttkUtils::CellVertexFromPoints(vtkDataSet *const dataSet,
329 vtkPoints *const points) {
330
331 if(dataSet == nullptr || points == nullptr) {
332 return 0;
333 }
334
335 if(!dataSet->IsA("vtkUnstructuredGrid") && !dataSet->IsA("vtkPolyData")) {
336 return 0;
337 }
338
339 const size_t nPoints = points->GetNumberOfPoints();
340 if(nPoints == 0) {
341 return 0;
342 }
343
344 vtkNew<vtkIdTypeArray> offsets{};
345 offsets->SetNumberOfTuples(nPoints + 1);
346 auto offsetsData = ttkUtils::GetPointer<vtkIdType>(offsets);
347 for(size_t i = 0; i <= nPoints; i++)
348 offsetsData[i] = i;
349
350 vtkNew<vtkIdTypeArray> connectivity{};
351 connectivity->SetNumberOfTuples(nPoints);
352 auto connectivityData = ttkUtils::GetPointer<vtkIdType>(connectivity);
353 for(size_t i = 0; i < nPoints; i++)
354 connectivityData[i] = i;
355
356 vtkNew<vtkCellArray> cells{};
357 cells->SetData(offsets, connectivity);
358
359 if(dataSet->IsA("vtkUnstructuredGrid")) {
360 const auto vtu{vtkUnstructuredGrid::SafeDownCast(dataSet)};
361 if(vtu != nullptr) {
362 vtu->SetPoints(points);
363 vtu->SetCells(VTK_VERTEX, cells);
364 }
365 } else if(dataSet->IsA("vtkPolyData")) {
366 const auto vtp{vtkPolyData::SafeDownCast(dataSet)};
367 if(vtp != nullptr) {
368 vtp->SetPoints(points);
369 vtp->SetVerts(cells);
370 }
371 }
372
373 return 1;
374}
static int replaceVariables(const std::string &iString, vtkFieldData *fieldData, std::string &oString, std::string &errorMsg)
Definition ttkUtils.cpp:73
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static int stringListToDoubleVector(const std::string &iString, std::vector< double > &v)
Definition ttkUtils.cpp:133
static vtkSmartPointer< vtkAbstractArray > SliceArray(vtkAbstractArray *array, vtkIdType idx)
Definition ttkUtils.cpp:243
static void FillCellArrayFromSingle(vtkIdType const *cells, vtkIdType ncells, vtkCellArray *cellArray)
Definition ttkUtils.cpp:294
static void * WritePointer(vtkDataArray *array, vtkIdType start, vtkIdType numValues)
Definition ttkUtils.cpp:268
static int replaceVariable(const std::string &iString, vtkFieldData *fieldData, std::string &oString, std::string &errorMsg)
Definition ttkUtils.cpp:18
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
Definition ttkUtils.cpp:328
static vtkSmartPointer< vtkAbstractArray > csvToVtkArray(const std::string &line)
Definition ttkUtils.cpp:152
static vtkSmartPointer< vtkDoubleArray > csvToDoubleArray(const std::string &line)
Definition ttkUtils.cpp:198
static void FillCellArrayFromDual(vtkIdType const *cells_co, vtkIdType const *cells_off, vtkIdType ncells, vtkCellArray *cellArray)
Definition ttkUtils.cpp:311
static void SetVoidArray(vtkDataArray *array, void *data, vtkIdType size, int save)
Definition ttkUtils.cpp:280
static void * WriteVoidPointer(vtkDataArray *array, vtkIdType start, vtkIdType numValues)
Definition ttkUtils.cpp:254
static int stringListToVector(const std::string &iString, std::vector< std::string > &v)
Definition ttkUtils.cpp:115
std::string to_string(__int128)
Definition ripserpy.cpp:99