TTK
Loading...
Searching...
No Matches
ttkCinemaImagingVTK.cpp
Go to the documentation of this file.
2
3#include <ttkCinemaImaging.h>
4
5#include <vtkDataSetSurfaceFilter.h>
6#include <vtkSmartPointer.h>
7
8#include <vtkMultiBlockDataSet.h>
9#include <vtkPointSet.h>
10
11#include <vtkActor.h>
12#include <vtkCamera.h>
13#include <vtkCameraPass.h>
14#include <vtkOpenGLRenderer.h>
15#include <vtkPolyDataMapper.h>
16#include <vtkRenderPassCollection.h>
17#include <vtkRenderWindow.h>
18#include <vtkRenderer.h>
19#include <vtkSequencePass.h>
20#include <vtkValuePass.h>
21#include <vtkWindowToImageFilter.h>
22
23#include <vtkDoubleArray.h>
24#include <vtkFloatArray.h>
25#include <vtkSignedCharArray.h>
26
27#include <ttkUtils.h>
28#include <vtkCellData.h>
29#include <vtkPointData.h>
30
35
37 vtkPointSet *object,
38 vtkCamera *camera) const {
40 mapper->SetInputDataObject(object);
41
42 auto actor = vtkSmartPointer<vtkActor>::New();
43 actor->SetMapper(mapper);
44
45 renderer->SetBackground(0.0, 0.0, 0.0);
46 renderer->GradientBackgroundOff();
47 renderer->AddActor(actor);
48 renderer->SetActiveCamera(camera);
49
50 return 1;
51}
52
53int ttk::ttkCinemaImagingVTK::setupWindow(vtkRenderWindow *window,
54 vtkRenderer *renderer,
55 const double resolution[2]) const {
56 window->SetSize(resolution[0], resolution[1]);
57 window->SetMultiSamples(0); // Disable AA
58 window->OffScreenRenderingOn();
59 window->AddRenderer(renderer);
60
61 return 1;
62};
63
65 vtkPointSet *object,
66 int fieldType,
67 vtkRenderPassCollection *valuePassCollection,
68 std::vector<std::string> &valuePassNames) const {
69 auto fd = fieldType == 0 ? static_cast<vtkFieldData *>(object->GetPointData())
70 : static_cast<vtkFieldData *>(object->GetCellData());
71
72 for(size_t i = 0, j = fd->GetNumberOfArrays(); i < j; i++) {
73 auto array = fd->GetArray(i);
74 if(!array)
75 continue;
76
77 std::string const name(array->GetName());
78
79 double minmax[2];
80 array->GetRange(minmax);
81
82 size_t const nComponents = array->GetNumberOfComponents();
83 for(size_t c = 0; c < nComponents; c++) {
84 auto valuePass = vtkSmartPointer<vtkValuePass>::New();
85 valuePass->SetInputArrayToProcess(fieldType == 0
86 ? VTK_SCALAR_MODE_USE_POINT_FIELD_DATA
87 : VTK_SCALAR_MODE_USE_CELL_FIELD_DATA,
88 name.data());
89 valuePass->SetInputComponentToProcess(c);
90
91 valuePassCollection->AddItem(valuePass);
92 valuePassNames.push_back(
93 nComponents == 1 ? name : name + "_" + std::to_string(c));
94 }
95 }
96
97 return 1;
98};
99
101 vtkMultiBlockDataSet *outputImages,
102
103 vtkPointSet *inputObject,
104 vtkPointSet *inputGrid) const {
105
106 auto inputObjectAsPD = vtkSmartPointer<vtkPolyData>::New();
107 if(inputObject->IsA("vtkPolyData")) {
108 inputObjectAsPD->ShallowCopy(inputObject);
109 } else {
111 surfaceFilter->SetInputDataObject(inputObject);
112 surfaceFilter->Update();
113 inputObjectAsPD->ShallowCopy(surfaceFilter->GetOutput());
114 }
115
116 // ---------------------------------------------------------------------------
117 // Prepare Field Data for Depth Values
118 // ---------------------------------------------------------------------------
119
120 // iterate over sampling locations
122 float *samplingPositions
123 = static_cast<float *>(ttkUtils::GetVoidPointer(inputGrid->GetPoints()));
124 int const nSamplingPositions = inputGrid->GetNumberOfPoints();
125 auto camParameters = inputGrid->GetPointData();
126 auto camUp = static_cast<double *>(
127 ttkUtils::GetVoidPointer(camParameters->GetArray("CamUp")));
128 auto camFocalPoint = static_cast<double *>(
129 ttkUtils::GetVoidPointer(camParameters->GetArray("CamFocalPoint")));
130 auto camHeight = static_cast<double *>(
131 ttkUtils::GetVoidPointer(camParameters->GetArray("CamHeight")));
132 auto camAngle = static_cast<double *>(
133 ttkUtils::GetVoidPointer(camParameters->GetArray("CamAngle")));
134 auto camNearFar = static_cast<double *>(
135 ttkUtils::GetVoidPointer(camParameters->GetArray("CamNearFar")));
136 auto resolution = static_cast<double *>(
137 ttkUtils::GetVoidPointer(camParameters->GetArray("Resolution")));
138 auto projectionMode = static_cast<double *>(
139 ttkUtils::GetVoidPointer(camParameters->GetArray("ProjectionMode")));
140
141 // Camera
142 auto camera = vtkSmartPointer<vtkCamera>::New();
143
144 // Depth Pass Elements
145 auto rendererDepth = vtkSmartPointer<vtkRenderer>::New();
146 this->setupRenderer(rendererDepth, inputObjectAsPD, camera);
147
148 auto windowDepth = vtkSmartPointer<vtkRenderWindow>::New();
149 this->setupWindow(windowDepth, rendererDepth, resolution);
150
151 auto windowDepthToImageFilter
153 windowDepthToImageFilter->SetInput(windowDepth);
154 windowDepthToImageFilter->SetInputBufferTypeToZBuffer();
155
156 // Value passes Elements
157 size_t nValuePasses = 0;
158
159 auto rendererScalars = vtkSmartPointer<vtkRenderer>::New();
160 this->setupRenderer(rendererScalars, inputObjectAsPD, camera);
161
162 auto windowScalars = vtkSmartPointer<vtkRenderWindow>::New();
163 this->setupWindow(windowScalars, rendererScalars, resolution);
164
165 if(windowScalars->SupportsOpenGL() == 0) {
166 // MS Windows in a VM does not support OpenGL
167 this->printErr("RenderWindow does not support OpenGL");
168 return 0;
169 }
170
171 auto valuePassCollection = vtkSmartPointer<vtkRenderPassCollection>::New();
172 std::vector<std::string> valuePassNames;
173 size_t firstValuePassIndex = 0;
174
175 auto preventVTKBug = [](vtkPointSet *object) {
176 auto pd = object->GetPointData();
177 auto cd = object->GetCellData();
178
179 if(pd->GetNumberOfArrays() < 1 && cd->GetNumberOfArrays() > 0) {
180 size_t const nP = object->GetNumberOfPoints();
181
183 fakeArray->SetName("Fake");
184 fakeArray->SetNumberOfComponents(1);
185 fakeArray->SetNumberOfTuples(nP);
186 auto fakeArrayData = (signed char *)fakeArray->GetVoidPointer(0);
187 for(size_t i = 0; i < nP; i++)
188 fakeArrayData[i] = 0;
189 pd->AddArray(fakeArray);
190 return 1;
191 }
192 return 0;
193 };
194
195 if(preventVTKBug(inputObjectAsPD)) {
196 firstValuePassIndex = 1;
197 };
198
199 this->addValuePass(inputObjectAsPD, 0, valuePassCollection, valuePassNames);
200 this->addValuePass(inputObjectAsPD, 1, valuePassCollection, valuePassNames);
201 nValuePasses = valuePassNames.size();
202
204 sequence->SetPasses(valuePassCollection);
205
206 auto cameraPass = vtkSmartPointer<vtkCameraPass>::New();
207 cameraPass->SetDelegatePass(sequence);
208
209 auto glRenderer = vtkOpenGLRenderer::SafeDownCast(rendererScalars);
210 glRenderer->SetPass(cameraPass);
211
212 // First pass to setup everything
213 windowScalars->Render();
214
215 for(int i = 0; i < nSamplingPositions; i++) {
216 ttk::Timer timer;
217 int const resX = resolution[i * 2];
218 int const resY = resolution[i * 2 + 1];
219
220 this->printMsg("Rendering Image ("
221 + std::string(projectionMode[i] ? "P" : "O") + "|"
222 + std::to_string(resX) + "x" + std::to_string(resY) + ")",
223 0, 0, this->threadNumber_, ttk::debug::LineMode::REPLACE);
224
225 double camPos[3]{samplingPositions[i * 3], samplingPositions[i * 3 + 1],
226 samplingPositions[i * 3 + 2]};
227
228 if(projectionMode[i] == 0) {
229 camera->SetParallelProjection(true);
230 camera->SetParallelScale(
231 camHeight[i]
232 * 0.5); // *0.5 to convert CamHeight to weird VTK convention
233 } else {
234 camera->SetParallelProjection(false);
235 camera->SetViewAngle(camAngle[i] * resolution[i * 2 + 1]
236 / resolution[i * 2]);
237 }
238
239 camera->SetPosition(camPos);
240 camera->SetViewUp(&camUp[i * 3]);
241 camera->SetFocalPoint(&camFocalPoint[i * 3]);
242 camera->SetClippingRange(&camNearFar[i * 2]);
243
244 // Initialize Output
245 auto outputImage = vtkSmartPointer<vtkImageData>::New();
246 outputImage->SetDimensions(resolution[i * 2], resolution[i * 2 + 1], 1);
247 outputImage->SetSpacing(1, 1, 1);
248 outputImage->SetOrigin(0, 0, 0);
249 outputImage->AllocateScalars(VTK_FLOAT, 1);
250
251 auto outputImagePD = outputImage->GetPointData();
252
253 // Initialize as depth image
254 windowDepthToImageFilter->Modified();
255 windowDepthToImageFilter->Update();
256 outputImage->DeepCopy(windowDepthToImageFilter->GetOutput());
257 outputImagePD->GetAbstractArray(0)->SetName("Depth");
258
259 ttkCinemaImaging::AddAllFieldDataArrays(inputGrid, outputImage, i);
260
261 // Render Scalar Images
262 if(nValuePasses > firstValuePassIndex) {
263 windowScalars->Render();
264
265 for(size_t j = firstValuePassIndex; j < nValuePasses; j++) {
266 auto valuePass
267 = vtkValuePass::SafeDownCast(valuePassCollection->GetItemAsObject(j));
268 auto newValueArray = vtkSmartPointer<vtkFloatArray>::New();
269 newValueArray->DeepCopy(
270 valuePass->GetFloatImageDataArray(rendererScalars));
271 newValueArray->SetName(valuePassNames[j].data());
272 outputImagePD->AddArray(newValueArray);
273 }
274 }
275
276 // Add to output collection
277 outputImages->SetBlock(i, outputImage);
278
279 this->printMsg("Rendering Image ("
280 + std::string(projectionMode[i] ? "P" : "O") + "|"
281 + std::to_string(resX) + "x" + std::to_string(resY) + ")",
282 1, timer.getElapsedTime(), this->threadNumber_);
283 }
285
286 return 1;
287};
static int AddAllFieldDataArrays(vtkPointSet *inputGrid, vtkImageData *image, int tupleIdx)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
double getElapsedTime()
Definition Timer.h:15
~ttkCinemaImagingVTK() override
int addValuePass(vtkPointSet *object, int fieldType, vtkRenderPassCollection *valuePassCollection, std::vector< std::string > &valuePassNames) const
int setupWindow(vtkRenderWindow *window, vtkRenderer *renderer, const double resolution[2]) const
int RenderVTKObject(vtkMultiBlockDataSet *outputImages, vtkPointSet *inputObject, vtkPointSet *inputGrid) const
int setupRenderer(vtkRenderer *renderer, vtkPointSet *object, vtkCamera *camera) const
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)