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