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