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