Salome HOME
update of messages
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
1 // Copyright (C) 2010-2014  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "SMDS_UnstructuredGrid.hxx"
21 #include "SMDS_Mesh.hxx"
22 #include "SMDS_MeshInfo.hxx"
23 #include "SMDS_Downward.hxx"
24 #include "SMDS_MeshVolume.hxx"
25
26 #include "utilities.h"
27
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 // OUV_PORTING_VTK6: seems to be useless
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       this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
289       newFaceLocations->Delete();
290       newFaces->Delete();
291     }
292   else
293   {
294     this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
295   }
296
297   newPoints->Delete();
298   newTypes->Delete();
299   newLocations->Delete();
300   newConnectivity->Delete();
301   this->BuildLinks();
302 }
303
304 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
305                                       int start, int end)
306 {
307   MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
308   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
309   void *source = this->Points->GetVoidPointer(3 * start);
310   int nbPoints = end - start;
311   if (nbPoints > 0)
312     {
313       memcpy(target, source, 3 * sizeof(double) * nbPoints);
314       for (int j = start; j < end; j++)
315         idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
316     }
317 }
318
319 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
320                                      std::vector<int>&     idCellsOldToNew,
321                                      std::vector<int>&     idNodesOldToNew,
322                                      vtkCellArray*         newConnectivity,
323                                      vtkIdTypeArray*       newLocations,
324                                      vtkIdType*            pointsCell,
325                                      int&                  alreadyCopied,
326                                      int                   start,
327                                      int                   end)
328 {
329   //MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
330   for (int j = start; j < end; j++)
331     {
332       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
333       idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
334       vtkIdType oldLoc = this->Locations->GetValue(j);
335       vtkIdType nbpts;
336       vtkIdType *oldPtsCell = 0;
337       this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
338       assert(nbpts < NBMAXNODESINCELL);
339       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
340       for (int l = 0; l < nbpts; l++)
341         {
342           int oldval = oldPtsCell[l];
343           pointsCell[l] = idNodesOldToNew[oldval];
344           //MESSAGE("   " << oldval << " " << pointsCell[l]);
345         }
346       /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
347       int newLoc = newConnectivity->GetInsertLocation(nbpts);
348       //MESSAGE(newcnt << " " << newLoc);
349       newLocations->SetValue(alreadyCopied, newLoc);
350       alreadyCopied++;
351     }
352 }
353
354 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
355 {
356   if((vtkCellId < 0) || (vtkCellId >= _cellIdToDownId.size()))
357     {
358       //MESSAGE("SMDS_UnstructuredGrid::CellIdToDownId structure not up to date: vtkCellId="
359       //    << vtkCellId << " max="<< _cellIdToDownId.size());
360       return -1;
361     }
362   return _cellIdToDownId[vtkCellId];
363 }
364
365 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
366 {
367   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
368   _cellIdToDownId[vtkCellId] = downId;
369 }
370
371 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
372 {
373   for (int i = 0; i < _downArray.size(); i++)
374     {
375       if (_downArray[i])
376         delete _downArray[i];
377       _downArray[i] = 0;
378     }
379   _cellIdToDownId.clear();
380 }
381
382 /*! Build downward connectivity: to do only when needed because heavy memory load.
383  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
384  *
385  */
386 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
387 {
388   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
389   // TODO calcul partiel sans edges
390
391   // --- erase previous data if any
392
393   this->CleanDownwardConnectivity();
394
395   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
396
397   _downArray.resize(VTK_MAXTYPE + 1, 0);
398
399   _downArray[VTK_LINE]                    = new SMDS_DownEdge(this);
400   _downArray[VTK_QUADRATIC_EDGE]          = new SMDS_DownQuadEdge(this);
401   _downArray[VTK_TRIANGLE]                = new SMDS_DownTriangle(this);
402   _downArray[VTK_QUADRATIC_TRIANGLE]      = new SMDS_DownQuadTriangle(this);
403   _downArray[VTK_BIQUADRATIC_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_BIQUADRATIC_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_BIQUADRATIC_TRIANGLE]   ->allocate(nbQuadTriaGuess);
464   _downArray[VTK_QUAD]                   ->allocate(nbLinQuadGuess);
465   _downArray[VTK_QUADRATIC_QUAD]         ->allocate(nbQuadQuadGuess);
466   _downArray[VTK_BIQUADRATIC_QUAD]       ->allocate(nbQuadQuadGuess);
467   _downArray[VTK_TETRA]                  ->allocate(nbLinTetra);
468   _downArray[VTK_QUADRATIC_TETRA]        ->allocate(nbQuadTetra);
469   _downArray[VTK_PYRAMID]                ->allocate(nbLinPyra);
470   _downArray[VTK_QUADRATIC_PYRAMID]      ->allocate(nbQuadPyra);
471   _downArray[VTK_WEDGE]                  ->allocate(nbLinPrism);
472   _downArray[VTK_QUADRATIC_WEDGE]        ->allocate(nbQuadPrism);
473   _downArray[VTK_HEXAHEDRON]             ->allocate(nbLinHexa);
474   _downArray[VTK_QUADRATIC_HEXAHEDRON]   ->allocate(nbQuadHexa);
475   _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
476   _downArray[VTK_HEXAGONAL_PRISM]        ->allocate(nbHexPrism);
477
478   // --- iteration on vtkUnstructuredGrid cells, only faces
479   //     for each vtk face:
480   //       create a downward face entry with its downward id.
481   //       compute vtk volumes, create downward volumes entry.
482   //       mark face in downward volumes
483   //       mark volumes in downward face
484
485   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
486   int cellSize = this->Types->GetNumberOfTuples();
487   _cellIdToDownId.resize(cellSize, -1);
488
489   for (int i = 0; i < cellSize; i++)
490     {
491       int vtkFaceType = this->GetCellType(i);
492       if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
493         {
494           int vtkFaceId = i;
495           //ASSERT(_downArray[vtkFaceType]);
496           int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
497           SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
498           downFace->setTempNodes(connFaceId, vtkFaceId);
499           int vols[2] = { -1, -1 };
500           int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
501           //MESSAGE("nbVolumes="<< nbVolumes);
502           for (int ivol = 0; ivol < nbVolumes; ivol++)
503             {
504               int vtkVolId = vols[ivol];
505               int vtkVolType = this->GetCellType(vtkVolId);
506               //ASSERT(_downArray[vtkVolType]);
507               int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
508               _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
509               _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
510               // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
511             }
512         }
513     }
514
515   // --- iteration on vtkUnstructuredGrid cells, only volumes
516   //     for each vtk volume:
517   //       create downward volumes entry if not already done
518   //       build a temporary list of faces described with their nodes
519   //       for each face
520   //         compute the vtk volumes containing this face
521   //         check if the face is already listed in the volumes (comparison of ordered list of nodes)
522   //         if not, create a downward face entry (resizing of structure required)
523   //         (the downward faces store a temporary list of nodes to ease the comparison)
524   //         create downward volumes entry if not already done
525   //         mark volumes in downward face
526   //         mark face in downward volumes
527
528   CHRONOSTOP(20);
529   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
530
531   for (int i = 0; i < cellSize; i++)
532     {
533       int vtkType = this->GetCellType(i);
534       if (SMDS_Downward::getCellDimension(vtkType) == 3)
535         {
536           //CHRONO(31);
537           int vtkVolId = i;
538           // MESSAGE("vtk volume " << vtkVolId);
539           //ASSERT(_downArray[vtkType]);
540           /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
541
542           // --- find all the faces of the volume, describe the faces by their nodes
543
544           SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
545           ListElemByNodesType facesWithNodes;
546           downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
547           // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
548           //CHRONOSTOP(31);
549           for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
550             {
551               // --- find the volumes containing the face
552
553               //CHRONO(32);
554               int vtkFaceType = facesWithNodes.elems[iface].vtkType;
555               SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
556               int vols[2] = { -1, -1 };
557               int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
558               int lg = facesWithNodes.elems[iface].nbNodes;
559               int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
560               // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
561
562               // --- check if face is registered in the volumes
563               //CHRONOSTOP(32);
564
565               //CHRONO(33);
566               int connFaceId = -1;
567               for (int ivol = 0; ivol < nbVolumes; ivol++)
568                 {
569                   int vtkVolId2 = vols[ivol];
570                   int vtkVolType = this->GetCellType(vtkVolId2);
571                   //ASSERT(_downArray[vtkVolType]);
572                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
573                   SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
574                   connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
575                   if (connFaceId >= 0)
576                     break; // --- face already created
577                 }//CHRONOSTOP(33);
578
579               // --- if face is not registered in the volumes, create face
580
581               //CHRONO(34);
582               if (connFaceId < 0)
583                 {
584                   connFaceId = _downArray[vtkFaceType]->addCell();
585                   downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
586                 }//CHRONOSTOP(34);
587
588               // --- mark volumes in downward face and mark face in downward volumes
589
590               //CHRONO(35);
591               for (int ivol = 0; ivol < nbVolumes; ivol++)
592                 {
593                   int vtkVolId2 = vols[ivol];
594                   int vtkVolType = this->GetCellType(vtkVolId2);
595                   //ASSERT(_downArray[vtkVolType]);
596                   int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
597                   _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
598                   _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
599                   // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
600                 }//CHRONOSTOP(35);
601             }
602         }
603     }
604
605   // --- iteration on vtkUnstructuredGrid cells, only edges
606   //     for each vtk edge:
607   //       create downward edge entry
608   //       store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
609   //       find downward faces
610   //       (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
611   //       mark edge in downward faces
612   //       mark faces in downward edge
613
614   CHRONOSTOP(21);
615   MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
616
617   for (int i = 0; i < cellSize; i++)
618     {
619       int vtkEdgeType = this->GetCellType(i);
620       if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
621         {
622           int vtkEdgeId = i;
623           //ASSERT(_downArray[vtkEdgeType]);
624           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
625           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
626           downEdge->setNodes(connEdgeId, vtkEdgeId);
627           vector<int> vtkIds;
628           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
629           int downFaces[1000];
630           unsigned char downTypes[1000];
631           int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
632           for (int n = 0; n < nbDownFaces; n++)
633             {
634               _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
635               _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
636             }
637         }
638     }
639
640   // --- iteration on downward faces (they are all listed now)
641   //     for each downward face:
642   //       build a temporary list of edges with their ordered list of nodes
643   //       for each edge:
644   //         find all the vtk cells containing this edge
645   //         then identify all the downward faces containing the edge, from the vtk cells
646   //         check if the edge is already listed in the faces (comparison of ordered list of nodes)
647   //         if not, create a downward edge entry with the node id's
648   //         mark edge in downward faces
649   //         mark downward faces in edge (size of list unknown, to be allocated)
650
651   CHRONOSTOP(22);CHRONO(23);
652
653   for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
654     {
655       if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
656         continue;
657
658       // --- find all the edges of the face, describe the edges by their nodes
659
660       SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
661       int maxId = downFace->getMaxId();
662       for (int faceId = 0; faceId < maxId; faceId++)
663         {
664           //CHRONO(40);
665           ListElemByNodesType edgesWithNodes;
666           downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
667           // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
668
669           //CHRONOSTOP(40);
670           for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
671             {
672
673               // --- check if the edge is already registered by exploration of the faces
674
675               //CHRONO(41);
676               vector<int> vtkIds;
677               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
678               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
679               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
680               int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
681               //CHRONOSTOP(41);CHRONO(42);
682               int downFaces[1000];
683               unsigned char downTypes[1000];
684               int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
685               //CHRONOSTOP(42);
686
687               //CHRONO(43);
688               int connEdgeId = -1;
689               for (int idf = 0; idf < nbDownFaces; idf++)
690                 {
691                   int faceId2 = downFaces[idf];
692                   int faceType = downTypes[idf];
693                   //ASSERT(_downArray[faceType]);
694                   SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
695                   connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
696                   if (connEdgeId >= 0)
697                     break; // --- edge already created
698                 }//CHRONOSTOP(43);
699
700               // --- if edge is not registered in the faces, create edge
701
702               if (connEdgeId < 0)
703                 {
704                   //CHRONO(44);
705                   connEdgeId = _downArray[vtkEdgeType]->addCell();
706                   downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
707                   //CHRONOSTOP(44);
708                 }
709
710               // --- mark faces in downward edge and mark edge in downward faces
711
712               //CHRONO(45);
713               for (int idf = 0; idf < nbDownFaces; idf++)
714                 {
715                   int faceId2 = downFaces[idf];
716                   int faceType = downTypes[idf];
717                   //ASSERT(_downArray[faceType]);
718                   //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
719                   _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
720                   _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
721                   // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
722                   //  " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
723                 }//CHRONOSTOP(45);
724             }
725         }
726     }
727
728   CHRONOSTOP(23);CHRONO(24);
729
730   // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
731   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
732
733   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
734     {
735       if (SMDS_Downward *down = _downArray[vtkType])
736         {
737           down->compactStorage();
738         }
739     }
740
741   // --- Statistics
742
743   for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
744     {
745       if (SMDS_Downward *down = _downArray[vtkType])
746         {
747           if (down->getMaxId())
748             {
749               MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
750                   << GuessSize[vtkType] << " real: " << down->getMaxId());
751             }
752         }
753     }CHRONOSTOP(24);CHRONOSTOP(2);
754   counters::stats();
755 }
756
757 /*! Get the neighbors of a cell.
758  * Only the neighbors having the dimension of the cell are taken into account
759  * (neighbors of a volume are the volumes sharing a face with this volume,
760  *  neighbors of a face are the faces sharing an edge with this face...).
761  * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
762  * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
763  * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
764  * @param vtkId the vtk id of the cell
765  * @return number of neighbors
766  */
767 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
768 {
769   int vtkType = this->GetCellType(vtkId);
770   int cellDim = SMDS_Downward::getCellDimension(vtkType);
771   if (cellDim <2)
772     return 0; // TODO voisins des edges = edges connectees
773   int cellId = this->_cellIdToDownId[vtkId];
774
775   int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
776   const int *downCells = _downArray[vtkType]->getDownCells(cellId);
777   const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
778
779   // --- iteration on faces of the 3D cell (or edges on the 2D cell).
780
781   int nb = 0;
782   for (int i = 0; i < nbCells; i++)
783     {
784       int downId = downCells[i];
785       int cellType = downTyp[i];
786       int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
787       const int *upCells = _downArray[cellType]->getUpCells(downId);
788       const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
789
790       // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
791       //    for a face, number of neighbors (connected faces) not known
792
793       for (int j = 0; j < nbUp; j++)
794         {
795           if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
796             continue;
797           int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
798           neighborsVtkIds[nb] = vtkNeighbor;
799           downIds[nb] = downId;
800           downTypes[nb] = cellType;
801           nb++;
802           if (nb >= NBMAXNEIGHBORS)
803             {
804               INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
805               return nb;
806             }
807         }
808       if (getSkin)
809         {
810           if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
811             {
812               neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
813               downIds[nb] = downId;
814               downTypes[nb] = cellType;
815               nb++;
816               if (nb >= NBMAXNEIGHBORS)
817                 {
818                   INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
819                   return nb;
820                 }
821             }
822         }
823     }
824   return nb;
825 }
826
827 /*! get the volumes containing a face or an edge of the grid
828  * The edge or face belongs to the vtkUnstructuredGrid
829  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
830  * @param vtkId vtk id of the face or edge
831  */
832 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
833 {
834   int vtkType = this->GetCellType(vtkId);
835   int dim = SMDS_Downward::getCellDimension(vtkType);
836   int nbFaces = 0;
837   unsigned char cellTypes[1000];
838   int downCellId[1000];
839   if (dim == 1)
840     {
841       int downId = this->CellIdToDownId(vtkId);
842       if (downId < 0)
843         {
844           MESSAGE("Downward structure not up to date: new edge not taken into account");
845           return 0;
846         }
847       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
848       const int *upCells = _downArray[vtkType]->getUpCells(downId);
849       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
850       for (int i=0; i< nbFaces; i++)
851         {
852           cellTypes[i] = upTypes[i];
853           downCellId[i] = upCells[i];
854         }
855     }
856   else if (dim == 2)
857     {
858       nbFaces = 1;
859       cellTypes[0] = this->GetCellType(vtkId);
860       int downId = this->CellIdToDownId(vtkId);
861       if (downId < 0)
862         {
863           MESSAGE("Downward structure not up to date: new face not taken into account");
864           return 0;
865         }
866       downCellId[0] = downId;
867     }
868
869   int nbvol =0;
870   for (int i=0; i<nbFaces; i++)
871     {
872       int vtkTypeFace = cellTypes[i];
873       int downId = downCellId[i];
874       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
875       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
876       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
877        for (int j=0; j<nv; j++)
878         {
879           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
880           if (vtkVolId >= 0)
881             volVtkIds[nbvol++] = vtkVolId;
882         }
883     }
884   return nbvol;
885 }
886
887 /*! get the volumes containing a face or an edge of the downward structure
888  * The edge or face does not necessary belong to the vtkUnstructuredGrid
889  * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
890  * @param downId id in the downward structure
891  * @param downType type of cell
892  */
893 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
894 {
895   int vtkType = downType;
896   int dim = SMDS_Downward::getCellDimension(vtkType);
897   int nbFaces = 0;
898   unsigned char cellTypes[1000];
899   int downCellId[1000];
900   if (dim == 1)
901     {
902       nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
903       const int *upCells = _downArray[vtkType]->getUpCells(downId);
904       const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
905       for (int i=0; i< nbFaces; i++)
906         {
907           cellTypes[i] = upTypes[i];
908           downCellId[i] = upCells[i];
909         }
910     }
911   else if (dim == 2)
912     {
913       nbFaces = 1;
914       cellTypes[0] = vtkType;
915       downCellId[0] = downId;
916     }
917
918   int nbvol =0;
919   for (int i=0; i<nbFaces; i++)
920     {
921       int vtkTypeFace = cellTypes[i];
922       int downId = downCellId[i];
923       int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
924       const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
925       const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
926        for (int j=0; j<nv; j++)
927         {
928           int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
929           if (vtkVolId >= 0)
930             volVtkIds[nbvol++] = vtkVolId;
931         }
932     }
933   return nbvol;
934 }
935
936 /*! get the node id's of a cell.
937  * The cell is defined by it's downward connectivity id and type.
938  * @param nodeSet set of of vtk node id's to fill.
939  * @param downId downward connectivity id of the cell.
940  * @param downType type of cell.
941  */
942 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
943 {
944   _downArray[downType]->getNodeIds(downId, nodeSet);
945 }
946
947 /*! change some nodes in cell without modifying type or internal connectivity.
948  * Nodes inverse connectivity is maintained up to date.
949  * @param vtkVolId vtk id of the cell
950  * @param localClonedNodeIds map old node id to new node id.
951  */
952 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
953 {
954   vtkIdType npts = 0;
955   vtkIdType *pts; // will refer to the point id's of the face
956   this->GetCellPoints(vtkVolId, npts, pts);
957   for (int i = 0; i < npts; i++)
958     {
959       if (localClonedNodeIds.count(pts[i]))
960         {
961           vtkIdType oldpt = pts[i];
962           pts[i] = localClonedNodeIds[oldpt];
963           //MESSAGE(oldpt << " --> " << pts[i]);
964           //this->RemoveReferenceToCell(oldpt, vtkVolId);
965           //this->AddReferenceToCell(pts[i], vtkVolId);
966         }
967     }
968 }
969
970 /*! reorder the nodes of a face
971  * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
972  * @param orderedNodes list of nodes to reorder (in out)
973  * @return size of the list
974  */
975 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
976 {
977   int vtkType = this->GetCellType(vtkVolId);
978   dim = SMDS_Downward::getCellDimension(vtkType);
979   if (dim == 3)
980     {
981       SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
982       int downVolId = this->_cellIdToDownId[vtkVolId];
983       downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
984     }
985   // else nothing to do;
986   return orderedNodes.size();
987 }
988
989 void SMDS_UnstructuredGrid::BuildLinks()
990 {
991   // Remove the old links if they are already built
992   if (this->Links)
993     {
994     this->Links->UnRegister(this);
995     }
996
997   this->Links = SMDS_CellLinks::New();
998   this->Links->Allocate(this->GetNumberOfPoints());
999   this->Links->Register(this);
1000   this->Links->BuildLinks(this, this->Connectivity);
1001   this->Links->Delete();
1002 }
1003
1004 /*! Create a volume (prism or hexahedron) by duplication of a face.
1005  * Designed for use in creation of flat elements separating volume domains.
1006  * A face separating two domains is shared by two volume cells.
1007  * All the nodes are already created (for the two faces).
1008  * Each original Node is associated to corresponding nodes in the domains.
1009  * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1010  * In that case, even some of the nodes to use for the original face may be changed.
1011  * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1012  * @param domain1: domain of the original face
1013  * @param domain2: domain of the duplicated face
1014  * @param originalNodes: the vtk node ids of the original face
1015  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1016  * @return ok if success.
1017  */
1018 SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1019                                                   int domain1,
1020                                                   int domain2,
1021                                                   std::set<int>& originalNodes,
1022                                                   std::map<int, std::map<int, int> >& nodeDomains,
1023                                                   std::map<int, std::map<long, int> >& nodeQuadDomains)
1024 {
1025   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1026   vector<vtkIdType> orderedOriginals;
1027   orderedOriginals.clear();
1028   set<int>::const_iterator it = originalNodes.begin();
1029   for (; it != originalNodes.end(); ++it)
1030     orderedOriginals.push_back(*it);
1031
1032   int dim = 0;
1033   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1034   vector<vtkIdType> orderedNodes;
1035
1036   bool isQuadratic = false;
1037   switch (orderedOriginals.size())
1038   {
1039     case 3:
1040       if (dim == 2)
1041         isQuadratic = true;
1042       break;
1043     case 6:
1044     case 8:
1045       isQuadratic = true;
1046       break;
1047     default:
1048       isQuadratic = false;
1049       break;
1050   }
1051
1052   if (isQuadratic)
1053     {
1054       long dom1 = domain1;
1055       long dom2 = domain2;
1056       long dom1_2; // for nodeQuadDomains
1057       if (domain1 < domain2)
1058         dom1_2 = dom1 + INT_MAX * dom2;
1059       else
1060         dom1_2 = dom2 + INT_MAX * dom1;
1061       //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1062       int ima = orderedOriginals.size();
1063       int mid = orderedOriginals.size() / 2;
1064       //cerr << "ima=" << ima << " mid=" << mid << endl;
1065       for (int i = 0; i < mid; i++)
1066         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1067       for (int i = 0; i < mid; i++)
1068         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1069       for (int i = mid; i < ima; i++)
1070         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1071       for (int i = mid; i < ima; i++)
1072         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1073       for (int i = 0; i < mid; i++)
1074         {
1075           int oldId = orderedOriginals[i];
1076           int newId;
1077           if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1078             newId = nodeQuadDomains[oldId][dom1_2];
1079           else
1080             {
1081               double *coords = this->GetPoint(oldId);
1082               SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1083               newId = newNode->getVtkId();
1084               std::map<long, int> emptyMap;
1085               nodeQuadDomains[oldId] = emptyMap;
1086               nodeQuadDomains[oldId][dom1_2] = newId;
1087             }
1088           orderedNodes.push_back(newId);
1089         }
1090     }
1091   else
1092     {
1093       for (int i = 0; i < nbNodes; i++)
1094         orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1095       if (dim == 3)
1096         for (int i = 0; i < nbNodes; i++)
1097           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1098       else
1099         for (int i = nbNodes-1; i >= 0; i--)
1100           orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1101
1102     }
1103
1104   if (dim == 3)
1105     {
1106       SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1107       return vol;
1108     }
1109   else if (dim == 2)
1110     {
1111       SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1112       return face;
1113     }
1114
1115   // TODO update sub-shape list of elements and nodes
1116   return 0;
1117 }
1118
1119 //================================================================================
1120 /*!
1121  * \brief Allocates data array for ball diameters
1122  *  \param MaxVtkID - max ID of a ball element
1123  */
1124 //================================================================================
1125
1126 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1127 {
1128   SetBallDiameter( MaxVtkID, 0 );
1129 }
1130
1131 //================================================================================
1132 /*!
1133  * \brief Sets diameter of a ball element
1134  *  \param vtkID - vtk id of the ball element
1135  *  \param diameter - diameter of the ball element
1136  */
1137 //================================================================================
1138
1139 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1140 {
1141   vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1142   if ( !array )
1143   {
1144     array = vtkDoubleArray::New();
1145     array->SetNumberOfComponents(1);
1146     vtkDataSet::CellData->SetScalars( array );
1147   }
1148   array->InsertValue( vtkID, diameter );
1149 }
1150
1151 //================================================================================
1152 /*!
1153  * \brief Returns diameter of a ball element
1154  *  \param vtkID - vtk id of the ball element
1155  */
1156 //================================================================================
1157
1158 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1159 {
1160   if ( vtkDataSet::CellData )
1161     return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
1162   return 0;
1163 }
1164