Salome HOME
Update of CheckDone
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
1 // Copyright (C) 2010-2024  CEA, EDF, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "SMDS_UnstructuredGrid.hxx"
21 #include "SMDS_Mesh.hxx"
22 #include "SMDS_MeshInfo.hxx"
23 #include "SMDS_Downward.hxx"
24 #include "SMDS_MeshVolume.hxx"
25
26 #include "utilities.h"
27 #include "chrono.hxx"
28
29 #include <vtkCellArray.h>
30 #include <vtkCellData.h>
31 #include <vtkCellLinks.h>
32 #include <vtkDoubleArray.h>
33 #include <vtkIdTypeArray.h>
34 #include <vtkUnsignedCharArray.h>
35
36 #include <list>
37 #include <climits>
38
39 SMDS_CellLinks* SMDS_CellLinks::New()
40 {
41   return new SMDS_CellLinks();
42 }
43
44 void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
45 {
46   if ( vtkID > this->MaxId )
47   {
48     this->MaxId = vtkID;
49     if ( vtkID >= this->Size )
50       vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
51   }
52 }
53
54 void SMDS_CellLinks::BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types)
55 {
56   // build links taking into account removed cells
57
58   vtkIdType numPts = data->GetNumberOfPoints();
59   vtkIdType j, cellId = 0;
60   unsigned short *linkLoc;
61   vtkIdType npts=0;
62   vtkIdType const *pts(nullptr);
63   vtkIdType loc = Connectivity->GetTraversalLocation();
64
65   // traverse data to determine number of uses of each point
66   cellId = 0;
67   for (Connectivity->InitTraversal();
68        Connectivity->GetNextCell(npts,pts); cellId++)
69   {
70     if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
71       for (j=0; j < npts; j++)
72       {
73         this->IncrementLinkCount(pts[j]);
74       }
75   }
76
77   // now allocate storage for the links
78   this->AllocateLinks(numPts);
79   this->MaxId = numPts - 1;
80
81   // fill out lists with references to cells
82   linkLoc = new unsigned short[numPts];
83   memset(linkLoc, 0, numPts*sizeof(unsigned short));
84
85   cellId = 0;
86   for (Connectivity->InitTraversal();
87        Connectivity->GetNextCell(npts,pts); cellId++)
88   {
89     if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
90       for (j=0; j < npts; j++)
91       {
92         this->InsertCellReference(pts[j], (linkLoc[pts[j]])++, cellId);
93       }
94   }
95   delete [] linkLoc;
96   Connectivity->SetTraversalLocation(loc);
97 }
98
99 SMDS_CellLinks::SMDS_CellLinks() :
100   vtkCellLinks()
101 {
102 }
103
104 SMDS_CellLinks::~SMDS_CellLinks()
105 {
106 }
107
108 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
109 {
110   return new SMDS_UnstructuredGrid();
111 }
112
113 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
114   vtkUnstructuredGrid()
115 {
116   _cellIdToDownId.clear();
117   _downTypes.clear();
118   _downArray.clear();
119   _mesh = 0;
120 }
121
122 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
123 {
124 }
125
126 vtkMTimeType SMDS_UnstructuredGrid::GetMTime()
127 {
128   vtkMTimeType mtime = vtkUnstructuredGrid::GetMTime();
129   return mtime;
130 }
131
132 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
133 {
134   // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
135   return this->Points;
136 }
137
138 vtkIdType SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
139 {
140   if ( !this->Links ) // don't create Links until they are needed
141   {
142     return this->InsertNextCell(type, npts, pts);
143   }
144
145   if ( type != VTK_POLYHEDRON )
146     return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
147
148   // --- type = VTK_POLYHEDRON
149   vtkIdType cellid = this->InsertNextCell(type, npts, pts);
150
151   std::set<vtkIdType> setOfNodes;
152   setOfNodes.clear();
153   int nbfaces = npts;
154   int i = 0;
155   for (int nf = 0; nf < nbfaces; nf++)
156   {
157     int nbnodes = pts[i];
158     i++;
159     for (int k = 0; k < nbnodes; k++)
160     {
161       if ( setOfNodes.insert( pts[i] ).second )
162       {
163         (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->ResizeCellList( pts[i], 1 );
164         (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->AddCellReference( cellid, pts[i] );
165       }
166       i++;
167     }
168   }
169
170   return cellid;
171 }
172
173 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
174 {
175   _mesh = mesh;
176 }
177
178 void SMDS_UnstructuredGrid::compactGrid(std::vector<smIdType>& idNodesOldToNew, smIdType newNodeSize,
179                                         std::vector<smIdType>& idCellsNewToOld, smIdType newCellSize)
180 {
181   this->DeleteLinks();
182
183   // IDs of VTK nodes always correspond to SMDS IDs but there can be "holes" in SMDS numeration.
184   // We compact only if there were holes
185
186   vtkIdType oldNodeSize = this->GetNumberOfPoints();
187   bool updateNodes = ( oldNodeSize > newNodeSize );
188   if ( true /*updateNodes*/ )
189   {
190     // 21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion
191     // Use double type for storing coordinates of nodes instead float.
192     vtkPoints *newPoints = vtkPoints::New();
193     newPoints->SetDataType( VTK_DOUBLE );
194     newPoints->SetNumberOfPoints( FromSmIdType<vtkIdType>(newNodeSize) );
195
196     vtkIdType i = 0, alreadyCopied = 0;
197     while ( i < oldNodeSize )
198     {
199       // skip a hole if any
200       while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
201         ++i;
202       vtkIdType startBloc = i;
203       // look for a block end
204       while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
205         ++i;
206       vtkIdType endBloc = i;
207       copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
208     }
209     this->SetPoints(newPoints);
210     newPoints->Delete();
211   }
212   else
213   {
214     this->Points->Squeeze();
215     this->Points->Modified();
216   }
217
218   // Compact cells if VTK IDs do not correspond to SMDS IDs or nodes compacted
219
220   vtkIdType oldCellSize = this->Types->GetNumberOfTuples();
221   bool      updateCells = ( updateNodes || newCellSize != oldCellSize );
222   for ( vtkIdType newID = 0, nbIDs = idCellsNewToOld.size(); newID < nbIDs &&  !updateCells; ++newID )
223     updateCells = ( idCellsNewToOld[ newID ] != newID );
224
225   if ( false /*!updateCells*/ ) // no holes in elements
226   {
227     this->Connectivity->Squeeze();
228     this->CellLocations->Squeeze();
229     this->Types->Squeeze();
230     if ( this->FaceLocations )
231     {
232       this->FaceLocations->Squeeze();
233       this->Faces->Squeeze();
234     }
235     this->Connectivity->Modified();
236     return;
237   }
238
239   if ((vtkIdType) idNodesOldToNew.size() < oldNodeSize )
240   {
241     idNodesOldToNew.reserve( oldNodeSize );
242     for ( vtkIdType i = idNodesOldToNew.size(); i < oldNodeSize; ++i )
243       idNodesOldToNew.push_back( i );
244   }
245
246   // --- create new compacted Connectivity, Locations and Types
247
248   vtkIdType newConnectivitySize = this->Connectivity->GetNumberOfConnectivityEntries();
249   if ( newCellSize != oldCellSize )
250     for ( vtkIdType i = 0; i < oldCellSize - 1; ++i )
251       if ( this->Types->GetValue( i ) == VTK_EMPTY_CELL )
252         newConnectivitySize -= this->Connectivity->GetCellSize( i );
253
254   vtkCellArray *newConnectivity = vtkCellArray::New();
255   newConnectivity->Initialize();
256   newConnectivity->Allocate( newConnectivitySize );
257
258   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
259   newTypes->Initialize();
260   newTypes->SetNumberOfValues(FromSmIdType<vtkIdType>(newCellSize));
261
262   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
263   newLocations->Initialize();
264   newLocations->SetNumberOfValues(FromSmIdType<vtkIdType>(newCellSize));
265
266   std::vector< vtkIdType > pointsCell(1024); // --- points id to fill a new cell
267
268   copyBloc(newTypes, idCellsNewToOld, idNodesOldToNew,
269            newConnectivity, newLocations, pointsCell );
270
271   if (vtkDoubleArray* diameters =
272       vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
273   {
274     vtkDoubleArray* newDiameters = vtkDoubleArray::New();
275     newDiameters->SetNumberOfComponents(1);
276     for ( vtkIdType newCellID = 0; newCellID < newCellSize; newCellID++ )
277     {
278       if ( newTypes->GetValue( newCellID ) == VTK_POLY_VERTEX )
279       {
280         vtkIdType oldCellID = idCellsNewToOld[ newCellID ];
281         newDiameters->InsertValue( newCellID, diameters->GetValue( oldCellID ));
282       }
283       vtkDataSet::CellData->SetScalars( newDiameters );
284     }
285   }
286
287   if (this->FaceLocations)
288   {
289     vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
290     newFaceLocations->Initialize();
291     newFaceLocations->Allocate(newTypes->GetSize());
292     vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
293     newFaces->Initialize();
294     newFaces->Allocate(this->Faces->GetSize());
295     for ( vtkIdType newCellID = 0; newCellID < newCellSize; newCellID++ )
296     {
297       if ( newTypes->GetValue( newCellID ) == VTK_POLYHEDRON )
298       {
299         smIdType oldCellId = idCellsNewToOld[ newCellID ];
300         newFaceLocations->InsertNextValue( newFaces->GetMaxId()+1 );
301         smIdType oldFaceLoc = this->FaceLocations->GetValue( FromSmIdType<vtkIdType>(oldCellId) );
302         smIdType nCellFaces = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
303         newFaces->InsertNextValue( FromSmIdType<vtkIdType>(nCellFaces) );
304         for ( int n = 0; n < nCellFaces; n++ )
305         {
306           int nptsInFace = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
307           newFaces->InsertNextValue( nptsInFace );
308           for ( int k = 0; k < nptsInFace; k++ )
309           {
310             vtkIdType oldpt = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
311             newFaces->InsertNextValue( idNodesOldToNew[ oldpt ]);
312           }
313         }
314       }
315       else
316       {
317         newFaceLocations->InsertNextValue(-1);
318       }
319     }
320     newFaceLocations->Squeeze();
321     newFaces->Squeeze();
322     this->SetCells( newTypes, newLocations, newConnectivity, newFaceLocations, newFaces );
323     this->CellLocations = newLocations;
324     newFaceLocations->Delete();
325     newFaces->Delete();
326   }
327   else
328   {
329     this->SetCells( newTypes, newLocations, newConnectivity, FaceLocations, Faces );
330     this->CellLocations = newLocations;
331   }
332
333   newTypes->Delete();
334   newLocations->Delete();
335   newConnectivity->Delete();
336 }
337
338 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *             newPoints,
339                                       std::vector<smIdType>& /*idNodesOldToNew*/,
340                                       vtkIdType&              alreadyCopied,
341                                       vtkIdType               start,
342                                       vtkIdType               end)
343 {
344   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
345   void *source = this->Points->GetVoidPointer(3 * start);
346   vtkIdType nbPoints = end - start;
347   if (nbPoints > 0)
348   {
349     memcpy(target, source, 3 * sizeof(double) * nbPoints);
350     alreadyCopied += nbPoints;
351   }
352 }
353
354 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *  newTypes,
355                                      const std::vector<smIdType>& idCellsNewToOld,
356                                      const std::vector<smIdType>& idNodesOldToNew,
357                                      vtkCellArray*           newConnectivity,
358                                      vtkIdTypeArray*         newLocations,
359                                      std::vector<vtkIdType>& pointsCell)
360 {
361   for ( size_t iNew = 0; iNew < idCellsNewToOld.size(); iNew++ )
362   {
363     vtkIdType iOld = idCellsNewToOld[ iNew ];
364     newTypes->SetValue( iNew, this->Types->GetValue( iOld ));
365
366     vtkIdType oldLoc = ((vtkIdTypeArray *)(this->Connectivity->GetOffsetsArray()))->GetValue( iOld );
367     vtkIdType nbpts;
368     vtkIdType const *oldPtsCell(nullptr);
369     this->Connectivity->GetCell( oldLoc+iOld, nbpts, oldPtsCell );
370     if ((vtkIdType) pointsCell.size() < nbpts )
371       pointsCell.resize( nbpts );
372     for ( int l = 0; l < nbpts; l++ )
373     {
374       vtkIdType oldval = oldPtsCell[l];
375       pointsCell[l] = idNodesOldToNew[oldval];
376     }
377     /*vtkIdType newcnt = */newConnectivity->InsertNextCell( nbpts, pointsCell.data() );
378     vtkIdType newLoc = newConnectivity->GetInsertLocation( nbpts );
379     newLocations->SetValue( iNew, newLoc );
380   }
381 }
382
383 int SMDS_UnstructuredGrid::CellIdToDownId(vtkIdType vtkCellId)
384 {
385   if ((vtkCellId < 0) || (vtkCellId >= (vtkIdType)_cellIdToDownId.size()))
386   {
387     return -1;
388   }
389   return _cellIdToDownId[vtkCellId];
390 }
391
392 void SMDS_UnstructuredGrid::setCellIdToDownId(vtkIdType vtkCellId, int downId)
393 {
394   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
395   _cellIdToDownId[vtkCellId] = downId;
396 }
397
398 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
399 {
400   for (size_t i = 0; i < _downArray.size(); i++)
401   {
402     if (_downArray[i])
403       delete _downArray[i];
404     _downArray[i] = 0;
405   }
406   _cellIdToDownId.clear();
407 }
408
409 /*! Build downward connectivity: to do only when needed because heavy memory load.
410  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
411  *
412  */
413 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool /*withEdges*/)
414 {
415   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
416   // TODO calcul partiel sans edges
417
418   // --- erase previous data if any
419
420   this->CleanDownwardConnectivity();
421
422   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
423
424   _downArray.resize(VTK_MAXTYPE + 1, 0);
425
426   _downArray[VTK_LINE]                    = new SMDS_DownEdge(this);
427   _downArray[VTK_QUADRATIC_EDGE]          = new SMDS_DownQuadEdge(this);
428   _downArray[VTK_TRIANGLE]                = new SMDS_DownTriangle(this);
429   _downArray[VTK_QUADRATIC_TRIANGLE]      = new SMDS_DownQuadTriangle(this);
430   _downArray[VTK_BIQUADRATIC_TRIANGLE]    = new SMDS_DownQuadTriangle(this);
431   _downArray[VTK_QUAD]                    = new SMDS_DownQuadrangle(this);
432   _downArray[VTK_QUADRATIC_QUAD]          = new SMDS_DownQuadQuadrangle(this);
433   _downArray[VTK_BIQUADRATIC_QUAD]        = new SMDS_DownQuadQuadrangle(this);
434   _downArray[VTK_TETRA]                   = new SMDS_DownTetra(this);
435   _downArray[VTK_QUADRATIC_TETRA]         = new SMDS_DownQuadTetra(this);
436   _downArray[VTK_PYRAMID]                 = new SMDS_DownPyramid(this);
437   _downArray[VTK_QUADRATIC_PYRAMID]       = new SMDS_DownQuadPyramid(this);
438   _downArray[VTK_WEDGE]                   = new SMDS_DownPenta(this);
439   _downArray[VTK_QUADRATIC_WEDGE]         = new SMDS_DownQuadPenta(this);
440   _downArray[VTK_HEXAHEDRON]              = new SMDS_DownHexa(this);
441   _downArray[VTK_QUADRATIC_HEXAHEDRON]    = new SMDS_DownQuadHexa(this);
442   _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
443   _downArray[VTK_HEXAGONAL_PRISM]         = new SMDS_DownPenta(this);
444
445   // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
446
447   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
448
449   int nbLinTetra  = FromSmIdType<int>(meshInfo.NbTetras  (ORDER_LINEAR));
450   int nbQuadTetra = FromSmIdType<int>(meshInfo.NbTetras  (ORDER_QUADRATIC));
451   int nbLinPyra   = FromSmIdType<int>(meshInfo.NbPyramids(ORDER_LINEAR));
452   int nbQuadPyra  = FromSmIdType<int>(meshInfo.NbPyramids(ORDER_QUADRATIC));
453   int nbLinPrism  = FromSmIdType<int>(meshInfo.NbPrisms  (ORDER_LINEAR));
454   int nbQuadPrism = FromSmIdType<int>(meshInfo.NbPrisms  (ORDER_QUADRATIC));
455   int nbLinHexa   = FromSmIdType<int>(meshInfo.NbHexas   (ORDER_LINEAR));
456   int nbQuadHexa  = FromSmIdType<int>(meshInfo.NbHexas   (ORDER_QUADRATIC));
457   int nbHexPrism  = FromSmIdType<int>(meshInfo.NbHexPrisms());
458
459   int nbLineGuess     = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
460   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
461   int nbLinTriaGuess  = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
462   int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
463   int nbLinQuadGuess  = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
464   int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
465
466   int GuessSize[VTK_MAXTYPE];
467   GuessSize[VTK_LINE]                    = nbLineGuess;
468   GuessSize[VTK_QUADRATIC_EDGE]          = nbQuadEdgeGuess;
469   GuessSize[VTK_TRIANGLE]                = nbLinTriaGuess;
470   GuessSize[VTK_QUADRATIC_TRIANGLE]      = nbQuadTriaGuess;
471   GuessSize[VTK_BIQUADRATIC_TRIANGLE]    = nbQuadTriaGuess;
472   GuessSize[VTK_QUAD]                    = nbLinQuadGuess;
473   GuessSize[VTK_QUADRATIC_QUAD]          = nbQuadQuadGuess;
474   GuessSize[VTK_BIQUADRATIC_QUAD]        = nbQuadQuadGuess;
475   GuessSize[VTK_TETRA]                   = nbLinTetra;
476   GuessSize[VTK_QUADRATIC_TETRA]         = nbQuadTetra;
477   GuessSize[VTK_PYRAMID]                 = nbLinPyra;
478   GuessSize[VTK_QUADRATIC_PYRAMID]       = nbQuadPyra;
479   GuessSize[VTK_WEDGE]                   = nbLinPrism;
480   GuessSize[VTK_QUADRATIC_WEDGE]         = nbQuadPrism;
481   GuessSize[VTK_HEXAHEDRON]              = nbLinHexa;
482   GuessSize[VTK_QUADRATIC_HEXAHEDRON]    = nbQuadHexa;
483   GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
484   GuessSize[VTK_HEXAGONAL_PRISM]         = nbHexPrism;
485   (void)GuessSize; // unused in Release mode
486
487   _downArray[VTK_LINE]                   ->allocate(nbLineGuess);
488   _downArray[VTK_QUADRATIC_EDGE]         ->allocate(nbQuadEdgeGuess);
489   _downArray[VTK_TRIANGLE]               ->allocate(nbLinTriaGuess);
490   _downArray[VTK_QUADRATIC_TRIANGLE]     ->allocate(nbQuadTriaGuess);
491   _downArray[VTK_BIQUADRATIC_TRIANGLE]   ->allocate(nbQuadTriaGuess);
492   _downArray[VTK_QUAD]                   ->allocate(nbLinQuadGuess);
493   _downArray[VTK_QUADRATIC_QUAD]         ->allocate(nbQuadQuadGuess);
494   _downArray[VTK_BIQUADRATIC_QUAD]       ->allocate(nbQuadQuadGuess);
495   _downArray[VTK_TETRA]                  ->allocate(nbLinTetra);
496   _downArray[VTK_QUADRATIC_TETRA]        ->allocate(nbQuadTetra);
497   _downArray[VTK_PYRAMID]                ->allocate(nbLinPyra);
498   _downArray[VTK_QUADRATIC_PYRAMID]      ->allocate(nbQuadPyra);
499   _downArray[VTK_WEDGE]                  ->allocate(nbLinPrism);
500   _downArray[VTK_QUADRATIC_WEDGE]        ->allocate(nbQuadPrism);
501   _downArray[VTK_HEXAHEDRON]             ->allocate(nbLinHexa);
502   _downArray[VTK_QUADRATIC_HEXAHEDRON]   ->allocate(nbQuadHexa);
503   _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
504   _downArray[VTK_HEXAGONAL_PRISM]        ->allocate(nbHexPrism);
505
506   // --- iteration on vtkUnstructuredGrid cells, only faces
507   //     for each vtk face:
508   //       create a downward face entry with its downward id.
509   //       compute vtk volumes, create downward volumes entry.
510   //       mark face in downward volumes
511   //       mark volumes in downward face
512
513   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
514   int cellSize = this->Types->GetNumberOfTuples();
515   _cellIdToDownId.resize(cellSize, -1);
516
517   for (int i = 0; i < cellSize; i++)
518     {
519       int vtkFaceType = this->GetCellType(i);
520       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
521         {
522           int vtkFaceId = i;
523           //ASSERT(_downArray[vtkFaceType]);
524           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
525           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
526           downFace->setTempNodes(connFaceId, vtkFaceId);
527           int vols[2] = { -1, -1 };
528           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
529           //MESSAGE("nbVolumes="<< nbVolumes);
530           for (int ivol = 0; ivol < nbVolumes; ivol++)
531             {
532               int vtkVolId = vols[ivol];
533               int vtkVolType = this->GetCellType(vtkVolId);
534               //ASSERT(_downArray[vtkVolType]);
535               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
536               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
537               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
538               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
539             }
540         }
541     }
542
543   // --- iteration on vtkUnstructuredGrid cells, only volumes
544   //     for each vtk volume:
545   //       create downward volumes entry if not already done
546   //       build a temporary list of faces described with their nodes
547   //       for each face
548   //         compute the vtk volumes containing this face
549   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
550   //         if not, create a downward face entry (resizing of structure required)
551   //         (the downward faces store a temporary list of nodes to ease the comparison)
552   //         create downward volumes entry if not already done
553   //         mark volumes in downward face
554   //         mark face in downward volumes
555
556   CHRONOSTOP(20);
557   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
558
559   for (int i = 0; i < cellSize; i++)
560     {
561       int vtkType = this->GetCellType(i);
562       if (SMDS_Downward::getCellDimension(vtkType) == 3)
563         {
564           //CHRONO(31);
565           int vtkVolId = i;
566           // MESSAGE("vtk volume " << vtkVolId);
567           //ASSERT(_downArray[vtkType]);
568           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
569
570           // --- find all the faces of the volume, describe the faces by their nodes
571
572           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
573           ListElemByNodesType facesWithNodes;
574           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
575           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
576           //CHRONOSTOP(31);
577           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
578             {
579               // --- find the volumes containing the face
580
581               //CHRONO(32);
582               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
583               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
584               int vols[2] = { -1, -1 };
585               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
586               int lg = facesWithNodes.elems[iface].nbNodes;
587               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
588               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
589
590               // --- check if face is registered in the volumes
591               //CHRONOSTOP(32);
592
593               //CHRONO(33);
594               int connFaceId = -1;
595               for (int ivol = 0; ivol < nbVolumes; ivol++)
596                 {
597                   int vtkVolId2 = vols[ivol];
598                   int vtkVolType = this->GetCellType(vtkVolId2);
599                   //ASSERT(_downArray[vtkVolType]);
600                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
601                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
602                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
603                   if (connFaceId >= 0)
604                     break; // --- face already created
605                 }//CHRONOSTOP(33);
606
607               // --- if face is not registered in the volumes, create face
608
609               //CHRONO(34);
610               if (connFaceId < 0)
611                 {
612                   connFaceId = _downArray[vtkFaceType]->addCell();
613                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
614                 }//CHRONOSTOP(34);
615
616               // --- mark volumes in downward face and mark face in downward volumes
617
618               //CHRONO(35);
619               for (int ivol = 0; ivol < nbVolumes; ivol++)
620                 {
621                   int vtkVolId2 = vols[ivol];
622                   int vtkVolType = this->GetCellType(vtkVolId2);
623                   //ASSERT(_downArray[vtkVolType]);
624                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
625                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
626                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
627                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
628                 }//CHRONOSTOP(35);
629             }
630         }
631     }
632
633   // --- iteration on vtkUnstructuredGrid cells, only edges
634   //     for each vtk edge:
635   //       create downward edge entry
636   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
637   //       find downward faces
638   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
639   //       mark edge in downward faces
640   //       mark faces in downward edge
641
642   CHRONOSTOP(21);
643   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
644
645   for (int i = 0; i < cellSize; i++)
646     {
647       int vtkEdgeType = this->GetCellType(i);
648       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
649         {
650           int vtkEdgeId = i;
651           //ASSERT(_downArray[vtkEdgeType]);
652           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
653           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
654           downEdge->setNodes(connEdgeId, vtkEdgeId);
655           std::vector<int> vtkIds;
656           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
657           int downFaces[1000];
658           unsigned char downTypes[1000];
659           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
660           for (int n = 0; n < nbDownFaces; n++)
661             {
662               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
663               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
664             }
665         }
666     }
667
668   // --- iteration on downward faces (they are all listed now)
669   //     for each downward face:
670   //       build a temporary list of edges with their ordered list of nodes
671   //       for each edge:
672   //         find all the vtk cells containing this edge
673   //         then identify all the downward faces containing the edge, from the vtk cells
674   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
675   //         if not, create a downward edge entry with the node id's
676   //         mark edge in downward faces
677   //         mark downward faces in edge (size of list unknown, to be allocated)
678
679   CHRONOSTOP(22);CHRONO(23);
680
681   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
682     {
683       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
684         continue;
685
686       // --- find all the edges of the face, describe the edges by their nodes
687
688       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
689       int maxId = downFace->getMaxId();
690       for (int faceId = 0; faceId < maxId; faceId++)
691         {
692           //CHRONO(40);
693           ListElemByNodesType edgesWithNodes;
694           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
695           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
696
697           //CHRONOSTOP(40);
698           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
699             {
700
701               // --- check if the edge is already registered by exploration of the faces
702
703               //CHRONO(41);
704               std::vector<int> vtkIds;
705               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
706               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
707               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
708               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
709               //CHRONOSTOP(41);CHRONO(42);
710               int downFaces[1000];
711               unsigned char downTypes[1000];
712               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
713               //CHRONOSTOP(42);
714
715               //CHRONO(43);
716               int connEdgeId = -1;
717               for (int idf = 0; idf < nbDownFaces; idf++)
718                 {
719                   int faceId2 = downFaces[idf];
720                   int faceType = downTypes[idf];
721                   //ASSERT(_downArray[faceType]);
722                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
723                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
724                   if (connEdgeId >= 0)
725                     break; // --- edge already created
726                 }//CHRONOSTOP(43);
727
728               // --- if edge is not registered in the faces, create edge
729
730               if (connEdgeId < 0)
731                 {
732                   //CHRONO(44);
733                   connEdgeId = _downArray[vtkEdgeType]->addCell();
734                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
735                   //CHRONOSTOP(44);
736                 }
737
738               // --- mark faces in downward edge and mark edge in downward faces
739
740               //CHRONO(45);
741               for (int idf = 0; idf < nbDownFaces; idf++)
742                 {
743                   int faceId2 = downFaces[idf];
744                   int faceType = downTypes[idf];
745                   //ASSERT(_downArray[faceType]);
746                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
747                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
748                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
749                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
750                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
751                 }//CHRONOSTOP(45);
752             }
753         }
754     }
755
756   CHRONOSTOP(23);CHRONO(24);
757
758   // compact downward connectivity structure: adjust downward arrays size, replace std::vector<vector int>> by a single std::vector<int>
759   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
760
761   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
762     {
763       if (SMDS_Downward *down = _downArray[vtkType])
764         {
765           down->compactStorage();
766         }
767     }
768
769   // --- Statistics
770
771   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
772     {
773       if (SMDS_Downward *down = _downArray[vtkType])
774         {
775           if (down->getMaxId())
776             {
777               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
778                   << GuessSize[vtkType] << " real: " << down->getMaxId());
779             }
780         }
781     }CHRONOSTOP(24);CHRONOSTOP(2);
782   counters::stats();
783 }
784
785 /*! Get the neighbors of a cell.
786  * Only the neighbors having the dimension of the cell are taken into account
787  * (neighbors of a volume are the volumes sharing a face with this volume,
788  *  neighbors of a face are the faces sharing an edge with this face...).
789  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
790  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
791  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
792  * @param vtkId the vtk id of the cell
793  * @return number of neighbors
794  */
795 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
796 {
797   int vtkType = this->GetCellType(vtkId);
798   int cellDim = SMDS_Downward::getCellDimension(vtkType);
799   if (cellDim <2)
800     return 0; // TODO voisins des edges = edges connectees
801   int cellId = this->_cellIdToDownId[vtkId];
802
803   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
804   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
805   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
806
807   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
808
809   int nb = 0;
810   for (int i = 0; i < nbCells; i++)
811     {
812       int downId = downCells[i];
813       int cellType = downTyp[i];
814       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
815       const int *upCells = _downArray[cellType]->getUpCells(downId);
816       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
817
818       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
819       //    for a face, number of neighbors (connected faces) not known
820
821       for (int j = 0; j < nbUp; j++)
822         {
823           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
824             continue;
825           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
826           neighborsVtkIds[nb] = vtkNeighbor;
827           downIds[nb] = downId;
828           downTypes[nb] = cellType;
829           nb++;
830           if (nb >= NBMAXNEIGHBORS)
831             {
832               INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
833               return nb;
834             }
835         }
836       if (getSkin)
837         {
838           if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
839             {
840               neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
841               downIds[nb] = downId;
842               downTypes[nb] = cellType;
843               nb++;
844               if (nb >= NBMAXNEIGHBORS)
845                 {
846                   INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
847                   return nb;
848                 }
849             }
850         }
851     }
852   return nb;
853 }
854
855 /*! get the volumes containing a face or an edge of the grid
856  * The edge or face belongs to the vtkUnstructuredGrid
857  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
858  * @param vtkId vtk id of the face or edge
859  */
860 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
861 {
862   int vtkType = this->GetCellType(vtkId);
863   int dim = SMDS_Downward::getCellDimension(vtkType);
864   int nbFaces = 0;
865   unsigned char cellTypes[1000];
866   int downCellId[1000];
867   if (dim == 1)
868     {
869       int downId = this->CellIdToDownId(vtkId);
870       if (downId < 0)
871         {
872           MESSAGE("Downward structure not up to date: new edge not taken into account");
873           return 0;
874         }
875       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
876       const int *upCells = _downArray[vtkType]->getUpCells(downId);
877       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
878       for (int i=0; i< nbFaces; i++)
879         {
880           cellTypes[i] = upTypes[i];
881           downCellId[i] = upCells[i];
882         }
883     }
884   else if (dim == 2)
885     {
886       nbFaces = 1;
887       cellTypes[0] = this->GetCellType(vtkId);
888       int downId = this->CellIdToDownId(vtkId);
889       if (downId < 0)
890         {
891           MESSAGE("Downward structure not up to date: new face not taken into account");
892           return 0;
893         }
894       downCellId[0] = downId;
895     }
896
897   int nbvol =0;
898   for (int i=0; i<nbFaces; i++)
899     {
900       int vtkTypeFace = cellTypes[i];
901       int downId = downCellId[i];
902       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
903       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
904       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
905        for (int j=0; j<nv; j++)
906         {
907           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
908           if (vtkVolId >= 0)
909             volVtkIds[nbvol++] = vtkVolId;
910         }
911     }
912   return nbvol;
913 }
914
915 /*! get the volumes containing a face or an edge of the downward structure
916  * The edge or face does not necessary belong to the vtkUnstructuredGrid
917  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
918  * @param downId id in the downward structure
919  * @param downType type of cell
920  */
921 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
922 {
923   int vtkType = downType;
924   int dim = SMDS_Downward::getCellDimension(vtkType);
925   int nbFaces = 0;
926   unsigned char cellTypes[1000];
927   int downCellId[1000];
928   if (dim == 1)
929     {
930       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
931       const int *upCells = _downArray[vtkType]->getUpCells(downId);
932       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
933       for (int i=0; i< nbFaces; i++)
934         {
935           cellTypes[i] = upTypes[i];
936           downCellId[i] = upCells[i];
937         }
938     }
939   else if (dim == 2)
940     {
941       nbFaces = 1;
942       cellTypes[0] = vtkType;
943       downCellId[0] = downId;
944     }
945
946   int nbvol =0;
947   for (int i=0; i<nbFaces; i++)
948     {
949       int vtkTypeFace = cellTypes[i];
950       int downId = downCellId[i];
951       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
952       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
953       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
954        for (int j=0; j<nv; j++)
955         {
956           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
957           if (vtkVolId >= 0)
958             volVtkIds[nbvol++] = vtkVolId;
959         }
960     }
961   return nbvol;
962 }
963
964 /*! get the node id's of a cell.
965  * The cell is defined by it's downward connectivity id and type.
966  * @param nodeSet set of of vtk node id's to fill.
967  * @param downId downward connectivity id of the cell.
968  * @param downType type of cell.
969  */
970 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
971 {
972   _downArray[downType]->getNodeIds(downId, nodeSet);
973 }
974
975 /*! change some nodes in cell without modifying type or internal connectivity.
976  * Nodes inverse connectivity is maintained up to date.
977  * @param vtkVolId vtk id of the cell
978  * @param localClonedNodeIds map old node id to new node id.
979  */
980 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
981 {
982   vtkIdType npts = 0;
983   vtkIdType const *tmp(nullptr); // will refer to the point id's of the face
984   vtkIdType *pts;                // will refer to the point id's of the face
985   this->GetCellPoints(vtkVolId, npts, tmp);
986   pts = const_cast< vtkIdType*>( tmp );
987   for (int i = 0; i < npts; i++)
988     {
989       if (localClonedNodeIds.count(pts[i]))
990         {
991           vtkIdType oldpt = pts[i];
992           pts[i] = localClonedNodeIds[oldpt];
993           //MESSAGE(oldpt << " --> " << pts[i]);
994           //this->RemoveReferenceToCell(oldpt, vtkVolId);
995           //this->AddReferenceToCell(pts[i], vtkVolId);
996         }
997     }
998 }
999
1000 /*! reorder the nodes of a face
1001  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
1002  * @param orderedNodes list of nodes to reorder (in out)
1003  * @return size of the list
1004  */
1005 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
1006 {
1007   int vtkType = this->GetCellType( vtkVolId );
1008   dim = SMDS_Downward::getCellDimension( vtkType );
1009   if (dim == 3)
1010   {
1011     SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
1012     int        downVolId = this->_cellIdToDownId[ vtkVolId ];
1013     downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
1014   }
1015   // else nothing to do;
1016   return orderedNodes.size();
1017 }
1018
1019 void SMDS_UnstructuredGrid::BuildLinks()
1020 {
1021   // Remove the old links if they are already built
1022   if (this->Links)
1023   {
1024     this->Links->UnRegister(this);
1025   }
1026
1027   SMDS_CellLinks* links;
1028   this->Links = links = SMDS_CellLinks::New();
1029   (static_cast< vtkCellLinks *>(this->Links.GetPointer()))->Allocate(this->GetNumberOfPoints());
1030   this->Links->Register(this);
1031   links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
1032   this->Links->Delete();
1033 }
1034
1035 void SMDS_UnstructuredGrid::DeleteLinks()
1036 {
1037   // Remove the old links if they are already built
1038   if (this->Links)
1039   {
1040     this->Links->UnRegister(this);
1041     this->Links = NULL;
1042   }
1043 }
1044 SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
1045 {
1046   if ( !this->Links )
1047     BuildLinks();
1048   return static_cast< SMDS_CellLinks* >( this->Links.GetPointer() );
1049 }
1050
1051 /*! Create a volume (prism or hexahedron) by duplication of a face.
1052  * Designed for use in creation of flat elements separating volume domains.
1053  * A face separating two domains is shared by two volume cells.
1054  * All the nodes are already created (for the two faces).
1055  * Each original Node is associated to corresponding nodes in the domains.
1056  * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1057  * In that case, even some of the nodes to use for the original face may be changed.
1058  * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1059  * @param domain1: domain of the original face
1060  * @param domain2: domain of the duplicated face
1061  * @param originalNodes: the vtk node ids of the original face
1062  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1063  * @return ok if success.
1064  */
1065 SMDS_MeshCell*
1066 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1067                                              int domain1,
1068                                              int domain2,
1069                                              std::set<int>& originalNodes,
1070                                              std::map<int, std::map<int, int> >& nodeDomains,
1071                                              std::map<int, std::map<long, int> >& nodeQuadDomains)
1072 {
1073   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1074   std::vector<vtkIdType> orderedOriginals( originalNodes.begin(), originalNodes.end() );
1075
1076   int dim = 0;
1077   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1078   std::vector<vtkIdType> orderedNodes;
1079
1080   bool isQuadratic = false;
1081   switch (orderedOriginals.size())
1082   {
1083   case 3:
1084     if (dim == 2)
1085       isQuadratic = true;
1086     break;
1087   case 6:
1088   case 8:
1089     isQuadratic = true;
1090     break;
1091   default:
1092     isQuadratic = false;
1093     break;
1094   }
1095
1096   if (isQuadratic)
1097   {
1098     long dom1 = domain1;
1099     long dom2 = domain2;
1100     long dom1_2; // for nodeQuadDomains
1101     if (domain1 < domain2)
1102       dom1_2 = dom1 + INT_MAX * dom2;
1103     else
1104       dom1_2 = dom2 + INT_MAX * dom1;
1105     //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1106     int ima = orderedOriginals.size();
1107     int mid = orderedOriginals.size() / 2;
1108     //cerr << "ima=" << ima << " mid=" << mid << endl;
1109     for (int i = 0; i < mid; i++)
1110       orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1111     for (int i = 0; i < mid; i++)
1112       orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1113     for (int i = mid; i < ima; i++)
1114       orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1115     for (int i = mid; i < ima; i++)
1116       orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1117     for (int i = 0; i < mid; i++)
1118     {
1119       int oldId = orderedOriginals[i];
1120       int newId;
1121       if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1122         newId = nodeQuadDomains[oldId][dom1_2];
1123       else
1124       {
1125         double *coords = this->GetPoint(oldId);
1126         SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1127         newId = newNode->GetVtkID();
1128         if (! nodeQuadDomains.count(oldId))
1129         {
1130           std::map<long, int> emptyMap;
1131           nodeQuadDomains[oldId] = emptyMap;
1132         }
1133         nodeQuadDomains[oldId][dom1_2] = newId;
1134       }
1135       orderedNodes.push_back(newId);
1136     }
1137   }
1138   else
1139   {
1140     for (int i = 0; i < nbNodes; i++)
1141       orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1142     if (dim == 3)
1143       for (int i = 0; i < nbNodes; i++)
1144         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1145     else
1146       for (int i = nbNodes-1; i >= 0; i--)
1147         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1148
1149   }
1150
1151   if (dim == 3)
1152   {
1153     SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1154     return vol;
1155   }
1156   else if (dim == 2)
1157   {
1158     // bos #24368
1159     // orient face by the original one, as getOrderedNodesOfFace() not implemented for faces
1160     const SMDS_MeshElement* origFace = _mesh->FindElementVtk( vtkVolId );
1161     int   i0 = origFace->GetNodeIndex( _mesh->FindNodeVtk( orderedNodes[0] ));
1162     int   i1 = origFace->GetNodeIndex( _mesh->FindNodeVtk( orderedNodes[1] ));
1163     int diff = i0 - i1;
1164     // order of nodes must be reverse in face and origFace
1165     bool oriOk = ( diff == 1 ) || ( diff == -3 );
1166     if ( !oriOk )
1167     {
1168       SMDSAbs_EntityType type = isQuadratic ? SMDSEntity_Quad_Quadrangle : SMDSEntity_Quadrangle;
1169       const std::vector<int>& interlace = SMDS_MeshCell::reverseSmdsOrder( type );
1170       SMDS_MeshCell::applyInterlace( interlace, orderedNodes );
1171     }
1172     SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1173     return face;
1174   }
1175
1176   // TODO update sub-shape list of elements and nodes
1177   return 0;
1178 }
1179
1180 //================================================================================
1181 /*!
1182  * \brief Allocates data array for ball diameters
1183  *  \param MaxVtkID - max ID of a ball element
1184  */
1185 //================================================================================
1186
1187 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1188 {
1189   SetBallDiameter( MaxVtkID, 0 );
1190 }
1191
1192 //================================================================================
1193 /*!
1194  * \brief Sets diameter of a ball element
1195  *  \param vtkID - vtk id of the ball element
1196  *  \param diameter - diameter of the ball element
1197  */
1198 //================================================================================
1199
1200 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1201 {
1202   vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1203   if ( !array )
1204   {
1205     array = vtkDoubleArray::New();
1206     array->SetNumberOfComponents(1);
1207     vtkDataSet::CellData->SetScalars( array );
1208   }
1209   array->InsertValue( vtkID, diameter );
1210 }
1211
1212 //================================================================================
1213 /*!
1214  * \brief Returns diameter of a ball element
1215  *  \param vtkID - vtk id of the ball element
1216  */
1217 //================================================================================
1218
1219 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1220 {
1221   if ( vtkDataSet::CellData )
1222     return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
1223   return 0;
1224 }