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