Salome HOME
23514: EDF 16031 - SMESH freezes
[modules/smesh.git] / src / SMDS / SMDS_Mesh.cxx
index 6582e23f156cf269ea2239d2131d406fec87b328..40caa2c8854d51f04cbb99b91ba3f5943a34f8cd 100644 (file)
@@ -56,7 +56,7 @@
 #include <iterator>
 using namespace std;
 
-#ifndef WIN32
+#if !defined WIN32 && !defined __APPLE__
 #include <sys/sysinfo.h>
 #endif
 
@@ -77,7 +77,7 @@ int SMDS_Mesh::chunkSize = 1024;
 
 int SMDS_Mesh::CheckMemory(const bool doNotRaise) throw (std::bad_alloc)
 {
-#ifndef WIN32
+#if !defined WIN32 && !defined __APPLE__
   struct sysinfo si;
   int err = sysinfo( &si );
   if ( err )
@@ -130,7 +130,6 @@ SMDS_Mesh::SMDS_Mesh():
   myNodeIDFactory(new SMDS_MeshNodeIDFactory()),
   myElementIDFactory(new SMDS_MeshElementIDFactory()),
   myModified(false), myModifTime(0), myCompactTime(0),
-  myNodeMin(0), myNodeMax(0),
   myHasConstructionEdges(false), myHasConstructionFaces(false),
   myHasInverseElements(true),
   xmin(0), xmax(0), ymin(0), ymax(0), zmin(0), zmax(0)
@@ -160,8 +159,19 @@ SMDS_Mesh::SMDS_Mesh():
   points->SetNumberOfPoints(0 /*SMDS_Mesh::chunkSize*/);
   myGrid->SetPoints( points );
   points->Delete();
-  myGrid->BuildLinks();
+  //myGrid->BuildLinks();
   this->Modified();
+
+  // initialize static maps in SMDS_MeshCell, to be thread-safe
+  if ( myMeshId == 0 )
+  {
+    SMDS_MeshCell::toVtkType( SMDSEntity_Node );
+    SMDS_MeshCell::toVtkOrder( SMDSEntity_Node );
+    SMDS_MeshCell::reverseSmdsOrder( SMDSEntity_Node );
+    SMDS_MeshCell::interlacedSmdsOrder( SMDSEntity_Node );
+    SMDS_MeshCell::toSmdsType( VTK_VERTEX );
+    SMDS_MeshCell::fromVtkOrder( SMDSEntity_Node );
+  }
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -373,7 +383,7 @@ SMDS_MeshEdge* SMDS_Mesh::AddEdgeWithID(const SMDS_MeshNode * n1,
   if ( !n1 || !n2 ) return 0;
   SMDS_MeshEdge * edge = 0;
 
-  // --- retreive nodes ID
+  // --- retrieve nodes ID
   vector<vtkIdType> nodeIds;
   nodeIds.clear();
   nodeIds.push_back(n1->getVtkId());
@@ -1546,9 +1556,15 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeFromVtkIdsWithID(const std::vector<vtkIdTyp
     case VTK_QUADRATIC_WEDGE:
       myInfo.myNbQuadPrisms++;
       break;
+    case VTK_BIQUADRATIC_QUADRATIC_WEDGE:
+      myInfo.myNbBiQuadPrisms++;
+      break;
     case VTK_QUADRATIC_HEXAHEDRON:
       myInfo.myNbQuadHexas++;
       break;
+    case VTK_TRIQUADRATIC_HEXAHEDRON:
+      myInfo.myNbTriQuadHexas++;
+      break;
 //#ifdef VTK_HAVE_POLYHEDRON
     case VTK_POLYHEDRON:
       myInfo.myNbPolyhedrons++;
@@ -1680,10 +1696,10 @@ const SMDS_MeshNode * SMDS_Mesh::FindNodeVtk(int vtkId) const
   return (const SMDS_MeshNode *)myNodes[vtkId+1];
 }
 
-///////////////////////////////////////////////////////////////////////////////
-///Create a triangle and add it to the current mesh. This method do not bind an
-///ID to the create triangle.
-///////////////////////////////////////////////////////////////////////////////
+//////////////////////////////////////////////////////////////////////////////
+///Create a triangle and add it to the current mesh. This method does not bind
+///an ID to the create triangle.
+//////////////////////////////////////////////////////////////////////////////
 SMDS_MeshFace * SMDS_Mesh::createTriangle(const SMDS_MeshNode * node1,
                                           const SMDS_MeshNode * node2,
                                           const SMDS_MeshNode * node3,
@@ -1730,10 +1746,10 @@ SMDS_MeshFace * SMDS_Mesh::createTriangle(const SMDS_MeshNode * node1,
   }
 }
 
-///////////////////////////////////////////////////////////////////////////////
-///Create a quadrangle and add it to the current mesh. This methode do not bind
-///a ID to the create triangle.
-///////////////////////////////////////////////////////////////////////////////
+////////////////////////////////////////////////////////////////////////////////
+///Create a quadrangle and add it to the current mesh. This method does not bind
+///an ID to the create triangle.
+////////////////////////////////////////////////////////////////////////////////
 SMDS_MeshFace * SMDS_Mesh::createQuadrangle(const SMDS_MeshNode * node1,
                                             const SMDS_MeshNode * node2,
                                             const SMDS_MeshNode * node3,
@@ -2448,6 +2464,49 @@ const SMDS_MeshElement* SMDS_Mesh::FindElement (const vector<const SMDS_MeshNode
   return NULL;
 }
 
+//================================================================================
+/*!
+ * \brief Return elements including all given nodes
+ *  \param [in] nodes - nodes to find elements around
+ *  \param [out] foundElems - the found elements
+ *  \param [in] type - type of elements to find
+ *  \return int - a number of found elements
+ */
+//================================================================================
+
+int SMDS_Mesh::GetElementsByNodes(const std::vector<const SMDS_MeshNode *>& nodes,
+                                  std::vector<const SMDS_MeshElement *>&    foundElems,
+                                  const SMDSAbs_ElementType                 type)
+{
+  // chose a node with minimal number of inverse elements
+  const SMDS_MeshNode* n0 = nodes[0];
+  int minNbInverse = n0 ? n0->NbInverseElements( type ) : 1000;
+  for ( size_t i = 1; i < nodes.size(); ++i )
+    if ( nodes[i] && nodes[i]->NbInverseElements( type ) < minNbInverse )
+    {
+      n0 = nodes[i];
+      minNbInverse = n0->NbInverseElements( type );
+    }
+
+  foundElems.clear();
+  if ( n0 )
+  {
+    foundElems.reserve( minNbInverse );
+    SMDS_ElemIteratorPtr eIt = n0->GetInverseElementIterator( type );
+    while ( eIt->more() )
+    {
+      const SMDS_MeshElement* e = eIt->next();
+      bool includeAll = true;
+      for ( size_t i = 0; i < nodes.size() &&  includeAll; ++i )
+        if ( nodes[i] != n0 && e->GetNodeIndex( nodes[i] ) < 0 )
+          includeAll = false;
+      if ( includeAll )
+        foundElems.push_back( e );
+    }
+  }
+  return foundElems.size();
+}
+
 //=======================================================================
 //function : DumpNodes
 //purpose  :
@@ -2560,6 +2619,13 @@ int SMDS_Mesh::NbNodes() const
   return myInfo.NbNodes();
 }
 
+///////////////////////////////////////////////////////////////////////////////
+/// Return the number of elements
+///////////////////////////////////////////////////////////////////////////////
+int SMDS_Mesh::NbElements() const
+{
+  return myInfo.NbElements();
+}
 ///////////////////////////////////////////////////////////////////////////////
 /// Return the number of 0D elements
 ///////////////////////////////////////////////////////////////////////////////
@@ -2661,52 +2727,26 @@ SMDS_Mesh::~SMDS_Mesh()
 void SMDS_Mesh::Clear()
 {
   if (myParent!=NULL)
-    {
+  {
     SMDS_ElemIteratorPtr eIt = elementsIterator();
     while ( eIt->more() )
-      {
-        const SMDS_MeshElement *elem = eIt->next();
-        myElementIDFactory->ReleaseID(elem->GetID(), elem->getVtkId());
-      }
+    {
+      const SMDS_MeshElement *elem = eIt->next();
+      myElementIDFactory->ReleaseID(elem->GetID(), elem->getVtkId());
+    }
     SMDS_NodeIteratorPtr itn = nodesIterator();
     while (itn->more())
-      {
-        const SMDS_MeshNode *node = itn->next();
-        myNodeIDFactory->ReleaseID(node->GetID(), node->getVtkId());
-      }
+    {
+      const SMDS_MeshNode *node = itn->next();
+      myNodeIDFactory->ReleaseID(node->GetID(), node->getVtkId());
     }
+  }
   else
-    {
+  {
     myNodeIDFactory->Clear();
     myElementIDFactory->Clear();
-    }
+  }
 
-  // SMDS_ElemIteratorPtr itv = elementsIterator();
-  // while (itv->more())
-  //   {
-  //     SMDS_MeshElement* elem = (SMDS_MeshElement*)(itv->next());
-  //     SMDSAbs_ElementType aType = elem->GetType();
-  //     switch (aType)
-  //     {
-  //       case SMDSAbs_0DElement:
-  //         delete elem;
-  //         break;
-  //       case SMDSAbs_Edge:
-  //          myEdgePool->destroy(static_cast<SMDS_VtkEdge*>(elem));
-  //         break;
-  //       case SMDSAbs_Face:
-  //         myFacePool->destroy(static_cast<SMDS_VtkFace*>(elem));
-  //         break;
-  //       case SMDSAbs_Volume:
-  //         myVolumePool->destroy(static_cast<SMDS_VtkVolume*>(elem));
-  //         break;
-  //       case SMDSAbs_Ball:
-  //         myBallPool->destroy(static_cast<SMDS_BallElement*>(elem));
-  //         break;
-  //       default:
-  //         break;
-  //     }
-  //   }
   myVolumePool->clear();
   myFacePool->clear();
   myEdgePool->clear();
@@ -2717,11 +2757,11 @@ void SMDS_Mesh::Clear()
 
   SMDS_NodeIteratorPtr itn = nodesIterator();
   while (itn->more())
-    {
-      SMDS_MeshNode *node = (SMDS_MeshNode*)(itn->next());
-      node->SetPosition(SMDS_SpacePosition::originSpacePosition());
-      //myNodePool->destroy(node);
-    }
+  {
+    SMDS_MeshNode *node = (SMDS_MeshNode*)(itn->next());
+    node->SetPosition(SMDS_SpacePosition::originSpacePosition());
+    //myNodePool->destroy(node);
+  }
   myNodePool->clear();
   clearVector( myNodes );
 
@@ -2743,10 +2783,10 @@ void SMDS_Mesh::Clear()
   // rnv: to fix bug "21125: EDF 1233 SMESH: Degrardation of precision in a test case for quadratic conversion"
   // using double type for storing coordinates of nodes instead float.
   points->SetDataType(VTK_DOUBLE);
-  points->SetNumberOfPoints(0 /*SMDS_Mesh::chunkSize*/);
+  points->SetNumberOfPoints( 0 );
   myGrid->SetPoints( points );
   points->Delete();
-  myGrid->BuildLinks();
+  myGrid->DeleteLinks();
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2756,7 +2796,7 @@ void SMDS_Mesh::Clear()
 ///////////////////////////////////////////////////////////////////////////////
 bool SMDS_Mesh::hasConstructionEdges()
 {
-        return myHasConstructionEdges;
+  return myHasConstructionEdges;
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2768,7 +2808,7 @@ bool SMDS_Mesh::hasConstructionEdges()
 ///////////////////////////////////////////////////////////////////////////////
 bool SMDS_Mesh::hasConstructionFaces()
 {
-        return myHasConstructionFaces;
+  return myHasConstructionFaces;
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2777,7 +2817,7 @@ bool SMDS_Mesh::hasConstructionFaces()
 ///////////////////////////////////////////////////////////////////////////////
 bool SMDS_Mesh::hasInverseElements()
 {
-        return myHasInverseElements;
+  return myHasInverseElements;
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2786,7 +2826,7 @@ bool SMDS_Mesh::hasInverseElements()
 ///////////////////////////////////////////////////////////////////////////////
 void SMDS_Mesh::setConstructionEdges(bool b)
 {
-        myHasConstructionEdges=b;
+  myHasConstructionEdges=b;
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2795,7 +2835,7 @@ void SMDS_Mesh::setConstructionEdges(bool b)
 ///////////////////////////////////////////////////////////////////////////////
 void SMDS_Mesh::setConstructionFaces(bool b)
 {
-         myHasConstructionFaces=b;
+  myHasConstructionFaces=b;
 }
 
 ///////////////////////////////////////////////////////////////////////////////
@@ -2828,7 +2868,7 @@ namespace {
     IdSortedIterator(const SMDS_MeshElementIDFactory& fact,
                      const SMDSAbs_ElementType        type, // SMDSAbs_All NOT allowed!!! 
                      const int                        totalNb)
-      :myIDFact( fact ),
+      :myIDFact( const_cast<SMDS_MeshElementIDFactory&>(fact) ),
        myID(1), myMaxID( myIDFact.GetMaxID() ),myNbFound(0), myTotalNb( totalNb ),
        myType( type ),
        myElem(0)
@@ -3242,7 +3282,7 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
           if (const SMDS_VtkEdge* vtkElem = dynamic_cast<const SMDS_VtkEdge*>(*it))
             myEdgePool->destroy((SMDS_VtkEdge*) vtkElem);
           else {
-            ((SMDS_MeshElement*) *it)->init( -1, -1, -1 ); // avoid reuse
+            ((SMDS_MeshElement*) *it)->init( 0, -1, -1 ); // avoid reuse
             delete (*it);
           }
           break;
@@ -3257,7 +3297,7 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
           if (const SMDS_VtkFace* vtkElem = dynamic_cast<const SMDS_VtkFace*>(*it))
             myFacePool->destroy((SMDS_VtkFace*) vtkElem);
           else {
-            ((SMDS_MeshElement*) *it)->init( -1, -1, -1 ); // avoid reuse
+            ((SMDS_MeshElement*) *it)->init( 0, -1, -1 ); // avoid reuse
             delete (*it);
           }
           break;
@@ -3272,7 +3312,7 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
           if (const SMDS_VtkVolume* vtkElem = dynamic_cast<const SMDS_VtkVolume*>(*it))
             myVolumePool->destroy((SMDS_VtkVolume*) vtkElem);
           else {
-            ((SMDS_MeshElement*) *it)->init( -1, -1, -1 ); // avoid reuse
+            ((SMDS_MeshElement*) *it)->init( 0, -1, -1 ); // avoid reuse
             delete (*it);
           }
           break;
@@ -3287,7 +3327,7 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
           if (const SMDS_BallElement* vtkElem = dynamic_cast<const SMDS_BallElement*>(*it))
             myBallPool->destroy(const_cast<SMDS_BallElement*>( vtkElem ));
           else {
-            ((SMDS_MeshElement*) *it)->init( -1, -1, -1 ); // avoid reuse
+            ((SMDS_MeshElement*) *it)->init( 0, -1, -1 ); // avoid reuse
             delete (*it);
           }
           break;
@@ -3338,22 +3378,24 @@ void SMDS_Mesh::RemoveElement(const SMDS_MeshElement *        elem,
 void SMDS_Mesh::RemoveFreeElement(const SMDS_MeshElement * elem)
 {
   int elemId = elem->GetID();
-  int vtkId = elem->getVtkId();
+  int  vtkId = elem->getVtkId();
   SMDSAbs_ElementType aType = elem->GetType();
-  SMDS_MeshElement* todest = (SMDS_MeshElement*)(elem);
-  if (aType == SMDSAbs_Node) {
+  SMDS_MeshElement*  todest = (SMDS_MeshElement*)(elem);
+  if ( aType == SMDSAbs_Node )
+  {
     // only free node can be removed by this method
     const SMDS_MeshNode* n = static_cast<SMDS_MeshNode*>(todest);
-    SMDS_ElemIteratorPtr itFe = n->GetInverseElementIterator();
-    if (!itFe->more()) { // free node
+    if ( n->NbInverseElements() == 0 ) { // free node
       myNodes[elemId] = 0;
       myInfo.myNbNodes--;
       ((SMDS_MeshNode*) n)->SetPosition(SMDS_SpacePosition::originSpacePosition());
-      ((SMDS_MeshNode*) n)->SMDS_MeshElement::init( -1, -1, -1 ); // avoid reuse
+      ((SMDS_MeshNode*) n)->SMDS_MeshElement::init( 0, -1, -1 ); // avoid reuse
       myNodePool->destroy(static_cast<SMDS_MeshNode*>(todest));
       myNodeIDFactory->ReleaseID(elemId, vtkId);
     }
-  } else {
+  }
+  else
+  {
     if (hasConstructionEdges() || hasConstructionFaces())
       // this methods is only for meshes without descendants
       return;
@@ -3403,7 +3445,7 @@ void SMDS_Mesh::RemoveFreeElement(const SMDS_MeshElement * elem)
     // --- to do: keep vtkid in a list of reusable cells
 
     if ( elem )
-      ((SMDS_MeshElement*) elem)->init( -1, -1, -1 ); // avoid reuse
+      ((SMDS_MeshElement*) elem)->init( 0, -1, -1 ); // avoid reuse
   }
 }
 
@@ -3413,7 +3455,7 @@ void SMDS_Mesh::RemoveFreeElement(const SMDS_MeshElement * elem)
  */
 bool SMDS_Mesh::Contains (const SMDS_MeshElement* elem) const
 {
-  // we should not imply on validity of *elem, so iterate on containers
+  // we should not rely on validity of *elem, so iterate on containers
   // of all types in the hope of finding <elem> somewhere there
   SMDS_NodeIteratorPtr itn = nodesIterator();
   while (itn->more())
@@ -3433,7 +3475,7 @@ bool SMDS_Mesh::Contains (const SMDS_MeshElement* elem) const
 
 int SMDS_Mesh::MaxNodeID() const
 {
-  return myNodeMax;
+  return myNodeIDFactory->GetMaxID();
 }
 
 //=======================================================================
@@ -3443,7 +3485,7 @@ int SMDS_Mesh::MaxNodeID() const
 
 int SMDS_Mesh::MinNodeID() const
 {
-  return myNodeMin;
+  return myNodeIDFactory->GetMinID();
 }
 
 //=======================================================================
@@ -4174,7 +4216,7 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
 
 //=======================================================================
 //function : AddVolume
-//purpose  :
+//purpose  : 2d order Pentahedron (prism) with 15 nodes
 //=======================================================================
 SMDS_MeshVolume* SMDS_Mesh::AddVolume(const SMDS_MeshNode * n1,
                                       const SMDS_MeshNode * n2,
@@ -4202,7 +4244,7 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolume(const SMDS_MeshNode * n1,
 
 //=======================================================================
 //function : AddVolumeWithID
-//purpose  :
+//purpose  : 2d order Pentahedron (prism) with 15 nodes
 //=======================================================================
 SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(int n1, int n2, int n3,
                                             int n4, int n5, int n6,
@@ -4231,7 +4273,7 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(int n1, int n2, int n3,
 
 //=======================================================================
 //function : AddVolumeWithID
-//purpose  : 2d order Pentahedron with 15 nodes
+//purpose  : 2d order Pentahedron (prism) with 15 nodes
 //=======================================================================
 SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
                                             const SMDS_MeshNode * n2,
@@ -4298,6 +4340,149 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
   return volvtk;
 }
 
+//=======================================================================
+//function : AddVolume
+//purpose  : 2d order Pentahedron (prism) with 18 nodes
+//=======================================================================
+SMDS_MeshVolume* SMDS_Mesh::AddVolume(const SMDS_MeshNode * n1,
+                                      const SMDS_MeshNode * n2,
+                                      const SMDS_MeshNode * n3,
+                                      const SMDS_MeshNode * n4,
+                                      const SMDS_MeshNode * n5,
+                                      const SMDS_MeshNode * n6,
+                                      const SMDS_MeshNode * n12,
+                                      const SMDS_MeshNode * n23,
+                                      const SMDS_MeshNode * n31,
+                                      const SMDS_MeshNode * n45,
+                                      const SMDS_MeshNode * n56,
+                                      const SMDS_MeshNode * n64,
+                                      const SMDS_MeshNode * n14,
+                                      const SMDS_MeshNode * n25,
+                                      const SMDS_MeshNode * n36,
+                                      const SMDS_MeshNode * n1245,
+                                      const SMDS_MeshNode * n2356,
+                                      const SMDS_MeshNode * n1346)
+{
+  //MESSAGE("AddVolume penta18");
+  int ID = myElementIDFactory->GetFreeID();
+  SMDS_MeshVolume * v =
+    SMDS_Mesh::AddVolumeWithID(n1, n2, n3, n4, n5, n6, n12, n23, n31,
+                               n45, n56, n64, n14, n25, n36, n1245, n2356, n1346, ID);
+  if(v==NULL) myElementIDFactory->ReleaseID(ID);
+  return v;
+}
+
+//=======================================================================
+//function : AddVolumeWithID
+//purpose  : 2d order Pentahedron (prism) with 18 nodes
+//=======================================================================
+SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(int n1, int n2, int n3,
+                                            int n4, int n5, int n6,
+                                            int n12,int n23,int n31,
+                                            int n45,int n56,int n64,
+                                            int n14,int n25,int n36,
+                                            int n1245, int n2356, int n1346, int ID)
+{
+  //MESSAGE("AddVolumeWithID penta18 " << ID);
+  return SMDS_Mesh::AddVolumeWithID
+    ((SMDS_MeshNode*) myNodeIDFactory->MeshElement(n1) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n2) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n3) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n4) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n5) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n6) ,
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n12),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n23),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n31),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n45),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n56),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n64),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n14),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n25),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n36),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n1245),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n2356),
+     (SMDS_MeshNode*) myNodeIDFactory->MeshElement(n1346),
+     ID);
+}
+
+//=======================================================================
+//function : AddVolumeWithID
+//purpose  : 2d order Pentahedron (prism) with 18 nodes
+//=======================================================================
+SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
+                                            const SMDS_MeshNode * n2,
+                                            const SMDS_MeshNode * n3,
+                                            const SMDS_MeshNode * n4,
+                                            const SMDS_MeshNode * n5,
+                                            const SMDS_MeshNode * n6,
+                                            const SMDS_MeshNode * n12,
+                                            const SMDS_MeshNode * n23,
+                                            const SMDS_MeshNode * n31,
+                                            const SMDS_MeshNode * n45,
+                                            const SMDS_MeshNode * n56,
+                                            const SMDS_MeshNode * n64,
+                                            const SMDS_MeshNode * n14,
+                                            const SMDS_MeshNode * n25,
+                                            const SMDS_MeshNode * n36,
+                                            const SMDS_MeshNode * n1245,
+                                            const SMDS_MeshNode * n2356,
+                                            const SMDS_MeshNode * n1346,
+                                            int ID)
+{
+  //MESSAGE("AddVolumeWithID penta18 "<< ID);
+  if (!n1 || !n2 || !n3 || !n4 || !n5 || !n6 || !n12 || !n23 ||
+      !n31 || !n45 || !n56 || !n64 || !n14 || !n25 || !n36 || !n1245 || !n2356 || !n1346)
+    return 0;
+  if(hasConstructionFaces()) {
+    // creation quadratic faces - not implemented
+    return 0;
+  }
+  // --- retrieve nodes ID
+  myNodeIds.resize(18);
+  myNodeIds[0] = n1->getVtkId();
+  myNodeIds[1] = n2->getVtkId();
+  myNodeIds[2] = n3->getVtkId();
+
+  myNodeIds[3] = n4->getVtkId();
+  myNodeIds[4] = n5->getVtkId();
+  myNodeIds[5] = n6->getVtkId();
+
+  myNodeIds[6] = n12->getVtkId();
+  myNodeIds[7] = n23->getVtkId();
+  myNodeIds[8] = n31->getVtkId();
+
+  myNodeIds[9] = n45->getVtkId();
+  myNodeIds[10] = n56->getVtkId();
+  myNodeIds[11] = n64->getVtkId();
+
+  myNodeIds[12] = n14->getVtkId();
+  myNodeIds[13] = n25->getVtkId();
+  myNodeIds[14] = n36->getVtkId();
+
+  myNodeIds[15] = n1245->getVtkId();
+  myNodeIds[16] = n2356->getVtkId();
+  myNodeIds[17] = n1346->getVtkId();
+
+  SMDS_VtkVolume *volvtk = myVolumePool->getNew();
+  volvtk->init(myNodeIds, this);
+  if (!this->registerElement(ID,volvtk))
+  {
+    this->myGrid->GetCellTypesArray()->SetValue(volvtk->getVtkId(), VTK_EMPTY_CELL);
+    myVolumePool->destroy(volvtk);
+    return 0;
+  }
+  adjustmyCellsCapacity(ID);
+  myCells[ID] = volvtk;
+  myInfo.myNbBiQuadPrisms++;
+
+  //  if (!registerElement(ID, volvtk)) {
+  //    RemoveElement(volvtk, false);
+  //    volvtk = NULL;
+  //  }
+  return volvtk;
+}
+
 
 //=======================================================================
 //function : AddVolume
@@ -4621,44 +4806,6 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
   return volvtk;
 }
 
-
-void SMDS_Mesh::updateNodeMinMax()
-{
-  myNodeMin = 0;
-  if (myNodes.size() == 0)
-  {
-    myNodeMax=0;
-    return;
-  }
-  while ( !myNodes[myNodeMin] && myNodeMin < (int)myNodes.size() )
-    myNodeMin++;
-  myNodeMax=myNodes.size()-1;
-  while (!myNodes[myNodeMax] && (myNodeMin>=0))
-    myNodeMin--;
-}
-
-void SMDS_Mesh::incrementNodesCapacity(int nbNodes)
-{
-//  int val = myCellIdSmdsToVtk.size();
-//  MESSAGE(" ------------------- resize myCellIdSmdsToVtk " << val << " --> " << val + nbNodes);
-//  myCellIdSmdsToVtk.resize(val + nbNodes, -1); // fill new elements with -1
-  int val = myNodes.size();
-  myNodes.resize(val +nbNodes, 0);
-}
-
-void SMDS_Mesh::incrementCellsCapacity(int nbCells)
-{
-  int val = myCellIdVtkToSmds.size();
-  myCellIdVtkToSmds.resize(val + nbCells, -1); // fill new elements with -1
-  val = myCells.size();
-  myNodes.resize(val +nbCells, 0);
-}
-
-void SMDS_Mesh::adjustStructure()
-{
-  myGrid->GetPoints()->GetData()->SetNumberOfTuples(myNodeIDFactory->GetMaxID());
-}
-
 void SMDS_Mesh::dumpGrid(string ficdump)
 {
   //  vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
@@ -4691,7 +4838,7 @@ void SMDS_Mesh::dumpGrid(string ficdump)
     ficcon << endl;
   }
   ficcon << "-------------------------------- connectivity " <<  nbPoints << endl;
-  vtkCellLinks *links = myGrid->GetCellLinks();
+  vtkCellLinks *links = myGrid->GetLinks();
   for (int i=0; i<nbPoints; i++)
   {
     int ncells = links->GetNcells(i);
@@ -4700,8 +4847,8 @@ void SMDS_Mesh::dumpGrid(string ficdump)
     for (int j=0; j<ncells; j++)
     {
       ficcon << " " << cells[j];
-        }
-        ficcon << endl;
+    }
+    ficcon << endl;
   }
   ficcon.close();
 
@@ -4709,7 +4856,7 @@ void SMDS_Mesh::dumpGrid(string ficdump)
 
 void SMDS_Mesh::compactMesh()
 {
-  MESSAGE("SMDS_Mesh::compactMesh do nothing!");
+  this->myCompactTime = this->myModifTime;
 }
 
 int SMDS_Mesh::fromVtkToSmds(int vtkid)
@@ -4719,28 +4866,28 @@ int SMDS_Mesh::fromVtkToSmds(int vtkid)
   throw SALOME_Exception(LOCALIZED ("vtk id out of bounds"));
 }
 
-void SMDS_Mesh::updateBoundingBox()
-{
-  xmin = 0; xmax = 0;
-  ymin = 0; ymax = 0;
-  zmin = 0; zmax = 0;
-  vtkPoints *points = myGrid->GetPoints();
-  int myNodesSize = this->myNodes.size();
-  for (int i = 0; i < myNodesSize; i++)
-    {
-      if (SMDS_MeshNode *n = myNodes[i])
-        {
-          double coords[3];
-          points->GetPoint(n->myVtkID, coords);
-          if (coords[0] < xmin) xmin = coords[0];
-          else if (coords[0] > xmax) xmax = coords[0];
-          if (coords[1] < ymin) ymin = coords[1];
-          else if (coords[1] > ymax) ymax = coords[1];
-          if (coords[2] < zmin) zmin = coords[2];
-          else if (coords[2] > zmax) zmax = coords[2];
-        }
-    }
-}
+// void SMDS_Mesh::updateBoundingBox()
+// {
+//   xmin = 0; xmax = 0;
+//   ymin = 0; ymax = 0;
+//   zmin = 0; zmax = 0;
+//   vtkPoints *points = myGrid->GetPoints();
+//   int myNodesSize = this->myNodes.size();
+//   for (int i = 0; i < myNodesSize; i++)
+//     {
+//       if (SMDS_MeshNode *n = myNodes[i])
+//         {
+//           double coords[3];
+//           points->GetPoint(n->myVtkID, coords);
+//           if (coords[0] < xmin) xmin = coords[0];
+//           else if (coords[0] > xmax) xmax = coords[0];
+//           if (coords[1] < ymin) ymin = coords[1];
+//           else if (coords[1] > ymax) ymax = coords[1];
+//           if (coords[2] < zmin) zmin = coords[2];
+//           else if (coords[2] > zmax) zmax = coords[2];
+//         }
+//     }
+// }
 
 double SMDS_Mesh::getMaxDim()
 {
@@ -4762,17 +4909,12 @@ void SMDS_Mesh::Modified()
 }
 
 //! get last modification timeStamp
-unsigned long SMDS_Mesh::GetMTime() const
+vtkMTimeType SMDS_Mesh::GetMTime() const
 {
   return this->myModifTime;
 }
 
 bool SMDS_Mesh::isCompacted()
 {
-  if (this->myModifTime > this->myCompactTime)
-  {
-    this->myCompactTime = this->myModifTime;
-    return false;
-  }
-  return true;
+  return this->myCompactTime == this->myModifTime;
 }