]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
PR: some debug with SMESH_TEST grid from branch V6_1_BR
authorprascle <prascle>
Fri, 10 Sep 2010 14:58:23 +0000 (14:58 +0000)
committerprascle <prascle>
Fri, 10 Sep 2010 14:58:23 +0000 (14:58 +0000)
14 files changed:
src/OBJECT/SMESH_Actor.cxx
src/OBJECT/SMESH_ActorUtils.cxx
src/OBJECT/SMESH_Object.cxx
src/SMDS/SMDS_Downward.cxx
src/SMDS/SMDS_Mesh.cxx
src/SMDS/SMDS_MeshElementIDFactory.cxx
src/SMDS/SMDS_MeshNode.cxx
src/SMDS/SMDS_UnstructuredGrid.cxx
src/SMDS/SMDS_UnstructuredGrid.hxx
src/SMDS/SMDS_VtkCellIterator.cxx
src/SMDS/SMDS_VtkVolume.cxx
src/SMDS/SMDS_VtkVolume.hxx
src/SMESH/SMESH_MeshEditor.cxx
src/SMESH_I/SMESH_MeshEditor_i.cxx

index cd9e0691967992a98cc2b113fd863498a0c67211..9424ac0dfca3da477d02ca3f2a78b7a487fec129 100644 (file)
@@ -30,6 +30,7 @@
 #include "SMESH_DeviceActor.h"
 #include "SMESH_ObjectDef.h"
 #include "SMESH_ControlsDef.hxx"
+#include "SMDS_UnstructuredGrid.hxx"
 #include "VTKViewer_ExtractUnstructuredGrid.h"
 #include "VTKViewer_FramedTextActor.h"
 #include "SALOME_InteractiveObject.hxx"
@@ -82,7 +83,7 @@
 #ifdef _DEBUG_
 static int MYDEBUG = 1;
 #else
-static int MYDEBUG = 0;
+static int MYDEBUG = 1;
 #endif
 
 static int aLineWidthInc = 2;
@@ -197,6 +198,10 @@ SMESH_ActorDef::SMESH_ActorDef()
   aFilter->RegisterCellsWithType(VTK_QUADRATIC_HEXAHEDRON);
   aFilter->RegisterCellsWithType(VTK_QUADRATIC_WEDGE);
   aFilter->RegisterCellsWithType(VTK_CONVEX_POINT_SET);
+#ifdef VTK_HAVE_POLYHEDRON
+  MESSAGE("RegisterCellsWithType(VTK_POLYHEDRON)");
+  aFilter->RegisterCellsWithType(VTK_POLYHEDRON);
+#endif
 
   //Definition 1D device of the actor
   //---------------------------------
@@ -1335,6 +1340,9 @@ void SMESH_ActorDef::SetEntityMode(unsigned int theMode)
     aFilter->RegisterCellsWithType(VTK_QUADRATIC_HEXAHEDRON);
     aFilter->RegisterCellsWithType(VTK_QUADRATIC_WEDGE);
     aFilter->RegisterCellsWithType(VTK_CONVEX_POINT_SET);
+#ifdef VTK_HAVE_POLYHEDRON
+    aFilter->RegisterCellsWithType(VTK_POLYHEDRON);
+#endif
     
     aHightFilter->RegisterCellsWithType(VTK_TETRA);
     aHightFilter->RegisterCellsWithType(VTK_VOXEL);
@@ -1345,6 +1353,9 @@ void SMESH_ActorDef::SetEntityMode(unsigned int theMode)
     aHightFilter->RegisterCellsWithType(VTK_QUADRATIC_HEXAHEDRON);
     aHightFilter->RegisterCellsWithType(VTK_QUADRATIC_WEDGE);
     aHightFilter->RegisterCellsWithType(VTK_CONVEX_POINT_SET);
+#ifdef VTK_HAVE_POLYHEDRON
+    aHightFilter->RegisterCellsWithType(VTK_POLYHEDRON);
+#endif
   }
   aFilter->Update();
   if (MYDEBUG) MESSAGE(aFilter->GetOutput()->GetNumberOfCells());
index 6808e145e1c668b6aecff46a237e1835cf05fd7b..899e413ec23266f1896593cd895ae906e5ebde98 100644 (file)
@@ -29,6 +29,7 @@
 #include "utilities.h"
 
 #include <vtkUnstructuredGrid.h>
+#include <vtkXMLUnstructuredGridWriter.h>
 #include <vtkUnstructuredGridWriter.h>
 
 #ifdef _DEBUG_
@@ -73,9 +74,10 @@ namespace SMESH
   WriteUnstructuredGrid(vtkUnstructuredGrid* theGrid, 
                         const char* theFileName)
   {
-    vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
+    vtkXMLUnstructuredGridWriter* aWriter = vtkXMLUnstructuredGridWriter::New();
     aWriter->SetFileName(theFileName);
     aWriter->SetInput(theGrid);
+    aWriter->SetDataModeToAscii();
     if(theGrid->GetNumberOfCells()){
       aWriter->Write();
     }
index 169d2ed07b0648e36a9dfeea24bdf9627b709d37..a2141b54e14841d33a9ebc63e65999a7d838a6af 100644 (file)
@@ -65,7 +65,7 @@ using namespace std;
 
 #ifdef _DEBUG_
 static int MYDEBUG = 1;
-static int MYDEBUGWITHFILES = 0;
+static int MYDEBUGWITHFILES = 1;
 #else
 static int MYDEBUG = 0;
 static int MYDEBUGWITHFILES = 0;
@@ -268,7 +268,7 @@ void SMESH_VisualObjDef::buildPrs(bool buildGrid)
        //MESSAGE(myGrid->GetReferenceCount());
        //MESSAGE( "Update - myGrid->GetNumberOfCells() = "<<myGrid->GetNumberOfCells() );
        //MESSAGE( "Update - myGrid->GetNumberOfPoints() = "<<myGrid->GetNumberOfPoints() );
-       if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"/tmp/buildPrs" );
+       if( MYDEBUGWITHFILES ) SMESH::WriteUnstructuredGrid( myGrid,"buildPrs.vtu" );
   }
 }
 
index 6875aebd5bf5e3404cb75c8bfa1287d7b6aeea08..f735c820aba60c33f90d211572fe426eac10d06a 100644 (file)
@@ -46,7 +46,7 @@ SMDS_Downward::SMDS_Downward(SMDS_UnstructuredGrid *grid, int nbDownCells) :
   this->_cellTypes.clear();
   if (_cellDimension.empty())
     {
-      _cellDimension.resize(VTK_QUADRATIC_PYRAMID + 1, 0);
+      _cellDimension.resize(VTK_MAXTYPE + 1, 0);
       _cellDimension[VTK_LINE] = 1;
       _cellDimension[VTK_QUADRATIC_EDGE] = 1;
       _cellDimension[VTK_TRIANGLE] = 2;
index 79f41cbbcaca9109e6aecc1856eae2e1c8174fc0..ce6145a50414fefb0e09488023d498d5ab11bf63 100644 (file)
@@ -1172,7 +1172,7 @@ SMDS_MeshVolume * SMDS_Mesh::AddPolyhedralVolumeWithID
 SMDS_MeshVolume* SMDS_Mesh::AddPolyhedralVolumeWithID
                             (vector<const SMDS_MeshNode*> nodes,
                              vector<int>                  quantities,
-                             const int                         ID)
+                             const int                    ID)
 {
   SMDS_MeshVolume* volume;
   //if ( myVolumes.Extent() % CHECKMEMORY_INTERVAL == 0 ) CheckMemory();
@@ -1183,9 +1183,21 @@ SMDS_MeshVolume* SMDS_Mesh::AddPolyhedralVolumeWithID
     MESSAGE("Error : Not implemented");
     return NULL;
   } else {
+#ifdef VTK_HAVE_POLYHEDRON
+    vector<vtkIdType> nodeIds;
+    nodeIds.clear();
+    vector<const SMDS_MeshNode*>::iterator it = nodes.begin();
+    for ( ; it != nodes.end(); ++it)
+      nodeIds.push_back((*it)->getId());
+
+    SMDS_VtkVolume *volvtk = myVolumePool->getNew();
+    volvtk->initPoly(nodeIds, quantities, this);
+    volume = volvtk;
+#else
     for ( int i = 0; i < nodes.size(); ++i )
       if ( !nodes[ i ] ) return 0;
     volume = new SMDS_PolyhedralVolumeOfNodes(nodes, quantities);
+#endif
     adjustmyCellsCapacity(ID);
     myCells[ID] = volume;
     myInfo.myNbPolyhedrons++;
@@ -1261,7 +1273,9 @@ const SMDS_MeshNode * SMDS_Mesh::FindNode(int ID) const
 {
   if (ID < 0 || ID >= myNodes.size())
   {
+    MESSAGE("------------------------------------------------------------------------- ");
     MESSAGE("----------------------------------- bad ID " << ID << " " << myNodes.size());
+    MESSAGE("------------------------------------------------------------------------- ");
     return 0;
   }
   return (const SMDS_MeshNode *)myNodes[ID];
@@ -1986,8 +2000,11 @@ const SMDS_MeshElement* SMDS_Mesh::FindElement(int IDelem) const
 {
   if ((IDelem < 0) || IDelem >= myCells.size())
   {
+    MESSAGE("--------------------------------------------------------------------------------- ");
     MESSAGE("----------------------------------- bad IDelem " << IDelem << " " << myCells.size());
-    assert(0);
+    MESSAGE("--------------------------------------------------------------------------------- ");
+    // TODO raise an exception
+    //assert(0);
     return 0;
   }
   return myCells[IDelem];
index cee7145626d7f9727a719c647bb48716dfe8d92e..901e6a0445a5a5604aa498f5d671fe1afd1171a4 100644 (file)
@@ -35,7 +35,7 @@
 
 #include "utilities.h"
 
-#include <vtkUnstructuredGrid.h>
+#include <SMDS_UnstructuredGrid.hxx>
 #include <vtkCellType.h>
 
 using namespace std;
@@ -69,7 +69,11 @@ SMDS_MeshElementIDFactory::SMDS_MeshElementIDFactory():
     myVtkCellTypes[SMDSEntity_Quad_Hexa]       = VTK_QUADRATIC_HEXAHEDRON;
     myVtkCellTypes[SMDSEntity_Penta]           = VTK_WEDGE;
     myVtkCellTypes[SMDSEntity_Quad_Penta]      = VTK_QUADRATIC_WEDGE;
+#ifdef VTK_HAVE_POLYHEDRON
+    myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_POLYHEDRON;
+#else
     myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_CONVEX_POINT_SET;
+#endif
     myVtkCellTypes[SMDSEntity_Quad_Polyhedra]  = VTK_CONVEX_POINT_SET;
 }
 
index c774deb8c36b39c53496c5dc527f945a77db158b..227c36732ef5e86f3d742bafcea22f67e1ce3c10 100644 (file)
@@ -134,43 +134,59 @@ const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
  */
 //=======================================================================
 
-class SMDS_MeshNode_MyInvIterator:public SMDS_ElemIterator
+class SMDS_MeshNode_MyInvIterator: public SMDS_ElemIterator
 {
 private:
   SMDS_Mesh* myMesh;
   vtkIdType* myCells;
-  int  myNcells;
-  SMDSAbs_ElementType                                 myType;
-  int  iter;
+  int myNcells;
+  SMDSAbs_ElementType myType;
+  int iter;
+  vector<vtkIdType> cellList;
 
- public:
-  SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh,
-                              vtkIdType* cells,
-                              int ncells,
-                              SMDSAbs_ElementType type):
+public:
+  SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
   {
-     //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
+    //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
+    if (type == SMDSAbs_All)
+      return;
+    cellList.clear();
+    for (int i = 0; i < ncells; i++)
+      {
+        int vtkId = cells[i];
+        int smdsId = myMesh->fromVtkToSmds(vtkId);
+        const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
+        if (elem->GetType() == type)
+          {
+            //MESSAGE("Add element vtkId " << vtkId << " " << elem->GetType())
+            cellList.push_back(vtkId);
+          }
+      }
+    myCells = &cellList[0];
+    myNcells = cellList.size();
+    //MESSAGE("myNcells="<<myNcells);
   }
 
   bool more()
   {
-      return (iter< myNcells);
+    //MESSAGE("iter " << iter << " ncells " << myNcells);
+    return (iter < myNcells);
   }
 
   const SMDS_MeshElement* next()
   {
-      int vtkId = myCells[iter];
-      int smdsId = myMesh->fromVtkToSmds(vtkId);
-      const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
-      if (!elem)
+    int vtkId = myCells[iter];
+    int smdsId = myMesh->fromVtkToSmds(vtkId);
+    const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
+    if (!elem)
       {
-          assert(0);
-          throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
+        assert(0);
+        throw SALOME_Exception("SMDS_MeshNode_MyInvIterator problem Null element");
       }
-      //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << (elem!=0));
-      iter++;
-      return elem;
+    //MESSAGE("vtkId " << vtkId << " smdsId " << smdsId << " " << elem->GetType());
+    iter++;
+    return elem;
   }
 };
 
index a543900a944722bd416f40ceba41ca945de03acd..a62267de92c5b00e51ab584bec94d086adcfdc75 100644 (file)
@@ -53,6 +53,51 @@ void 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)
+{
+  if (type != VTK_POLYHEDRON)
+    return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
+
+  // --- type = VTK_POLYHEDRON
+  MESSAGE("InsertNextLinkedCell VTK_POLYHEDRON");
+  int cellid = this->InsertNextCell(type, npts, pts);
+
+  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++)
+        {
+          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;
@@ -61,6 +106,8 @@ void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
                                         std::vector<int>& idCellsOldToNew, int newCellSize)
 {
+  // TODO utiliser mieux vtk pour faire plus simple (plus couteux ?)
+
   MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);CHRONO(1);
   int startHole = 0;
   int endHole = 0;
@@ -238,7 +285,29 @@ void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int n
       this->SetPoints(newPoints);
       MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
     }
+#ifdef VTK_HAVE_POLYHEDRON
+  // TODO compact faces for Polyhedrons
+  // refaire completement faces et faceLocation
+  // pour chaque cell, recup oldCellId, oldFacesId, recopie dans newFaces de la faceStream
+  // en changeant les numeros de noeuds
+//  vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
+//  newFaceLocations->Initialize();
+//  vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
+//  newFaces->Initialize();
+//  newFaceLocations->DeepCopy(this->FaceLocations);
+//  newFaces->DeepCopy(this->Faces);
+//  this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
+//  newFaceLocations->Delete();
+//  newFaces->Delete();
+  if (this->FaceLocations) this->FaceLocations->Register(this);
+  if (this->Faces) this->Faces->Register(this);
+  this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
+#else
   this->SetCells(newTypes, newLocations, newConnectivity);
+#endif
+  newTypes->Delete();
+  newLocations->Delete();
+  newConnectivity->Delete();
   this->BuildLinks();
 }
 
@@ -319,7 +388,7 @@ void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
 
   // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
 
-  _downArray.resize(VTK_QUADRATIC_PYRAMID + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
+  _downArray.resize(VTK_MAXTYPE + 1, 0); // --- max. type value = VTK_QUADRATIC_PYRAMID
 
   _downArray[VTK_LINE] = new SMDS_DownEdge(this);
   _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
index e8f871f17a2b44fbe4fa1613d7ce9ea77306b173..214e65310f78a5bc843d1e109a0004419507da9b 100644 (file)
@@ -8,12 +8,19 @@
 #ifndef _SMDS_UNSTRUCTUREDGRID_HXX
 #define        _SMDS_UNSTRUCTUREDGRID_HXX
 
+#include <vtkUnstructuredGrid.h>
+#include "chrono.hxx"
+
 #include <vector>
 #include <set>
 #include <map>
 
-#include <vtkUnstructuredGrid.h>
-#include "chrono.hxx"
+#define VTK_HAVE_POLYHEDRON
+#ifdef VTK_HAVE_POLYHEDRON
+  #define VTK_MAXTYPE VTK_POLYHEDRON
+#else
+  #define VTK_MAXTYPE VTK_QUADRATIC_PYRAMID
+#endif
 
 #define NBMAXNEIGHBORS 10
 
@@ -30,6 +37,11 @@ public:
   virtual unsigned long GetMTime();
   virtual void Update();
   virtual void UpdateInformation();
+  virtual vtkPoints *GetPoints();
+
+#ifdef VTK_HAVE_POLYHEDRON
+  int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
+#endif
 
   int CellIdToDownId(int vtkCellId);
   void setCellIdToDownId(int vtkCellId, int downId);
index d888a58bae699268ba3523c8a69c6e99678514df..e2188f36bb8c4ad64276fae42f61f2bbd994ce34 100644 (file)
@@ -1,7 +1,8 @@
 #include "SMDS_VtkCellIterator.hxx"
+#include "utilities.h"
 
-SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId,
-                                           SMDSAbs_EntityType aType) :
+
+SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId, SMDSAbs_EntityType aType) :
   _mesh(mesh), _cellId(vtkCellId), _index(0), _type(aType)
 {
   vtkUnstructuredGrid* grid = _mesh->getGrid();
@@ -59,6 +60,7 @@ SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId,
       }
     case SMDSEntity_Quad_Hexa:
       {
+        MESSAGE("SMDS_VtkCellIterator Quad_Hexa");
         this->exchange(1, 3);
         this->exchange(5, 7);
         this->exchange(8, 11);
@@ -68,6 +70,11 @@ SMDS_VtkCellIterator::SMDS_VtkCellIterator(SMDS_Mesh* mesh, int vtkCellId,
         this->exchange(17, 19);
         break;
       }
+    case SMDSEntity_Polyhedra:
+      {
+        MESSAGE("SMDS_VtkCellIterator Polyhedra");
+        break;
+      }
     default:
       break;
   }
index f333a8267ee795abd6a9acdc4a4b605d7079b25f..d2a15d5e70e3d754642e52194b380b5a663cc5c4 100644 (file)
@@ -58,8 +58,30 @@ void SMDS_VtkVolume::init(std::vector<vtkIdType> nodeIds, SMDS_Mesh* mesh)
   myVtkID = grid->InsertNextLinkedCell(aType, nodeIds.size(), &nodeIds[0]);
 }
 
-bool SMDS_VtkVolume::ChangeNodes(const SMDS_MeshNode* nodes[],
-                                 const int nbNodes)
+#ifdef VTK_HAVE_POLYHEDRON
+void SMDS_VtkVolume::initPoly(std::vector<vtkIdType> nodeIds, std::vector<int> nbNodesPerFace, SMDS_Mesh* mesh)
+{
+  MESSAGE("SMDS_VtkVolume::initPoly");
+  SMDS_UnstructuredGrid* grid = mesh->getGrid();
+  vector<vtkIdType> ptIds;
+  ptIds.clear();
+  vtkIdType nbFaces = nbNodesPerFace.size();
+  int k = 0;
+  for (int i = 0; i < nbFaces; i++)
+    {
+      int nf = nbNodesPerFace[i];
+      ptIds.push_back(nf);
+      for (int n = 0; n < nf; n++)
+        {
+          ptIds.push_back(nodeIds[k]);
+          k++;
+        }
+    }
+  myVtkID = grid->InsertNextLinkedCell(VTK_POLYHEDRON, nbFaces, &ptIds[0]);
+}
+#endif
+
+bool SMDS_VtkVolume::ChangeNodes(const SMDS_MeshNode* nodes[], const int nbNodes)
 {
   // TODO utilise dans SMDS_Mesh
   return true;
@@ -76,6 +98,7 @@ void SMDS_VtkVolume::Print(ostream & OS) const
 
 int SMDS_VtkVolume::NbFaces() const
 {
+  // TODO polyedres
   switch (NbNodes())
   {
     case 4:
@@ -106,6 +129,7 @@ int SMDS_VtkVolume::NbNodes() const
 
 int SMDS_VtkVolume::NbEdges() const
 {
+  // TODO polyedres
   switch (NbNodes())
   {
     case 4:
@@ -132,11 +156,7 @@ SMDS_ElemIteratorPtr SMDS_VtkVolume::elementsIterator(SMDSAbs_ElementType type)
   switch (type)
   {
     case SMDSAbs_Node:
-      return SMDS_ElemIteratorPtr(
-                                  new SMDS_VtkCellIterator(
-                                                           SMDS_Mesh::_meshList[myMeshId],
-                                                           myVtkID,
-                                                           GetEntityType()));
+      return SMDS_ElemIteratorPtr(new SMDS_VtkCellIterator(SMDS_Mesh::_meshList[myMeshId], myVtkID, GetEntityType()));
     default:
       MESSAGE("ERROR : Iterator not implemented")
       ;
@@ -161,6 +181,7 @@ const SMDS_MeshNode* SMDS_VtkVolume::GetNode(const int ind) const
 
 bool SMDS_VtkVolume::IsQuadratic() const
 {
+  // TODO polyedres
   if (this->NbNodes() > 9)
     return true;
   else
@@ -169,35 +190,43 @@ bool SMDS_VtkVolume::IsQuadratic() const
 
 SMDSAbs_EntityType SMDS_VtkVolume::GetEntityType() const
 {
+  // TODO see SMDS_MeshElementIDFactory::GetVtkCellType
+  vtkIdType aVtkType = this->GetVtkType();
+
   SMDSAbs_EntityType aType = SMDSEntity_Tetra;
-  switch (NbNodes())
+  switch (aVtkType)
   {
-    case 4:
+    case VTK_TETRA:
       aType = SMDSEntity_Tetra;
       break;
-    case 5:
+    case VTK_PYRAMID:
       aType = SMDSEntity_Pyramid;
       break;
-    case 6:
+    case VTK_WEDGE:
       aType = SMDSEntity_Penta;
       break;
-    case 8:
+    case VTK_HEXAHEDRON:
       aType = SMDSEntity_Hexa;
       break;
-    case 10:
+    case VTK_QUADRATIC_TETRA:
       aType = SMDSEntity_Quad_Tetra;
       break;
-    case 13:
+    case VTK_QUADRATIC_PYRAMID:
       aType = SMDSEntity_Quad_Pyramid;
       break;
-    case 15:
+    case VTK_QUADRATIC_WEDGE:
       aType = SMDSEntity_Quad_Penta;
       break;
-    case 20:
+    case VTK_QUADRATIC_HEXAHEDRON:
       aType = SMDSEntity_Quad_Hexa;
       break;
+#ifdef VTK_HAVE_POLYHEDRON
+    case VTK_POLYHEDRON:
+      aType = SMDSEntity_Polyhedra;
+      break;
+#endif
     default:
-      aType = SMDSEntity_Hexa;
+      aType = SMDSEntity_Polyhedra;
       break;
   }
   return aType;
@@ -205,36 +234,7 @@ SMDSAbs_EntityType SMDS_VtkVolume::GetEntityType() const
 
 vtkIdType SMDS_VtkVolume::GetVtkType() const
 {
-  vtkIdType aType = VTK_TETRA;
-  switch (NbNodes())
-  {
-    case 4:
-      aType = VTK_TETRA;
-      break;
-    case 5:
-      aType = VTK_PYRAMID;
-      break;
-    case 6:
-      aType = VTK_WEDGE;
-      break;
-    case 8:
-      aType = VTK_HEXAHEDRON;
-      break;
-    case 10:
-      aType = VTK_QUADRATIC_TETRA;
-      break;
-    case 13:
-      aType = VTK_QUADRATIC_PYRAMID;
-      break;
-    case 15:
-      aType = VTK_QUADRATIC_WEDGE;
-      break;
-    case 20:
-      aType = VTK_QUADRATIC_HEXAHEDRON;
-      break;
-    default:
-      aType = VTK_HEXAHEDRON;
-      break;
-  }
+  vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
+  vtkIdType aType = grid->GetCellType(myVtkID);
   return aType;
 }
index 9ee9f81dfb1165c83dc1f6a82a54b4239d63d7b3..54d01b605d97ed664b52957cd8388d603395b2d6 100644 (file)
@@ -4,7 +4,7 @@
 #include "SMESH_SMDS.hxx"
 
 #include "SMDS_MeshVolume.hxx"
-#include <vtkUnstructuredGrid.h>
+#include "SMDS_UnstructuredGrid.hxx"
 #include <vector>
 
 class SMDS_EXPORT SMDS_VtkVolume: public SMDS_MeshVolume
@@ -14,6 +14,9 @@ public:
   SMDS_VtkVolume(std::vector<vtkIdType> nodeIds, SMDS_Mesh* mesh);
   ~SMDS_VtkVolume();
   void init(std::vector<vtkIdType> nodeIds, SMDS_Mesh* mesh);
+#ifdef VTK_HAVE_POLYHEDRON
+  void initPoly(std::vector<vtkIdType> nodeIds, std::vector<int> nbNodesPerFace, SMDS_Mesh* mesh);
+#endif
   bool ChangeNodes(const SMDS_MeshNode* nodes[], const int nbNodes);
 
   void Print(std::ostream & OS) const;
index edbe6da4d150a4bbd87aeb5cac2b7328b18dc801..5b4bfbaf52954b919b9b9d123ba362a2247f7c69 100644 (file)
@@ -7770,8 +7770,10 @@ SMESH_MeshEditor::FindFaceInSet(const SMDS_MeshNode*    n1,
   const SMDS_MeshElement* face = 0;
 
   SMDS_ElemIteratorPtr invElemIt = n1->GetInverseElementIterator(SMDSAbs_Face);
+  //MESSAGE("n1->GetInverseElementIterator(SMDSAbs_Face) " << invElemIt);
   while ( invElemIt->more() && !face ) // loop on inverse faces of n1
   {
+    //MESSAGE("in while ( invElemIt->more() && !face )");
     const SMDS_MeshElement* elem = invElemIt->next();
     if (avoidSet.count( elem ))
       continue;
index 17ba61a84767104617f24e9bcf264666102669e6..054b424f7cc57b87430f809c842bdd0f6ec29bba 100644 (file)
@@ -558,7 +558,11 @@ CORBA::Long SMESH_MeshEditor_i::AddPolyhedralVolume (const SMESH::long_array & I
   int NbNodes = IDsOfNodes.length();
   std::vector<const SMDS_MeshNode*> n (NbNodes);
   for (int i = 0; i < NbNodes; i++)
-    n[i] = GetMeshDS()->FindNode(IDsOfNodes[i]);
+    {
+      const SMDS_MeshNode* aNode = GetMeshDS()->FindNode(IDsOfNodes[i]);
+      if (!aNode) return 0;
+      n[i] = aNode;
+    }
 
   int NbFaces = Quantities.length();
   std::vector<int> q (NbFaces);