]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
PR: debug compact vtkUnstructuredGrid in progress
authorprascle <prascle>
Thu, 15 Apr 2010 15:52:18 +0000 (15:52 +0000)
committerprascle <prascle>
Thu, 15 Apr 2010 15:52:18 +0000 (15:52 +0000)
14 files changed:
src/SMDS/SMDS_Mesh.cxx
src/SMDS/SMDS_Mesh.hxx
src/SMDS/SMDS_MeshElement.hxx
src/SMDS/SMDS_MeshIDFactory.cxx
src/SMDS/SMDS_MeshIDFactory.hxx
src/SMDS/SMDS_MeshNode.cxx
src/SMDS/SMDS_UnstructuredGrid.cxx
src/SMDS/SMDS_UnstructuredGrid.hxx
src/SMDS/SMDS_VtkCellIterator.cxx
src/SMESH/SMESH_Gen.cxx
src/SMESHDS/SMESHDS_Mesh.cxx
src/SMESHDS/SMESHDS_Mesh.hxx
src/SMESHDS/SMESHDS_SubMesh.cxx
src/SMESHDS/SMESHDS_SubMesh.hxx

index c6f18287e5ee81d0915ecfd591394cbf464d1f29..1f704a681d20b72e46de6438efc7fcfbe99d5aed 100644 (file)
@@ -123,7 +123,7 @@ SMDS_Mesh::SMDS_Mesh()
          myElementIDFactory(new SMDS_MeshElementIDFactory()),
          myHasConstructionEdges(false), myHasConstructionFaces(false),
          myHasInverseElements(true),
-         myNodeMin(0), myNodeMax(0), myCellLinksSize(0),
+         myNodeMin(0), myNodeMax(0),
          myNodePool(0), myEdgePool(0), myFacePool(0), myVolumePool(0)
 {
   myMeshId = _meshList.size();         // --- index of the mesh to push back in the vector
@@ -796,11 +796,11 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
     vector<vtkIdType> nodeIds;
     nodeIds.clear();
     nodeIds.push_back(n1->getId());
-    nodeIds.push_back(n3->getId());
     nodeIds.push_back(n2->getId());
+    nodeIds.push_back(n3->getId());
     nodeIds.push_back(n4->getId());
-    nodeIds.push_back(n6->getId());
     nodeIds.push_back(n5->getId());
+    nodeIds.push_back(n6->getId());
 
     SMDS_VtkVolume *volvtk = myVolumePool->getNew();
     volvtk->init(nodeIds, this);
@@ -1979,6 +1979,7 @@ const SMDS_MeshElement* SMDS_Mesh::FindElement(int IDelem) const
   if ((IDelem < 0) || IDelem >= myCells.size())
   {
     MESSAGE("----------------------------------- bad IDelem " << IDelem << " " << myCells.size());
+    assert(0);
     return 0;
   }
   return myCells[IDelem];
@@ -2687,7 +2688,7 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
                               list<const SMDS_MeshElement *>& removedNodes,
                               bool                            removenodes)
 {
-  MESSAGE("SMDS_Mesh::RemoveElement " << elem->GetID() << " " << removenodes);
+  //MESSAGE("SMDS_Mesh::RemoveElement " << elem->GetID() << " " << removenodes);
   // get finite elements built on elem
   set<const SMDS_MeshElement*> * s1;
   if (elem->GetType() == SMDSAbs_0DElement ||
@@ -3530,24 +3531,24 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
   vector<vtkIdType> nodeIds;
   nodeIds.clear();
   nodeIds.push_back(n1->getId());
-  nodeIds.push_back(n3->getId());
   nodeIds.push_back(n2->getId());
+  nodeIds.push_back(n3->getId());
 
   nodeIds.push_back(n4->getId());
-  nodeIds.push_back(n6->getId());
   nodeIds.push_back(n5->getId());
+  nodeIds.push_back(n6->getId());
 
-  nodeIds.push_back(n31->getId());
-  nodeIds.push_back(n23->getId());
   nodeIds.push_back(n12->getId());
+  nodeIds.push_back(n23->getId());
+  nodeIds.push_back(n31->getId());
 
-  nodeIds.push_back(n64->getId());
-  nodeIds.push_back(n56->getId());
   nodeIds.push_back(n45->getId());
+  nodeIds.push_back(n56->getId());
+  nodeIds.push_back(n64->getId());
 
   nodeIds.push_back(n14->getId());
-  nodeIds.push_back(n36->getId());
   nodeIds.push_back(n25->getId());
+  nodeIds.push_back(n36->getId());
 
   SMDS_VtkVolume *volvtk = myVolumePool->getNew();
   volvtk->init(nodeIds, this);
@@ -3747,14 +3748,14 @@ void SMDS_Mesh::adjustStructure()
 void SMDS_Mesh::dumpGrid(string ficdump)
 {
        MESSAGE("SMDS_Mesh::dumpGrid " << ficdump);
-  vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
-  aWriter->SetFileName(ficdump.c_str());
-  aWriter->SetInput(myGrid);
-  if(myGrid->GetNumberOfCells())
-  {
-    aWriter->Write();
-  }
-  aWriter->Delete();
+//  vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
+//  aWriter->SetFileName(ficdump.c_str());
+//  aWriter->SetInput(myGrid);
+//  if(myGrid->GetNumberOfCells())
+//  {
+//    aWriter->Write();
+//  }
+//  aWriter->Delete();
   ficdump = ficdump + "_connectivity";
   ofstream ficcon(ficdump.c_str(), ios::out);
   int nbPoints = myGrid->GetNumberOfPoints();
@@ -3767,6 +3768,8 @@ void SMDS_Mesh::dumpGrid(string ficdump)
   ficcon << "-------------------------------- cells " <<  nbCells << endl;
   for (int i=0; i<nbCells; i++)
   {
+//     MESSAGE(i << " " << myGrid->GetCell(i));
+//     MESSAGE("  " << myGrid->GetCell(i)->GetCellType());
        ficcon << i << " - " << myGrid->GetCell(i)->GetCellType() << " -";
        int nbptcell = myGrid->GetCell(i)->GetNumberOfPoints();
        vtkIdList *listid = myGrid->GetCell(i)->GetPointIds();
index f32a366fecf8143e8052735815c675edd0c6d668..357ee8eebcd792389336e5a9acd083c17a7b37e5 100644 (file)
@@ -43,6 +43,7 @@
 #include "SMDS_VtkFace.hxx"
 #include "SMDS_VtkVolume.hxx"
 #include "ObjectPool.hxx"
+#include "SMDS_UnstructuredGrid.hxx"
 
 #include <boost/shared_ptr.hpp>
 #include <set>
@@ -53,8 +54,6 @@
 #include "Utils_SALOME_Exception.hxx"
 #define MYASSERT(val) if (!(val)) throw SALOME_Exception(LOCALIZED("assertion not verified"));
 
-class vtkUnstructuredGrid;
-
 class SMDS_EXPORT SMDS_Mesh:public SMDS_MeshObject{
 public:
   friend class SMDS_MeshElementIDFactory;
@@ -570,12 +569,10 @@ public:
   void incrementCellsCapacity(int nbCells);
   void adjustStructure();
   void dumpGrid(string ficdump="dumpGrid");
-  
-  int myCellLinksSize;
 
   static int chunkSize;
 
-private:
+protected:
   SMDS_Mesh(SMDS_Mesh * parent);
 
   SMDS_MeshFace * createTriangle(const SMDS_MeshNode * node1,
@@ -617,7 +614,7 @@ private:
   int myMeshId;
 
   //! actual nodes coordinates, cells definition and reverse connectivity are stored in a vtkUnstructuredGrid
-  vtkUnstructuredGrid*      myGrid;
+  SMDS_UnstructuredGrid*      myGrid;
 
   //! Small objects like SMDS_MeshNode are allocated by chunks to limit memory costs of new
   ObjectPool<SMDS_MeshNode>* myNodePool;
index 9fea64ceafcba2c3e282a5c80dd210a95783f51f..939eabd5cadb1df58072fb89a7755872b84f5edc 100644 (file)
@@ -81,6 +81,7 @@ public:
   friend SMDS_EXPORT std::ostream & operator <<(std::ostream & OS, const SMDS_MeshElement *);
   friend SMDS_EXPORT bool SMDS_MeshElementIDFactory::BindID(int ID,SMDS_MeshElement* elem);
   friend class SMDS_Mesh;
+  friend class SMESHDS_Mesh;
 
   // ===========================
   //  Access to nodes by index
@@ -134,6 +135,7 @@ public:
   inline void setIdInShape(int id) { myIdInShape = id; };
 
 protected:
+  inline void setId(int id) {myID = id; };
   SMDS_MeshElement(int ID=-1);
   SMDS_MeshElement(int id, ShortType meshId, ShortType shapeId=-1);
   virtual void Print(std::ostream & OS) const;
index 7465a1a8d2addd358e85d5aec5fddd5e7d0d621f..8bde46aa5fba2b0c9a5063292d0b1dfae0bb9653 100644 (file)
@@ -90,16 +90,23 @@ void SMDS_MeshIDFactory::ReleaseID(const int ID)
 
 void SMDS_MeshIDFactory::Clear()
 {
-  myMaxID = -1;
-  myPoolOfID.clear();
+       myMaxID = -1;
+       myPoolOfID.clear();
 }
 
-  void SMDS_MeshIDFactory::SetMesh(SMDS_Mesh *mesh)
-  {
-      myMesh = mesh;
-  }
+void SMDS_MeshIDFactory::SetMesh(SMDS_Mesh *mesh)
+{
+       myMesh = mesh;
+}
 
-  SMDS_Mesh* SMDS_MeshIDFactory::GetMesh()
-  {
-      return myMesh;
-  }
+SMDS_Mesh* SMDS_MeshIDFactory::GetMesh()
+{
+       return myMesh;
+}
+
+void SMDS_MeshIDFactory::emptyPool(int maxId)
+{
+       MESSAGE("SMDS_MeshIDFactory::emptyPool " << myMaxID << " --> " << maxId-1);
+       myMaxID = maxId-1;
+       myPoolOfID.clear();
+}
index 09ed074b91c45ae1d5743abec785c25461d0ac80..9e433c9518684e03b4508a57a29d59e2cb1419f6 100644 (file)
@@ -42,6 +42,8 @@ public:
 
   void SetMesh(SMDS_Mesh *mesh);
   SMDS_Mesh* GetMesh();
+  inline bool isPoolIdEmpty() { return myPoolOfID.empty(); };
+  void emptyPool(int maxId);
 protected:
   SMDS_MeshIDFactory();
   int myMaxID;
index 888d8967bd51323a6340bac216866fe3f596158f..f1c100046be18ea73676daea1ece4c404c82f953 100644 (file)
@@ -71,42 +71,10 @@ void SMDS_MeshNode::init(int id, int meshId, int shapeId, double x, double y, do
   SMDS_Mesh* mesh = SMDS_Mesh::_meshList[myMeshId];
   vtkUnstructuredGrid * grid = mesh->getGrid();
   vtkPoints *points = grid->GetPoints();
-  //int nbp = points->GetNumberOfPoints();
   points->InsertPoint(myID, x, y, z);
-  if (myID >= mesh->myCellLinksSize)
-  {
-      //points->SetNumberOfPoints(myID+SMDS_Mesh::chunkSize);
-      vtkCellLinks *cellLinks = grid->GetCellLinks();
-
-//      int imax = cellLinks->Size;
-//      for (int i =0; i<imax; i++)
-//      {
-//        vtkCellLinks::Link &ilink = cellLinks->GetLink(i);
-//        int ncells = ilink.ncells;
-//        int *cells = ilink.cells;
-//        MESSAGE("NODE " << i << " " << cellLinks << " " << cells << " " << ncells);
-//        for (int j=0; j< ncells; j++)
-//          MESSAGE("             " << j << " " << cells[j]);
-//      }
-
-      cellLinks->Resize(myID+SMDS_Mesh::chunkSize);
-
-//      cellLinks = grid->GetCellLinks();
-//      imax = cellLinks->Size;
-//      for (int i =0; i<imax; i++)
-//      {
-//        vtkCellLinks::Link &ilink = cellLinks->GetLink(i);
-//        int ncells = ilink.ncells;
-//        int *cells = ilink.cells;
-//        MESSAGE("NODE " << i << " " << cellLinks << " " << cells << " " << ncells);
-//        for (int j=0; j< ncells; j++)
-//          MESSAGE("             " << j << " " << cells[j]);
-//      }
-
-      mesh->myCellLinksSize = cellLinks->Size;
-      //MESSAGE(" -------------------------------------- resize CellLinks " << myID << " --> " << mesh->myCellLinksSize);
-  }
-  //setXYZ(x, y, z);
+  vtkCellLinks *cellLinks = grid->GetCellLinks();
+  if (myID >=cellLinks->Size)
+         cellLinks->Resize(myID+SMDS_Mesh::chunkSize);
 }
 
 SMDS_MeshNode::~SMDS_MeshNode()
index 1c219a27d61e1a077bb6941aeae9b65320462411..7fe994f2108b438560de2eebeb76e82b139f7cf7 100644 (file)
@@ -3,75 +3,18 @@
 #include "SMDS_UnstructuredGrid.hxx"
 #include "utilities.h"
 
-using namespace std;
+#include <vtkCellArray.h>
+#include <vtkIdTypeArray.h>
+#include <vtkUnsignedCharArray.h>
+
+#include <list>
 
-//vtkCellLinks::Link* SMDS_CellLinks::AdjustSize(vtkIdType sz)
-//{
-//  vtkIdType i;
-//  vtkCellLinks::Link *newArray;
-//  vtkIdType newSize = sz;
-//  vtkCellLinks::Link linkInit = {0,NULL};
-//
-//  newArray = new vtkCellLinks::Link[newSize];
-//
-//  for (i=0; i<sz && i<this->Size; i++)
-//    {
-//    newArray[i] = this->Array[i];
-//    }
-//
-//  for (i=this->Size; i < newSize ; i++)
-//    {
-//    newArray[i] = linkInit;
-//    }
-//
-//  this->Size = newSize;
-//  delete [] this->Array;
-//  this->Array = newArray;
-//
-//  return this->Array;
-//}
-//
-//SMDS_CellLinks* SMDS_CellLinks::New()
-//{
-//  return new SMDS_CellLinks();
-//}
-//
-//SMDS_CellLinks::SMDS_CellLinks() : vtkCellLinks()
-//{
-//}
-//
-//SMDS_CellLinks::~SMDS_CellLinks()
-//{
-//}
-//
-//
-///*! initialize an SMDS_CellLinks instance instead of a vtkCellLinks instance
-// *
-// */
-//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());
-//  this->Links->Register(this);
-//  this->Links->BuildLinks(this, this->Connectivity);
-//  this->Links->Delete();
-//}
-//
-//SMDS_CellLinks* SMDS_UnstructuredGrid::GetCellLinks()
-//{
-//  return static_cast<SMDS_CellLinks*>(this->Links);
-//}
+using namespace std;
 
 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
 {
        MESSAGE("SMDS_UnstructuredGrid::New");
-  return new SMDS_UnstructuredGrid();
+       return new SMDS_UnstructuredGrid();
 }
 
 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() : vtkUnstructuredGrid()
@@ -95,3 +38,189 @@ void SMDS_UnstructuredGrid::UpdateInformation()
        MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
        return vtkUnstructuredGrid::UpdateInformation();
 }
+
+void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
+                                                                               std::vector<int>& idCellsOldToNew, int newCellSize)
+{
+       MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);
+
+       int startHole = 0;
+       int endHole = 0;
+       int startBloc = 0;
+       int endBloc = 0;
+       int alreadyCopied = 0;
+       typedef enum {lookHoleStart, lookHoleEnd, lookBlocEnd} enumState;
+       enumState compactState = lookHoleStart;
+
+       // --- if newNodeSize, create a new compacted vtkPoints
+
+       vtkPoints *newPoints = 0;
+       if (newNodeSize)
+       {
+               MESSAGE("-------------- compactGrid, newNodeSize");
+               newPoints = vtkPoints::New();
+               newPoints->Initialize();
+               newPoints->Allocate(newNodeSize);
+               newPoints->SetNumberOfPoints(newNodeSize);
+               int oldNodeSize = idNodesOldToNew.size();
+
+               for (int i=0; i< oldNodeSize; i++)
+               {
+                       switch(compactState)
+                       {
+                       case lookHoleStart:
+                               if (idNodesOldToNew[i] < 0)
+                               {
+                                       startHole = i;
+                                       compactState = lookHoleEnd;
+                               }
+                               break;
+                       case lookHoleEnd:
+                               if (idNodesOldToNew[i] >= 0)
+                               {
+                                       endHole = i;
+                                       startBloc = i;
+                                       compactState = lookBlocEnd;
+                               }
+                               break;
+                       case lookBlocEnd:
+                               if (idNodesOldToNew[i] < 0) endBloc = i; // see nbPoints below
+                               else if (i == (oldNodeSize-1)) endBloc = i+1;
+                               if (endBloc)
+                               {
+                                       MESSAGE("-------------- newNodeSize, endbloc");
+                                       void *target = newPoints->GetVoidPointer(3*alreadyCopied);
+                                       void *source = this->Points->GetVoidPointer(3*startBloc);
+                                       int nbPoints = endBloc - startBloc;
+                                       memcpy(target, source, 3*sizeof(float)*nbPoints);
+                                       for (int j=startBloc; j<endBloc; j++)
+                                               idNodesOldToNew[j] = alreadyCopied++;
+                                       compactState = lookHoleStart;
+                                       startHole = i;
+                                       endHole = 0;
+                                       startBloc = 0;
+                                       endBloc = 0;
+                               }
+                               break;
+                       }
+               }
+               if (!alreadyCopied) // no hole, but shorter, no need to modify idNodesOldToNew
+               {
+                       MESSAGE("------------- newNodeSize, shorter")
+                       void *target = newPoints->GetVoidPointer(0);
+                       void *source = this->Points->GetVoidPointer(0);
+                       int nbPoints = newNodeSize;
+                       memcpy(target, source, 3*sizeof(float)*nbPoints);
+               }
+       }
+
+       // --- create new compacted Connectivity, Locations and Types
+
+    int oldCellSize = this->Types->GetNumberOfTuples();
+
+       vtkCellArray *newConnectivity = vtkCellArray::New();
+       newConnectivity->Initialize();
+       int oldCellDataSize = this->Connectivity->GetData()->GetSize();
+       newConnectivity->Allocate(oldCellDataSize);
+       MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
+
+       vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
+       newTypes->Initialize();
+       //newTypes->Allocate(oldCellSize);
+       newTypes->SetNumberOfValues(newCellSize);
+
+       vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
+       newLocations->Initialize();
+       //newLocations->Allocate(oldCellSize);
+       newLocations->SetNumberOfValues(newCellSize);
+
+       startHole = 0;
+       endHole = 0;
+       startBloc = 0;
+       endBloc = 0;
+       alreadyCopied = 0;
+       compactState = lookHoleStart;
+
+       vtkIdType tmpid[50];
+       vtkIdType *pointsCell =&tmpid[0]; // --- points id to fill a new cell
+
+    for (int i=0; i<oldCellSize; i++)
+    {
+               switch(compactState)
+               {
+               case lookHoleStart:
+                       if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
+                       {
+                               startHole = i;
+                               compactState = lookHoleEnd;
+                       }
+                       break;
+               case lookHoleEnd:
+                       if (this->Types->GetValue(i) != VTK_EMPTY_CELL)
+                       {
+                               endHole = i;
+                               startBloc = i;
+                               compactState = lookBlocEnd;
+                       }
+                       break;
+               case lookBlocEnd:
+                       if (this->Types->GetValue(i) == VTK_EMPTY_CELL) endBloc =i;
+                       else if (i == (oldCellSize-1)) endBloc = i+1;
+                       if (endBloc)
+                       {
+                       MESSAGE(" -------- newCellSize, endBloc");
+                               for (int j=startBloc; j<endBloc; j++)
+                               {
+                                       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
+                                       idCellsOldToNew[j] = alreadyCopied;
+                                       int oldLoc = this->Locations->GetValue(j);
+                                       int nbpts;
+                                       this->Connectivity->GetCell(oldLoc, nbpts, pointsCell);
+                                       for (int l=0; l<nbpts; l++)
+                                               pointsCell[l] = idNodesOldToNew[pointsCell[l]];
+                                       int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
+                                       int newLoc = newConnectivity->GetInsertLocation(nbpts);
+                                       newLocations->SetValue(alreadyCopied, newLoc);
+                                       alreadyCopied++;
+                               }
+                       }
+                       break;
+               }
+    }
+    if (!alreadyCopied) // no hole, but shorter
+    {
+       MESSAGE(" -------- newCellSize, shorter");
+               for (int j=0; j<oldCellSize; j++)
+               {
+                       newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
+                       idCellsOldToNew[j] = alreadyCopied;
+                       int oldLoc = this->Locations->GetValue(j);
+                       int nbpts;
+                       this->Connectivity->GetCell(oldLoc, nbpts, pointsCell);
+                       //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
+                       for (int l=0; l<nbpts; l++)
+                       {
+                               int oldval = pointsCell[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++;
+               }
+    }
+
+    newConnectivity->Squeeze();
+    //newTypes->Squeeze();
+    //newLocations->Squeeze();
+
+    if (newNodeSize)
+    {
+       MESSAGE("------- newNodeSize, setPoints");
+       this->SetPoints(newPoints);
+    }
+    this->SetCells(newTypes, newLocations, newConnectivity);
+    this->BuildLinks();
+}
index ce4fed74c0af77162fb38a086cdd5cb95c182a39..e2bc01d0d0c16b8cba92b0c940ad0ec5dfde297b 100644 (file)
@@ -8,39 +8,24 @@
 #ifndef _SMDS_UNSTRUCTUREDGRID_HXX
 #define        _SMDS_UNSTRUCTUREDGRID_HXX
 
+#include <vector>
+
 #include <vtkUnstructuredGrid.h>
-//#include <vtkCellLinks.h>
-
-//class SMDS_CellLinks: public vtkCellLinks
-//{
-//public:
-//    Link *AdjustSize(vtkIdType sz);
-//
-//    //virtual void Delete();
-//    static SMDS_CellLinks* New();
-//protected:
-//    SMDS_CellLinks();
-//    ~SMDS_CellLinks();
-//};
 
 class SMDS_UnstructuredGrid: public vtkUnstructuredGrid
 {
 public:
-//    void BuildLinks(); // initialise un SMDS_CellLinks;
-//    SMDS_CellLinks* GetCellLinks();
-//
-//    vtkIdType GetCellArraySize() { return (this->Connectivity ? this->Connectivity->GetSize() : 0); };
+       void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
+                                        std::vector<int>& idCellsOldToNew, int newCellSize);
 
        virtual unsigned long GetMTime();
        virtual void UpdateInformation();
 
-    //virtual void Delete();
     static SMDS_UnstructuredGrid* New();
 protected:
     SMDS_UnstructuredGrid();
     ~SMDS_UnstructuredGrid();
 };
 
-
 #endif /* _SMDS_UNSTRUCTUREDGRID_HXX */
 
index 8d6fe266ff64c242dea9026478bb483b465e7f34..1d0138f1be067e767fd7c908a3b2061ce38e50df 100644 (file)
@@ -22,8 +22,8 @@ SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId, SMDSA
     }
     case SMDSEntity_Penta:
     {
-      this->exchange(1, 2);
-      this->exchange(4, 5);
+      //this->exchange(1, 2);
+      //this->exchange(4, 5);
       break;
     }
     case SMDSEntity_Hexa:
@@ -49,11 +49,11 @@ SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId, SMDSA
     }
     case SMDSEntity_Quad_Penta:
     {
-      this->exchange(1, 2);
-      this->exchange(4, 5);
-      this->exchange(6, 8);
-      this->exchange(9, 11);
-      this->exchange(13, 14);
+      //this->exchange(1, 2);
+      //this->exchange(4, 5);
+      //this->exchange(6, 8);
+      //this->exchange(9, 11);
+      //this->exchange(13, 14);
       break;
     }
     case SMDSEntity_Quad_Hexa:
index 0cbac04e861b3810b88992836397e468be0e04dc..a86d65b3538add601f70c9a96f776ef887b37806 100644 (file)
@@ -296,6 +296,8 @@ bool SMESH_Gen::Compute(SMESH_Mesh &          aMesh,
 
   SMESHDS_Mesh *myMesh = aMesh.GetMeshDS();
   myMesh->adjustStructure();
+  myMesh->compactMesh();
+  //myMesh->adjustStructure();
   list<int> listind = myMesh->SubMeshIndices();
   list<int>::iterator it = listind.begin();
   int total = 0;
@@ -304,9 +306,9 @@ bool SMESH_Gen::Compute(SMESH_Mesh &          aMesh,
       ::SMESHDS_SubMesh *subMesh = myMesh->MeshElements(*it);
       total +=  subMesh->getSize();
     }
-  cerr << "total elements and nodes in submesh sets:" << total << endl;
-  cerr << "Number of node objects " << SMDS_MeshNode::nbNodes << endl;
-  cerr << "Number of cell objects " << SMDS_MeshCell::nbCells << endl;
+  MESSAGE("total elements and nodes in submesh sets:" << total);
+  MESSAGE("Number of node objects " << SMDS_MeshNode::nbNodes);
+  MESSAGE("Number of cell objects " << SMDS_MeshCell::nbCells);
   //myMesh->dumpGrid();
   return ret;
 }
index 34118b82907afab02291f4d9651c74a35fa22e27..d3bcea07b7ffde59e53674165f04c721343fbd80 100644 (file)
@@ -1775,4 +1775,92 @@ SMDS_MeshVolume* SMESHDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
                          ID);
 }
 
+void SMESHDS_Mesh::compactMesh()
+{
+       int newNodeSize = 0;
+       int nbNodes = myNodes.size();
+       vector<int> idNodesOldToNew;
+       idNodesOldToNew.clear();
+       idNodesOldToNew.resize(nbNodes, -1); // all unused id will be -1
+
+       bool areNodesModified = ! myNodeIDFactory->isPoolIdEmpty();
+       MESSAGE("------------------------------------------------- SMESHDS_Mesh::compactMesh " << areNodesModified);
+       if (areNodesModified)
+       {
+               for (int i=0; i<nbNodes; i++)
+               {
+                       if (myNodes[i])
+                               {
+                                       idNodesOldToNew[i] = i;  // all valid id are >= 0
+                                       newNodeSize++;
+                               }
+               }
+       }
+       else
+       {
+               for (int i=0; i<nbNodes; i++)
+                       idNodesOldToNew[i] = i;
+       }
+
+       int newCellSize = 0;
+       int nbCells = myCells.size();
+       vector<int> idCellsOldToNew;
+       idCellsOldToNew.clear();
+       idCellsOldToNew.resize(nbCells, -1);             // all unused id will be -1
+
+       for (int i=0; i<nbCells; i++)
+       {
+               if (myCells[i])
+                       {
+                               idCellsOldToNew[i] = myVtkIndex[i];  // valid vtk indexes are > = 0
+                               newCellSize++;
+                       }
+       }
+       myGrid->compactGrid(idNodesOldToNew, newNodeSize, idCellsOldToNew, newCellSize);
+
+       // --- SMDS_MeshNode and myNodes (id in SMDS and in VTK are the same), myNodeIdFactory
+
+       if (areNodesModified)
+       {
+               MESSAGE("-------------- modify myNodes");
+               for (int i=0; i<nbNodes; i++)
+               {
+                       if (myNodes[i])
+                       {
+                               int newid = idNodesOldToNew[i];
+                               if (newid != i)
+                               {
+                                       MESSAGE(i << " --> " << newid);
+                                       myNodes[i]->setId(newid);
+                                       ASSERT(!myNodes[newid]);
+                                       myNodes[newid] = myNodes[i];
+                               }
+                       }
+               }
+               this->myNodeIDFactory->emptyPool(newNodeSize);
+       }
+
+       // --- SMDS_MeshCell, myIDElements and myVtkIndex (myCells and myElementIdFactory are not compacted)
+
+       for (int oldVtkId=0; oldVtkId<nbCells; oldVtkId++)
+       {
+               int smdsId = this->myVtkIndex[oldVtkId];
+               if (smdsId >=0)
+               {
+                       int newVtkId = idCellsOldToNew[oldVtkId];
+                       myCells[smdsId]->setVtkId(newVtkId);
+                       myIDElements[smdsId] = newVtkId;
+                       myVtkIndex[newVtkId] = smdsId;
+               }
+       }
+
+       // ---TODO: myNodes, myElements in submeshes
+
+//    map<int,SMESHDS_SubMesh*>::iterator it = myShapeIndexToSubMesh.begin();
+//    for(; it != myShapeIndexToSubMesh.end(); ++it)
+//    {
+//     (*it).second->compactList(idNodesOldToNew, newNodeSize, idCellsOldToNew, newCellSize);
+//    }
+
+}
 
index 116deb607081194435817d4c2b8cfb289ac614f0..fa8d97d2040802d2f30e38555e7960bb6dafc755 100644 (file)
@@ -438,6 +438,8 @@ public:
 
   bool IsGroupOfSubShapes (const TopoDS_Shape& aSubShape) const;
 
+  void compactMesh();
+
   ~SMESHDS_Mesh();
   
 private:
index f06398de151fbfb4df7115723400c92ab44afd3b..02141893d96c91b030cf51ccebda7b5c0da86860 100644 (file)
@@ -412,8 +412,12 @@ int SMESHDS_SubMesh::getSize()
 {
   int c = NbNodes();
   int d = NbElements();
-  cerr << "SMESHDS_SubMesh::NbNodes " << c << endl;
-  cerr << "SMESHDS_SubMesh::NbElements " << d << endl;
+  //cerr << "SMESHDS_SubMesh::NbNodes " << c << endl;
+  //cerr << "SMESHDS_SubMesh::NbElements " << d << endl;
   return c+d;
 }
 
+void SMESHDS_SubMesh::compactList()
+{
+       // todo : compact vector of nodes and elements
+}
index ca5e3098d0375ff380914763f25fdfdbc53d4c0a..84b93d62d1ecb4355301de0d9643ce09e0cac5e5 100644 (file)
@@ -66,6 +66,7 @@ class SMESHDS_EXPORT SMESHDS_SubMesh
   // clear the contents
   void Clear();
   int getSize();
+  void compactList();
 
  private: