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