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