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