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