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