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