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