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