Salome HOME
Porting to ParaView 5.8
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
1 // Copyright (C) 2010-2019  CEA/DEN, EDF R&D, 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 int 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   int 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<int>& idNodesOldToNew, int newNodeSize,
179                                         std::vector<int>& idCellsNewToOld, int 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   int 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( newNodeSize );
195
196     int 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       int startBloc = i;
203       // look for a block end
204       while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
205         ++i;
206       int 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   int  oldCellSize = this->Types->GetNumberOfTuples();
221   bool updateCells = ( updateNodes || newCellSize != oldCellSize );
222   for ( int 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 ((int) idNodesOldToNew.size() < oldNodeSize )
240   {
241     idNodesOldToNew.reserve( oldNodeSize );
242     for ( int i = idNodesOldToNew.size(); i < oldNodeSize; ++i )
243       idNodesOldToNew.push_back( i );
244   }
245
246   // --- create new compacted Connectivity, Locations and Types
247
248   int newConnectivitySize = this->Connectivity->GetNumberOfConnectivityEntries();
249   if ( newCellSize != oldCellSize )
250     for ( int 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(newCellSize);
261
262   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
263   newLocations->Initialize();
264   newLocations->SetNumberOfValues(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 ( int newCellID = 0; newCellID < newCellSize; newCellID++ )
277     {
278       if ( newTypes->GetValue( newCellID ) == VTK_POLY_VERTEX )
279       {
280         int 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 ( int newCellID = 0; newCellID < newCellSize; newCellID++ )
296     {
297       if ( newTypes->GetValue( newCellID ) == VTK_POLYHEDRON )
298       {
299         int oldCellId = idCellsNewToOld[ newCellID ];
300         newFaceLocations->InsertNextValue( newFaces->GetMaxId()+1 );
301         int oldFaceLoc = this->FaceLocations->GetValue( oldCellId );
302         int nCellFaces = this->Faces->GetValue( oldFaceLoc++ );
303         newFaces->InsertNextValue( nCellFaces );
304         for ( int n = 0; n < nCellFaces; n++ )
305         {
306           int nptsInFace = this->Faces->GetValue( oldFaceLoc++ );
307           newFaces->InsertNextValue( nptsInFace );
308           for ( int k = 0; k < nptsInFace; k++ )
309           {
310             int oldpt = this->Faces->GetValue( 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<int>& idNodesOldToNew,
340                                       int&              alreadyCopied,
341                                       int               start,
342                                       int               end)
343 {
344   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
345   void *source = this->Points->GetVoidPointer(3 * start);
346   int 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<int>& idCellsNewToOld,
356                                      const std::vector<int>& 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     int 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       int oldval = oldPtsCell[l];
375       pointsCell[l] = idNodesOldToNew[oldval];
376     }
377     /*int newcnt = */newConnectivity->InsertNextCell( nbpts, pointsCell.data() );
378     int newLoc = newConnectivity->GetInsertLocation( nbpts );
379     newLocations->SetValue( iNew, newLoc );
380   }
381 }
382
383 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
384 {
385   if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
386   {
387     return -1;
388   }
389   return _cellIdToDownId[vtkCellId];
390 }
391
392 void SMDS_UnstructuredGrid::setCellIdToDownId(int 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  = meshInfo.NbTetras  (ORDER_LINEAR);
450   int nbQuadTetra = meshInfo.NbTetras  (ORDER_QUADRATIC);
451   int nbLinPyra   = meshInfo.NbPyramids(ORDER_LINEAR);
452   int nbQuadPyra  = meshInfo.NbPyramids(ORDER_QUADRATIC);
453   int nbLinPrism  = meshInfo.NbPrisms  (ORDER_LINEAR);
454   int nbQuadPrism = meshInfo.NbPrisms  (ORDER_QUADRATIC);
455   int nbLinHexa   = meshInfo.NbHexas   (ORDER_LINEAR);
456   int nbQuadHexa  = meshInfo.NbHexas   (ORDER_QUADRATIC);
457   int nbHexPrism  = 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
486   _downArray[VTK_LINE]                   ->allocate(nbLineGuess);
487   _downArray[VTK_QUADRATIC_EDGE]         ->allocate(nbQuadEdgeGuess);
488   _downArray[VTK_TRIANGLE]               ->allocate(nbLinTriaGuess);
489   _downArray[VTK_QUADRATIC_TRIANGLE]     ->allocate(nbQuadTriaGuess);
490   _downArray[VTK_BIQUADRATIC_TRIANGLE]   ->allocate(nbQuadTriaGuess);
491   _downArray[VTK_QUAD]                   ->allocate(nbLinQuadGuess);
492   _downArray[VTK_QUADRATIC_QUAD]         ->allocate(nbQuadQuadGuess);
493   _downArray[VTK_BIQUADRATIC_QUAD]       ->allocate(nbQuadQuadGuess);
494   _downArray[VTK_TETRA]                  ->allocate(nbLinTetra);
495   _downArray[VTK_QUADRATIC_TETRA]        ->allocate(nbQuadTetra);
496   _downArray[VTK_PYRAMID]                ->allocate(nbLinPyra);
497   _downArray[VTK_QUADRATIC_PYRAMID]      ->allocate(nbQuadPyra);
498   _downArray[VTK_WEDGE]                  ->allocate(nbLinPrism);
499   _downArray[VTK_QUADRATIC_WEDGE]        ->allocate(nbQuadPrism);
500   _downArray[VTK_HEXAHEDRON]             ->allocate(nbLinHexa);
501   _downArray[VTK_QUADRATIC_HEXAHEDRON]   ->allocate(nbQuadHexa);
502   _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
503   _downArray[VTK_HEXAGONAL_PRISM]        ->allocate(nbHexPrism);
504
505   // --- iteration on vtkUnstructuredGrid cells, only faces
506   //     for each vtk face:
507   //       create a downward face entry with its downward id.
508   //       compute vtk volumes, create downward volumes entry.
509   //       mark face in downward volumes
510   //       mark volumes in downward face
511
512   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
513   int cellSize = this->Types->GetNumberOfTuples();
514   _cellIdToDownId.resize(cellSize, -1);
515
516   for (int i = 0; i < cellSize; i++)
517     {
518       int vtkFaceType = this->GetCellType(i);
519       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
520         {
521           int vtkFaceId = i;
522           //ASSERT(_downArray[vtkFaceType]);
523           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
524           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
525           downFace->setTempNodes(connFaceId, vtkFaceId);
526           int vols[2] = { -1, -1 };
527           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
528           //MESSAGE("nbVolumes="<< nbVolumes);
529           for (int ivol = 0; ivol < nbVolumes; ivol++)
530             {
531               int vtkVolId = vols[ivol];
532               int vtkVolType = this->GetCellType(vtkVolId);
533               //ASSERT(_downArray[vtkVolType]);
534               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
535               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
536               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
537               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
538             }
539         }
540     }
541
542   // --- iteration on vtkUnstructuredGrid cells, only volumes
543   //     for each vtk volume:
544   //       create downward volumes entry if not already done
545   //       build a temporary list of faces described with their nodes
546   //       for each face
547   //         compute the vtk volumes containing this face
548   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
549   //         if not, create a downward face entry (resizing of structure required)
550   //         (the downward faces store a temporary list of nodes to ease the comparison)
551   //         create downward volumes entry if not already done
552   //         mark volumes in downward face
553   //         mark face in downward volumes
554
555   CHRONOSTOP(20);
556   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
557
558   for (int i = 0; i < cellSize; i++)
559     {
560       int vtkType = this->GetCellType(i);
561       if (SMDS_Downward::getCellDimension(vtkType) == 3)
562         {
563           //CHRONO(31);
564           int vtkVolId = i;
565           // MESSAGE("vtk volume " << vtkVolId);
566           //ASSERT(_downArray[vtkType]);
567           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
568
569           // --- find all the faces of the volume, describe the faces by their nodes
570
571           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
572           ListElemByNodesType facesWithNodes;
573           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
574           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
575           //CHRONOSTOP(31);
576           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
577             {
578               // --- find the volumes containing the face
579
580               //CHRONO(32);
581               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
582               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
583               int vols[2] = { -1, -1 };
584               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
585               int lg = facesWithNodes.elems[iface].nbNodes;
586               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
587               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
588
589               // --- check if face is registered in the volumes
590               //CHRONOSTOP(32);
591
592               //CHRONO(33);
593               int connFaceId = -1;
594               for (int ivol = 0; ivol < nbVolumes; ivol++)
595                 {
596                   int vtkVolId2 = vols[ivol];
597                   int vtkVolType = this->GetCellType(vtkVolId2);
598                   //ASSERT(_downArray[vtkVolType]);
599                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
600                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
601                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
602                   if (connFaceId >= 0)
603                     break; // --- face already created
604                 }//CHRONOSTOP(33);
605
606               // --- if face is not registered in the volumes, create face
607
608               //CHRONO(34);
609               if (connFaceId < 0)
610                 {
611                   connFaceId = _downArray[vtkFaceType]->addCell();
612                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
613                 }//CHRONOSTOP(34);
614
615               // --- mark volumes in downward face and mark face in downward volumes
616
617               //CHRONO(35);
618               for (int ivol = 0; ivol < nbVolumes; ivol++)
619                 {
620                   int vtkVolId2 = vols[ivol];
621                   int vtkVolType = this->GetCellType(vtkVolId2);
622                   //ASSERT(_downArray[vtkVolType]);
623                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
624                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
625                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
626                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
627                 }//CHRONOSTOP(35);
628             }
629         }
630     }
631
632   // --- iteration on vtkUnstructuredGrid cells, only edges
633   //     for each vtk edge:
634   //       create downward edge entry
635   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
636   //       find downward faces
637   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
638   //       mark edge in downward faces
639   //       mark faces in downward edge
640
641   CHRONOSTOP(21);
642   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
643
644   for (int i = 0; i < cellSize; i++)
645     {
646       int vtkEdgeType = this->GetCellType(i);
647       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
648         {
649           int vtkEdgeId = i;
650           //ASSERT(_downArray[vtkEdgeType]);
651           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
652           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
653           downEdge->setNodes(connEdgeId, vtkEdgeId);
654           std::vector<int> vtkIds;
655           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
656           int downFaces[1000];
657           unsigned char downTypes[1000];
658           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
659           for (int n = 0; n < nbDownFaces; n++)
660             {
661               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
662               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
663             }
664         }
665     }
666
667   // --- iteration on downward faces (they are all listed now)
668   //     for each downward face:
669   //       build a temporary list of edges with their ordered list of nodes
670   //       for each edge:
671   //         find all the vtk cells containing this edge
672   //         then identify all the downward faces containing the edge, from the vtk cells
673   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
674   //         if not, create a downward edge entry with the node id's
675   //         mark edge in downward faces
676   //         mark downward faces in edge (size of list unknown, to be allocated)
677
678   CHRONOSTOP(22);CHRONO(23);
679
680   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
681     {
682       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
683         continue;
684
685       // --- find all the edges of the face, describe the edges by their nodes
686
687       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
688       int maxId = downFace->getMaxId();
689       for (int faceId = 0; faceId < maxId; faceId++)
690         {
691           //CHRONO(40);
692           ListElemByNodesType edgesWithNodes;
693           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
694           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
695
696           //CHRONOSTOP(40);
697           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
698             {
699
700               // --- check if the edge is already registered by exploration of the faces
701
702               //CHRONO(41);
703               std::vector<int> vtkIds;
704               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
705               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
706               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
707               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
708               //CHRONOSTOP(41);CHRONO(42);
709               int downFaces[1000];
710               unsigned char downTypes[1000];
711               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
712               //CHRONOSTOP(42);
713
714               //CHRONO(43);
715               int connEdgeId = -1;
716               for (int idf = 0; idf < nbDownFaces; idf++)
717                 {
718                   int faceId2 = downFaces[idf];
719                   int faceType = downTypes[idf];
720                   //ASSERT(_downArray[faceType]);
721                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
722                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
723                   if (connEdgeId >= 0)
724                     break; // --- edge already created
725                 }//CHRONOSTOP(43);
726
727               // --- if edge is not registered in the faces, create edge
728
729               if (connEdgeId < 0)
730                 {
731                   //CHRONO(44);
732                   connEdgeId = _downArray[vtkEdgeType]->addCell();
733                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
734                   //CHRONOSTOP(44);
735                 }
736
737               // --- mark faces in downward edge and mark edge in downward faces
738
739               //CHRONO(45);
740               for (int idf = 0; idf < nbDownFaces; idf++)
741                 {
742                   int faceId2 = downFaces[idf];
743                   int faceType = downTypes[idf];
744                   //ASSERT(_downArray[faceType]);
745                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
746                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
747                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
748                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
749                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
750                 }//CHRONOSTOP(45);
751             }
752         }
753     }
754
755   CHRONOSTOP(23);CHRONO(24);
756
757   // compact downward connectivity structure: adjust downward arrays size, replace std::vector<vector int>> by a single std::vector<int>
758   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
759
760   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
761     {
762       if (SMDS_Downward *down = _downArray[vtkType])
763         {
764           down->compactStorage();
765         }
766     }
767
768   // --- Statistics
769
770   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
771     {
772       if (SMDS_Downward *down = _downArray[vtkType])
773         {
774           if (down->getMaxId())
775             {
776               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
777                   << GuessSize[vtkType] << " real: " << down->getMaxId());
778             }
779         }
780     }CHRONOSTOP(24);CHRONOSTOP(2);
781   counters::stats();
782 }
783
784 /*! Get the neighbors of a cell.
785  * Only the neighbors having the dimension of the cell are taken into account
786  * (neighbors of a volume are the volumes sharing a face with this volume,
787  *  neighbors of a face are the faces sharing an edge with this face...).
788  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
789  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
790  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
791  * @param vtkId the vtk id of the cell
792  * @return number of neighbors
793  */
794 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
795 {
796   int vtkType = this->GetCellType(vtkId);
797   int cellDim = SMDS_Downward::getCellDimension(vtkType);
798   if (cellDim <2)
799     return 0; // TODO voisins des edges = edges connectees
800   int cellId = this->_cellIdToDownId[vtkId];
801
802   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
803   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
804   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
805
806   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
807
808   int nb = 0;
809   for (int i = 0; i < nbCells; i++)
810     {
811       int downId = downCells[i];
812       int cellType = downTyp[i];
813       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
814       const int *upCells = _downArray[cellType]->getUpCells(downId);
815       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
816
817       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
818       //    for a face, number of neighbors (connected faces) not known
819
820       for (int j = 0; j < nbUp; j++)
821         {
822           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
823             continue;
824           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
825           neighborsVtkIds[nb] = vtkNeighbor;
826           downIds[nb] = downId;
827           downTypes[nb] = cellType;
828           nb++;
829           if (nb >= NBMAXNEIGHBORS)
830             {
831               INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
832               return nb;
833             }
834         }
835       if (getSkin)
836         {
837           if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
838             {
839               neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
840               downIds[nb] = downId;
841               downTypes[nb] = cellType;
842               nb++;
843               if (nb >= NBMAXNEIGHBORS)
844                 {
845                   INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
846                   return nb;
847                 }
848             }
849         }
850     }
851   return nb;
852 }
853
854 /*! get the volumes containing a face or an edge of the grid
855  * The edge or face belongs to the vtkUnstructuredGrid
856  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
857  * @param vtkId vtk id of the face or edge
858  */
859 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
860 {
861   int vtkType = this->GetCellType(vtkId);
862   int dim = SMDS_Downward::getCellDimension(vtkType);
863   int nbFaces = 0;
864   unsigned char cellTypes[1000];
865   int downCellId[1000];
866   if (dim == 1)
867     {
868       int downId = this->CellIdToDownId(vtkId);
869       if (downId < 0)
870         {
871           MESSAGE("Downward structure not up to date: new edge not taken into account");
872           return 0;
873         }
874       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
875       const int *upCells = _downArray[vtkType]->getUpCells(downId);
876       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
877       for (int i=0; i< nbFaces; i++)
878         {
879           cellTypes[i] = upTypes[i];
880           downCellId[i] = upCells[i];
881         }
882     }
883   else if (dim == 2)
884     {
885       nbFaces = 1;
886       cellTypes[0] = this->GetCellType(vtkId);
887       int downId = this->CellIdToDownId(vtkId);
888       if (downId < 0)
889         {
890           MESSAGE("Downward structure not up to date: new face not taken into account");
891           return 0;
892         }
893       downCellId[0] = downId;
894     }
895
896   int nbvol =0;
897   for (int i=0; i<nbFaces; i++)
898     {
899       int vtkTypeFace = cellTypes[i];
900       int downId = downCellId[i];
901       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
902       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
903       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
904        for (int j=0; j<nv; j++)
905         {
906           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
907           if (vtkVolId >= 0)
908             volVtkIds[nbvol++] = vtkVolId;
909         }
910     }
911   return nbvol;
912 }
913
914 /*! get the volumes containing a face or an edge of the downward structure
915  * The edge or face does not necessary belong to the vtkUnstructuredGrid
916  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
917  * @param downId id in the downward structure
918  * @param downType type of cell
919  */
920 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
921 {
922   int vtkType = downType;
923   int dim = SMDS_Downward::getCellDimension(vtkType);
924   int nbFaces = 0;
925   unsigned char cellTypes[1000];
926   int downCellId[1000];
927   if (dim == 1)
928     {
929       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
930       const int *upCells = _downArray[vtkType]->getUpCells(downId);
931       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
932       for (int i=0; i< nbFaces; i++)
933         {
934           cellTypes[i] = upTypes[i];
935           downCellId[i] = upCells[i];
936         }
937     }
938   else if (dim == 2)
939     {
940       nbFaces = 1;
941       cellTypes[0] = vtkType;
942       downCellId[0] = downId;
943     }
944
945   int nbvol =0;
946   for (int i=0; i<nbFaces; i++)
947     {
948       int vtkTypeFace = cellTypes[i];
949       int downId = downCellId[i];
950       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
951       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
952       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
953        for (int j=0; j<nv; j++)
954         {
955           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
956           if (vtkVolId >= 0)
957             volVtkIds[nbvol++] = vtkVolId;
958         }
959     }
960   return nbvol;
961 }
962
963 /*! get the node id's of a cell.
964  * The cell is defined by it's downward connectivity id and type.
965  * @param nodeSet set of of vtk node id's to fill.
966  * @param downId downward connectivity id of the cell.
967  * @param downType type of cell.
968  */
969 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
970 {
971   _downArray[downType]->getNodeIds(downId, nodeSet);
972 }
973
974 /*! change some nodes in cell without modifying type or internal connectivity.
975  * Nodes inverse connectivity is maintained up to date.
976  * @param vtkVolId vtk id of the cell
977  * @param localClonedNodeIds map old node id to new node id.
978  */
979 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
980 {
981   vtkIdType npts = 0;
982   vtkIdType const *tmp(nullptr); // will refer to the point id's of the face
983   vtkIdType *pts;                // will refer to the point id's of the face
984   this->GetCellPoints(vtkVolId, npts, tmp);
985   pts = const_cast< vtkIdType*>( tmp );
986   for (int i = 0; i < npts; i++)
987     {
988       if (localClonedNodeIds.count(pts[i]))
989         {
990           vtkIdType oldpt = pts[i];
991           pts[i] = localClonedNodeIds[oldpt];
992           //MESSAGE(oldpt << " --> " << pts[i]);
993           //this->RemoveReferenceToCell(oldpt, vtkVolId);
994           //this->AddReferenceToCell(pts[i], vtkVolId);
995         }
996     }
997 }
998
999 /*! reorder the nodes of a face
1000  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
1001  * @param orderedNodes list of nodes to reorder (in out)
1002  * @return size of the list
1003  */
1004 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
1005 {
1006   int vtkType = this->GetCellType(vtkVolId);
1007   dim = SMDS_Downward::getCellDimension(vtkType);
1008   if (dim == 3)
1009     {
1010       SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
1011       int downVolId = this->_cellIdToDownId[vtkVolId];
1012       downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
1013     }
1014   // else nothing to do;
1015   return orderedNodes.size();
1016 }
1017
1018 void SMDS_UnstructuredGrid::BuildLinks()
1019 {
1020   // Remove the old links if they are already built
1021   if (this->Links)
1022   {
1023     this->Links->UnRegister(this);
1024   }
1025
1026   SMDS_CellLinks* links;
1027   this->Links = links = SMDS_CellLinks::New();
1028   (static_cast< vtkCellLinks *>(this->Links.GetPointer()))->Allocate(this->GetNumberOfPoints());
1029   this->Links->Register(this);
1030   links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
1031   this->Links->Delete();
1032 }
1033
1034 void SMDS_UnstructuredGrid::DeleteLinks()
1035 {
1036   // Remove the old links if they are already built
1037   if (this->Links)
1038   {
1039     this->Links->UnRegister(this);
1040     this->Links = NULL;
1041   }
1042 }
1043 SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
1044 {
1045   if ( !this->Links )
1046     BuildLinks();
1047   return static_cast< SMDS_CellLinks* >( this->Links.GetPointer() );
1048 }
1049
1050 /*! Create a volume (prism or hexahedron) by duplication of a face.
1051  * Designed for use in creation of flat elements separating volume domains.
1052  * A face separating two domains is shared by two volume cells.
1053  * All the nodes are already created (for the two faces).
1054  * Each original Node is associated to corresponding nodes in the domains.
1055  * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1056  * In that case, even some of the nodes to use for the original face may be changed.
1057  * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1058  * @param domain1: domain of the original face
1059  * @param domain2: domain of the duplicated face
1060  * @param originalNodes: the vtk node ids of the original face
1061  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1062  * @return ok if success.
1063  */
1064 SMDS_MeshCell*
1065 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1066                                              int domain1,
1067                                              int domain2,
1068                                              std::set<int>& originalNodes,
1069                                              std::map<int, std::map<int, int> >& nodeDomains,
1070                                              std::map<int, std::map<long, int> >& nodeQuadDomains)
1071 {
1072   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1073   std::vector<vtkIdType> orderedOriginals( originalNodes.begin(), originalNodes.end() );
1074
1075   int dim = 0;
1076   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1077   std::vector<vtkIdType> orderedNodes;
1078
1079   bool isQuadratic = false;
1080   switch (orderedOriginals.size())
1081   {
1082   case 3:
1083     if (dim == 2)
1084       isQuadratic = true;
1085     break;
1086   case 6:
1087   case 8:
1088     isQuadratic = true;
1089     break;
1090   default:
1091     isQuadratic = false;
1092       break;
1093   }
1094
1095   if (isQuadratic)
1096     {
1097       long dom1 = domain1;
1098       long dom2 = domain2;
1099       long dom1_2; // for nodeQuadDomains
1100       if (domain1 < domain2)
1101         dom1_2 = dom1 + INT_MAX * dom2;
1102       else
1103         dom1_2 = dom2 + INT_MAX * dom1;
1104       //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1105       int ima = orderedOriginals.size();
1106       int mid = orderedOriginals.size() / 2;
1107       //cerr << "ima=" << ima << " mid=" << mid << endl;
1108       for (int i = 0; i < mid; i++)
1109         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1110       for (int i = 0; i < mid; i++)
1111         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1112       for (int i = mid; i < ima; i++)
1113         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1114       for (int i = mid; i < ima; i++)
1115         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1116       for (int i = 0; i < mid; i++)
1117         {
1118           int oldId = orderedOriginals[i];
1119           int newId;
1120           if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1121             newId = nodeQuadDomains[oldId][dom1_2];
1122           else
1123             {
1124               double *coords = this->GetPoint(oldId);
1125               SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1126               newId = newNode->GetVtkID();
1127               if (! nodeQuadDomains.count(oldId))
1128                 {
1129                   std::map<long, int> emptyMap;
1130                   nodeQuadDomains[oldId] = emptyMap;
1131                 }
1132               nodeQuadDomains[oldId][dom1_2] = newId;
1133             }
1134           orderedNodes.push_back(newId);
1135         }
1136     }
1137   else
1138     {
1139       for (int i = 0; i < nbNodes; i++)
1140         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1141       if (dim == 3)
1142         for (int i = 0; i < nbNodes; i++)
1143           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1144       else
1145         for (int i = nbNodes-1; i >= 0; i--)
1146           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1147
1148     }
1149
1150   if (dim == 3)
1151     {
1152       SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1153       return vol;
1154     }
1155   else if (dim == 2)
1156     {
1157       SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1158       return face;
1159     }
1160
1161   // TODO update sub-shape list of elements and nodes
1162   return 0;
1163 }
1164
1165 //================================================================================
1166 /*!
1167  * \brief Allocates data array for ball diameters
1168  *  \param MaxVtkID - max ID of a ball element
1169  */
1170 //================================================================================
1171
1172 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1173 {
1174   SetBallDiameter( MaxVtkID, 0 );
1175 }
1176
1177 //================================================================================
1178 /*!
1179  * \brief Sets diameter of a ball element
1180  *  \param vtkID - vtk id of the ball element
1181  *  \param diameter - diameter of the ball element
1182  */
1183 //================================================================================
1184
1185 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1186 {
1187   vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1188   if ( !array )
1189   {
1190     array = vtkDoubleArray::New();
1191     array->SetNumberOfComponents(1);
1192     vtkDataSet::CellData->SetScalars( array );
1193   }
1194   array->InsertValue( vtkID, diameter );
1195 }
1196
1197 //================================================================================
1198 /*!
1199  * \brief Returns diameter of a ball element
1200  *  \param vtkID - vtk id of the ball element
1201  */
1202 //================================================================================
1203
1204 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1205 {
1206   if ( vtkDataSet::CellData )
1207     return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
1208   return 0;
1209 }