Salome HOME
bos #24368 EDF 23667 - Duplicates nodes
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
index 52916ec724c12d491ed92982b26e50e30fae467b..124d069f1b3e407335b4a4c78fbb26f31fd3e37b 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2010-2015  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2010-2021  CEA/DEN, EDF R&D, OPEN CASCADE
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
 #include <list>
 #include <climits>
 
-using namespace std;
-
 SMDS_CellLinks* SMDS_CellLinks::New()
 {
-  MESSAGE("SMDS_CellLinks::New");
   return new SMDS_CellLinks();
 }
 
@@ -49,11 +46,56 @@ void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
   if ( vtkID > this->MaxId )
   {
     this->MaxId = vtkID;
-    if ( vtkID >= this->Size ) 
+    if ( vtkID >= this->Size )
       vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
   }
 }
 
+void SMDS_CellLinks::BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types)
+{
+  // build links taking into account removed cells
+
+  vtkIdType numPts = data->GetNumberOfPoints();
+  vtkIdType j, cellId = 0;
+  unsigned short *linkLoc;
+  vtkIdType npts=0;
+  vtkIdType const *pts(nullptr);
+  vtkIdType loc = Connectivity->GetTraversalLocation();
+
+  // traverse data to determine number of uses of each point
+  cellId = 0;
+  for (Connectivity->InitTraversal();
+       Connectivity->GetNextCell(npts,pts); cellId++)
+  {
+    if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
+      for (j=0; j < npts; j++)
+      {
+        this->IncrementLinkCount(pts[j]);
+      }
+  }
+
+  // now allocate storage for the links
+  this->AllocateLinks(numPts);
+  this->MaxId = numPts - 1;
+
+  // fill out lists with references to cells
+  linkLoc = new unsigned short[numPts];
+  memset(linkLoc, 0, numPts*sizeof(unsigned short));
+
+  cellId = 0;
+  for (Connectivity->InitTraversal();
+       Connectivity->GetNextCell(npts,pts); cellId++)
+  {
+    if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
+      for (j=0; j < npts; j++)
+      {
+        this->InsertCellReference(pts[j], (linkLoc[pts[j]])++, cellId);
+      }
+  }
+  delete [] linkLoc;
+  Connectivity->SetTraversalLocation(loc);
+}
+
 SMDS_CellLinks::SMDS_CellLinks() :
   vtkCellLinks()
 {
@@ -65,7 +107,6 @@ SMDS_CellLinks::~SMDS_CellLinks()
 
 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
 {
-  MESSAGE("SMDS_UnstructuredGrid::New");
   return new SMDS_UnstructuredGrid();
 }
 
@@ -82,288 +123,273 @@ SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
 {
 }
 
-unsigned long SMDS_UnstructuredGrid::GetMTime()
+vtkMTimeType SMDS_UnstructuredGrid::GetMTime()
 {
-  unsigned long mtime = vtkUnstructuredGrid::GetMTime();
-  MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
+  vtkMTimeType mtime = vtkUnstructuredGrid::GetMTime();
   return mtime;
 }
-// OUV_PORTING_VTK6: seems to be useless
-/*
-void SMDS_UnstructuredGrid::Update()
-{
-  MESSAGE("SMDS_UnstructuredGrid::Update");
-  return vtkUnstructuredGrid::Update();
-}
 
-void SMDS_UnstructuredGrid::UpdateInformation()
-{
-  MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
-  return vtkUnstructuredGrid::UpdateInformation();
-}
-*/
 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
 {
   // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
-  //MESSAGE("*********************** SMDS_UnstructuredGrid::GetPoints " << this->Points << " " << vtkUnstructuredGrid::GetPoints());
   return this->Points;
 }
 
-//#ifdef VTK_HAVE_POLYHEDRON
-int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
+vtkIdType SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
 {
-  if (type != VTK_POLYHEDRON)
+  if ( !this->Links ) // don't create Links until they are needed
+  {
+    return this->InsertNextCell(type, npts, pts);
+  }
+
+  if ( type != VTK_POLYHEDRON )
     return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
 
   // --- type = VTK_POLYHEDRON
-  //MESSAGE("InsertNextLinkedCell VTK_POLYHEDRON");
-  int cellid = this->InsertNextCell(type, npts, pts);
+  vtkIdType cellid = this->InsertNextCell(type, npts, pts);
 
-  set<vtkIdType> setOfNodes;
+  std::set<vtkIdType> setOfNodes;
   setOfNodes.clear();
   int nbfaces = npts;
   int i = 0;
   for (int nf = 0; nf < nbfaces; nf++)
+  {
+    int nbnodes = pts[i];
+    i++;
+    for (int k = 0; k < nbnodes; k++)
     {
-      int nbnodes = pts[i];
+      if ( setOfNodes.insert( pts[i] ).second )
+      {
+        (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->ResizeCellList( pts[i], 1 );
+        (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->AddCellReference( cellid, pts[i] );
+      }
       i++;
-      for (int k = 0; k < nbnodes; k++)
-        {
-          //MESSAGE(" cell " << cellid << " face " << nf << " node " << pts[i]);
-          setOfNodes.insert(pts[i]);
-          i++;
-        }
-    }
-
-  set<vtkIdType>::iterator it = setOfNodes.begin();
-  for (; it != setOfNodes.end(); ++it)
-    {
-      //MESSAGE("reverse link for node " << *it << " cell " << cellid);
-      this->Links->ResizeCellList(*it, 1);
-      this->Links->AddCellReference(cellid, *it);
     }
+  }
 
   return cellid;
 }
-//#endif
 
 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
 {
   _mesh = mesh;
 }
 
-void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
-                                        std::vector<int>& idCellsOldToNew, int newCellSize)
+void SMDS_UnstructuredGrid::compactGrid(std::vector<smIdType>& idNodesOldToNew, smIdType newNodeSize,
+                                        std::vector<smIdType>& idCellsNewToOld, smIdType newCellSize)
 {
-  MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);//CHRONO(1);
-  int alreadyCopied = 0;
+  this->DeleteLinks();
 
-  // --- if newNodeSize, create a new compacted vtkPoints
+  // IDs of VTK nodes always correspond to SMDS IDs but there can be "holes" in SMDS numeration.
+  // We compact only if there were holes
 
-  vtkPoints *newPoints = vtkPoints::New();
-  newPoints->SetDataType(VTK_DOUBLE);
-  newPoints->SetNumberOfPoints(newNodeSize);
-  if (newNodeSize)
+  vtkIdType oldNodeSize = this->GetNumberOfPoints();
+  bool updateNodes = ( oldNodeSize > newNodeSize );
+  if ( true /*updateNodes*/ )
+  {
+    // 21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion
+    // Use double type for storing coordinates of nodes instead float.
+    vtkPoints *newPoints = vtkPoints::New();
+    newPoints->SetDataType( VTK_DOUBLE );
+    newPoints->SetNumberOfPoints( FromSmIdType<vtkIdType>(newNodeSize) );
+
+    vtkIdType i = 0, alreadyCopied = 0;
+    while ( i < oldNodeSize )
     {
-      MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
-      // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
-      // using double type for storing coordinates of nodes instead float.
-      int oldNodeSize = idNodesOldToNew.size();
+      // skip a hole if any
+      while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
+        ++i;
+      vtkIdType startBloc = i;
+      // look for a block end
+      while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
+        ++i;
+      vtkIdType endBloc = i;
+      copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
+    }
+    this->SetPoints(newPoints);
+    newPoints->Delete();
+  }
+  else
+  {
+    this->Points->Squeeze();
+    this->Points->Modified();
+  }
 
-      int i = 0;
-      while ( i < oldNodeSize )
-      {
-        // skip a hole if any
-        while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
-          ++i;
-        int startBloc = i;
-        // look for a block end
-        while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
-          ++i;
-        int endBloc = i;
-        copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
-      }
-      newPoints->Squeeze();
+  // Compact cells if VTK IDs do not correspond to SMDS IDs or nodes compacted
+
+  vtkIdType oldCellSize = this->Types->GetNumberOfTuples();
+  bool      updateCells = ( updateNodes || newCellSize != oldCellSize );
+  for ( vtkIdType newID = 0, nbIDs = idCellsNewToOld.size(); newID < nbIDs &&  !updateCells; ++newID )
+    updateCells = ( idCellsNewToOld[ newID ] != newID );
+
+  if ( false /*!updateCells*/ ) // no holes in elements
+  {
+    this->Connectivity->Squeeze();
+    this->CellLocations->Squeeze();
+    this->Types->Squeeze();
+    if ( this->FaceLocations )
+    {
+      this->FaceLocations->Squeeze();
+      this->Faces->Squeeze();
     }
+    this->Connectivity->Modified();
+    return;
+  }
+
+  if ((vtkIdType) idNodesOldToNew.size() < oldNodeSize )
+  {
+    idNodesOldToNew.reserve( oldNodeSize );
+    for ( vtkIdType i = idNodesOldToNew.size(); i < oldNodeSize; ++i )
+      idNodesOldToNew.push_back( i );
+  }
 
   // --- create new compacted Connectivity, Locations and Types
 
-  int oldCellSize = this->Types->GetNumberOfTuples();
+  vtkIdType newConnectivitySize = this->Connectivity->GetNumberOfConnectivityEntries();
+  if ( newCellSize != oldCellSize )
+    for ( vtkIdType i = 0; i < oldCellSize - 1; ++i )
+      if ( this->Types->GetValue( i ) == VTK_EMPTY_CELL )
+        newConnectivitySize -= this->Connectivity->GetCellSize( i );
 
   vtkCellArray *newConnectivity = vtkCellArray::New();
   newConnectivity->Initialize();
-  int oldCellDataSize = this->Connectivity->GetData()->GetSize();
-  newConnectivity->Allocate(oldCellDataSize);
-  MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
+  newConnectivity->Allocate( newConnectivitySize );
 
   vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
   newTypes->Initialize();
-  newTypes->SetNumberOfValues(newCellSize);
+  newTypes->SetNumberOfValues(FromSmIdType<vtkIdType>(newCellSize));
 
   vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
   newLocations->Initialize();
-  newLocations->SetNumberOfValues(newCellSize);
+  newLocations->SetNumberOfValues(FromSmIdType<vtkIdType>(newCellSize));
 
-  // TODO some polyhedron may be huge (only in some tests)
-  vtkIdType tmpid[NBMAXNODESINCELL];
-  vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
+  std::vector< vtkIdType > pointsCell(1024); // --- points id to fill a new cell
 
-  alreadyCopied = 0;
-  int i = 0;
-  while ( i < oldCellSize )
-  {
-    // skip a hole if any
-    while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
-      ++i;
-    int startBloc = i;
-    // look for a block end
-    while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
-      ++i;
-    int endBloc = i;
-    if ( endBloc > startBloc )
-      copyBloc(newTypes,
-               idCellsOldToNew, idNodesOldToNew,
-               newConnectivity, newLocations,
-               pointsCell, alreadyCopied,
-               startBloc, endBloc);
-  }
-  newConnectivity->Squeeze();
-
-  if (1/*newNodeSize*/)
-    {
-      MESSAGE("------- newNodeSize, setPoints");
-      this->SetPoints(newPoints);
-      MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
-    }
+  copyBloc(newTypes, idCellsNewToOld, idNodesOldToNew,
+           newConnectivity, newLocations, pointsCell );
 
   if (vtkDoubleArray* diameters =
       vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
   {
-    for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
+    vtkDoubleArray* newDiameters = vtkDoubleArray::New();
+    newDiameters->SetNumberOfComponents(1);
+    for ( vtkIdType newCellID = 0; newCellID < newCellSize; newCellID++ )
     {
-      if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
-        continue;
-      int newCellId = idCellsOldToNew[ oldCellID ];
-      if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
-        diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
+      if ( newTypes->GetValue( newCellID ) == VTK_POLY_VERTEX )
+      {
+        vtkIdType oldCellID = idCellsNewToOld[ newCellID ];
+        newDiameters->InsertValue( newCellID, diameters->GetValue( oldCellID ));
+      }
+      vtkDataSet::CellData->SetScalars( newDiameters );
     }
   }
 
   if (this->FaceLocations)
+  {
+    vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
+    newFaceLocations->Initialize();
+    newFaceLocations->Allocate(newTypes->GetSize());
+    vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
+    newFaces->Initialize();
+    newFaces->Allocate(this->Faces->GetSize());
+    for ( vtkIdType newCellID = 0; newCellID < newCellSize; newCellID++ )
     {
-      vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
-      newFaceLocations->Initialize();
-      newFaceLocations->Allocate(newTypes->GetSize());
-      vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
-      newFaces->Initialize();
-      newFaces->Allocate(this->Faces->GetSize());
-      for (int i = 0; i < oldCellSize; i++)
+      if ( newTypes->GetValue( newCellID ) == VTK_POLYHEDRON )
+      {
+        smIdType oldCellId = idCellsNewToOld[ newCellID ];
+        newFaceLocations->InsertNextValue( newFaces->GetMaxId()+1 );
+        smIdType oldFaceLoc = this->FaceLocations->GetValue( FromSmIdType<vtkIdType>(oldCellId) );
+        smIdType nCellFaces = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
+        newFaces->InsertNextValue( FromSmIdType<vtkIdType>(nCellFaces) );
+        for ( int n = 0; n < nCellFaces; n++ )
         {
-          if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
-            continue;
-          int newCellId = idCellsOldToNew[i];
-          if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
-            {
-               newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
-               int oldFaceLoc = this->FaceLocations->GetValue(i);
-               int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
-               newFaces->InsertNextValue(nCellFaces);
-               for (int n=0; n<nCellFaces; n++)
-                 {
-                   int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
-                   newFaces->InsertNextValue(nptsInFace);
-                   for (int k=0; k<nptsInFace; k++)
-                     {
-                       int oldpt = this->Faces->GetValue(oldFaceLoc++);
-                       newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
-                     }
-                 }
-            }
-          else
-            {
-               newFaceLocations->InsertNextValue(-1);
-            }
+          int nptsInFace = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
+          newFaces->InsertNextValue( nptsInFace );
+          for ( int k = 0; k < nptsInFace; k++ )
+          {
+            vtkIdType oldpt = this->Faces->GetValue( FromSmIdType<vtkIdType>(oldFaceLoc++) );
+            newFaces->InsertNextValue( idNodesOldToNew[ oldpt ]);
+          }
         }
-      newFaceLocations->Squeeze();
-      newFaces->Squeeze();
-      this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
-      newFaceLocations->Delete();
-      newFaces->Delete();
+      }
+      else
+      {
+        newFaceLocations->InsertNextValue(-1);
+      }
     }
+    newFaceLocations->Squeeze();
+    newFaces->Squeeze();
+    this->SetCells( newTypes, newLocations, newConnectivity, newFaceLocations, newFaces );
+    this->CellLocations = newLocations;
+    newFaceLocations->Delete();
+    newFaces->Delete();
+  }
   else
   {
-    this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
+    this->SetCells( newTypes, newLocations, newConnectivity, FaceLocations, Faces );
+    this->CellLocations = newLocations;
   }
 
-  newPoints->Delete();
   newTypes->Delete();
   newLocations->Delete();
   newConnectivity->Delete();
-  this->BuildLinks();
 }
 
-void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
-                                      int start, int end)
+void SMDS_UnstructuredGrid::copyNodes(vtkPoints *             newPoints,
+                                      std::vector<smIdType>& /*idNodesOldToNew*/,
+                                      vtkIdType&              alreadyCopied,
+                                      vtkIdType               start,
+                                      vtkIdType               end)
 {
-  MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
   void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
   void *source = this->Points->GetVoidPointer(3 * start);
-  int nbPoints = end - start;
+  vtkIdType nbPoints = end - start;
   if (nbPoints > 0)
-    {
-      memcpy(target, source, 3 * sizeof(double) * nbPoints);
-      for (int j = start; j < end; j++)
-        idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
-    }
+  {
+    memcpy(target, source, 3 * sizeof(double) * nbPoints);
+    alreadyCopied += nbPoints;
+  }
 }
 
-void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
-                                     std::vector<int>&     idCellsOldToNew,
-                                     std::vector<int>&     idNodesOldToNew,
-                                     vtkCellArray*         newConnectivity,
-                                     vtkIdTypeArray*       newLocations,
-                                     vtkIdType*            pointsCell,
-                                     int&                  alreadyCopied,
-                                     int                   start,
-                                     int                   end)
+void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *  newTypes,
+                                     const std::vector<smIdType>& idCellsNewToOld,
+                                     const std::vector<smIdType>& idNodesOldToNew,
+                                     vtkCellArray*           newConnectivity,
+                                     vtkIdTypeArray*         newLocations,
+                                     std::vector<vtkIdType>& pointsCell)
 {
-  //MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
-  for (int j = start; j < end; j++)
+  for ( size_t iNew = 0; iNew < idCellsNewToOld.size(); iNew++ )
+  {
+    vtkIdType iOld = idCellsNewToOld[ iNew ];
+    newTypes->SetValue( iNew, this->Types->GetValue( iOld ));
+
+    vtkIdType oldLoc = ((vtkIdTypeArray *)(this->Connectivity->GetOffsetsArray()))->GetValue( iOld );
+    vtkIdType nbpts;
+    vtkIdType const *oldPtsCell(nullptr);
+    this->Connectivity->GetCell( oldLoc+iOld, nbpts, oldPtsCell );
+    if ((vtkIdType) pointsCell.size() < nbpts )
+      pointsCell.resize( nbpts );
+    for ( int l = 0; l < nbpts; l++ )
     {
-      newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
-      idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
-      vtkIdType oldLoc = this->Locations->GetValue(j);
-      vtkIdType nbpts;
-      vtkIdType *oldPtsCell = 0;
-      this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
-      assert(nbpts < NBMAXNODESINCELL);
-      //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
-      for (int l = 0; l < nbpts; l++)
-        {
-          int oldval = oldPtsCell[l];
-          pointsCell[l] = idNodesOldToNew[oldval];
-          //MESSAGE("   " << oldval << " " << pointsCell[l]);
-        }
-      /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
-      int newLoc = newConnectivity->GetInsertLocation(nbpts);
-      //MESSAGE(newcnt << " " << newLoc);
-      newLocations->SetValue(alreadyCopied, newLoc);
-      alreadyCopied++;
+      vtkIdType oldval = oldPtsCell[l];
+      pointsCell[l] = idNodesOldToNew[oldval];
     }
+    /*vtkIdType newcnt = */newConnectivity->InsertNextCell( nbpts, pointsCell.data() );
+    vtkIdType newLoc = newConnectivity->GetInsertLocation( nbpts );
+    newLocations->SetValue( iNew, newLoc );
+  }
 }
 
-int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
+int SMDS_UnstructuredGrid::CellIdToDownId(vtkIdType vtkCellId)
 {
-  if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
+  if ((vtkCellId < 0) || (vtkCellId >= (vtkIdType)_cellIdToDownId.size()))
   {
-    //MESSAGE("SMDS_UnstructuredGrid::CellIdToDownId structure not up to date: vtkCellId="
-    //    << vtkCellId << " max="<< _cellIdToDownId.size());
     return -1;
   }
   return _cellIdToDownId[vtkCellId];
 }
 
-void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
+void SMDS_UnstructuredGrid::setCellIdToDownId(vtkIdType vtkCellId, int downId)
 {
   // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
   _cellIdToDownId[vtkCellId] = downId;
@@ -384,7 +410,7 @@ void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
  *  Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
  *
  */
-void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
+void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool /*withEdges*/)
 {
   MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
   // TODO calcul partiel sans edges
@@ -420,15 +446,15 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
 
   const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
 
-  int nbLinTetra  = meshInfo.NbTetras  (ORDER_LINEAR);
-  int nbQuadTetra = meshInfo.NbTetras  (ORDER_QUADRATIC);
-  int nbLinPyra   = meshInfo.NbPyramids(ORDER_LINEAR);
-  int nbQuadPyra  = meshInfo.NbPyramids(ORDER_QUADRATIC);
-  int nbLinPrism  = meshInfo.NbPrisms  (ORDER_LINEAR);
-  int nbQuadPrism = meshInfo.NbPrisms  (ORDER_QUADRATIC);
-  int nbLinHexa   = meshInfo.NbHexas   (ORDER_LINEAR);
-  int nbQuadHexa  = meshInfo.NbHexas   (ORDER_QUADRATIC);
-  int nbHexPrism  = meshInfo.NbHexPrisms();
+  int nbLinTetra  = FromSmIdType<int>(meshInfo.NbTetras  (ORDER_LINEAR));
+  int nbQuadTetra = FromSmIdType<int>(meshInfo.NbTetras  (ORDER_QUADRATIC));
+  int nbLinPyra   = FromSmIdType<int>(meshInfo.NbPyramids(ORDER_LINEAR));
+  int nbQuadPyra  = FromSmIdType<int>(meshInfo.NbPyramids(ORDER_QUADRATIC));
+  int nbLinPrism  = FromSmIdType<int>(meshInfo.NbPrisms  (ORDER_LINEAR));
+  int nbQuadPrism = FromSmIdType<int>(meshInfo.NbPrisms  (ORDER_QUADRATIC));
+  int nbLinHexa   = FromSmIdType<int>(meshInfo.NbHexas   (ORDER_LINEAR));
+  int nbQuadHexa  = FromSmIdType<int>(meshInfo.NbHexas   (ORDER_QUADRATIC));
+  int nbHexPrism  = FromSmIdType<int>(meshInfo.NbHexPrisms());
 
   int nbLineGuess     = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
   int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
@@ -456,6 +482,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
   GuessSize[VTK_QUADRATIC_HEXAHEDRON]    = nbQuadHexa;
   GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
   GuessSize[VTK_HEXAGONAL_PRISM]         = nbHexPrism;
+  (void)GuessSize; // unused in Release mode
 
   _downArray[VTK_LINE]                   ->allocate(nbLineGuess);
   _downArray[VTK_QUADRATIC_EDGE]         ->allocate(nbQuadEdgeGuess);
@@ -625,7 +652,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
           int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
           SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
           downEdge->setNodes(connEdgeId, vtkEdgeId);
-          vector<int> vtkIds;
+          std::vector<int> vtkIds;
           int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
           int downFaces[1000];
           unsigned char downTypes[1000];
@@ -674,7 +701,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
               // --- check if the edge is already registered by exploration of the faces
 
               //CHRONO(41);
-              vector<int> vtkIds;
+              std::vector<int> vtkIds;
               unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
               int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
               SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
@@ -728,7 +755,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
 
   CHRONOSTOP(23);CHRONO(24);
 
-  // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
+  // compact downward connectivity structure: adjust downward arrays size, replace std::vector<vector int>> by a single std::vector<int>
   // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
 
   for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
@@ -953,8 +980,10 @@ void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsig
 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
 {
   vtkIdType npts = 0;
-  vtkIdType *pts; // will refer to the point id's of the face
-  this->GetCellPoints(vtkVolId, npts, pts);
+  vtkIdType const *tmp(nullptr); // will refer to the point id's of the face
+  vtkIdType *pts;                // will refer to the point id's of the face
+  this->GetCellPoints(vtkVolId, npts, tmp);
+  pts = const_cast< vtkIdType*>( tmp );
   for (int i = 0; i < npts; i++)
     {
       if (localClonedNodeIds.count(pts[i]))
@@ -975,14 +1004,14 @@ void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> loc
  */
 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
 {
-  int vtkType = this->GetCellType(vtkVolId);
-  dim = SMDS_Downward::getCellDimension(vtkType);
+  int vtkType = this->GetCellType( vtkVolId );
+  dim = SMDS_Downward::getCellDimension( vtkType );
   if (dim == 3)
-    {
-      SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
-      int downVolId = this->_cellIdToDownId[vtkVolId];
-      downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
-    }
+  {
+    SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
+    int        downVolId = this->_cellIdToDownId[ vtkVolId ];
+    downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
+  }
   // else nothing to do;
   return orderedNodes.size();
 }
@@ -991,17 +1020,34 @@ void SMDS_UnstructuredGrid::BuildLinks()
 {
   // Remove the old links if they are already built
   if (this->Links)
-    {
+  {
     this->Links->UnRegister(this);
-    }
+  }
 
-  this->Links = SMDS_CellLinks::New();
-  this->Links->Allocate(this->GetNumberOfPoints());
+  SMDS_CellLinks* links;
+  this->Links = links = SMDS_CellLinks::New();
+  (static_cast< vtkCellLinks *>(this->Links.GetPointer()))->Allocate(this->GetNumberOfPoints());
   this->Links->Register(this);
-  this->Links->BuildLinks(this, this->Connectivity);
+  links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
   this->Links->Delete();
 }
 
+void SMDS_UnstructuredGrid::DeleteLinks()
+{
+  // Remove the old links if they are already built
+  if (this->Links)
+  {
+    this->Links->UnRegister(this);
+    this->Links = NULL;
+  }
+}
+SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
+{
+  if ( !this->Links )
+    BuildLinks();
+  return static_cast< SMDS_CellLinks* >( this->Links.GetPointer() );
+}
+
 /*! Create a volume (prism or hexahedron) by duplication of a face.
  * Designed for use in creation of flat elements separating volume domains.
  * A face separating two domains is shared by two volume cells.
@@ -1016,105 +1062,116 @@ void SMDS_UnstructuredGrid::BuildLinks()
  * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
  * @return ok if success.
  */
-SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
-                                                  int domain1,
-                                                  int domain2,
-                                                  std::set<int>& originalNodes,
-                                                  std::map<int, std::map<int, int> >& nodeDomains,
-                                                  std::map<int, std::map<long, int> >& nodeQuadDomains)
+SMDS_MeshCell*
+SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
+                                             int domain1,
+                                             int domain2,
+                                             std::set<int>& originalNodes,
+                                             std::map<int, std::map<int, int> >& nodeDomains,
+                                             std::map<int, std::map<long, int> >& nodeQuadDomains)
 {
   //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
-  vector<vtkIdType> orderedOriginals;
-  orderedOriginals.clear();
-  set<int>::const_iterator it = originalNodes.begin();
-  for (; it != originalNodes.end(); ++it)
-    orderedOriginals.push_back(*it);
+  std::vector<vtkIdType> orderedOriginals( originalNodes.begin(), originalNodes.end() );
 
   int dim = 0;
   int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
-  vector<vtkIdType> orderedNodes;
+  std::vector<vtkIdType> orderedNodes;
 
   bool isQuadratic = false;
   switch (orderedOriginals.size())
   {
-    case 3:
-      if (dim == 2)
-        isQuadratic = true;
-      break;
-    case 6:
-    case 8:
+  case 3:
+    if (dim == 2)
       isQuadratic = true;
-      break;
-    default:
-      isQuadratic = false;
-      break;
+    break;
+  case 6:
+  case 8:
+    isQuadratic = true;
+    break;
+  default:
+    isQuadratic = false;
+    break;
   }
 
   if (isQuadratic)
+  {
+    long dom1 = domain1;
+    long dom2 = domain2;
+    long dom1_2; // for nodeQuadDomains
+    if (domain1 < domain2)
+      dom1_2 = dom1 + INT_MAX * dom2;
+    else
+      dom1_2 = dom2 + INT_MAX * dom1;
+    //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
+    int ima = orderedOriginals.size();
+    int mid = orderedOriginals.size() / 2;
+    //cerr << "ima=" << ima << " mid=" << mid << endl;
+    for (int i = 0; i < mid; i++)
+      orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+    for (int i = 0; i < mid; i++)
+      orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+    for (int i = mid; i < ima; i++)
+      orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+    for (int i = mid; i < ima; i++)
+      orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+    for (int i = 0; i < mid; i++)
     {
-      long dom1 = domain1;
-      long dom2 = domain2;
-      long dom1_2; // for nodeQuadDomains
-      if (domain1 < domain2)
-        dom1_2 = dom1 + INT_MAX * dom2;
+      int oldId = orderedOriginals[i];
+      int newId;
+      if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
+        newId = nodeQuadDomains[oldId][dom1_2];
       else
-        dom1_2 = dom2 + INT_MAX * dom1;
-      //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
-      int ima = orderedOriginals.size();
-      int mid = orderedOriginals.size() / 2;
-      //cerr << "ima=" << ima << " mid=" << mid << endl;
-      for (int i = 0; i < mid; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-      for (int i = 0; i < mid; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
-      for (int i = mid; i < ima; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-      for (int i = mid; i < ima; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
-      for (int i = 0; i < mid; i++)
+      {
+        double *coords = this->GetPoint(oldId);
+        SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
+        newId = newNode->GetVtkID();
+        if (! nodeQuadDomains.count(oldId))
         {
-          int oldId = orderedOriginals[i];
-          int newId;
-          if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
-            newId = nodeQuadDomains[oldId][dom1_2];
-          else
-            {
-              double *coords = this->GetPoint(oldId);
-              SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
-              newId = newNode->getVtkId();
-              if (! nodeQuadDomains.count(oldId))
-                {
-                  std::map<long, int> emptyMap;
-                  nodeQuadDomains[oldId] = emptyMap;
-                }
-              nodeQuadDomains[oldId][dom1_2] = newId;
-            }
-          orderedNodes.push_back(newId);
+          std::map<long, int> emptyMap;
+          nodeQuadDomains[oldId] = emptyMap;
         }
+        nodeQuadDomains[oldId][dom1_2] = newId;
+      }
+      orderedNodes.push_back(newId);
     }
+  }
   else
-    {
+  {
+    for (int i = 0; i < nbNodes; i++)
+      orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
+    if (dim == 3)
       for (int i = 0; i < nbNodes; i++)
-        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
-      if (dim == 3)
-        for (int i = 0; i < nbNodes; i++)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
-      else
-        for (int i = nbNodes-1; i >= 0; i--)
-          orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
+    else
+      for (int i = nbNodes-1; i >= 0; i--)
+        orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
 
-    }
+  }
 
   if (dim == 3)
-    {
-      SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
-      return vol;
-    }
+  {
+    SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
+    return vol;
+  }
   else if (dim == 2)
+  {
+    // bos #24368
+    // orient face by the original one, as getOrderedNodesOfFace() not implemented for faces
+    const SMDS_MeshElement* origFace = _mesh->FindElementVtk( vtkVolId );
+    int   i0 = origFace->GetNodeIndex( _mesh->FindNodeVtk( orderedNodes[0] ));
+    int   i1 = origFace->GetNodeIndex( _mesh->FindNodeVtk( orderedNodes[1] ));
+    int diff = i0 - i1;
+    // order of nodes must be reverse in face and origFace
+    bool oriOk = ( diff == 1 ) || ( diff == -3 );
+    if ( !oriOk )
     {
-      SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
-      return face;
+      SMDSAbs_EntityType type = isQuadratic ? SMDSEntity_Quad_Quadrangle : SMDSEntity_Quadrangle;
+      const std::vector<int>& interlace = SMDS_MeshCell::reverseSmdsOrder( type );
+      SMDS_MeshCell::applyInterlace( interlace, orderedNodes );
     }
+    SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
+    return face;
+  }
 
   // TODO update sub-shape list of elements and nodes
   return 0;
@@ -1165,4 +1222,3 @@ double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
     return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );
   return 0;
 }
-