TTK
Loading...
Searching...
No Matches
ttkCinemaImaging.cpp
Go to the documentation of this file.
1#include <ttkCinemaImaging.h>
2
3#include <vtkInformation.h>
4
5#include <vtkImageData.h>
6#include <vtkMultiBlockDataSet.h>
7#include <vtkPolyData.h>
8#include <vtkUnstructuredGrid.h>
9
10#include <vtkCellData.h>
11#include <vtkPointData.h>
12
13#include <ttkUtils.h>
14#include <vtkDoubleArray.h>
15
18#include <ttkCinemaImagingVTK.h>
19
20#include <CinemaImaging.h>
21
23
25 this->setDebugMsgPrefix("CinemaImaging");
26 this->SetNumberOfInputPorts(2);
27 this->SetNumberOfOutputPorts(1);
28};
29
31;
32
33int ttkCinemaImaging::FillInputPortInformation(int port, vtkInformation *info) {
34 if(port == 0) {
35 info->Remove(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE());
36 info->Append(
37 vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
38 info->Append(
39 vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
40 info->Append(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPolyData");
41 } else if(port == 1)
42 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPointSet");
43 else
44 return 0;
45 return 1;
46};
47
49 vtkInformation *info) {
50 if(port == 0) {
51 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
52 return 1;
53 }
54 return 0;
55};
56
57int ttkCinemaImaging::RequestData(vtkInformation *ttkNotUsed(request),
58 vtkInformationVector **inputVector,
59 vtkInformationVector *outputVector) {
60
61 auto inputObject = vtkDataObject::GetData(inputVector[0]);
62 auto inputGrid = vtkPointSet::GetData(inputVector[1]);
63 auto outputCollection = vtkMultiBlockDataSet::GetData(outputVector);
64 if(!inputObject || !inputGrid) {
65 this->printErr("Unsupported input object types.");
66 return 0;
67 }
68
69 auto inputObjectAsMB = vtkSmartPointer<vtkMultiBlockDataSet>::New();
70 if(inputObject->IsA("vtkMultiBlockDataSet"))
71 inputObjectAsMB->ShallowCopy(inputObject);
72 else
73 inputObjectAsMB->SetBlock(0, inputObject);
74 const size_t nInputObjects = inputObjectAsMB->GetNumberOfBlocks();
75
76 // get default parameters
77 std::vector<double> defaultFocalPoint{
78 this->FocalPoint[0], this->FocalPoint[1], this->FocalPoint[2]};
79 std::vector<double> defaultNearFar{this->NearFar[0], this->NearFar[1]};
80 double defaultHeight = this->Height;
81 double const defaultAngle = this->Angle;
82
83 if(this->AutoFocalPoint || this->AutoNearFar || this->AutoHeight) {
84
85 double objectBounds[6];
86 inputObjectAsMB->GetBounds(objectBounds);
87
88 auto norm = [](const double v[3]) {
89 return std::sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
90 };
91
92 const double d[3]{objectBounds[1] - objectBounds[0],
93 objectBounds[3] - objectBounds[2],
94 objectBounds[5] - objectBounds[4]};
95 const double objectDiameter = norm(d);
96
97 const double c[3]{objectBounds[0] + 0.5 * d[0],
98 objectBounds[2] + 0.5 * d[1],
99 objectBounds[4] + 0.5 * d[2]};
100
101 if(this->AutoFocalPoint) {
102 defaultFocalPoint[0] = c[0];
103 defaultFocalPoint[1] = c[1];
104 defaultFocalPoint[2] = c[2];
105 }
106
107 if(this->AutoNearFar) {
108 double gridBounds[6];
109 inputGrid->GetBounds(gridBounds);
110
111 const double l[3]{
112 gridBounds[0] - c[0],
113 gridBounds[2] - c[1],
114 gridBounds[4] - c[2],
115 };
116 const double ld = norm(l);
117
118 defaultNearFar[0] = std::max(0.01, ld - objectDiameter * 0.5);
119 defaultNearFar[1] = ld + objectDiameter * 0.5;
120 }
121
122 if(this->AutoHeight) {
123 defaultHeight = objectDiameter;
124 }
125 }
126
127 // ensure grid contains all camera parameters as field data
128 auto aInputGrid
129 = vtkSmartPointer<vtkPointSet>::Take(inputGrid->NewInstance());
130 aInputGrid->ShallowCopy(inputGrid);
131 int const n = aInputGrid->GetNumberOfPoints();
132 auto aInputGridPD = aInputGrid->GetPointData();
133
135 aInputGridPD, "Resolution", n,
136 {(double)this->Resolution[0], (double)this->Resolution[1]});
138 aInputGridPD, "ProjectionMode", n, {(double)this->ProjectionMode});
139
141 aInputGridPD, "CamNearFar", n, defaultNearFar);
143 aInputGridPD, "CamHeight", n, {defaultHeight});
144 ttkCinemaImaging::EnsureGridData(aInputGridPD, "CamAngle", n, {defaultAngle});
145 ttkCinemaImaging::EnsureGridData(aInputGridPD, "CamUp", n, {0, 1, 0});
146
147 const bool gridHasDirection = aInputGridPD->HasArray("CamDirection");
148 const bool gridHasFocalPoint = aInputGridPD->HasArray("CamFocalPoint");
149 if(gridHasFocalPoint || !gridHasDirection) {
151 aInputGridPD, "CamFocalPoint", n, defaultFocalPoint);
153 } else {
155 aInputGridPD, "CamDirection", n, {0, 0, -1});
156 }
157
158 // print parameters
159 {
160 this->printMsg(
161 {{"#Objects", std::to_string(nInputObjects)},
162 {"Backend", this->Backend == 0 ? std::string("VTK_OPENGL")
163 : this->Backend == 1 ? std::string("EMBREE")
164 : std::string("NATIVE")},
165 {"Resolution", std::to_string(this->Resolution[0]) + " x "
166 + std::to_string(this->Resolution[1])},
167 {"Projection", this->ProjectionMode ? "Perspective" : "Orthographic"},
168 {"CamNearFar", std::to_string(defaultNearFar[0]) + " / "
169 + std::to_string(defaultNearFar[1])},
170 {"CamFocalPoint", gridHasDirection || gridHasFocalPoint
171 ? "Grid Data"
172 : "(" + std::to_string(defaultFocalPoint[0]) + ", "
173 + std::to_string(defaultFocalPoint[1]) + ", "
174 + std::to_string(defaultFocalPoint[2]) + ")"},
175 {std::string(this->ProjectionMode == 0 ? "CamHeight" : "CamAngle"),
176 std::to_string(this->ProjectionMode == 0 ? defaultHeight
177 : defaultAngle)}});
179 }
180
182 for(size_t b = 0; b < nInputObjects; b++) {
183 auto outputImages = vtkSmartPointer<vtkMultiBlockDataSet>::New();
184 auto inputPointSet
185 = vtkPointSet::SafeDownCast(inputObjectAsMB->GetBlock(b));
186 if(!inputPointSet->IsA("vtkUnstructuredGrid")
187 && !inputPointSet->IsA("vtkPolyData")) {
188 this->printErr("Unsupported input object type.");
189 return 0;
190 }
191 this->RequestDataSingle(outputImages,
192
193 inputPointSet, aInputGrid, defaultFocalPoint,
194 defaultNearFar, defaultHeight, defaultAngle);
195 outputAsMB->SetBlock(b, outputImages);
196 }
197
198 if(inputObject->IsA("vtkMultiBlockDataSet"))
199 outputCollection->ShallowCopy(outputAsMB);
200 else
201 outputCollection->ShallowCopy(outputAsMB->GetBlock(0));
202
203 return 1;
204}
205
207 vtkMultiBlockDataSet *outputImages,
208
209 vtkPointSet *inputObject,
210 vtkPointSet *inputGrid,
211 const std::vector<double> &ttkNotUsed(defaultFocalPoint),
212 const std::vector<double> &ttkNotUsed(defaultNearFar),
213 const double ttkNotUsed(defaultHeight),
214 const double ttkNotUsed(defaultAngle)) {
215 ttk::Timer const globalTimer;
216
217 auto cells = ttkCinemaImaging::GetCells(inputObject);
218 if(!cells)
219 return 0;
220
221 size_t const nTriangles = cells->GetNumberOfCells();
222 // make sure that cells consists only of triangles
223 {
224 auto offsets = static_cast<vtkIdType *>(
225 ttkUtils::GetVoidPointer(cells->GetOffsetsArray()));
226 vtkIdType j = 0;
227 for(size_t i = 0; i < nTriangles; i++, j += 3) {
228 if(j != offsets[i]) {
229 this->printErr("Connectivity list has to contain only triangles.");
230 return 0;
231 }
232 }
233 }
234
235 int status = 0;
236 if(this->Backend == 0) {
238 renderer.setDebugLevel(this->debugLevel_);
239 renderer.setThreadNumber(this->threadNumber_);
240 status = renderer.RenderVTKObject(outputImages, inputObject, inputGrid);
241 } else if(this->Backend == 1) {
243 renderer.setDebugLevel(this->debugLevel_);
244 renderer.setThreadNumber(this->threadNumber_);
245 status = renderer.RenderVTKObject(outputImages, inputObject, inputGrid);
246 } else {
248 renderer.setDebugLevel(this->debugLevel_);
249 renderer.setThreadNumber(this->threadNumber_);
250 status = renderer.RenderVTKObject(outputImages, inputObject, inputGrid);
251 }
252
253 if(!status)
254 return 0;
255
256 return 1;
257}
258
259vtkCellArray *ttkCinemaImaging::GetCells(vtkPointSet *pointSet) {
260 switch(pointSet->GetDataObjectType()) {
261 case VTK_UNSTRUCTURED_GRID:
262 return static_cast<vtkUnstructuredGrid *>(pointSet)->GetCells();
263 case VTK_POLY_DATA:
264 return static_cast<vtkPolyData *>(pointSet)->GetPolys();
265 }
266 return nullptr;
267};
268
270 vtkDataArray *array,
271 int tupleIdx,
272 const std::string &name) {
273 if(!array)
274 return 0;
275
276 size_t const nComponents = array->GetNumberOfComponents();
277
278 auto newArray = vtkSmartPointer<vtkDoubleArray>::New();
279 newArray->SetName(name.empty() ? array->GetName() : name.data());
280 newArray->SetNumberOfComponents(nComponents);
281 newArray->SetNumberOfTuples(1);
282
283 if(newArray->GetDataType() == array->GetDataType()) {
284 newArray->SetTuple(0, tupleIdx, array);
285 } else {
286 for(size_t i = 0; i < nComponents; i++)
287 newArray->SetValue(
288 i, array->GetVariantValue(tupleIdx * nComponents + i).ToDouble());
289 }
290
291 fd->AddArray(newArray);
292
293 return 1;
294};
295
297 vtkImageData *image,
298 int tupleIdx) {
299 auto imageFD = image->GetFieldData();
300
301 auto inputGridPD = inputGrid->GetPointData();
302 for(int i = 0; i < inputGridPD->GetNumberOfArrays(); i++) {
304 imageFD, inputGridPD->GetArray(i), tupleIdx);
305 }
306
308 imageFD, inputGrid->GetPoints()->GetData(), tupleIdx, "CamPosition");
309
310 return 1;
311};
312
314 auto pos
315 = static_cast<float *>(ttkUtils::GetVoidPointer(inputGrid->GetPoints()));
316 auto focal = static_cast<double *>(ttkUtils::GetVoidPointer(
317 inputGrid->GetPointData()->GetArray("CamFocalPoint")));
318
319 int const nTuples = inputGrid->GetNumberOfPoints();
320
321 auto newArray = vtkSmartPointer<vtkDoubleArray>::New();
322 newArray->SetName("CamDirection");
323 newArray->SetNumberOfComponents(3);
324 newArray->SetNumberOfTuples(nTuples);
325
326 auto dir = static_cast<double *>(ttkUtils::GetVoidPointer(newArray));
327
328 for(int i = 0, j = nTuples * 3; i < j; i++)
329 dir[i] = focal[i] - pos[i];
330
331 inputGrid->GetPointData()->AddArray(newArray);
332
333 return 1;
334};
335
337 const std::string &name,
338 int nTuples,
339 const std::vector<double> &defaultValues) {
340 auto array = vtkDoubleArray::SafeDownCast(fd->GetArray(name.data()));
341
342 if(!array) {
343 int const nComponents = defaultValues.size();
344
345 auto newArray = vtkSmartPointer<vtkDoubleArray>::New();
346 newArray->SetName(name.data());
347 newArray->SetNumberOfComponents(nComponents);
348 newArray->SetNumberOfTuples(nTuples);
349
350 auto data = static_cast<double *>(ttkUtils::GetVoidPointer(newArray));
351 for(int i = 0, q = 0; i < nTuples; i++)
352 for(int j = 0; j < nComponents; j++)
353 data[q++] = defaultValues[j];
354
355 fd->AddArray(newArray);
356 }
357
358 return 1;
359};
360
361int ttkCinemaImaging::Normalize(vtkDataArray *depthArray,
362 const double nearFar[2]) {
363 if(!depthArray->IsA("vtkFloatArray")
364 || depthArray->GetNumberOfComponents() != 1)
365 return 0;
366
367 if(nearFar[0] == 0.0 && nearFar[1] == 0.0)
368 return 1;
369
370 const size_t nValues = depthArray->GetNumberOfTuples();
371 auto data = static_cast<float *>(ttkUtils::GetVoidPointer(depthArray));
372
373 const float near = (float)nearFar[0];
374 const float far = (float)nearFar[1];
375 const float delta = far - near;
376
377 for(size_t i = 0; i < nValues; i++) {
378 if(std::isnan(data[i])) {
379 data[i] = 1.0;
380 } else {
381 data[i] = (data[i] - near) / delta;
382 if(data[i] > 1.0)
383 data[i] = 1.0;
384 else if(data[i] < 0.0)
385 data[i] = 0.0;
386 }
387 }
388
389 return 1;
390};
391
393 vtkImageData *outputImage,
394
395 vtkPointSet *inputObject,
396 const ttk::CinemaImaging *renderer,
397 const unsigned int *primitiveIdArray,
398 const float *barycentricCoordinates,
399 const vtkIdType *inputObjectConnectivityList) {
400
401 auto inputObjectPD = inputObject->GetPointData();
402 auto inputObjectCD = inputObject->GetCellData();
403 auto outputImagePD = outputImage->GetPointData();
404 int dim[3];
405 outputImage->GetDimensions(dim);
406 size_t const nPixels = dim[0] * dim[1];
407
408 const size_t nInputObjectPDArrays = inputObjectPD->GetNumberOfArrays();
409 const size_t nInputObjectCDArrays = inputObjectCD->GetNumberOfArrays();
410
411 int status = 0;
412
413 // Map Point Data
414 for(size_t j = 0; j < nInputObjectPDArrays; j++) {
415 auto inputArray = inputObjectPD->GetArray(j);
416 auto outputArray
417 = vtkSmartPointer<vtkDataArray>::Take(inputArray->NewInstance());
418 outputArray->SetName(inputArray->GetName());
419 outputArray->SetNumberOfComponents(inputArray->GetNumberOfComponents());
420 outputArray->SetNumberOfTuples(nPixels);
421
422 outputImagePD->AddArray(outputArray);
423
424 switch(outputArray->GetDataType()) {
425 vtkTemplateMacro(status = renderer->interpolateArray(
426 (VTK_TT *)ttkUtils::GetVoidPointer(outputArray),
427
428 primitiveIdArray, barycentricCoordinates,
429 inputObjectConnectivityList,
430
431 (const VTK_TT *)ttkUtils::GetVoidPointer(inputArray),
432 nPixels, inputArray->GetNumberOfComponents()));
433 }
434
435 if(!status)
436 return 0;
437 }
438
439 // Map Cell Data
440 for(size_t j = 0; j < nInputObjectCDArrays; j++) {
441 auto inputArray = inputObjectCD->GetArray(j);
442 auto outputArray
443 = vtkSmartPointer<vtkDataArray>::Take(inputArray->NewInstance());
444 outputArray->SetName(inputArray->GetName());
445 outputArray->SetNumberOfComponents(inputArray->GetNumberOfComponents());
446 outputArray->SetNumberOfTuples(nPixels);
447
448 outputImagePD->AddArray(outputArray);
449
450 switch(outputArray->GetDataType()) {
451 vtkTemplateMacro(status = renderer->lookupArray(
452 (VTK_TT *)ttkUtils::GetVoidPointer(outputArray),
453
454 primitiveIdArray,
455 (const VTK_TT *)ttkUtils::GetVoidPointer(inputArray),
456 nPixels, inputArray->GetNumberOfComponents()));
457 }
458
459 if(!status)
460 return 0;
461 }
462
463 return 1;
464};
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
TTK VTK-filter that generates images of a vtkDataSet.
static int AddFieldDataArray(vtkFieldData *fd, vtkDataArray *array, int tupleIdx, const std::string &name="")
static int EnsureGridData(vtkPointData *fd, const std::string &name, int nTuples, const std::vector< double > &Values)
static int AddAllFieldDataArrays(vtkPointSet *inputGrid, vtkImageData *image, int tupleIdx)
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
static int MapPointAndCellData(vtkImageData *outputImage, vtkPointSet *inputObject, const ttk::CinemaImaging *renderer, const unsigned int *primitiveIdArray, const float *barycentricCoordinates, const vtkIdType *inputObjectConnectivityList)
~ttkCinemaImaging() override
int RequestDataSingle(vtkMultiBlockDataSet *collection, vtkPointSet *object, vtkPointSet *grid, const std::vector< double > &defaultFocal, const std::vector< double > &defaultNearFar, const double defaultHeight, const double defaultAngle)
static vtkCellArray * GetCells(vtkPointSet *pointSet)
static int Normalize(vtkDataArray *depthArray, const double nearFar[2])
static int ComputeDirFromFocalPoint(vtkPointSet *inputGrid)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
virtual int setThreadNumber(const int threadNumber)
Definition BaseClass.h:80
TTK modules that generates images of a dataset.
int interpolateArray(DT *outputArray, const unsigned int *primitiveIds, const float *barycentricCoordinates, const IT *connectivityList, const DT *inputArray, const size_t &nTuples, const size_t &nComponents=1, const DT &missingValue=std::numeric_limits< DT >::has_quiet_NaN ? std::numeric_limits< DT >::quiet_NaN() :std::numeric_limits< DT >::max()) const
int lookupArray(DT *outputArray, const unsigned int *primitiveIds, const DT *inputArray, const size_t &nTuples, const size_t &nComponents=1, const DT &missingValue=std::numeric_limits< DT >::has_quiet_NaN ? std::numeric_limits< DT >::quiet_NaN() :std::numeric_limits< DT >::max()) const
int debugLevel_
Definition Debug.h:379
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
int RenderVTKObject(vtkMultiBlockDataSet *outputImages, vtkPointSet *inputObject, vtkPointSet *inputGrid) const
int RenderVTKObject(vtkMultiBlockDataSet *outputImages, vtkPointSet *inputObject, vtkPointSet *inputGrid) const
int RenderVTKObject(vtkMultiBlockDataSet *outputImages, vtkPointSet *inputObject, vtkPointSet *inputGrid) const
vtkStandardNewMacro(ttkCinemaImaging)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)