Salome HOME
7d139de0727e0000a820b41267b38454aae04a6b
[modules/gui.git] / src / VTKViewer / VTKViewer_GeometryFilter.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  File   : VTKViewer_GeometryFilter.cxx
24 //  Author : Michael ZORIN
25 //  Module : SALOME
26 //
27 #include "VTKViewer_GeometryFilter.h"
28 #include "VTKViewer_ConvexTool.h"
29 #include "VTKViewer_ArcBuilder.h"
30
31 #include <vtkSmartPointer.h>
32 #include <vtkCellArray.h>
33 #include <vtkCellData.h>
34 #include <vtkGenericCell.h>
35 #include <vtkHexagonalPrism.h>
36 #include <vtkHexahedron.h>
37 #include <vtkInformation.h>
38 #include <vtkInformationVector.h>
39 #include <vtkMergePoints.h>
40 #include <vtkObjectFactory.h>
41 #include <vtkPointData.h>
42 #include <vtkPolyData.h>
43 #include <vtkPolygon.h>
44 #include <vtkPyramid.h>
45 #include <vtkStructuredGrid.h>
46 #include <vtkTetra.h>
47 #include <vtkUnsignedCharArray.h>
48 #include <vtkUnstructuredGrid.h>
49 #include <vtkVoxel.h>
50 #include <vtkWedge.h>
51
52 #include <algorithm>
53 #include <iterator>
54 #include <vector>
55 #include <map>
56 #include <set>
57
58 #include "utilities.h"
59
60 #if defined __GNUC__
61   #if __GNUC__ == 2
62     #define __GNUC_2__
63   #endif
64 #endif
65
66 #define VTK_XVERSION (VTK_MAJOR_VERSION*10000+VTK_MINOR_VERSION*100+VTK_BUILD_VERSION)
67
68 //#define __MYDEBUG__
69 //#define USE_ROBUST_TRIANGULATION
70
71 ///////////////////////////////////////////////////////////////////////////////////////////////
72 // VSR 26/10/2012: fix of regression (issue 21924) - increased memory consumption
73 // for displaying of 3d elements, introduced by fix for issue 20314.
74 // ...
75 // The macro SHOW_COINCIDING_3D_PAL20314, when defined, allows correct visualization of
76 // coincident 3d elements but causes substantial increasing of memory consumption, as all 3d 
77 // elements are always shown, even if they are totally covered by surrounding faces.
78 // If this macro is not defined (commented), the behaviour is defined by another macro -
79 // SHOW_COINCIDING_3D_PAL21924, as follows:
80 // - If SHOW_COINCIDING_3D_PAL21924 is defined, an alternative solution for computing 
81 //   visibility of 3d elements is used; this solution allows to fix problem with visibility
82 //   of coinciding 3d elements in most cases (though some cases might not work), while not
83 //   causing significant increasing of memory consumption.
84 // - If SHOW_COINCIDING_3D_PAL21924 is not defined (commented), coinciding 3d elements are 
85 //   not shown at all (this corresponds to the state before issue 20314 fixing).
86 ///////////////////////////////////////////////////////////////////////////////////////////////
87 //#define SHOW_COINCIDING_3D_PAL20314
88 #ifndef SHOW_COINCIDING_3D_PAL20314
89 #define SHOW_COINCIDING_3D_PAL21924
90 #endif
91 ///////////////////////////////////////////////////////////////////////////////////////////////
92
93 vtkCxxRevisionMacro(VTKViewer_GeometryFilter, "$Revision$");
94 vtkStandardNewMacro(VTKViewer_GeometryFilter);
95
96 VTKViewer_GeometryFilter
97 ::VTKViewer_GeometryFilter():
98   myShowInside(0),
99   myStoreMapping(0),
100   myAppendCoincident3D(0),
101   myIsWireframeMode(0),
102   myIsBuildArc(false),
103   myMaxArcAngle(2)
104 {}
105
106
107 VTKViewer_GeometryFilter
108 ::~VTKViewer_GeometryFilter()
109 {}
110
111
112 int
113 VTKViewer_GeometryFilter
114 ::RequestData(
115   vtkInformation *request,
116   vtkInformationVector **inputVector,
117   vtkInformationVector *outputVector)
118 {
119   // get the info objects
120   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
121   vtkInformation *outInfo = outputVector->GetInformationObject(0);
122
123   // get the input and ouptut
124   vtkDataSet *input = vtkDataSet::SafeDownCast(
125     inInfo->Get(vtkDataObject::DATA_OBJECT()));
126   vtkPolyData *output = vtkPolyData::SafeDownCast(
127     outInfo->Get(vtkDataObject::DATA_OBJECT()));
128
129   vtkIdType numCells=input->GetNumberOfCells();
130
131   if (numCells == 0)
132     {
133       return 0;
134     }
135
136   if (input->GetDataObjectType() == VTK_UNSTRUCTURED_GRID){
137     return this->UnstructuredGridExecute(input, output, outInfo);
138   }else
139     return Superclass::RequestData(request,inputVector,outputVector);
140
141   return 1;
142 }
143
144
145 int
146 VTKViewer_GeometryFilter
147 ::UnstructuredGridExecute(
148                           vtkDataSet *dataSetInput,
149                           vtkPolyData *output,
150                           vtkInformation *outInfo)
151 {
152
153   vtkUnstructuredGrid *input= (vtkUnstructuredGrid *)dataSetInput;
154   vtkCellArray *Connectivity = input->GetCells();
155   // Check input
156   if ( Connectivity == NULL )
157     {
158     vtkDebugMacro(<<"Nothing to extract");
159     return 0;
160     }
161
162   vtkIdType cellId;
163   int i;
164   int allVisible;
165   vtkIdType npts = 0;
166   vtkIdType *pts = 0;
167   vtkPoints *p = input->GetPoints();
168   vtkIdType numCells=input->GetNumberOfCells();
169   vtkPointData *pd = input->GetPointData();
170   vtkCellData *cd = input->GetCellData();
171   vtkPointData *outputPD = output->GetPointData();
172
173   VTKViewer_OrderedTriangulator anOrderedTriangulator;
174   VTKViewer_DelaunayTriangulator aDelaunayTriangulator;
175
176   vtkCellData *outputCD = output->GetCellData();
177   vtkGenericCell *cell = vtkGenericCell::New();
178
179   vtkIdList *cellIds = vtkIdList::New();
180   vtkIdList *faceIds = vtkIdList::New();
181   vtkIdList *cellIdsTmp = vtkIdList::New();
182   vtkIdList *faceIdsTmp = vtkIdList::New();
183
184   char *cellVis;
185   vtkIdType newCellId;
186   int faceId, *faceVerts, numFacePts;
187   vtkFloatingPointType *x;
188   vtkIdType PixelConvert[4];
189   // Change the type from int to vtkIdType in order to avoid compilation errors while using VTK
190   // from ParaView-3.4.0 compiled on 64-bit Debian platform with VTK_USE_64BIT_IDS = ON
191   vtkIdType aNewPts[VTK_CELL_SIZE];
192   // ghost cell stuff
193   unsigned char  updateLevel = (unsigned char)(output->GetUpdateGhostLevel());
194   unsigned char  *cellGhostLevels = 0;
195
196   PixelConvert[0] = 0;
197   PixelConvert[1] = 1;
198   PixelConvert[2] = 3;
199   PixelConvert[3] = 2;
200
201   vtkDebugMacro(<<"Executing geometry filter for unstructured grid input");
202
203   vtkDataArray* temp = 0;
204   if (cd)
205     {
206     temp = cd->GetArray("vtkGhostLevels");
207     }
208   if ( (!temp) || (temp->GetDataType() != VTK_UNSIGNED_CHAR)
209     || (temp->GetNumberOfComponents() != 1))
210     {
211     vtkDebugMacro("No appropriate ghost levels field available.");
212     }
213   else
214     {
215     cellGhostLevels = ((vtkUnsignedCharArray*)temp)->GetPointer(0);
216     }
217
218   // Determine nature of what we have to do
219   if ( (!this->CellClipping) && (!this->PointClipping) &&
220        (!this->ExtentClipping) )
221     {
222     allVisible = 1;
223     cellVis = NULL;
224     }
225   else
226     {
227     allVisible = 0;
228     cellVis = new char[numCells];
229     }
230
231   // Issue 0020115: [CEA 308] Quadratic elements visualization
232   // Fix of remark described in note 0005222 - SIGSEGV
233   vtkPoints* outputPoints = vtkPoints::New();
234   outputPoints->DeepCopy(input->GetPoints());
235   output->SetPoints(outputPoints);
236   outputPoints->Delete();
237
238   outputPD->PassData(pd);
239
240   outputCD->CopyAllocate(cd,numCells,numCells/2);
241
242   output->Allocate(numCells/4+1,numCells);
243
244   // Loop over the cells determining what's visible
245   if (!allVisible)
246     {
247     for (cellId=0, Connectivity->InitTraversal();
248          Connectivity->GetNextCell(npts,pts);
249          cellId++)
250       {
251       cellVis[cellId] = 1;
252       if ( ( this->CellClipping && cellId < this->CellMinimum ) ||
253            cellId > this->CellMaximum )
254         {
255         cellVis[cellId] = 0;
256         }
257       else
258         {
259         for (i=0; i < npts; i++)
260           {
261           x = p->GetPoint(pts[i]);
262           if ( ( ( ( this->PointClipping && (pts[i] < this->PointMinimum ) ) ||
263                                              pts[i] > this->PointMaximum) ) ||
264                ( this->ExtentClipping &&
265                 ( x[0] < this->Extent[0] || x[0] > this->Extent[1] ||
266                   x[1] < this->Extent[2] || x[1] > this->Extent[3] ||
267                   x[2] < this->Extent[4] || x[2] > this->Extent[5] )) )
268             {
269             cellVis[cellId] = 0;
270             break;
271             }//point/extent clipping
272           }//for each point
273         }//if point clipping needs checking
274       }//for all cells
275     }//if not all visible
276
277   // Loop over all cells now that visibility is known
278   // (Have to compute visibility first for 3D cell boundarys)
279   int progressInterval = numCells/20 + 1;
280   TMapOfVectorId aDimension2VTK2ObjIds;
281   if(myStoreMapping){
282     myVTK2ObjIds.clear();
283     myVTK2ObjIds.reserve(numCells);
284   }
285   for (cellId=0, Connectivity->InitTraversal();
286        Connectivity->GetNextCell(npts,pts);
287        cellId++)
288     {
289     //Progress and abort method support
290     if ( !(cellId % progressInterval) )
291       {
292       vtkDebugMacro(<<"Process cell #" << cellId);
293       this->UpdateProgress ((float)cellId/numCells);
294       }
295
296     // Handle ghost cells here.  Another option was used cellVis array.
297     if (cellGhostLevels && cellGhostLevels[cellId] > updateLevel)
298       { // Do not create surfaces in outer ghost cells.
299       continue;
300       }
301
302     if (allVisible || cellVis[cellId])  //now if visible extract geometry
303       {
304       //special code for nonlinear cells - rarely occurs, so right now it
305       //is slow.
306       vtkIdType aCellType = input->GetCellType(cellId);
307       switch (aCellType)
308         {
309         case VTK_EMPTY_CELL:
310           break;
311
312         case VTK_VERTEX:
313         case VTK_POLY_VERTEX:
314           newCellId = output->InsertNextCell(aCellType,npts,pts);
315           if(myStoreMapping){
316             InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
317           }
318           outputCD->CopyData(cd,cellId,newCellId);
319           break;
320
321         case VTK_LINE:
322         case VTK_POLY_LINE:
323           newCellId = output->InsertNextCell(aCellType,npts,pts);
324           if(myStoreMapping)
325             InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
326           outputCD->CopyData(cd,cellId,newCellId);
327           break;
328
329         case VTK_TRIANGLE:
330         case VTK_QUAD:
331         case VTK_POLYGON:
332           newCellId = output->InsertNextCell(aCellType,npts,pts);
333           if(myStoreMapping)
334             InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
335           outputCD->CopyData(cd,cellId,newCellId);
336           break;
337
338         case VTK_TRIANGLE_STRIP:
339           newCellId = output->InsertNextCell(aCellType,npts,pts);
340           if(myStoreMapping)
341             InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
342           outputCD->CopyData(cd,cellId,newCellId);
343           break;
344
345         case VTK_PIXEL:
346           newCellId = output->InsertNextCell(aCellType,npts,pts);
347           if(myStoreMapping)
348             InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
349           outputCD->CopyData(cd,cellId,newCellId);
350           break;
351
352         case VTK_CONVEX_POINT_SET: {
353           bool anIsOk = anOrderedTriangulator.Execute(input,
354                                                       cd,
355                                                       cellId,
356                                                       myShowInside,
357                                                       allVisible,
358                                                       GetAppendCoincident3D(),
359                                                       cellVis,
360                                                       output,
361                                                       outputCD,
362                                                       myStoreMapping,
363                                                       myVTK2ObjIds,
364                                                       aDimension2VTK2ObjIds,
365                                                       true);
366           if(!anIsOk)
367             aDelaunayTriangulator.Execute(input,
368                                           cd,
369                                           cellId,
370                                           myShowInside,
371                                           allVisible,
372                                           GetAppendCoincident3D(),
373                                           cellVis,
374                                           output,
375                                           outputCD,
376                                           myStoreMapping,
377                                           myVTK2ObjIds,
378                                           aDimension2VTK2ObjIds,
379                                           false);
380
381           break;
382         }
383         case VTK_TETRA: {
384 #ifdef SHOW_COINCIDING_3D_PAL21924
385           faceIdsTmp->Reset();
386           for (int ai=0; ai<npts; ai++)
387             faceIdsTmp->InsertNextId(pts[ai]);
388           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
389 #endif
390           for (faceId = 0; faceId < 4; faceId++)
391             {
392             faceIds->Reset();
393             faceVerts = vtkTetra::GetFaceArray(faceId);
394             faceIds->InsertNextId(pts[faceVerts[0]]);
395             faceIds->InsertNextId(pts[faceVerts[1]]);
396             faceIds->InsertNextId(pts[faceVerts[2]]);
397             aCellType = VTK_TRIANGLE;
398             numFacePts = 3;
399             input->GetCellNeighbors(cellId, faceIds, cellIds);
400 #ifdef SHOW_COINCIDING_3D_PAL21924
401             int nbNeighbors = 0;
402             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
403               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
404             }
405             bool process = nbNeighbors <= 0;
406 #else
407             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
408 #endif
409             if ( process || myShowInside ||
410                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
411               {
412               for ( i=0; i < numFacePts; i++)
413                 aNewPts[i] = pts[faceVerts[i]];
414               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
415               if(myStoreMapping)
416                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
417               outputCD->CopyData(cd,cellId,newCellId);
418               }
419             }
420           break;
421         }
422         case VTK_VOXEL: {
423 #ifdef SHOW_COINCIDING_3D_PAL21924
424           faceIdsTmp->Reset();
425           for (int ai=0; ai<npts; ai++)
426             faceIdsTmp->InsertNextId(pts[ai]);
427           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
428 #endif
429           for (faceId = 0; faceId < 6; faceId++)
430             {
431             faceIds->Reset();
432             faceVerts = vtkVoxel::GetFaceArray(faceId);
433             faceIds->InsertNextId(pts[faceVerts[0]]);
434             faceIds->InsertNextId(pts[faceVerts[1]]);
435             faceIds->InsertNextId(pts[faceVerts[2]]);
436             faceIds->InsertNextId(pts[faceVerts[3]]);
437             aCellType = VTK_QUAD;
438             numFacePts = 4;
439             input->GetCellNeighbors(cellId, faceIds, cellIds);
440 #ifdef SHOW_COINCIDING_3D_PAL21924
441             int nbNeighbors = 0;
442             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
443               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
444             }
445             bool process = nbNeighbors <= 0;
446 #else
447             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
448 #endif
449             if ( process || myShowInside ||
450                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
451               {
452               for ( i=0; i < numFacePts; i++)
453                 aNewPts[i] = pts[faceVerts[PixelConvert[i]]];
454               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
455               if(myStoreMapping)
456                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
457               outputCD->CopyData(cd,cellId,newCellId);
458               }
459             }
460           break;
461         }
462         case VTK_HEXAHEDRON: {
463 #ifdef SHOW_COINCIDING_3D_PAL21924
464           faceIdsTmp->Reset();
465           for (int ai=0; ai<npts; ai++)
466             faceIdsTmp->InsertNextId(pts[ai]);
467           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
468 #endif
469           for (faceId = 0; faceId < 6; faceId++)
470             {
471             faceIds->Reset();
472             faceVerts = vtkHexahedron::GetFaceArray(faceId);
473             faceIds->InsertNextId(pts[faceVerts[0]]);
474             faceIds->InsertNextId(pts[faceVerts[1]]);
475             faceIds->InsertNextId(pts[faceVerts[2]]);
476             faceIds->InsertNextId(pts[faceVerts[3]]);
477             aCellType = VTK_QUAD;
478             numFacePts = 4;
479             input->GetCellNeighbors(cellId, faceIds, cellIds);
480 #ifdef SHOW_COINCIDING_3D_PAL21924
481             int nbNeighbors = 0;
482             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
483               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
484             }
485             bool process = nbNeighbors <= 0;
486 #else
487             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
488 #endif
489             if ( process || myShowInside ||
490                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
491               {
492               for ( i=0; i < numFacePts; i++)
493                 aNewPts[i] = pts[faceVerts[i]];
494               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
495               if(myStoreMapping)
496                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
497               outputCD->CopyData(cd,cellId,newCellId);
498               }
499             }
500           break;
501         }
502         case VTK_WEDGE: {
503 #ifdef SHOW_COINCIDING_3D_PAL21924
504           faceIdsTmp->Reset();
505           for (int ai=0; ai<npts; ai++)
506             faceIdsTmp->InsertNextId(pts[ai]);
507           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
508 #endif
509           for (faceId = 0; faceId < 5; faceId++)
510             {
511             faceIds->Reset();
512             faceVerts = vtkWedge::GetFaceArray(faceId);
513             faceIds->InsertNextId(pts[faceVerts[0]]);
514             faceIds->InsertNextId(pts[faceVerts[1]]);
515             faceIds->InsertNextId(pts[faceVerts[2]]);
516             aCellType = VTK_TRIANGLE;
517             numFacePts = 3;
518             if (faceVerts[3] >= 0)
519               {
520               faceIds->InsertNextId(pts[faceVerts[3]]);
521               aCellType = VTK_QUAD;
522               numFacePts = 4;
523               }
524
525             input->GetCellNeighbors(cellId, faceIds, cellIds);
526 #ifdef SHOW_COINCIDING_3D_PAL21924
527             int nbNeighbors = 0;
528             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
529               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
530             }
531             bool process = nbNeighbors <= 0;
532 #else
533             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
534 #endif
535             if ( process || myShowInside ||
536                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
537               {
538               for ( i=0; i < numFacePts; i++)
539                 aNewPts[i] = pts[faceVerts[i]];
540               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
541               if(myStoreMapping)
542                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
543               outputCD->CopyData(cd,cellId,newCellId);
544               }
545             }
546           break;
547         }
548         case VTK_HEXAGONAL_PRISM: {
549 #ifdef SHOW_COINCIDING_3D_PAL21924
550           faceIdsTmp->Reset();
551           for (int ai=0; ai<npts; ai++)
552             faceIdsTmp->InsertNextId(pts[ai]);
553           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
554 #endif
555           for (faceId = 0; faceId < 8; faceId++)
556           {
557             faceVerts = vtkHexagonalPrism::GetFaceArray(faceId);
558             faceIds->Reset();
559             faceIds->InsertNextId(pts[faceVerts[0]]);
560             faceIds->InsertNextId(pts[faceVerts[1]]);
561             faceIds->InsertNextId(pts[faceVerts[2]]);
562             faceIds->InsertNextId(pts[faceVerts[3]]);
563             aCellType = VTK_QUAD;
564             numFacePts = 4;
565             if (faceVerts[5] >= 0)
566             {
567               faceIds->InsertNextId(pts[faceVerts[4]]);
568               faceIds->InsertNextId(pts[faceVerts[5]]);
569               aCellType = VTK_POLYGON;
570               numFacePts = 6;
571             }
572             input->GetCellNeighbors(cellId, faceIds, cellIds);
573 #ifdef SHOW_COINCIDING_3D_PAL21924
574             int nbNeighbors = 0;
575             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
576               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
577             }
578             bool process = nbNeighbors <= 0;
579 #else
580             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
581 #endif
582             if ( process || myShowInside ||
583                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
584             {
585               for ( i=0; i < numFacePts; i++)
586                 aNewPts[i] = pts[faceVerts[i]];
587               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
588               if(myStoreMapping)
589                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
590               outputCD->CopyData(cd,cellId,newCellId);
591             }
592             }
593           break;
594         }
595         case VTK_PYRAMID: {
596 #ifdef SHOW_COINCIDING_3D_PAL21924
597           faceIdsTmp->Reset();
598           for (int ai=0; ai<npts; ai++)
599             faceIdsTmp->InsertNextId(pts[ai]);
600           input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
601 #endif
602           for (faceId = 0; faceId < 5; faceId++)
603             {
604             faceIds->Reset();
605             faceVerts = vtkPyramid::GetFaceArray(faceId);
606             faceIds->InsertNextId(pts[faceVerts[0]]);
607             faceIds->InsertNextId(pts[faceVerts[1]]);
608             faceIds->InsertNextId(pts[faceVerts[2]]);
609             aCellType = VTK_TRIANGLE;
610             numFacePts = 3;
611             if (faceVerts[3] >= 0)
612               {
613               faceIds->InsertNextId(pts[faceVerts[3]]);
614               aCellType = VTK_QUAD;
615               numFacePts = 4;
616               }
617             input->GetCellNeighbors(cellId, faceIds, cellIds);
618 #ifdef SHOW_COINCIDING_3D_PAL21924
619             int nbNeighbors = 0;
620             for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
621               if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
622             }
623             bool process = nbNeighbors <= 0;
624 #else
625             bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
626 #endif
627             if ( process || myShowInside ||
628                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
629               {
630               for ( i=0; i < numFacePts; i++)
631                 aNewPts[i] = pts[faceVerts[i]];
632               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
633               if(myStoreMapping)
634                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
635               outputCD->CopyData(cd,cellId,newCellId);
636               }
637             }
638           break;
639         }
640
641 #if VTK_XVERSION > 50700
642         case VTK_POLYHEDRON:
643           {
644             //MESSAGE("VTK_POLYHEDRON geometry filter");
645             vtkIdType nFaces = 0;
646             vtkIdType* ptIds = 0;
647             int idp = 0;
648             input->GetFaceStream(cellId, nFaces, ptIds);
649 #ifdef SHOW_COINCIDING_3D_PAL21924
650             faceIdsTmp->Reset();
651             for (int ai=0; ai<npts; ai++)
652               faceIdsTmp->InsertNextId(pts[ai]);
653             input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
654 #endif
655             for (faceId = 0; faceId < nFaces; faceId++)
656               {
657                 faceIds->Reset();
658                 numFacePts = ptIds[idp];
659                 //MESSAGE("numFacePts="<< numFacePts);
660                 int pt0 = ++idp;
661                 for (i = 0; i < numFacePts; i++)
662                   {
663                     //MESSAGE("ptIds[" << idp + i << "]=" << ptIds[idp + i]);
664                     faceIds->InsertNextId(ptIds[idp + i]);
665                   }
666                 idp += numFacePts;
667                 switch (numFacePts)
668                   {
669                   case 3:
670                     aCellType = VTK_TRIANGLE;
671                     break;
672                   case 4:
673                     aCellType = VTK_QUAD;
674                     break;
675                   default:
676                     aCellType = VTK_POLYGON;
677                     break;
678                   }
679                 // TODO understand and fix display of several polyhedrons                
680                 input->GetCellNeighbors(cellId, faceIds, cellIds);
681 #ifdef SHOW_COINCIDING_3D_PAL21924
682                 int nbNeighbors = 0;
683                 for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
684                   if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
685                 }
686                 bool process = nbNeighbors <= 0;
687 #else
688                 bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
689 #endif
690                 if (process || myShowInside
691                     || (!allVisible && !cellVis[cellIds->GetId(0)]))
692                   {
693                     for (i = 0; i < numFacePts; i++)
694                       aNewPts[i] = ptIds[pt0 + i];
695                     newCellId = output->InsertNextCell(aCellType, numFacePts, aNewPts);
696                     if (myStoreMapping)
697                       InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
698                     outputCD->CopyData(cd, cellId, newCellId);
699                   }
700               }
701             break;
702           }
703 #endif
704         //Quadratic cells
705         case VTK_QUADRATIC_EDGE:
706         case VTK_QUADRATIC_TRIANGLE:
707         case VTK_QUADRATIC_QUAD:
708         case VTK_BIQUADRATIC_QUAD:
709         case VTK_QUADRATIC_TETRA:
710         case VTK_QUADRATIC_HEXAHEDRON:
711         case VTK_TRIQUADRATIC_HEXAHEDRON:
712         case VTK_QUADRATIC_WEDGE:
713         case VTK_QUADRATIC_PYRAMID:
714           if(!myIsWireframeMode) {
715             input->GetCell(cellId,cell);
716             vtkIdList *lpts = vtkIdList::New();
717             vtkPoints *coords = vtkPoints::New();
718             vtkIdList *cellIds = vtkIdList::New();
719             vtkIdType newCellId;
720
721             if ( cell->GetCellDimension() == 1 ) {
722               vtkIdType arcResult = -1;
723               if(myIsBuildArc) {
724                 arcResult = Build1DArc(cellId, input, output, pts, myMaxArcAngle);
725                 newCellId = arcResult;
726               }
727
728               if(!myIsBuildArc || arcResult == -1 ) {
729                 aCellType = VTK_LINE;
730                 numFacePts = 2;
731                 cell->Triangulate(0,lpts,coords);
732                 for (i=0; i < lpts->GetNumberOfIds(); i+=2) {
733                   aNewPts[0] = lpts->GetId(i);
734                   aNewPts[1] = lpts->GetId(i+1);
735                   newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
736                   if(myStoreMapping)
737                     InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
738                   outputCD->CopyData(cd,cellId,newCellId);
739                 }
740               }
741               else {
742                 if(myStoreMapping)
743                   InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
744                 outputCD->CopyData(cd,cellId,newCellId);
745               }
746             }
747             else if ( cell->GetCellDimension() == 2 ) {
748               if(!myIsBuildArc) {
749                 aCellType = VTK_TRIANGLE;
750                 numFacePts = 3;
751                 cell->Triangulate(0,lpts,coords);
752                 for (i=0; i < lpts->GetNumberOfIds(); i+=3) {
753                   aNewPts[0] = lpts->GetId(i);
754                   aNewPts[1] = lpts->GetId(i+1);
755                   aNewPts[2] = lpts->GetId(i+2);
756                   newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
757                   if(myStoreMapping)
758                     InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
759                   outputCD->CopyData(cd,cellId,newCellId);
760                 }
761               }
762               else{
763                 BuildArcedPolygon(cellId,input,output,aDimension2VTK2ObjIds,true);
764               }
765             }
766             else //3D nonlinear cell
767             {
768 #ifdef SHOW_COINCIDING_3D_PAL21924
769               faceIdsTmp->Reset();
770               int npts1 = 0;
771               switch (aCellType ){
772               case VTK_QUADRATIC_TETRA:         npts1 = 4; break;
773               case VTK_QUADRATIC_HEXAHEDRON:    npts1 = 8; break;
774               case VTK_TRIQUADRATIC_HEXAHEDRON: npts1 = 8; break;
775               case VTK_QUADRATIC_WEDGE:         npts1 = 6; break;
776               case VTK_QUADRATIC_PYRAMID:       npts1 = 5; break;
777               }
778               if ( npts1 > 0 ) {
779                 for (int ai=0; ai<npts; ai++)
780                   faceIdsTmp->InsertNextId(pts[ai]);
781                 input->GetCellNeighbors(cellId, faceIdsTmp, cellIdsTmp);
782               }
783 #endif
784               aCellType = VTK_TRIANGLE;
785               numFacePts = 3;
786               for (int j=0; j < cell->GetNumberOfFaces(); j++){
787                 vtkCell *face = cell->GetFace(j);
788                 input->GetCellNeighbors(cellId, face->PointIds, cellIds);
789 #ifdef SHOW_COINCIDING_3D_PAL21924
790                 int nbNeighbors = 0;
791                 for(int ai=0;ai<cellIds->GetNumberOfIds();ai++) {
792                   if (cellIdsTmp->IsId(cellIds->GetId(ai)) == -1) nbNeighbors++;
793                 }
794                 bool process = nbNeighbors <= 0;
795 #else
796                 bool process = cellIds->GetNumberOfIds() <= 0 || GetAppendCoincident3D();
797 #endif
798                 if ( process || myShowInside ) {
799                   face->Triangulate(0,lpts,coords);
800                   for (i=0; i < lpts->GetNumberOfIds(); i+=3) {
801                     aNewPts[0] = lpts->GetId(i);
802                     aNewPts[1] = lpts->GetId(i+1);
803                     aNewPts[2] = lpts->GetId(i+2);
804                     newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
805                     if(myStoreMapping)
806                       InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
807                     outputCD->CopyData(cd,cellId,newCellId);
808                   }
809                 }
810               }
811             } //3d nonlinear cell
812             cellIds->Delete();
813             coords->Delete();
814             lpts->Delete();
815             break;
816           } else { // wireframe
817             switch(aCellType) {
818             case VTK_QUADRATIC_EDGE: {
819               vtkIdType arcResult =-1;
820               if(myIsBuildArc) {
821                arcResult = Build1DArc(cellId, input, output, pts,myMaxArcAngle);
822                newCellId = arcResult;
823               }
824               if(!myIsBuildArc || arcResult == -1) {
825                 aCellType = VTK_POLY_LINE;
826                 numFacePts = 3;
827
828                 aNewPts[0] = pts[0];
829                 aNewPts[2] = pts[1];
830                 aNewPts[1] = pts[2];
831
832                 newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
833               }
834
835               if(myStoreMapping)
836                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
837
838               outputCD->CopyData(cd,cellId,newCellId);
839               break;
840             }
841             case VTK_QUADRATIC_TRIANGLE: {
842               if(!myIsBuildArc) {
843                 aCellType = VTK_POLYGON;
844                 numFacePts = 6;
845
846                 aNewPts[0] = pts[0];
847                 aNewPts[1] = pts[3];
848                 aNewPts[2] = pts[1];
849                 aNewPts[3] = pts[4];
850                 aNewPts[4] = pts[2];
851                 aNewPts[5] = pts[5];
852
853                 newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
854                 if(myStoreMapping)
855                   InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
856
857                 outputCD->CopyData(cd,cellId,newCellId);
858               }
859               else
860                 BuildArcedPolygon(cellId,input,output,aDimension2VTK2ObjIds);
861               break;
862             }
863             case VTK_BIQUADRATIC_QUAD:
864             case VTK_QUADRATIC_QUAD: {
865               if(!myIsBuildArc) {
866                 aCellType = VTK_POLYGON;
867                 numFacePts = 8;
868
869                 aNewPts[0] = pts[0];
870                 aNewPts[1] = pts[4];
871                 aNewPts[2] = pts[1];
872                 aNewPts[3] = pts[5];
873                 aNewPts[4] = pts[2];
874                 aNewPts[5] = pts[6];
875                 aNewPts[6] = pts[3];
876                 aNewPts[7] = pts[7];
877
878                 newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
879                 if(myStoreMapping)
880                   InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
881
882                 outputCD->CopyData(cd,cellId,newCellId);
883               }
884               else
885                 BuildArcedPolygon(cellId,input,output,aDimension2VTK2ObjIds);
886               break;
887             }
888             case VTK_QUADRATIC_TETRA: {
889               aCellType = VTK_POLYGON;
890               numFacePts = 6;
891
892               //---------------------------------------------------------------
893               aNewPts[0] = pts[0];
894               aNewPts[1] = pts[4];
895               aNewPts[2] = pts[1];
896               aNewPts[3] = pts[5];
897               aNewPts[4] = pts[2];
898               aNewPts[5] = pts[6];
899
900               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
901               if(myStoreMapping)
902                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
903
904               outputCD->CopyData(cd,cellId,newCellId);
905
906               //---------------------------------------------------------------
907               aNewPts[0] = pts[0];
908               aNewPts[1] = pts[7];
909               aNewPts[2] = pts[3];
910               aNewPts[3] = pts[8];
911               aNewPts[4] = pts[1];
912               aNewPts[5] = pts[4];
913
914               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
915               if(myStoreMapping)
916                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
917
918               outputCD->CopyData(cd,cellId,newCellId);
919
920               //---------------------------------------------------------------
921               aNewPts[0] = pts[1];
922               aNewPts[1] = pts[8];
923               aNewPts[2] = pts[3];
924               aNewPts[3] = pts[9];
925               aNewPts[4] = pts[2];
926               aNewPts[5] = pts[5];
927
928               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
929               if(myStoreMapping)
930                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
931
932               outputCD->CopyData(cd,cellId,newCellId);
933
934               //---------------------------------------------------------------
935               aNewPts[0] = pts[2];
936               aNewPts[1] = pts[9];
937               aNewPts[2] = pts[3];
938               aNewPts[3] = pts[7];
939               aNewPts[4] = pts[0];
940               aNewPts[5] = pts[6];
941
942               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
943               if(myStoreMapping)
944                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
945
946               outputCD->CopyData(cd,cellId,newCellId);
947
948               break;
949             }
950             case VTK_QUADRATIC_WEDGE: {
951               aCellType = VTK_POLYGON;
952               numFacePts = 6;
953               //---------------------------------------------------------------
954               //Face 1
955               aNewPts[0] = pts[0];
956               aNewPts[1] = pts[6];
957               aNewPts[2] = pts[1];
958               aNewPts[3] = pts[7];
959               aNewPts[4] = pts[2];
960               aNewPts[5] = pts[8];
961               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
962               if(myStoreMapping)
963                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
964               outputCD->CopyData(cd,cellId,newCellId);
965
966               //---------------------------------------------------------------
967               //Face 2
968               aNewPts[0] = pts[3];
969               aNewPts[1] = pts[9];
970               aNewPts[2] = pts[4];
971               aNewPts[3] = pts[10];
972               aNewPts[4] = pts[5];
973               aNewPts[5] = pts[11];
974               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
975               if(myStoreMapping)
976                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
977               outputCD->CopyData(cd,cellId,newCellId);
978
979               //---------------------------------------------------------------
980               //Face 3
981               numFacePts = 8;
982               aNewPts[0] = pts[0];
983               aNewPts[1] = pts[8];
984               aNewPts[2] = pts[2];
985               aNewPts[3] = pts[14];
986               aNewPts[4] = pts[5];
987               aNewPts[5] = pts[11];
988               aNewPts[6] = pts[3];
989               aNewPts[7] = pts[12];
990               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
991               if(myStoreMapping)
992                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
993               outputCD->CopyData(cd,cellId,newCellId);
994
995               //---------------------------------------------------------------
996               //Face 4
997               aNewPts[0] = pts[1];
998               aNewPts[1] = pts[13];
999               aNewPts[2] = pts[4];
1000               aNewPts[3] = pts[10];
1001               aNewPts[4] = pts[5];
1002               aNewPts[5] = pts[14];
1003               aNewPts[6] = pts[2];
1004               aNewPts[7] = pts[7];
1005               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1006               if(myStoreMapping)
1007                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1008               outputCD->CopyData(cd,cellId,newCellId);
1009
1010               //---------------------------------------------------------------
1011               //Face 5
1012               aNewPts[0] = pts[0];
1013               aNewPts[1] = pts[12];
1014               aNewPts[2] = pts[3];
1015               aNewPts[3] = pts[9];
1016               aNewPts[4] = pts[4];
1017               aNewPts[5] = pts[13];
1018               aNewPts[6] = pts[1];
1019               aNewPts[7] = pts[6];
1020               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1021               if(myStoreMapping)
1022                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1023               outputCD->CopyData(cd,cellId,newCellId);
1024               break;
1025             }
1026             case VTK_TRIQUADRATIC_HEXAHEDRON:
1027             case VTK_QUADRATIC_HEXAHEDRON: {
1028               aCellType = VTK_POLYGON;
1029               numFacePts = 8;
1030
1031               //---------------------------------------------------------------
1032               aNewPts[0] = pts[0];
1033               aNewPts[1] = pts[8];
1034               aNewPts[2] = pts[1];
1035               aNewPts[3] = pts[17];
1036               aNewPts[4] = pts[5];
1037               aNewPts[5] = pts[12];
1038               aNewPts[6] = pts[4];
1039               aNewPts[7] = pts[16];
1040
1041               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1042               if(myStoreMapping)
1043                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1044
1045               outputCD->CopyData(cd,cellId,newCellId);
1046
1047               //---------------------------------------------------------------
1048               aNewPts[0] = pts[1];
1049               aNewPts[1] = pts[9];
1050               aNewPts[2] = pts[2];
1051               aNewPts[3] = pts[18];
1052               aNewPts[4] = pts[6];
1053               aNewPts[5] = pts[13];
1054               aNewPts[6] = pts[5];
1055               aNewPts[7] = pts[17];
1056
1057               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1058               if(myStoreMapping)
1059                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1060
1061               outputCD->CopyData(cd,cellId,newCellId);
1062
1063               //---------------------------------------------------------------
1064               aNewPts[0] = pts[2];
1065               aNewPts[1] = pts[10];
1066               aNewPts[2] = pts[3];
1067               aNewPts[3] = pts[19];
1068               aNewPts[4] = pts[7];
1069               aNewPts[5] = pts[14];
1070               aNewPts[6] = pts[6];
1071               aNewPts[7] = pts[18];
1072
1073               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1074               if(myStoreMapping)
1075                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1076
1077               outputCD->CopyData(cd,cellId,newCellId);
1078
1079               //---------------------------------------------------------------
1080               aNewPts[0] = pts[3];
1081               aNewPts[1] = pts[11];
1082               aNewPts[2] = pts[0];
1083               aNewPts[3] = pts[16];
1084               aNewPts[4] = pts[4];
1085               aNewPts[5] = pts[15];
1086               aNewPts[6] = pts[7];
1087               aNewPts[7] = pts[19];
1088
1089               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1090               if(myStoreMapping)
1091                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1092
1093               outputCD->CopyData(cd,cellId,newCellId);
1094
1095               //---------------------------------------------------------------
1096               aNewPts[0] = pts[0];
1097               aNewPts[1] = pts[8];
1098               aNewPts[2] = pts[1];
1099               aNewPts[3] = pts[9];
1100               aNewPts[4] = pts[2];
1101               aNewPts[5] = pts[10];
1102               aNewPts[6] = pts[3];
1103               aNewPts[7] = pts[11];
1104
1105               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1106               if(myStoreMapping)
1107                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1108
1109               outputCD->CopyData(cd,cellId,newCellId);
1110
1111               //---------------------------------------------------------------
1112               aNewPts[0] = pts[4];
1113               aNewPts[1] = pts[12];
1114               aNewPts[2] = pts[5];
1115               aNewPts[3] = pts[13];
1116               aNewPts[4] = pts[6];
1117               aNewPts[5] = pts[14];
1118               aNewPts[6] = pts[7];
1119               aNewPts[7] = pts[15];
1120
1121               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1122               if(myStoreMapping)
1123                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1124
1125               outputCD->CopyData(cd,cellId,newCellId);
1126
1127               break;
1128             }
1129             case VTK_QUADRATIC_PYRAMID: {
1130               aCellType = VTK_POLYGON;
1131               numFacePts = 6;
1132
1133               //---------------------------------------------------------------
1134               aNewPts[0] = pts[0];
1135               aNewPts[1] = pts[8];
1136               aNewPts[2] = pts[3];
1137               aNewPts[3] = pts[12];
1138               aNewPts[4] = pts[4];
1139               aNewPts[5] = pts[9];
1140
1141               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1142               if(myStoreMapping)
1143                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1144
1145               outputCD->CopyData(cd,cellId,newCellId);
1146
1147               //---------------------------------------------------------------
1148               aNewPts[0] = pts[0];
1149               aNewPts[1] = pts[9];
1150               aNewPts[2] = pts[4];
1151               aNewPts[3] = pts[10];
1152               aNewPts[4] = pts[1];
1153               aNewPts[5] = pts[5];
1154
1155               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1156               if(myStoreMapping)
1157                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1158
1159               outputCD->CopyData(cd,cellId,newCellId);
1160
1161               //---------------------------------------------------------------
1162               aNewPts[0] = pts[1];
1163               aNewPts[1] = pts[10];
1164               aNewPts[2] = pts[4];
1165               aNewPts[3] = pts[11];
1166               aNewPts[4] = pts[2];
1167               aNewPts[5] = pts[6];
1168
1169               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1170               if(myStoreMapping)
1171                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1172
1173               outputCD->CopyData(cd,cellId,newCellId);
1174
1175               //---------------------------------------------------------------
1176               aNewPts[0] = pts[2];
1177               aNewPts[1] = pts[11];
1178               aNewPts[2] = pts[4];
1179               aNewPts[3] = pts[12];
1180               aNewPts[4] = pts[3];
1181               aNewPts[5] = pts[7];
1182
1183               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1184               if(myStoreMapping)
1185                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1186
1187               outputCD->CopyData(cd,cellId,newCellId);
1188
1189               //---------------------------------------------------------------
1190               numFacePts = 8;
1191               aNewPts[0] = pts[0];
1192               aNewPts[1] = pts[5];
1193               aNewPts[2] = pts[1];
1194               aNewPts[3] = pts[6];
1195               aNewPts[4] = pts[2];
1196               aNewPts[5] = pts[7];
1197               aNewPts[6] = pts[3];
1198               aNewPts[7] = pts[8];
1199
1200               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1201               if(myStoreMapping)
1202                 InsertId( cellId, aCellType, myVTK2ObjIds, aDimension2VTK2ObjIds );
1203
1204               outputCD->CopyData(cd,cellId,newCellId);
1205
1206               break;
1207             } // case VTK_QUADRATIC_PYRAMID:
1208             } // switch by type
1209           } // end WIREFRAME
1210           break;
1211         } //switch
1212       } //if visible
1213     } //for all cells
1214
1215   output->Squeeze();
1216
1217   vtkDebugMacro(<<"Extracted " << input->GetNumberOfPoints() << " points,"
1218                 << output->GetNumberOfCells() << " cells.");
1219
1220   cell->Delete();
1221
1222   cellIds->Delete();
1223   faceIds->Delete();
1224   cellIdsTmp->Delete();
1225   faceIdsTmp->Delete();
1226
1227   if ( cellVis )
1228   {
1229     delete [] cellVis;
1230   }
1231
1232   // to sort myVTK2ObjIds vector by cell dimension (ascending)
1233   if( myStoreMapping && !aDimension2VTK2ObjIds.empty() )
1234   {
1235     myVTK2ObjIds.clear();
1236     for( vtkIdType aDimension = 0; aDimension <= 2; aDimension++ )
1237     {
1238       TMapOfVectorId::const_iterator anIter = aDimension2VTK2ObjIds.find( aDimension );
1239       if( anIter != aDimension2VTK2ObjIds.end() )
1240       {
1241         const TVectorId& aCellIds = anIter->second;
1242         TVectorId::const_iterator anIdIter, anIdIterEnd = aCellIds.end();
1243         for( anIdIter = aCellIds.begin(); anIdIter != anIdIterEnd; anIdIter++ )
1244         {
1245           const vtkIdType aCellId = *anIdIter;
1246           myVTK2ObjIds.push_back( aCellId );
1247         }
1248       }
1249     }
1250   }
1251
1252   return 1;
1253 }
1254
1255 void
1256 VTKViewer_GeometryFilter
1257 ::InsertId( const vtkIdType theCellId,
1258             const vtkIdType theCellType,
1259             TVectorId& theVTK2ObjIds,
1260             TMapOfVectorId& theDimension2VTK2ObjIds )
1261 {
1262   theVTK2ObjIds.push_back( theCellId );
1263
1264   int aDimension = 0;
1265   switch( theCellType )
1266   {
1267     case VTK_VERTEX:
1268     case VTK_POLY_VERTEX:
1269       aDimension = 0;
1270       break;
1271     case VTK_LINE:
1272     case VTK_POLY_LINE:
1273       aDimension = 1;
1274       break;
1275     case VTK_TRIANGLE:
1276     case VTK_TRIANGLE_STRIP:
1277     case VTK_POLYGON:
1278     case VTK_PIXEL:
1279     case VTK_QUAD:
1280       aDimension = 2;
1281       break;
1282   }
1283
1284   TVectorId& aCellIds = theDimension2VTK2ObjIds[ aDimension ];
1285   aCellIds.push_back( theCellId );
1286 }
1287
1288 void
1289 VTKViewer_GeometryFilter
1290 ::SetInside(int theShowInside)
1291 {
1292   if(myShowInside == theShowInside)
1293     return;
1294
1295   myShowInside = theShowInside;
1296   this->Modified();
1297 }
1298
1299 int
1300 VTKViewer_GeometryFilter
1301 ::GetInside()
1302 {
1303   return myShowInside;
1304 }
1305
1306
1307 void
1308 VTKViewer_GeometryFilter
1309 ::SetWireframeMode(int theIsWireframeMode)
1310 {
1311   if(myIsWireframeMode == theIsWireframeMode)
1312     return;
1313
1314   myIsWireframeMode = theIsWireframeMode;
1315   this->Modified();
1316 }
1317
1318 int
1319 VTKViewer_GeometryFilter
1320 ::GetWireframeMode()
1321 {
1322   return myIsWireframeMode;
1323 }
1324
1325
1326 void
1327 VTKViewer_GeometryFilter
1328 ::SetStoreMapping(int theStoreMapping)
1329 {
1330   if(myStoreMapping == theStoreMapping)
1331     return;
1332
1333   myStoreMapping = theStoreMapping;
1334   this->Modified();
1335 }
1336
1337 int
1338 VTKViewer_GeometryFilter
1339 ::GetStoreMapping()
1340 {
1341   return myStoreMapping;
1342 }
1343
1344
1345 vtkIdType VTKViewer_GeometryFilter::GetElemObjId( int theVtkID )
1346 {
1347   if( theVtkID < 0 || theVtkID >= (int)myVTK2ObjIds.size() )
1348     return -1;
1349   return myVTK2ObjIds[theVtkID];
1350 }
1351
1352
1353 void VTKViewer_GeometryFilter::BuildArcedPolygon(vtkIdType cellId,
1354                                                  vtkUnstructuredGrid* input,
1355                                                  vtkPolyData *output,
1356                                                  TMapOfVectorId& theDimension2VTK2ObjIds,
1357                                                  bool triangulate){
1358   vtkIdType aCellType = VTK_POLYGON;
1359   vtkIdType *aNewPoints = NULL;
1360   vtkIdType aNbPoints = 0;
1361   vtkIdType newCellId;
1362
1363   //Input and output cell data
1364   vtkCellData *cd = input->GetCellData();
1365   vtkCellData *outputCD = output->GetCellData();
1366
1367   //Input and output scalars on point data
1368   vtkDataArray* inputScalars = input->GetPointData()->GetScalars();
1369   vtkDataArray* outputScalars = output->GetPointData()->GetScalars();
1370
1371   std::vector<vtkPoints*> aCollection;
1372   std::vector< std::vector<double> > aScalarCollection;
1373
1374   vtkCell* aCell = input->GetCell(cellId);
1375   switch(aCell->GetCellType()) {
1376     case VTK_QUADRATIC_TRIANGLE:
1377     {
1378       //Get All points from input cell
1379       Pnt P0 = CreatePnt( aCell, inputScalars, 0 );
1380       Pnt P1 = CreatePnt( aCell, inputScalars, 1 );
1381       Pnt P2 = CreatePnt( aCell, inputScalars, 2 );
1382       Pnt P3 = CreatePnt( aCell, inputScalars, 3 );
1383       Pnt P4 = CreatePnt( aCell, inputScalars, 4 );
1384       Pnt P5 = CreatePnt( aCell, inputScalars, 5 );
1385
1386       VTKViewer_ArcBuilder aBuilder1(P0,P3,P1,myMaxArcAngle); //Build arc using 0, 3 and 1 points
1387 #ifdef __MYDEBUG__
1388       cout << "Quadrangle arc 1 " << ( aBuilder1.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1389 #endif
1390
1391       VTKViewer_ArcBuilder aBuilder2(P1,P4,P2,myMaxArcAngle); //Build arc using 1, 4 and 2 points
1392 #ifdef __MYDEBUG__
1393       cout << "Quadrangle arc 2 " << ( aBuilder2.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1394 #endif
1395
1396       VTKViewer_ArcBuilder aBuilder3(P2,P5,P0,myMaxArcAngle); //Build arc using 2, 5 and 0 points
1397 #ifdef __MYDEBUG__
1398       cout << "Quadrangle arc 3 " << ( aBuilder3.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1399 #endif
1400
1401       aCollection.push_back(aBuilder1.GetPoints());
1402       aCollection.push_back(aBuilder2.GetPoints());
1403       aCollection.push_back(aBuilder3.GetPoints());
1404
1405       aScalarCollection.push_back(aBuilder1.GetScalarValues());
1406       aScalarCollection.push_back(aBuilder2.GetScalarValues());
1407       aScalarCollection.push_back(aBuilder3.GetScalarValues());
1408       break;
1409     }
1410     case VTK_BIQUADRATIC_QUAD:
1411     case VTK_QUADRATIC_QUAD:
1412     {
1413       //Get All points from input cell
1414       Pnt P0 = CreatePnt( aCell, inputScalars, 0 );
1415       Pnt P1 = CreatePnt( aCell, inputScalars, 1 );
1416       Pnt P2 = CreatePnt( aCell, inputScalars, 2 );
1417       Pnt P3 = CreatePnt( aCell, inputScalars, 3 );
1418       Pnt P4 = CreatePnt( aCell, inputScalars, 4 );
1419       Pnt P5 = CreatePnt( aCell, inputScalars, 5 );
1420       Pnt P6 = CreatePnt( aCell, inputScalars, 6 );
1421       Pnt P7 = CreatePnt( aCell, inputScalars, 7 );
1422
1423       VTKViewer_ArcBuilder aBuilder1(P0,P4,P1,myMaxArcAngle); //Build arc using 0, 4 and 1 points
1424 #ifdef __MYDEBUG__
1425       cout << "Quadrangle arc 1 " << ( aBuilder1.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1426 #endif
1427
1428       VTKViewer_ArcBuilder aBuilder2(P1,P5,P2,myMaxArcAngle); //Build arc using 1, 5 and 2 points
1429 #ifdef __MYDEBUG__
1430       cout << "Quadrangle arc 2 " << ( aBuilder2.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1431 #endif
1432
1433       VTKViewer_ArcBuilder aBuilder3(P2,P6,P3,myMaxArcAngle); //Build arc using 2, 6 and 3 points
1434 #ifdef __MYDEBUG__
1435       cout << "Quadrangle arc 3 " << ( aBuilder3.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1436 #endif
1437
1438       VTKViewer_ArcBuilder aBuilder4(P3,P7,P0,myMaxArcAngle); //Build arc using 3, 7 and 0 points
1439 #ifdef __MYDEBUG__
1440       cout << "Quadrangle arc 4 " << ( aBuilder4.GetStatus() == VTKViewer_ArcBuilder::Arc_Done ? "" : "NOT " ) << "done !!!" << endl;
1441 #endif
1442
1443       aCollection.push_back(aBuilder1.GetPoints());
1444       aCollection.push_back(aBuilder2.GetPoints());
1445       aCollection.push_back(aBuilder3.GetPoints());
1446       aCollection.push_back(aBuilder4.GetPoints());
1447
1448       aScalarCollection.push_back(aBuilder1.GetScalarValues());
1449       aScalarCollection.push_back(aBuilder2.GetScalarValues());
1450       aScalarCollection.push_back(aBuilder3.GetScalarValues());
1451       aScalarCollection.push_back(aBuilder4.GetScalarValues());
1452       break;
1453     }
1454     default: //Unsupported cell type
1455       return;
1456   }
1457
1458   if(triangulate){
1459     const vtkIdType numFacePts = 3;
1460     vtkIdList *pts = vtkIdList::New();
1461     vtkPoints *coords = vtkPoints::New();
1462     aCellType = VTK_TRIANGLE;
1463     vtkIdType aNewPts[numFacePts];
1464     vtkIdType aTriangleId;
1465
1466     vtkPolygon *aPlg = vtkPolygon::New();
1467     std::map<int, double> aPntId2ScalarValue;
1468     aNbPoints = MergevtkPoints(aCollection, aScalarCollection, aPlg->GetPoints(), aPntId2ScalarValue, aNewPoints);
1469     aPlg->GetPointIds()->SetNumberOfIds(aNbPoints);
1470
1471     for(vtkIdType i = 0; i < aNbPoints;i++) {
1472       aPlg->GetPointIds()->SetId(i, aNewPoints[i]);
1473     }
1474
1475     aPlg->Triangulate(0,pts,coords);
1476
1477     for (vtkIdType i=0; i < pts->GetNumberOfIds(); i+=3) {
1478       aNewPts[0] = output->GetPoints()->InsertNextPoint(coords->GetPoint(i));
1479       aNewPts[1] = output->GetPoints()->InsertNextPoint(coords->GetPoint(i+1));
1480       aNewPts[2] = output->GetPoints()->InsertNextPoint(coords->GetPoint(i+2));
1481
1482       if(outputScalars) {
1483         outputScalars->InsertNextTuple1(aPntId2ScalarValue[pts->GetId(i)]);
1484         outputScalars->InsertNextTuple1(aPntId2ScalarValue[pts->GetId(i+1)]);
1485         outputScalars->InsertNextTuple1(aPntId2ScalarValue[pts->GetId(i+2)]);
1486       }
1487
1488       aTriangleId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
1489
1490       if(myStoreMapping)
1491         InsertId( cellId, aCellType, myVTK2ObjIds, theDimension2VTK2ObjIds );
1492       outputCD->CopyData(cd,cellId,aTriangleId);
1493     }
1494     pts->Delete();
1495     coords->Delete();
1496     aPlg->Delete();
1497   }
1498   else {
1499     std::map<int, double> aPntId2ScalarValue;
1500     aNbPoints = MergevtkPoints(aCollection, aScalarCollection, output->GetPoints(), aPntId2ScalarValue, aNewPoints);
1501     if(outputScalars)
1502       for(vtkIdType i = 0; i < aNbPoints; i++)
1503         outputScalars->InsertNextTuple1(aPntId2ScalarValue[aNewPoints[i]]);
1504     newCellId = output->InsertNextCell(aCellType,aNbPoints,aNewPoints);
1505     outputCD->CopyData(cd,cellId,newCellId);
1506
1507     if(myStoreMapping)
1508       InsertId( cellId, aCellType, myVTK2ObjIds, theDimension2VTK2ObjIds );
1509   }
1510
1511   if (aNewPoints)
1512     delete [] aNewPoints;
1513 }
1514
1515
1516 void VTKViewer_GeometryFilter::SetQuadraticArcMode(bool theFlag)
1517 {
1518   if(myIsBuildArc != theFlag) {
1519     myIsBuildArc = theFlag;
1520     this->Modified();
1521   }
1522 }
1523 bool VTKViewer_GeometryFilter::GetQuadraticArcMode() const
1524 {
1525   return myIsBuildArc;
1526 }
1527
1528 void VTKViewer_GeometryFilter::SetQuadraticArcAngle(vtkFloatingPointType theMaxAngle)
1529 {
1530   if(myMaxArcAngle != theMaxAngle) {
1531     myMaxArcAngle = theMaxAngle;
1532     this->Modified();
1533   }
1534 }
1535
1536 vtkFloatingPointType VTKViewer_GeometryFilter:: GetQuadraticArcAngle() const
1537 {
1538   return myMaxArcAngle;
1539 }
1540
1541 int VTKViewer_GeometryFilter::GetAppendCoincident3D() const {
1542 // VSR 26/10/2012: see description of SHOW_COINCIDING_3D_PAL20314
1543 // in the top of this file
1544 #ifdef SHOW_COINCIDING_3D_PAL20314
1545   return myAppendCoincident3D;
1546 #else
1547   return false;
1548 #endif
1549 }
1550
1551 void VTKViewer_GeometryFilter::SetAppendCoincident3D(int theFlag) {
1552   if(myAppendCoincident3D != theFlag){
1553     myAppendCoincident3D = theFlag;
1554     this->Modified();
1555   }
1556 }