TTK
Loading...
Searching...
No Matches
ttkWRLExporter.cpp
Go to the documentation of this file.
1// VTK includes
2#include <vtkActor.h>
3#include <vtkActorCollection.h>
4#include <vtkAssemblyPath.h>
5#include <vtkCamera.h>
6#include <vtkCellArray.h>
7#include <vtkCompositeDataGeometryFilter.h>
8#include <vtkDataArray.h>
9#include <vtkDataObject.h>
10#include <vtkGeometryFilter.h>
11#include <vtkLightCollection.h>
12#include <vtkMapper.h>
13#include <vtkMath.h>
14#include <vtkObjectFactory.h>
15#include <vtkPointData.h>
16#include <vtkPoints.h>
17#include <vtkPolyDataMapper.h>
18#include <vtkRenderWindow.h>
19#include <vtkRenderer.h>
20#include <vtkRendererCollection.h>
21#include <vtkSmartPointer.h>
22#include <vtkTransform.h>
23#include <vtkUnsignedCharArray.h>
24#include <vtkVRMLExporter.h>
25#include <vtkVersionMacros.h>
26
27// base code includes
28#include <Debug.h>
29
30#include <ttkWRLExporter.h>
31
33
34TTKWRLEXPORTER_EXPORT vtkPolyData *ttkWRLExporterPolyData_ = nullptr;
35
36// Over-ride the appropriate functions of the vtkVRMLExporter class.
37TTKWRLEXPORTER_EXPORT void vtkVRMLExporter::WriteAnActor(vtkActor *anActor,
38 FILE *fp) {
39
40 ttk::Debug debugInfo;
41 debugInfo.setDebugMsgPrefix("WRLExporter");
42 debugInfo.printWrn("Using TTK fix for VRML export...");
43
45 vtkPointData *pntData;
46 vtkPoints *points;
47 vtkDataArray *normals = nullptr;
48 vtkDataArray *tcoords = nullptr;
49 int i, i1, i2;
50 double *tempd;
51 vtkCellArray *cells;
52 vtkIdType npts = 0;
53#ifdef VTK_CELL_ARRAY_V2
54 vtkIdType const *indx = nullptr;
55#else
56 vtkIdType *indx = nullptr;
57#endif
58 int pointDataWritten = 0;
59 vtkPolyDataMapper *pm;
60 vtkUnsignedCharArray *colors;
61 double *p;
62 unsigned char *c;
63 vtkTransform *trans;
64
65 // see if the actor has a mapper. it could be an assembly
66 if(anActor->GetMapper() == nullptr) {
67 return;
68 }
69
70 if(anActor->GetVisibility() == 0) {
71 return;
72 }
73
74 // first stuff out the transform
75 trans = vtkTransform::New();
76 trans->SetMatrix(anActor->vtkProp3D::GetMatrix());
77
78 fprintf(fp, " Transform {\n");
79 tempd = trans->GetPosition();
80 fprintf(fp, " translation %g %g %g\n", tempd[0], tempd[1], tempd[2]);
81 tempd = trans->GetOrientationWXYZ();
82 fprintf(fp, " rotation %g %g %g %g\n", tempd[1], tempd[2], tempd[3],
83 tempd[0] * vtkMath::Pi() / 180.0);
84 tempd = trans->GetScale();
85 fprintf(fp, " scale %g %g %g\n", tempd[0], tempd[1], tempd[2]);
86 fprintf(fp, " children [\n");
87 trans->Delete();
88
89 vtkDataObject *inputDO = anActor->GetMapper()->GetInputDataObject(0, 0);
90
91 if(inputDO == nullptr) {
92 return;
93 }
94
95 // we really want polydata
96 if(inputDO->IsA("vtkCompositeDataSet")) {
97 vtkCompositeDataGeometryFilter *gf = vtkCompositeDataGeometryFilter::New();
98 gf->SetInputConnection(anActor->GetMapper()->GetInputConnection(0, 0));
99 gf->Update();
100 pd = gf->GetOutput();
101 gf->Delete();
102 } else if(inputDO->GetDataObjectType() != VTK_POLY_DATA) {
103 vtkGeometryFilter *gf = vtkGeometryFilter::New();
104 gf->SetInputConnection(anActor->GetMapper()->GetInputConnection(0, 0));
105 gf->Update();
106 pd = gf->GetOutput();
107 gf->Delete();
108 } else {
109 anActor->GetMapper()->Update();
110 pd = static_cast<vtkPolyData *>(inputDO);
111 }
112
113 // BUG fix here
114 ttkWRLExporterPolyData_ = static_cast<vtkPolyData *>(pd);
115 // end of BUG fix here
116
117 pm = vtkPolyDataMapper::New();
118 pm->SetInputData(pd);
119 pm->SetScalarRange(anActor->GetMapper()->GetScalarRange());
120 pm->SetScalarVisibility(anActor->GetMapper()->GetScalarVisibility());
121 pm->SetLookupTable(anActor->GetMapper()->GetLookupTable());
122 pm->SetScalarMode(anActor->GetMapper()->GetScalarMode());
123
124 if(pm->GetScalarMode() == VTK_SCALAR_MODE_USE_POINT_FIELD_DATA
125 || pm->GetScalarMode() == VTK_SCALAR_MODE_USE_CELL_FIELD_DATA) {
126 if(anActor->GetMapper()->GetArrayAccessMode() == VTK_GET_ARRAY_BY_ID) {
127 pm->ColorByArrayComponent(anActor->GetMapper()->GetArrayId(),
128 anActor->GetMapper()->GetArrayComponent());
129 } else {
130 pm->ColorByArrayComponent(anActor->GetMapper()->GetArrayName(),
131 anActor->GetMapper()->GetArrayComponent());
132 }
133 }
134
135 points = pd->GetPoints();
136 pntData = pd->GetPointData();
137 normals = pntData->GetNormals();
138 tcoords = pntData->GetTCoords();
139 colors = pm->MapScalars(1.0);
140
141 // write out polys if any
142 if(pd->GetNumberOfPolys() > 0) {
143 WriteShapeBegin(anActor, fp, pd, pntData, colors);
144 fprintf(fp, " geometry IndexedFaceSet {\n");
145 // two sided lighting ? for now assume it is on
146 fprintf(fp, " solid FALSE\n");
147
148 if(!pointDataWritten) {
149
150#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 4, 20250513)
151 this->WritePointData(points, normals, tcoords, colors, false, fp);
152#else
153 this->WritePointData(points, normals, tcoords, colors, fp);
154#endif
155 pointDataWritten = 1;
156 } else {
157 fprintf(fp, " coord USE VTKcoordinates\n");
158
159 if(normals) {
160 fprintf(fp, " normal USE VTKnormals\n");
161 }
162
163 if(tcoords) {
164 fprintf(fp, " texCoord USE VTKtcoords\n");
165 }
166
167 if(colors) {
168 fprintf(fp, " color USE VTKcolors\n");
169 }
170 }
171
172 fprintf(fp, " coordIndex [\n");
173
174 cells = pd->GetPolys();
175
176 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
177 fprintf(fp, " ");
178
179 for(i = 0; i < npts; i++) {
180 // treating vtkIdType as int
181 fprintf(fp, "%i, ", static_cast<int>(indx[i]));
182 }
183
184 fprintf(fp, "-1,\n");
185 }
186
187 fprintf(fp, " ]\n");
188 fprintf(fp, " }\n");
189 WriteShapeEnd(fp);
190 }
191
192 // write out tstrips if any
193 if(pd->GetNumberOfStrips() > 0) {
194 WriteShapeBegin(anActor, fp, pd, pntData, colors);
195 fprintf(fp, " geometry IndexedFaceSet {\n");
196
197 if(!pointDataWritten) {
198
199#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 4, 20250513)
200 this->WritePointData(points, normals, tcoords, colors, false, fp);
201#else
202 this->WritePointData(points, normals, tcoords, colors, fp);
203#endif
204
205 pointDataWritten = 1;
206 } else {
207 fprintf(fp, " coord USE VTKcoordinates\n");
208
209 if(normals) {
210 fprintf(fp, " normal USE VTKnormals\n");
211 }
212
213 if(tcoords) {
214 fprintf(fp, " texCoord USE VTKtcoords\n");
215 }
216
217 if(colors) {
218 fprintf(fp, " color USE VTKcolors\n");
219 }
220 }
221
222 fprintf(fp, " coordIndex [\n");
223 cells = pd->GetStrips();
224
225 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
226 for(i = 2; i < npts; i++) {
227 if(i % 2) {
228 i1 = i - 1;
229 i2 = i - 2;
230 } else {
231 i1 = i - 2;
232 i2 = i - 1;
233 }
234
235 // treating vtkIdType as int
236 fprintf(fp, " %i, %i, %i, -1,\n",
237 static_cast<int>(indx[i1]), static_cast<int>(indx[i2]),
238 static_cast<int>(indx[i]));
239 }
240 }
241
242 fprintf(fp, " ]\n");
243 fprintf(fp, " }\n");
244 WriteShapeEnd(fp);
245 }
246
247 // write out lines if any
248 if(pd->GetNumberOfLines() > 0) {
249 WriteShapeBegin(anActor, fp, pd, pntData, colors);
250 fprintf(fp, " geometry IndexedLineSet {\n");
251
252 if(!pointDataWritten) {
253
254#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 4, 20250513)
255 this->WritePointData(points, nullptr, nullptr, colors, false, fp);
256#else
257 this->WritePointData(points, nullptr, nullptr, colors, fp);
258#endif
259 } else {
260 fprintf(fp, " coord USE VTKcoordinates\n");
261
262 if(colors) {
263 fprintf(fp, " color USE VTKcolors\n");
264 }
265 }
266
267 fprintf(fp, " coordIndex [\n");
268
269 cells = pd->GetLines();
270
271 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
272 fprintf(fp, " ");
273
274 for(i = 0; i < npts; i++) {
275 // treating vtkIdType as int
276 fprintf(fp, "%i, ", static_cast<int>(indx[i]));
277 }
278
279 fprintf(fp, "-1,\n");
280 }
281
282 fprintf(fp, " ]\n");
283 fprintf(fp, " }\n");
284 WriteShapeEnd(fp);
285 }
286
287 // write out verts if any
288 if(pd->GetNumberOfVerts() > 0) {
289 WriteShapeBegin(anActor, fp, pd, pntData, colors);
290 fprintf(fp, " geometry PointSet {\n");
291 cells = pd->GetVerts();
292 fprintf(fp, " coord Coordinate {");
293 fprintf(fp, " point [");
294
295 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
296 fprintf(fp, " ");
297
298 for(i = 0; i < npts; i++) {
299 p = points->GetPoint(indx[i]);
300 fprintf(fp, " %g %g %g,\n", p[0], p[1], p[2]);
301 }
302 }
303
304 fprintf(fp, " ]\n");
305 fprintf(fp, " }\n");
306
307 if(colors) {
308 fprintf(fp, " color Color {");
309 fprintf(fp, " color [");
310
311 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
312 fprintf(fp, " ");
313
314 for(i = 0; i < npts; i++) {
315 c = colors->GetPointer(4 * indx[i]);
316 fprintf(fp, " %g %g %g,\n", c[0] / 255.0, c[1] / 255.0,
317 c[2] / 255.0);
318 }
319 }
320
321 fprintf(fp, " ]\n");
322 fprintf(fp, " }\n");
323 }
324
325 fprintf(fp, " }\n");
326 WriteShapeEnd(fp);
327 }
328
329 fprintf(fp, " ]\n"); // close the original transforms children
330 fprintf(fp, " }\n"); // close the original transform
331
332 pm->Delete();
333}
334
335TTKWRLEXPORTER_EXPORT void vtkVRMLExporter::WriteData() {
336
337 vtkRenderer *ren;
338 vtkActorCollection *ac;
339 vtkActor *anActor, *aPart;
340 vtkLightCollection *lc;
341 vtkLight *aLight;
342 // // vtkCamera *cam;
343 // double *tempd;
344 FILE *fp;
345
346 // make sure the user specified a FileName or FilePointer
347 if(!this->FilePointer && (this->FileName == nullptr)) {
348 vtkErrorMacro(<< "Please specify FileName to use");
349 return;
350 }
351
352 // Always pick the first renderer
353 // first make sure there is only one renderer in this rendering window
354 // if (this->RenderWindow->GetRenderers()->GetNumberOfItems() > 1)
355 // {
356 // vtkErrorMacro(<< "VRML files only support one renderer per window.");
357 // return;
358 // }
359
360 // get the renderer
361 ren = this->RenderWindow->GetRenderers()->GetFirstRenderer();
362
363 // make sure it has at least one actor
364 if(ren->GetActors()->GetNumberOfItems() < 1) {
365 vtkErrorMacro(<< "no actors found for writing VRML file.");
366 return;
367 }
368
369 // try opening the files
370 if(!this->FilePointer) {
371 fp = fopen(this->FileName, "w");
372
373 if(!fp) {
374 vtkErrorMacro(<< "unable to open VRML file " << this->FileName);
375 return;
376 }
377 } else {
378 fp = this->FilePointer;
379 }
380
381 //
382 // Write header
383 //
384 vtkDebugMacro("Writing VRML file");
385 fprintf(fp, "#VRML V2.0 utf8\n");
386 fprintf(fp, "# VRML file written by the visualization toolkit\n\n");
387
388 // Start write the Background
389 double background[3];
390 ren->GetBackground(background);
391 fprintf(fp, " Background {\n ");
392 fprintf(fp, " skyColor [%f %f %f, ]\n", background[0], background[1],
393 background[2]);
394 fprintf(fp, " }\n ");
395 // End of Background
396
397 // BUG fix
398 // do the camera
399 // cam = ren->GetActiveCamera();
400 // fprintf(fp, " Viewpoint\n {\n fieldOfView %f\n",
401 // cam->GetViewAngle()*vtkMath::Pi() / 180.0);
402 // fprintf(fp, " position %f %f %f\n", cam->GetPosition()[0],
403 // cam->GetPosition()[1], cam->GetPosition()[2]);
404 // fprintf(fp, " description \"Default View\"\n");
405 // tempd = cam->GetOrientationWXYZ();
406 // fprintf(fp, " orientation %g %g %g %g\n }\n", tempd[1],
407 // tempd[2],
408 // tempd[3], tempd[0]*vtkMath::Pi() / 180.0);
409 //
410 // // do the lights first the ambient then the others
411 // fprintf(fp,
412 // " NavigationInfo {\n type [\"EXAMINE\",\"FLY\"]\n speed
413 // %f\n",
414 // this->Speed);
415 //
416 // if (ren->GetLights()->GetNumberOfItems() == 0){
417 // fprintf(fp, " headlight TRUE}\n\n");
418 // }
419 // else{
420 // fprintf(fp, " headlight FALSE}\n\n");
421 // }
422 // end of BUG fix
423
424 fprintf(
425 fp, "DirectionalLight { ambientIntensity 1 intensity 0 # ambient light\n");
426 fprintf(fp, " color %f %f %f }\n\n", ren->GetAmbient()[0],
427 ren->GetAmbient()[1], ren->GetAmbient()[2]);
428
429 // make sure we have a default light
430 // if we dont then use a headlight
431 lc = ren->GetLights();
432 vtkCollectionSimpleIterator lsit;
433
434 for(lc->InitTraversal(lsit); (aLight = lc->GetNextLight(lsit));) {
435 this->WriteALight(aLight, fp);
436 }
437
438 // do the actors now
439 ac = ren->GetActors();
440 vtkAssemblyPath *apath;
441 vtkCollectionSimpleIterator ait;
442
443 for(ac->InitTraversal(ait); (anActor = ac->GetNextActor(ait));) {
444 for(anActor->InitPathTraversal(); (apath = anActor->GetNextPath());) {
445 aPart = static_cast<vtkActor *>(apath->GetLastNode()->GetViewProp());
446 this->WriteAnActor(aPart, fp);
447 }
448 }
449
450 if(!this->FilePointer) {
451 fclose(fp);
452 }
453}
454
455TTKWRLEXPORTER_EXPORT void
456 vtkVRMLExporter::WritePointData(vtkPoints *points,
457 vtkDataArray *normals,
458 vtkDataArray *tcoords,
459 vtkUnsignedCharArray *colors,
460#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 4, 20250513)
461 bool vtkNotUsed(cellData),
462#endif
463 FILE *fp) {
464
465 double *p;
466 unsigned char *c;
467
468 // write out the points
469 fprintf(fp, " coord DEF VTKcoordinates Coordinate {\n");
470 fprintf(fp, " point [\n");
471 for(int i = 0; i < points->GetNumberOfPoints(); i++) {
472 p = points->GetPoint(i);
473 fprintf(fp, " %g %g %g,\n", p[0], p[1], p[2]);
474 }
475 fprintf(fp, " ]\n");
476 fprintf(fp, " }\n");
477
478 // write out the point data
479 if(normals) {
480
481 fprintf(fp, " normal DEF VTKnormals Normal {\n");
482 fprintf(fp, " vector [\n");
483 for(int i = 0; i < normals->GetNumberOfTuples(); i++) {
484 p = normals->GetTuple(i);
485 fprintf(fp, " %g %g %g,\n", p[0], p[1], p[2]);
486 }
487 fprintf(fp, " ]\n");
488 fprintf(fp, " }\n");
489 }
490
491 // write out the point data
492 if(tcoords) {
493 fprintf(fp, " texCoord DEF VTKtcoords TextureCoordinate {\n");
494 fprintf(fp, " point [\n");
495 for(int i = 0; i < tcoords->GetNumberOfTuples(); i++) {
496 p = tcoords->GetTuple(i);
497 fprintf(fp, " %g %g,\n", p[0], p[1]);
498 }
499 fprintf(fp, " ]\n");
500 fprintf(fp, " }\n");
501
502 // BUG fix here.
503 if(ttkWRLExporterPolyData_) {
504 fprintf(fp, " texCoordIndex[\n");
505 vtkCellArray *cells = ttkWRLExporterPolyData_->GetPolys();
506 vtkIdType npts = 0;
507#ifdef VTK_CELL_ARRAY_V2
508 vtkIdType const *indx = nullptr;
509#else
510 vtkIdType *indx = nullptr;
511#endif
512 for(cells->InitTraversal(); cells->GetNextCell(npts, indx);) {
513 fprintf(fp, " ");
514 for(int i = 0; i < npts; i++) {
515 fprintf(fp, "%i, ", static_cast<int>(indx[i]));
516 }
517 fprintf(fp, "-1,\n");
518 }
519 fprintf(fp, " ]\n");
520 }
521 // end of BUG fix here.
522 }
523
524 // write out the point data
525 if(colors) {
526 fprintf(fp, " color DEF VTKcolors Color {\n");
527 fprintf(fp, " color [\n");
528 for(int i = 0; i < colors->GetNumberOfTuples(); i++) {
529 c = colors->GetPointer(4 * i);
530 fprintf(
531 fp, " %g %g %g,\n", c[0] / 255.0, c[1] / 255.0, c[2] / 255.0);
532 }
533 fprintf(fp, " ]\n");
534 fprintf(fp, " }\n");
535 }
536}
537
Minimalist debugging class.
Definition Debug.h:88
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364