Salome HOME
0022098: EDF 2036 SMESH: Create groups from none conected parts of a mesh
[modules/smesh.git] / src / StdMeshers / StdMeshers_Penta_3D.cxx
index ad6ebc83198609b882539342061257d2b33e4f5d..3db3f2dd9e9a53ce62d9d9aece172c3c958a6957 100644 (file)
@@ -1,29 +1,29 @@
-//  SMESH StdMeshers_Penta_3D implementaion of SMESH idl descriptions
+// Copyright (C) 2007-2013  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+// Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
 //
-//  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
-//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
-// 
-//  This library is free software; you can redistribute it and/or 
-//  modify it under the terms of the GNU Lesser General Public 
-//  License as published by the Free Software Foundation; either 
-//  version 2.1 of the License. 
-// 
-//  This library is distributed in the hope that it will be useful, 
-//  but WITHOUT ANY WARRANTY; without even the implied warranty of 
-//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
-//  Lesser General Public License for more details. 
-// 
-//  You should have received a copy of the GNU Lesser General Public 
-//  License along with this library; if not, write to the Free Software 
-//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
-// 
-//  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
+// This library is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 2.1 of the License.
 //
+// This library is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+// Lesser General Public License for more details.
 //
+// You should have received a copy of the GNU Lesser General Public
+// License along with this library; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
 //
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+//  SMESH StdMeshers_Penta_3D implementaion of SMESH idl descriptions
 //  File   : StdMeshers_Penta_3D.cxx
 //  Module : SMESH
-
+//
 #include "StdMeshers_Penta_3D.hxx"
 
 #include "utilities.h"
 #include "SMDS_VolumeOfNodes.hxx"
 #include "SMDS_VolumeTool.hxx"
 #include "SMESHDS_SubMesh.hxx"
+#include "SMESH_Comment.hxx"
 #include "SMESH_Mesh.hxx"
+#include "SMESH_MeshAlgos.hxx"
+#include "SMESH_MesherHelper.hxx"
 #include "SMESH_subMesh.hxx"
-#include "SMESH_MeshEditor.hxx"
+#include "SMESH_subMeshEventListener.hxx"
 
 #include <BRep_Tool.hxx>
-#include <TopAbs_ShapeEnum.hxx>
 #include <TopExp.hxx>
+#include <TopExp_Explorer.hxx>
 #include <TopTools_IndexedDataMapOfShapeListOfShape.hxx>
 #include <TopTools_IndexedMapOfShape.hxx>
-#include <TopTools_IndexedMapOfShape.hxx>
 #include <TopTools_ListIteratorOfListOfShape.hxx>
 #include <TopTools_ListOfShape.hxx>
+#include <TopTools_SequenceOfShape.hxx>
+#include <TopTools_MapOfShape.hxx>
 #include <TopoDS.hxx>
 #include <TopoDS_Edge.hxx>
 #include <TopoDS_Shell.hxx>
@@ -60,57 +64,75 @@ using namespace std;
 typedef map < int, int, less<int> >::iterator   \
   StdMeshers_IteratorOfDataMapOfIntegerInteger;
 
+enum { NB_WALL_FACES = 4 };
+
 //=======================================================================
 //function : StdMeshers_Penta_3D
 //purpose  : 
 //=======================================================================
 StdMeshers_Penta_3D::StdMeshers_Penta_3D()
-: myErrorStatus(1)
+  : myErrorStatus(SMESH_ComputeError::New())
 {
   myTol3D=0.1;
   myWallNodesMaps.resize( SMESH_Block::NbFaces() );
   myShapeXYZ.resize( SMESH_Block::NbSubShapes() );
+  myTool = 0;
+}
+
+//=======================================================================
+//function : ~StdMeshers_Penta_3D
+//purpose  : 
+//=======================================================================
+
+StdMeshers_Penta_3D::~StdMeshers_Penta_3D()
+{
 }
+
 //=======================================================================
 //function : Compute
 //purpose  : 
 //=======================================================================
 bool StdMeshers_Penta_3D::Compute(SMESH_Mesh& aMesh, 
-                                 const TopoDS_Shape& aShape)
+                                  const TopoDS_Shape& aShape)
 {
   MESSAGE("StdMeshers_Penta_3D::Compute()");
   //
-  myErrorStatus=0;
-  //
   bool bOK=false;
   //
   myShape=aShape;
   SetMesh(aMesh);
   //
   CheckData();
-  if (myErrorStatus){
+  if (!myErrorStatus->IsOK()) {
     return bOK;
   }
+
   //
   MakeBlock();
-    if (myErrorStatus){
+  if (!myErrorStatus->IsOK()) {
     return bOK;
   }
   //
-  MakeNodes();
-  if (myErrorStatus){
+  ClearMeshOnFxy1();
+  if (!myErrorStatus->IsOK()) {
     return bOK;
   }
+
+  // now unnecessary faces removed, we can load medium nodes
+  SMESH_MesherHelper helper(aMesh);
+  myTool = &helper;
+  myCreateQuadratic = myTool->IsQuadraticSubMesh(aShape);
+
   //
-  MakeConnectingMap();
-  //
-  ClearMeshOnFxy1();
-  if (myErrorStatus) {
+  MakeNodes();
+  if (!myErrorStatus->IsOK()) {
     return bOK;
   }
   //
+  MakeConnectingMap();
+  //
   MakeMeshOnFxy1();
-  if (myErrorStatus) {
+  if (!myErrorStatus->IsOK()) {
     return bOK;
   }
   //
@@ -118,14 +140,13 @@ bool StdMeshers_Penta_3D::Compute(SMESH_Mesh& aMesh,
   //
   return !bOK;
 }
+
 //=======================================================================
 //function : MakeNodes
 //purpose  : 
 //=======================================================================
 void StdMeshers_Penta_3D::MakeNodes()
 {
-  myErrorStatus=0;
-  //
   const int aNbSIDs=9;
   int i, j, k, ij, iNbN, aNodeID, aSize, iErr;
   double aX, aY, aZ;
@@ -144,12 +165,24 @@ void StdMeshers_Penta_3D::MakeNodes()
   // 1.1 Horizontal size
   myJSize=0;
   for (i=0; i<aNbSIDs; ++i) {
-    const TopoDS_Shape& aS=myBlock.Shape(aSIDs[i]);
+    const TopoDS_Shape& aS = myBlock.Shape(aSIDs[i]);
     SMESH_subMesh *aSubMesh = pMesh->GetSubMeshContaining(aS);
     ASSERT(aSubMesh);
-    SMESHDS_SubMesh *aSM=aSubMesh->GetSubMeshDS();
-    iNbN=aSM->NbNodes();
-    myJSize+=iNbN;
+    SMESHDS_SubMesh *aSM = aSubMesh->GetSubMeshDS();
+    if(!myCreateQuadratic) {
+      iNbN = aSM->NbNodes();
+    }
+    else {
+      iNbN = 0;
+      SMDS_NodeIteratorPtr itn = aSM->GetNodes();
+      while(itn->more()) {
+        const SMDS_MeshNode* aNode = itn->next();
+        if(myTool->IsMedium(aNode))
+          continue;
+        iNbN++;
+      }
+    }
+    myJSize += iNbN;
   }
   //printf("***  Horizontal: number of nodes summary=%d\n", myJSize);
   //
@@ -159,9 +192,21 @@ void StdMeshers_Penta_3D::MakeNodes()
     const TopoDS_Shape& aS=myBlock.Shape(SMESH_Block::ID_E00z);
     SMESH_subMesh *aSubMesh = pMesh->GetSubMeshContaining(aS);
     ASSERT(aSubMesh);
-    SMESHDS_SubMesh *aSM=aSubMesh->GetSubMeshDS();
-    iNbN=aSM->NbNodes();
-    myISize+=iNbN;
+    SMESHDS_SubMesh *aSM = aSubMesh->GetSubMeshDS();
+    if(!myCreateQuadratic) {
+      iNbN = aSM->NbNodes();
+    }
+    else {
+      iNbN = 0;
+      SMDS_NodeIteratorPtr itn = aSM->GetNodes();
+      while(itn->more()) {
+        const SMDS_MeshNode* aNode = itn->next();
+        if(myTool->IsMedium(aNode))
+          continue;
+        iNbN++;
+      }
+    }
+    myISize += iNbN;
   }
   //printf("***  Vertical: number of nodes on edges and vertices=%d\n", myISize);
   //
@@ -177,20 +222,21 @@ void StdMeshers_Penta_3D::MakeNodes()
   // vertices
   for (k=0; k<aNbSIDs; ++k) {
     aSID=aSIDs[k];
-    const TopoDS_Shape& aS=myBlock.Shape(aSID);
-    SMDS_NodeIteratorPtr ite =pMesh->GetSubMeshContaining(aS)->GetSubMeshDS()->GetNodes();
+    const TopoDS_Shape& aS = myBlock.Shape(aSID);
+    SMDS_NodeIteratorPtr ite = pMesh->GetSubMeshContaining(aS)->GetSubMeshDS()->GetNodes();
     while(ite->more()) {
       const SMDS_MeshNode* aNode = ite->next();
+      if(myTool->IsMedium(aNode))
+        continue;
       aNodeID=aNode->GetID();
       //
       aTNode.SetNode(aNode);
       aTNode.SetShapeSupportID(aSID);
       aTNode.SetBaseNodeID(aNodeID);
       //
-      if ( SMESH_Block::IsEdgeID (aSID))
-      {
+      if ( SMESH_Block::IsEdgeID (aSID)) {
         const SMDS_EdgePosition* epos =
-          static_cast<const SMDS_EdgePosition*>(aNode->GetPosition().get());
+          static_cast<const SMDS_EdgePosition*>(aNode->GetPosition());
         myBlock.ComputeParameters( epos->GetUParameter(), aS, aCoords );
       }
       else {
@@ -200,11 +246,11 @@ void StdMeshers_Penta_3D::MakeNodes()
         aP3D.SetCoord(aX, aY, aZ);
         myBlock.ComputeParameters(aP3D, aS, aCoords);
       }
-      iErr=myBlock.ErrorStatus();
+      iErr = myBlock.ErrorStatus();
       if (iErr) {
         MESSAGE("StdMeshers_Penta_3D::MakeNodes()," <<
                 "SMESHBlock: ComputeParameters operation failed");
-        myErrorStatus=101; // SMESHBlock: ComputeParameters operation failed
+        myErrorStatus=myBlock.GetError();
         return;
       }
       aTNode.SetNormCoord(aCoords);
@@ -213,97 +259,35 @@ void StdMeshers_Penta_3D::MakeNodes()
       ++j;
     }
   }
-  /*
-  //DEB
-  {
-    int iShapeSupportID, iBaseNodeID;
-    //
-    //printf("\n\n*** Base Face\n");
-    i=0;
-    for (j=0; j<myJSize; ++j) {
-      ij=i*myJSize+j;
-      const StdMeshers_TNode& aTNode=myTNodes[ij];
-      iShapeSupportID=aTNode.ShapeSupportID();
-      iBaseNodeID=aTNode.BaseNodeID();
-      const gp_XYZ& aXYZ=aTNode.NormCoord();
-      printf("*** j:%d bID#%d iSS:%d { %lf %lf %lf }\n",
-            j,  iBaseNodeID, iShapeSupportID, aXYZ.X(),  aXYZ.Y(), aXYZ.Z());
-    }
-  }
-  */
-  //DEB
-  //return; //zz
-  //
-  // 3. Finding of Z-layers
-//   vector<double> aZL(myISize);
-//   vector<double>::iterator aItZL1, aItZL2 ;
-//   //
-//   const TopoDS_Shape& aE00z=myBlock.Shape(SMESH_Block::ID_E00z);
-//   SMDS_NodeIteratorPtr aItaE00z =
-//     pMesh->GetSubMeshContaining(aE00z)->GetSubMeshDS()->GetNodes();
-//   //
-//   aZL[0]=0.;
-//   i=1;
-//   while (aItaE00z->more()) {
-//     const SMDS_MeshNode* aNode=aItaE00z->next();
-//     const SMDS_EdgePosition* epos =
-//       static_cast<const SMDS_EdgePosition*>(aNode->GetPosition().get());
-//     myBlock.ComputeParameters( epos->GetUParameter(), aE00z, aCoords );
-//     iErr=myBlock.ErrorStatus();
-//     if (iErr) {
-//       MESSAGE("StdMeshers_Penta_3D::MakeNodes()," <<
-//               "SMESHBlock: ComputeParameters operation failed");
-//       myErrorStatus=101; // SMESHBlock: ComputeParameters operation failed
-//       return;
-//     }
-//     aZL[i]=aCoords.Z();
-//     ++i;
-//   }
-//   aZL[i]=1.;
-//   //
-//   aItZL1=aZL.begin();
-//   aItZL2=aZL.end();
-//   //
-//   // Sorting the layers
-//   sort(aItZL1, aItZL2);
-  //DEB
-  /*
-  printf("** \n\n Layers begin\n");
-  for(i=0, aItZL=aItZL1; aItZL!=aItZL2; ++aItZL, ++i) {
-    printf(" #%d : %lf\n", i, *aItZL);
-  } 
-  printf("** Layers end\n");
-  */
-  //DEB
-  //
-  //
 
   // 3.1 Fill maps of wall nodes
-  SMESH_Block::TShapeID wallFaceID[4] = {
+  SMESH_Block::TShapeID wallFaceID[ NB_WALL_FACES ] = {
     SMESH_Block::ID_Fx0z, SMESH_Block::ID_Fx1z,
     SMESH_Block::ID_F0yz, SMESH_Block::ID_F1yz
-    };
-  SMESH_Block::TShapeID baseEdgeID[4] = {
+  };
+  SMESH_Block::TShapeID baseEdgeID[ NB_WALL_FACES ] = {
     SMESH_Block::ID_Ex00, SMESH_Block::ID_Ex10,
     SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0
-    };
-  for ( i = 0; i < 4; ++i ) {
+  };
+  for ( i = 0; i < NB_WALL_FACES ; ++i ) {
     int fIndex = SMESH_Block::ShapeIndex( wallFaceID[ i ]);
     bool ok = LoadIJNodes (myWallNodesMaps[ fIndex ],
                            TopoDS::Face( myBlock.Shape( wallFaceID[ i ] )),
                            TopoDS::Edge( myBlock.Shape( baseEdgeID[ i ] )),
                            pMesh->GetMeshDS());
     if ( !ok ) {
-      myErrorStatus = i + 1;
-      MESSAGE(" Cant LoadIJNodes() from a wall face " << myErrorStatus );
+      myErrorStatus->myName = COMPERR_BAD_INPUT_MESH;
+      myErrorStatus->myComment = SMESH_Comment() <<
+        "Can't find regular quadrangle mesh on a side face #" <<
+        pMesh->GetMeshDS()->ShapeToIndex( myBlock.Shape( wallFaceID[ i ]));
       return;
     }
   }
 
   // 3.2 find node columns for vertical edges and edge IDs
-  vector<const SMDS_MeshNode*> * verticEdgeNodes[ 4 ];
-  SMESH_Block::TShapeID          verticEdgeID   [ 4 ];
-  for ( i = 0; i < 4; ++i ) { // 4 first base nodes are nodes on vertices
+  vector<const SMDS_MeshNode*> * verticEdgeNodes[ NB_WALL_FACES ];
+  SMESH_Block::TShapeID          verticEdgeID   [ NB_WALL_FACES ];
+  for ( i = 0; i < NB_WALL_FACES ; ++i ) { // 4 first base nodes are nodes on vertices
     // edge ID
     SMESH_Block::TShapeID eID, vID = aSIDs[ i ];
     ShapeSupportID(false, vID, eID);
@@ -311,7 +295,7 @@ void StdMeshers_Penta_3D::MakeNodes()
     // column nodes
     StdMeshers_TNode& aTNode = myTNodes[ i ];
     verticEdgeNodes[ i ] = 0;
-    for ( j = 0; j < 4; ++j ) { // loop on 4 wall faces
+    for ( j = 0; j < NB_WALL_FACES ; ++j ) { // loop on 4 wall faces
       int fIndex = SMESH_Block::ShapeIndex( wallFaceID[ j ]);
       StdMeshers_IJNodeMap & ijNodes= myWallNodesMaps[ fIndex ];
       if ( ijNodes.begin()->second[0] == aTNode.Node() )
@@ -325,8 +309,7 @@ void StdMeshers_Penta_3D::MakeNodes()
 
   // 3.3 set XYZ of vertices, and initialize of the rest
   SMESHDS_Mesh* aMesh = GetMesh()->GetMeshDS();
-  for ( int id = SMESH_Block::ID_V000; id < SMESH_Block::ID_Shell; ++id )
-  {
+  for ( int id = SMESH_Block::ID_V000; id < SMESH_Block::ID_Shell; ++id ) {
     if ( SMESH_Block::IsVertexID( id )) {
       TopoDS_Shape V = myBlock.Shape( id );
       SMESHDS_SubMesh* sm = aMesh->MeshElements( V );
@@ -344,23 +327,46 @@ void StdMeshers_Penta_3D::MakeNodes()
   SMESH_Block::TShapeID aSSID, aBNSSID;
   StdMeshers_TNode aTN;
   //
-  for (j=0; j<myJSize; ++j)
-  {
+
+  // create top face and find UV for it's corners
+  const TopoDS_Face& TopFace = TopoDS::Face(myBlock.Shape(SMESH_Block::ID_Fxy1));
+  SMESHDS_Mesh* meshDS = pMesh->GetMeshDS();
+  int topfaceID = meshDS->ShapeToIndex(TopFace);
+  const TopoDS_Vertex& v001 = TopoDS::Vertex(myBlock.Shape(SMESH_Block::ID_V001));
+  SMDS_NodeIteratorPtr itn = pMesh->GetSubMeshContaining(v001)->GetSubMeshDS()->GetNodes();
+  const SMDS_MeshNode* N = itn->next();
+  gp_XY UV001 = myTool->GetNodeUV(TopFace,N);
+  const TopoDS_Vertex& v101 = TopoDS::Vertex(myBlock.Shape(SMESH_Block::ID_V101));
+  itn = pMesh->GetSubMeshContaining(v101)->GetSubMeshDS()->GetNodes();
+  N = itn->next();
+  gp_XY UV101 = myTool->GetNodeUV(TopFace,N);
+  const TopoDS_Vertex& v011 = TopoDS::Vertex(myBlock.Shape(SMESH_Block::ID_V011));
+  itn = pMesh->GetSubMeshContaining(v011)->GetSubMeshDS()->GetNodes();
+  N = itn->next();
+  gp_XY UV011 = myTool->GetNodeUV(TopFace,N);
+  const TopoDS_Vertex& v111 = TopoDS::Vertex(myBlock.Shape(SMESH_Block::ID_V111));
+  itn = pMesh->GetSubMeshContaining(v111)->GetSubMeshDS()->GetNodes();
+  N = itn->next();
+  gp_XY UV111 = myTool->GetNodeUV(TopFace,N);
+
+  for (j=0; j<myJSize; ++j) { // loop on all nodes of the base face (ID_Fxy0)
     // base node info
-    const StdMeshers_TNode& aBN=myTNodes[j];
-    aBNSSID=(SMESH_Block::TShapeID)aBN.ShapeSupportID();
-    iBNID=aBN.BaseNodeID();
-    const gp_XYZ& aBNXYZ=aBN.NormCoord();
-    bool createNode = ( aBNSSID == SMESH_Block::ID_Fxy0 );
+    const StdMeshers_TNode& aBN = myTNodes[j];
+    aBNSSID = (SMESH_Block::TShapeID)aBN.ShapeSupportID();
+    iBNID = aBN.BaseNodeID();
+    const gp_XYZ& aBNXYZ = aBN.NormCoord();
+    bool createNode = ( aBNSSID == SMESH_Block::ID_Fxy0 ); // if base node is inside a bottom face
     //
     // set XYZ on horizontal edges and get node columns of faces:
     // 2 columns for each face, between which a base node is located
     vector<const SMDS_MeshNode*>* nColumns[8];
-    double ratio[4]; // base node position between columns [0.-1.]
-    if ( createNode )
-      for ( k = 0; k < 4; ++k )
+    double ratio[ NB_WALL_FACES ]; // base node position between columns [0.-1.]
+    if ( createNode ) {
+      for ( k = 0; k < NB_WALL_FACES ; ++k ) {
         ratio[ k ] = SetHorizEdgeXYZ (aBNXYZ, wallFaceID[ k ],
                                       nColumns[k*2], nColumns[k*2+1]);
+      }
+    }
     //
     // XYZ on the bottom and top faces
     const SMDS_MeshNode* n = aBN.Node();
@@ -368,20 +374,42 @@ void StdMeshers_Penta_3D::MakeNodes()
     myShapeXYZ[ SMESH_Block::ID_Fxy1 ].SetCoord( 0., 0., 0. );
     //
     // first create or find a top node, then the rest ones in a column
-    for (i=myISize-1; i>0; --i)
+    for (i=myISize-1; i>0; --i) // vertical loop, from top to bottom
     {
-      if ( createNode ) {
+      bIsUpperLayer = (i==(myISize-1));
+      gp_XY UV_Ex01, UV_Ex11, UV_E0y1, UV_E1y1;
+      if ( createNode ) // a base node is inside a top face
+      {
         // set XYZ on vertical edges and faces
-        for ( k = 0; k < 4; ++k ) {
+        for ( k = 0; k < NB_WALL_FACES ; ++k ) {
+          // XYZ on a vertical edge 
           const SMDS_MeshNode* n = (*verticEdgeNodes[ k ]) [ i ];
           myShapeXYZ[ verticEdgeID[ k ] ].SetCoord( n->X(), n->Y(), n->Z() );
-          //
+          // XYZ on a face (part 1 from one column)
           n = (*nColumns[k*2]) [ i ];
           gp_XYZ xyz( n->X(), n->Y(), n->Z() );
           myShapeXYZ[ wallFaceID[ k ]] = ( 1. - ratio[ k ]) * xyz;
+          gp_XY tmp1;
+          if( bIsUpperLayer ) {
+            tmp1 = myTool->GetNodeUV(TopFace,n);
+            tmp1 = ( 1. - ratio[ k ]) * tmp1;
+          }
+          // XYZ on a face (part 2 from other column)
           n = (*nColumns[k*2+1]) [ i ];
           xyz.SetCoord( n->X(), n->Y(), n->Z() );
           myShapeXYZ[ wallFaceID[ k ]] += ratio[ k ] * xyz;
+          if( bIsUpperLayer ) {
+            gp_XY tmp2 = myTool->GetNodeUV(TopFace,n);
+            tmp1 +=  ratio[ k ] * tmp2;
+            if( k==0 )
+              UV_Ex01 = tmp1;
+            else if( k==1 )
+              UV_Ex11 = tmp1;
+            else if( k==2 )
+              UV_E0y1 = tmp1;
+            else
+              UV_E1y1 = tmp1;
+          }
         }
       }
       // fill current node info
@@ -395,11 +423,10 @@ void StdMeshers_Penta_3D::MakeNodes()
       aCoords.SetCoord(aX, aY, aZ);
       //
       //   suporting shape ID
-      bIsUpperLayer=(i==(myISize-1));
       ShapeSupportID(bIsUpperLayer, aBNSSID, aSSID);
-      if (myErrorStatus) {
+      if (!myErrorStatus->IsOK()) {
         MESSAGE("StdMeshers_Penta_3D::MakeNodes() ");
-       return;
+        return;
       }
       //
       aTN.SetShapeSupportID(aSSID);
@@ -407,77 +434,74 @@ void StdMeshers_Penta_3D::MakeNodes()
       aTN.SetBaseNodeID(iBNID);
       //
       if (aSSID!=SMESH_Block::ID_NONE){
-       // try to find the node
-       const TopoDS_Shape& aS=myBlock.Shape((int)aSSID);
-       FindNodeOnShape(aS, aCoords, i, aTN);
+        // try to find the node
+        const TopoDS_Shape& aS=myBlock.Shape((int)aSSID);
+        FindNodeOnShape(aS, aCoords, i, aTN);
       }
       else{
-       // create node and get it id
-       CreateNode (bIsUpperLayer, aCoords, aTN);
+        // create node and get its id
+        CreateNode (bIsUpperLayer, aCoords, aTN);
         //
         if ( bIsUpperLayer ) {
           const SMDS_MeshNode* n = aTN.Node();
           myShapeXYZ[ SMESH_Block::ID_Fxy1 ].SetCoord( n->X(), n->Y(), n->Z() );
+          // set node on top face:
+          // find UV parameter for this node
+          //              UV_Ex11
+          //   UV011+-----+----------+UV111
+          //        |                |
+          //        |                |
+          // UV_E0y1+     +node      +UV_E1y1
+          //        |                |
+          //        |                |
+          //        |                |
+          //   UV001+-----+----------+UV101
+          //              UV_Ex01
+          gp_Pnt2d aP;
+          double u = aCoords.X(), v = aCoords.Y();
+          double u1 = ( 1. - u ), v1 = ( 1. - v );
+          aP.ChangeCoord()  = UV_Ex01 * v1;
+          aP.ChangeCoord() += UV_Ex11 * v;
+          aP.ChangeCoord() += UV_E0y1 * u1;
+          aP.ChangeCoord() += UV_E1y1 * u;
+          aP.ChangeCoord() -= UV001 * u1 * v1;
+          aP.ChangeCoord() -= UV101 * u  * v1;
+          aP.ChangeCoord() -= UV011 * u1 * v;
+          aP.ChangeCoord() -= UV111 * u  * v;
+          meshDS->SetNodeOnFace((SMDS_MeshNode*)n, topfaceID, aP.X(), aP.Y());
         }
       }
-      if (myErrorStatus) {
+      if (!myErrorStatus->IsOK()) {
         MESSAGE("StdMeshers_Penta_3D::MakeNodes() ");
-       return;
+        return;
       }
       //
       myTNodes[ij]=aTN;
     }
   }
-  //DEB
-  /*
-  {
-    int iSSID, iBNID, aID;
-    //
-    for (i=0; i<myISize; ++i) {
-      printf(" Layer# %d\n", i);
-      for (j=0; j<myJSize; ++j) {
-       ij=i*myJSize+j; 
-       const StdMeshers_TNode& aTN=myTNodes[ij];
-       //const StdMeshers_TNode& aTN=aTNodes[ij];
-       const gp_XYZ& aXYZ=aTN.NormCoord();
-       iSSID=aTN.ShapeSupportID();
-       iBNID=aTN.BaseNodeID();
-       //
-       const SMDS_MeshNode* aNode=aTN.Node();
-       aID=aNode->GetID(); 
-       aX=aNode->X();
-       aY=aNode->Y();
-       aZ=aNode->Z();
-       printf("*** j:%d BNID#%d iSSID:%d ID:%d { %lf %lf %lf },  { %lf %lf %lf }\n",
-              j,  iBNID, iSSID, aID, aXYZ.X(),  aXYZ.Y(), aXYZ.Z(), aX, aY, aZ);
-      }
-    }
-  }
-  */
-  //DEB t
 }
+
+
 //=======================================================================
 //function : FindNodeOnShape
 //purpose  : 
 //=======================================================================
+
 void StdMeshers_Penta_3D::FindNodeOnShape(const TopoDS_Shape& aS,
-                                         const gp_XYZ&       aParams,
+                                          const gp_XYZ&       aParams,
                                           const int           z,
-                                         StdMeshers_TNode&   aTN)
+                                          StdMeshers_TNode&   aTN)
 {
-  myErrorStatus=0;
-  //
   double aX, aY, aZ, aD, aTol2, minD;
   gp_Pnt aP1, aP2;
   //
-  SMESH_Mesh* pMesh=GetMesh();
-  aTol2=myTol3D*myTol3D;
+  SMESH_Mesh* pMesh = GetMesh();
+  aTol2 = myTol3D*myTol3D;
   minD = 1.e100;
-  SMDS_MeshNode* pNode=NULL;
+  SMDS_MeshNode* pNode = NULL;
   //
   if ( aS.ShapeType() == TopAbs_FACE ||
-       aS.ShapeType() == TopAbs_EDGE )
-  {
+       aS.ShapeType() == TopAbs_EDGE ) {
     // find a face ID to which aTN belongs to
     int faceID;
     if ( aS.ShapeType() == TopAbs_FACE )
@@ -492,13 +516,13 @@ void StdMeshers_Penta_3D::FindNodeOnShape(const TopoDS_Shape& aS,
     }
     ASSERT( SMESH_Block::IsFaceID( faceID ));
     int fIndex = SMESH_Block::ShapeIndex( faceID );
-    StdMeshers_IJNodeMap & ijNodes= myWallNodesMaps[ fIndex ];
+    StdMeshers_IJNodeMap & ijNodes = myWallNodesMaps[ fIndex ];
     // look for a base node in ijNodes
     const SMDS_MeshNode* baseNode = pMesh->GetMeshDS()->FindNode( aTN.BaseNodeID() );
     StdMeshers_IJNodeMap::const_iterator par_nVec = ijNodes.begin();
     for ( ; par_nVec != ijNodes.end(); par_nVec++ )
       if ( par_nVec->second[ 0 ] == baseNode ) {
-        pNode=(SMDS_MeshNode*)par_nVec->second.at( z );
+        pNode = (SMDS_MeshNode*)par_nVec->second.at( z );
         aTN.SetNode(pNode);
         return;
       }
@@ -510,6 +534,8 @@ void StdMeshers_Penta_3D::FindNodeOnShape(const TopoDS_Shape& aS,
     pMesh->GetSubMeshContaining(aS)->GetSubMeshDS()->GetNodes();
   while(ite->more()) {
     const SMDS_MeshNode* aNode = ite->next();
+    if(myTool->IsMedium(aNode))
+      continue;
     aX=aNode->X();
     aY=aNode->Y();
     aZ=aNode->Z();
@@ -532,6 +558,7 @@ void StdMeshers_Penta_3D::FindNodeOnShape(const TopoDS_Shape& aS,
   //myErrorStatus=11; // can not find the node;
 }
 
+
 //=======================================================================
 //function : SetHorizEdgeXYZ
 //purpose  : 
@@ -543,12 +570,15 @@ double StdMeshers_Penta_3D::SetHorizEdgeXYZ(const gp_XYZ&                  aBase
                                             vector<const SMDS_MeshNode*>*& aCol2)
 {
   // find base and top edges of the face
+  enum { BASE = 0, TOP };
   vector< int > edgeVec; // 0-base, 1-top
   SMESH_Block::GetFaceEdgesIDs( aFaceID, edgeVec );
   //
-  int coord = SMESH_Block::GetCoordIndOnEdge( edgeVec[ 0 ] );
+  int coord = SMESH_Block::GetCoordIndOnEdge( edgeVec[ BASE ] );
+  bool isForward = myBlock.IsForwadEdge( edgeVec[ BASE ] );
+
   double param = aBaseNodeParams.Coord( coord );
-  if ( !myBlock.IsForwadEdge( edgeVec[ 0 ] ))
+  if ( !isForward)
     param = 1. - param;
   //
   // look for columns around param
@@ -567,47 +597,61 @@ double StdMeshers_Penta_3D::SetHorizEdgeXYZ(const gp_XYZ&                  aBase
   aCol1 = & par_nVec_1->second;
   aCol2 = & par_nVec_2->second;
 
-  // base edge
-  const SMDS_MeshNode* n1 = aCol1->front();
-  const SMDS_MeshNode* n2 = aCol2->front();
-  gp_XYZ xyz1( n1->X(), n1->Y(), n1->Z() ), xyz2( n2->X(), n2->Y(), n2->Z() );
-  myShapeXYZ[ edgeVec[ 0 ] ] = ( 1. - r ) * xyz1 + r * xyz2;
-
   // top edge
-  n1 = aCol1->back();
-  n2 = aCol2->back();
-  xyz1.SetCoord( n1->X(), n1->Y(), n1->Z() );
-  xyz2.SetCoord( n2->X(), n2->Y(), n2->Z() );
-  myShapeXYZ[ edgeVec[ 1 ] ] = ( 1. - r ) * xyz1 + r * xyz2;
+  if (1) {
+    // this variant is better for cases with curved edges and
+    // different nodes distribution on top and base edges
+    const SMDS_MeshNode* n1 = aCol1->back();
+    const SMDS_MeshNode* n2 = aCol2->back();
+    gp_XYZ xyz1( n1->X(), n1->Y(), n1->Z() );
+    gp_XYZ xyz2( n2->X(), n2->Y(), n2->Z() );
+    myShapeXYZ[ edgeVec[ 1 ] ] = ( 1. - r ) * xyz1 + r * xyz2;
+  }
+  else {
+    // this variant is better for other cases
+    //   SMESH_MesherHelper helper( *GetMesh() );
+    //   const TopoDS_Edge & edge = TopoDS::Edge( myBlock.Shape( edgeVec[ TOP ]));
+    //   double u1 = helper.GetNodeU( edge, n1 );
+    //   double u2 = helper.GetNodeU( edge, n2 );
+    //   double u = ( 1. - r ) * u1 + r * u2;
+    //   gp_XYZ topNodeParams;
+    //   myBlock.Block().EdgeParameters( edgeVec[ TOP ], u, topNodeParams );
+    //   myBlock.Block().EdgePoint( edgeVec[ TOP ],
+    //                              topNodeParams,
+    //                              myShapeXYZ[ edgeVec[ TOP ]]);
+  }
 
+  // base edge
+  myBlock.Block().EdgePoint( edgeVec[ BASE ],
+                             aBaseNodeParams,
+                             myShapeXYZ[ edgeVec[ BASE ]]);
   return r;
 }
 
+
 //=======================================================================
 //function : MakeVolumeMesh
 //purpose  : 
 //=======================================================================
 void StdMeshers_Penta_3D::MakeVolumeMesh()
 {
-  myErrorStatus=0;
-  //
   int i, j, ij, ik, i1, i2, aSSID; 
   //
-  SMESH_Mesh*   pMesh =GetMesh();
-  SMESHDS_Mesh* meshDS=pMesh->GetMeshDS();
+  SMESH_Mesh*   pMesh = GetMesh();
+  SMESHDS_Mesh* meshDS = pMesh->GetMeshDS();
   //
   int shapeID = meshDS->ShapeToIndex( myShape );
   //
   // 1. Set Node In Volume
-  ik=myISize-1;
+  ik = myISize-1;
   for (i=1; i<ik; ++i){
     for (j=0; j<myJSize; ++j){
       ij=i*myJSize+j;
-      const StdMeshers_TNode& aTN=myTNodes[ij];
+      const StdMeshers_TNode& aTN = myTNodes[ij];
       aSSID=aTN.ShapeSupportID();
       if (aSSID==SMESH_Block::ID_NONE) {
-       SMDS_MeshNode* aNode=(SMDS_MeshNode*)aTN.Node();
-       meshDS->SetNodeInVolume(aNode, shapeID);
+        SMDS_MeshNode* aNode = (SMDS_MeshNode*)aTN.Node();
+        meshDS->SetNodeInVolume(aNode, shapeID);
       }
     }
   }
@@ -621,44 +665,44 @@ void StdMeshers_Penta_3D::MakeVolumeMesh()
   const TopoDS_Face& aFxy0=
     TopoDS::Face(myBlock.Shape(SMESH_Block::ID_Fxy0));
   SMESH_subMesh *aSubMesh0 = pMesh->GetSubMeshContaining(aFxy0);
-  SMESHDS_SubMesh *aSM0=aSubMesh0->GetSubMeshDS();
+  SMESHDS_SubMesh *aSM0 = aSubMesh0->GetSubMeshDS();
   //
-  itf=aSM0->GetElements();
+  itf = aSM0->GetElements();
   while(itf->more()) {
-    const SMDS_MeshElement* pE0=itf->next();
+    const SMDS_MeshElement* pE0 = itf->next();
     //
     int nbFaceNodes = pE0->NbNodes();
+    if(myCreateQuadratic)
+      nbFaceNodes = nbFaceNodes/2;
     if ( aN.size() < nbFaceNodes * 2 )
       aN.resize( nbFaceNodes * 2 );
     //
-    k=0;
-    aItNodes=pE0->nodesIterator();
-    while (aItNodes->more()) {
-      const SMDS_MeshElement* pNode=aItNodes->next();
-      aID0=pNode->GetID();
-      aJ[k]=GetIndexOnLayer(aID0);
-      if (myErrorStatus) {
+    for ( k=0; k<nbFaceNodes; ++k ) {
+      const SMDS_MeshNode* pNode = pE0->GetNode(k);
+//       if(myTool->IsMedium(pNode))
+//         continue;
+      aID0 = pNode->GetID();
+      aJ[k] = GetIndexOnLayer(aID0);
+      if (!myErrorStatus->IsOK()) {
         MESSAGE("StdMeshers_Penta_3D::MakeVolumeMesh");
-       return;
+        return;
       }
-      //
-      ++k;
     }
     //
     bool forward = true;
-    for (i=0; i<ik; ++i){
+    for (i=0; i<ik; ++i) {
       i1=i;
       i2=i+1;
       for(j=0; j<nbFaceNodes; ++j) {
-       ij=i1*myJSize+aJ[j];
-       const StdMeshers_TNode& aTN1=myTNodes[ij];
-       const SMDS_MeshNode* aN1=aTN1.Node();
-       aN[j]=aN1;
-       //
-       ij=i2*myJSize+aJ[j];
-       const StdMeshers_TNode& aTN2=myTNodes[ij];
-       const SMDS_MeshNode* aN2=aTN2.Node();
-       aN[j+nbFaceNodes]=aN2;
+        ij = i1*myJSize+aJ[j];
+        const StdMeshers_TNode& aTN1 = myTNodes[ij];
+        const SMDS_MeshNode* aN1 = aTN1.Node();
+        aN[j]=aN1;
+        //
+        ij=i2*myJSize+aJ[j];
+        const StdMeshers_TNode& aTN2 = myTNodes[ij];
+        const SMDS_MeshNode* aN2 = aTN2.Node();
+        aN[j+nbFaceNodes] = aN2;
       }
       // check if volume orientation will be ok
       if ( i == 0 ) {
@@ -685,20 +729,30 @@ void StdMeshers_Penta_3D::MakeVolumeMesh()
       SMDS_MeshVolume* aV = 0;
       switch ( nbFaceNodes ) {
       case 3:
-        if ( forward )
-          aV = meshDS->AddVolume(aN[0], aN[1], aN[2],
-                                 aN[3], aN[4], aN[5]);
-        else
-          aV = meshDS->AddVolume(aN[0], aN[2], aN[1],
-                                 aN[3], aN[5], aN[4]);
+        if ( forward ) {
+          //aV = meshDS->AddVolume(aN[0], aN[1], aN[2],
+          //                       aN[3], aN[4], aN[5]);
+          aV = myTool->AddVolume(aN[0], aN[1], aN[2], aN[3], aN[4], aN[5]);
+        }
+        else {
+          //aV = meshDS->AddVolume(aN[0], aN[2], aN[1],
+          //                       aN[3], aN[5], aN[4]);
+          aV = myTool->AddVolume(aN[0], aN[2], aN[1], aN[3], aN[5], aN[4]);
+        }
         break;
       case 4:
-        if ( forward )
-          aV = meshDS->AddVolume(aN[0], aN[1], aN[2], aN[3],
+        if ( forward ) {
+          //aV = meshDS->AddVolume(aN[0], aN[1], aN[2], aN[3],
+          //                       aN[4], aN[5], aN[6], aN[7]);
+          aV = myTool->AddVolume(aN[0], aN[1], aN[2], aN[3],
                                  aN[4], aN[5], aN[6], aN[7]);
-        else
-          aV = meshDS->AddVolume(aN[0], aN[3], aN[2], aN[1],
+        }
+        else {
+          //aV = meshDS->AddVolume(aN[0], aN[3], aN[2], aN[1],
+          //                       aN[4], aN[7], aN[6], aN[5]);
+          aV = myTool->AddVolume(aN[0], aN[3], aN[2], aN[1],
                                  aN[4], aN[7], aN[6], aN[5]);
+        }
         break;
       default:
         continue;
@@ -714,8 +768,6 @@ void StdMeshers_Penta_3D::MakeVolumeMesh()
 //=======================================================================
 void StdMeshers_Penta_3D::MakeMeshOnFxy1()
 {
-  myErrorStatus=0;
-  //
   int aID0, aJ, aLevel, ij, aNbNodes, k;
   //
   SMDS_NodeIteratorPtr itn;
@@ -727,89 +779,89 @@ void StdMeshers_Penta_3D::MakeMeshOnFxy1()
   const TopoDS_Face& aFxy1=
     TopoDS::Face(myBlock.Shape(SMESH_Block::ID_Fxy1));
   //
-  SMESH_Mesh* pMesh=GetMesh();
+  SMESH_Mesh* pMesh = GetMesh();
   SMESHDS_Mesh * meshDS = pMesh->GetMeshDS();
   //
+  SMESH_subMesh *aSubMesh1 = pMesh->GetSubMeshContaining(aFxy1);
   SMESH_subMesh *aSubMesh0 = pMesh->GetSubMeshContaining(aFxy0);
-  SMESHDS_SubMesh *aSM0=aSubMesh0->GetSubMeshDS();
+  SMESHDS_SubMesh *aSM0 = aSubMesh0->GetSubMeshDS();
   //
   // set nodes on aFxy1
-  aLevel=myISize-1;
-  itn=aSM0->GetNodes();
-  aNbNodes=aSM0->NbNodes();
+  aLevel = myISize-1;
+  itn = aSM0->GetNodes();
+  aNbNodes = aSM0->NbNodes();
   //printf("** aNbNodes=%d\n", aNbNodes);
-  while(itn->more()) {
-    const SMDS_MeshNode* aN0=itn->next();
-    aID0=aN0->GetID();
-    aJ=GetIndexOnLayer(aID0);
-    if (myErrorStatus) {
-      MESSAGE("StdMeshers_Penta_3D::MakeMeshOnFxy1() ");
-      return;
-    }
-    //
-    ij=aLevel*myJSize+aJ;
-    const StdMeshers_TNode& aTN1=myTNodes[ij];
-    SMDS_MeshNode* aN1=(SMDS_MeshNode*)aTN1.Node();
-    //
-    meshDS->SetNodeOnFace(aN1, aFxy1);
-  }
+  myTool->SetSubShape( aFxy1 ); // to set medium nodes to aFxy1
   //
   // set elements on aFxy1
   vector<const SMDS_MeshNode*> aNodes1;
   //
-  itf=aSM0->GetElements();
+  itf = aSM0->GetElements();
   while(itf->more()) {
-    const SMDS_MeshElement * pE0=itf->next();
-    aElementType=pE0->GetType();
+    const SMDS_MeshElement* pE0 = itf->next();
+    aElementType = pE0->GetType();
     if (!aElementType==SMDSAbs_Face) {
       continue;
     }
-    aNbNodes=pE0->NbNodes();
-//     if (aNbNodes!=3) {
-//       continue;
-//     }
+    aNbNodes = pE0->NbNodes();
+    if(myCreateQuadratic)
+      aNbNodes = aNbNodes/2;
     if ( aNodes1.size() < aNbNodes )
       aNodes1.resize( aNbNodes );
     //
-    k=aNbNodes-1; // reverse a face
-    aItNodes=pE0->nodesIterator();
+    k = aNbNodes-1; // reverse a face
+    aItNodes = pE0->nodesIterator();
     while (aItNodes->more()) {
-      const SMDS_MeshElement* pNode=aItNodes->next();
-      aID0=pNode->GetID();
-      aJ=GetIndexOnLayer(aID0);
-      if (myErrorStatus) {
+      const SMDS_MeshNode* pNode =
+        static_cast<const SMDS_MeshNode*> (aItNodes->next());
+      if(myTool->IsMedium(pNode))
+        continue;
+      aID0 = pNode->GetID();
+      aJ = GetIndexOnLayer(aID0);
+      if (!myErrorStatus->IsOK()) {
         MESSAGE("StdMeshers_Penta_3D::MakeMeshOnFxy1() ");
-       return;
+        return;
       }
       //
-      ij=aLevel*myJSize+aJ;
-      const StdMeshers_TNode& aTN1=myTNodes[ij];
-      const SMDS_MeshNode* aN1=aTN1.Node();
-      aNodes1[k]=aN1;
+      ij = aLevel*myJSize + aJ;
+      const StdMeshers_TNode& aTN1 = myTNodes[ij];
+      const SMDS_MeshNode* aN1 = aTN1.Node();
+      aNodes1[k] = aN1;
       --k;
     }
     SMDS_MeshFace * face = 0;
     switch ( aNbNodes ) {
     case 3:
-      face = meshDS->AddFace(aNodes1[0], aNodes1[1], aNodes1[2]);
+      face = myTool->AddFace(aNodes1[0], aNodes1[1], aNodes1[2]);
       break;
     case 4:
-      face = meshDS->AddFace(aNodes1[0], aNodes1[1], aNodes1[2], aNodes1[3]);
+      face = myTool->AddFace(aNodes1[0], aNodes1[1], aNodes1[2], aNodes1[3]);
       break;
     default:
       continue;
     }
     meshDS->SetMeshElementOnShape(face, aFxy1);
   }
+  myTool->SetSubShape( myShape );
+
+  // update compute state of top face submesh
+  aSubMesh1->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
+
+  // assure that mesh on the top face will be cleaned when it is cleaned
+  // on the bottom face
+  SMESH_subMesh* volSM = pMesh->GetSubMesh( myTool->GetSubShape() );
+  volSM->SetEventListener( new SMESH_subMeshEventListener(true, // deletable by SMESH_subMesh
+                                                          "StdMeshers_Penta_3D"),
+                           SMESH_subMeshEventListenerData::MakeData( aSubMesh1 ),
+                           aSubMesh0 ); // translate CLEAN event of aSubMesh0 to aSubMesh1
 }
+
 //=======================================================================
 //function : ClearMeshOnFxy1
 //purpose  : 
 //=======================================================================
 void StdMeshers_Penta_3D::ClearMeshOnFxy1()
 {
-  myErrorStatus=0;
-  //
   SMESH_subMesh* aSubMesh;
   SMESH_Mesh* pMesh=GetMesh();
   //
@@ -825,19 +877,19 @@ void StdMeshers_Penta_3D::ClearMeshOnFxy1()
 //=======================================================================
 int StdMeshers_Penta_3D::GetIndexOnLayer(const int aID)
 {
-  myErrorStatus=0;
-  //
   int j=-1;
   StdMeshers_IteratorOfDataMapOfIntegerInteger aMapIt;
   //
   aMapIt=myConnectingMap.find(aID);
   if (aMapIt==myConnectingMap.end()) {
-    myErrorStatus=200;
+    myErrorStatus->myName    = 200;
+    myErrorStatus->myComment = "Internal error of StdMeshers_Penta_3D";
     return j;
   }
   j=(*aMapIt).second;
   return j;
 }
+
 //=======================================================================
 //function : MakeConnectingMap
 //purpose  : 
@@ -852,17 +904,15 @@ void StdMeshers_Penta_3D::MakeConnectingMap()
     myConnectingMap[aBNID]=j;
   }
 }
+
 //=======================================================================
 //function : CreateNode
 //purpose  : 
 //=======================================================================
 void StdMeshers_Penta_3D::CreateNode(const bool bIsUpperLayer,
-                                    const gp_XYZ& aParams,
-                                    StdMeshers_TNode& aTN)
+                                     const gp_XYZ& aParams,
+                                     StdMeshers_TNode& aTN)
 {
-  myErrorStatus=0;
-  //
-  // int iErr;
   double aX, aY, aZ;
   //
   gp_Pnt aP;
@@ -870,17 +920,16 @@ void StdMeshers_Penta_3D::CreateNode(const bool bIsUpperLayer,
   SMDS_MeshNode* pNode=NULL; 
   aTN.SetNode(pNode);  
   //
-//   if (bIsUpperLayer) {
-//     // point on face Fxy1
-//     const TopoDS_Shape& aS=myBlock.Shape(SMESH_Block::ID_Fxy1);
-//     myBlock.Point(aParams, aS, aP);
-//   }
-//   else {
-//     // point inside solid
-//     myBlock.Point(aParams, aP);
-//   }
-  if (bIsUpperLayer)
-  {
+  //   if (bIsUpperLayer) {
+  //     // point on face Fxy1
+  //     const TopoDS_Shape& aS=myBlock.Shape(SMESH_Block::ID_Fxy1);
+  //     myBlock.Point(aParams, aS, aP);
+  //   }
+  //   else {
+  //     // point inside solid
+  //     myBlock.Point(aParams, aP);
+  //   }
+  if (bIsUpperLayer) {
     double u = aParams.X(), v = aParams.Y();
     double u1 = ( 1. - u ), v1 = ( 1. - v );
     aP.ChangeCoord()  = myShapeXYZ[ SMESH_Block::ID_Ex01 ] * v1;
@@ -893,67 +942,67 @@ void StdMeshers_Penta_3D::CreateNode(const bool bIsUpperLayer,
     aP.ChangeCoord() -= myShapeXYZ[ SMESH_Block::ID_V011 ] * u1 * v;
     aP.ChangeCoord() -= myShapeXYZ[ SMESH_Block::ID_V111 ] * u  * v;
   }
-  else
-  {
+  else {
     SMESH_Block::ShellPoint( aParams, myShapeXYZ, aP.ChangeCoord() );
   }
   //
-//   iErr=myBlock.ErrorStatus();
-//   if (iErr) {
-//     myErrorStatus=12; // can not find the node point;
-//     return;
-//   }
+  //   iErr=myBlock.ErrorStatus();
+  //   if (iErr) {
+  //     myErrorStatus=12; // can not find the node point;
+  //     return;
+  //   }
   //
   aX=aP.X(); aY=aP.Y(); aZ=aP.Z(); 
   //
-  SMESH_Mesh* pMesh=GetMesh();
-  SMESHDS_Mesh* pMeshDS=pMesh->GetMeshDS();
+  SMESH_Mesh* pMesh = GetMesh();
+  SMESHDS_Mesh* pMeshDS = pMesh->GetMeshDS();
   //
   pNode = pMeshDS->AddNode(aX, aY, aZ);
+
   aTN.SetNode(pNode);
 }
+
 //=======================================================================
 //function : ShapeSupportID
 //purpose  : 
 //=======================================================================
 void StdMeshers_Penta_3D::ShapeSupportID(const bool bIsUpperLayer,
-                                        const SMESH_Block::TShapeID aBNSSID,
-                                        SMESH_Block::TShapeID& aSSID)
+                                         const SMESH_Block::TShapeID aBNSSID,
+                                         SMESH_Block::TShapeID& aSSID)
 {
-  myErrorStatus=0;
-  //
   switch (aBNSSID) {
-    case SMESH_Block::ID_V000:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V001 : SMESH_Block::ID_E00z;
-      break;
-    case SMESH_Block::ID_V100:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V101 : SMESH_Block::ID_E10z;
-      break; 
-    case SMESH_Block::ID_V110:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V111 : SMESH_Block::ID_E11z;
-      break;
-    case SMESH_Block::ID_V010:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V011 : SMESH_Block::ID_E01z;
-      break;
-    case SMESH_Block::ID_Ex00:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_Ex01 : SMESH_Block::ID_Fx0z;
-      break;
-    case SMESH_Block::ID_Ex10:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_Ex11 : SMESH_Block::ID_Fx1z;
-      break; 
-    case SMESH_Block::ID_E0y0:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_E0y1 : SMESH_Block::ID_F0yz;
-      break; 
-    case SMESH_Block::ID_E1y0:
-      aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_E1y1 : SMESH_Block::ID_F1yz;
-      break; 
-    case SMESH_Block::ID_Fxy0:
-      aSSID=SMESH_Block::ID_NONE;//(bIsUpperLayer) ?  Shape_ID_Fxy1 : Shape_ID_NONE;
-      break;   
-    default:
-      aSSID=SMESH_Block::ID_NONE;
-      myErrorStatus=10; // Can not find supporting shape ID
-      break;
+  case SMESH_Block::ID_V000:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V001 : SMESH_Block::ID_E00z;
+    break;
+  case SMESH_Block::ID_V100:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V101 : SMESH_Block::ID_E10z;
+    break; 
+  case SMESH_Block::ID_V110:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V111 : SMESH_Block::ID_E11z;
+    break;
+  case SMESH_Block::ID_V010:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_V011 : SMESH_Block::ID_E01z;
+    break;
+  case SMESH_Block::ID_Ex00:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_Ex01 : SMESH_Block::ID_Fx0z;
+    break;
+  case SMESH_Block::ID_Ex10:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_Ex11 : SMESH_Block::ID_Fx1z;
+    break; 
+  case SMESH_Block::ID_E0y0:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_E0y1 : SMESH_Block::ID_F0yz;
+    break; 
+  case SMESH_Block::ID_E1y0:
+    aSSID=(bIsUpperLayer) ?  SMESH_Block::ID_E1y1 : SMESH_Block::ID_F1yz;
+    break; 
+  case SMESH_Block::ID_Fxy0:
+    aSSID=SMESH_Block::ID_NONE;//(bIsUpperLayer) ?  Shape_ID_Fxy1 : Shape_ID_NONE;
+    break;   
+  default:
+    aSSID=SMESH_Block::ID_NONE;
+    myErrorStatus->myName=10; // Can not find supporting shape ID
+    myErrorStatus->myComment = "Internal error of StdMeshers_Penta_3D";
+    break;
   }
   return;
 }
@@ -963,8 +1012,6 @@ void StdMeshers_Penta_3D::ShapeSupportID(const bool bIsUpperLayer,
 //=======================================================================
 void StdMeshers_Penta_3D::MakeBlock()
 {
-  myErrorStatus=0;
-  //
   bool bFound;
   int i, j, iNbEV, iNbE, iErr, iCnt, iNbNodes, iNbF;
   //
@@ -980,32 +1027,193 @@ void StdMeshers_Penta_3D::MakeBlock()
   SMDSAbs_ElementType aElementType;
   SMESH_Mesh* pMesh=GetMesh();
   //
-  iCnt=0;
-  iNbF=aM.Extent();
+  iCnt = 0;
+  iNbF = aM.Extent();
   for (i=1; i<=iNbF; ++i) {
-    const TopoDS_Shape& aF=aM(i);
+    const TopoDS_Shape& aF = aM(i);
     SMESH_subMesh *aSubMesh = pMesh->GetSubMeshContaining(aF);
     ASSERT(aSubMesh);
-    SMESHDS_SubMesh *aSM=aSubMesh->GetSubMeshDS();
-    SMDS_ElemIteratorPtr itf=aSM->GetElements();
+    SMESHDS_SubMesh *aSM = aSubMesh->GetSubMeshDS();
+    SMDS_ElemIteratorPtr itf = aSM->GetElements();
     while(itf->more()) {
-      const SMDS_MeshElement * pElement=itf->next();
-      aElementType=pElement->GetType();
+      const SMDS_MeshElement * pElement = itf->next();
+      aElementType = pElement->GetType();
       if (aElementType==SMDSAbs_Face) {
-       iNbNodes=pElement->NbNodes();
-       if (iNbNodes==3) {
-         aFTr=aF;
-         ++iCnt;
-         if (iCnt>1) {
-            MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-           myErrorStatus=5; // more than one face has triangulation
-           return;
-         }
-         break; // next face
-       }
+        iNbNodes = pElement->NbNodes();
+        if ( iNbNodes==3 || (pElement->IsQuadratic() && iNbNodes==6) ) {
+          aFTr = aF;
+          ++iCnt;
+          if (iCnt>1) {
+            // \begin{E.A.}
+            // The current algorithm fails if there is more that one
+            // face wich contains triangles ...
+            // In that case, replace return by break to try another
+            // method (coded in "if (iCnt != 1) { ... }")
+            //
+            // MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
+            // myErrorStatus=5; // more than one face has triangulation
+            // return;
+            break;
+            // \end{E.A.}
+          }
+          break; // next face
+        }
       }
     }
   }
+  //
+  // \begin{E.A.}
+  // The current algorithm fails if "iCnt != 1", the case "iCnt == 0"
+  // was not reached 'cause it was not called from Hexa_3D ... Now it
+  // can occurs and in my opinion, it is the most common case.
+  //
+  if (iCnt != 1) {
+    // The suggested algorithm is the following :
+    //
+    // o Check that nb_of_faces == 6 and nb_of_edges == 12
+    //   then the shape is tologically equivalent to a box
+    // o In a box, there are three set of four // edges ...
+    //   In the cascade notation, it seems to be the edges
+    //   numbered : 
+    //     - 1, 3, 5, 7
+    //     - 2, 4, 6, 8
+    //     - 9, 10, 11, 12
+    // o For each one of this set, check if the four edges
+    //   have the same number of element.
+    // o If so, check if the "corresponding" // faces contains
+    //   only quads. It's the faces numbered:
+    //     - 1, 2, 3, 4
+    //     - 1, 2, 5, 6
+    //     - 3, 4, 5, 6
+    // o If so, check if the opposite edges of each // faces
+    //   have the same number of elements. It is the edges
+    //   numbered :
+    //     - 2 and 4, 6 and 8, 9 and 10, 11 and 12
+    //     - 1 and 3, 5 and 7, 9 and 11, 10 and 12
+    //     - 1 and 5, 3 and 7, 4 and 8, 2 and 6
+    // o If so, check if the two other faces have the same
+    //   number of elements. It is the faces numbered:
+    //     - 5, 6
+    //     - 3, 4
+    //     - 1, 2
+    //   This test should be improved to test if the nodes
+    //   of the two faces are really "en face".
+    // o If so, one of the two faces is a candidate to an extrusion,
+    //   It is the faces numbered :
+    //     - 5
+    //     - 3
+    //     - 1
+    // o Finally, if there is only one candidate, let do the
+    //   extrusion job for the corresponding face
+    //
+    int isOK = 0;
+    //
+    int iNbF = aM.Extent();
+    if (iNbF == 6) {
+      //
+      int nb_f1 = pMesh->GetSubMeshContaining(aM(1))->GetSubMeshDS()->NbElements();
+      int nb_f2 = pMesh->GetSubMeshContaining(aM(2))->GetSubMeshDS()->NbElements();
+      int nb_f3 = pMesh->GetSubMeshContaining(aM(3))->GetSubMeshDS()->NbElements();
+      int nb_f4 = pMesh->GetSubMeshContaining(aM(4))->GetSubMeshDS()->NbElements();
+      int nb_f5 = pMesh->GetSubMeshContaining(aM(5))->GetSubMeshDS()->NbElements();
+      int nb_f6 = pMesh->GetSubMeshContaining(aM(6))->GetSubMeshDS()->NbElements();
+      //
+      int has_only_quad_f1 = 1;
+      int has_only_quad_f2 = 1;
+      int has_only_quad_f3 = 1;
+      int has_only_quad_f4 = 1;
+      int has_only_quad_f5 = 1;
+      int has_only_quad_f6 = 1;
+      //
+      for (i=1; i<=iNbF; ++i) {
+        int ok = 1;
+        const TopoDS_Shape& aF = aM(i);
+        SMESH_subMesh *aSubMesh = pMesh->GetSubMeshContaining(aF);
+        SMESHDS_SubMesh *aSM = aSubMesh->GetSubMeshDS();
+        SMDS_ElemIteratorPtr itf = aSM->GetElements();
+        while(itf->more()) {
+          const SMDS_MeshElement * pElement = itf->next();
+          aElementType = pElement->GetType();
+          if (aElementType==SMDSAbs_Face) {
+            iNbNodes = pElement->NbNodes();
+            if ( iNbNodes!=4 ) {
+              ok = 0;
+              break ;
+            }
+          }
+        }
+        if (i==1) has_only_quad_f1 = ok ;
+        if (i==2) has_only_quad_f2 = ok ;
+        if (i==3) has_only_quad_f3 = ok ;
+        if (i==4) has_only_quad_f4 = ok ;
+        if (i==5) has_only_quad_f5 = ok ;
+        if (i==6) has_only_quad_f6 = ok ;
+      }
+      //
+      TopTools_IndexedMapOfShape aE;
+      TopExp::MapShapes(myShape, TopAbs_EDGE, aE);
+      int iNbE = aE.Extent();
+      if (iNbE == 12) {
+        //
+        int nb_e01 = pMesh->GetSubMeshContaining(aE(1))->GetSubMeshDS()->NbElements();
+        int nb_e02 = pMesh->GetSubMeshContaining(aE(2))->GetSubMeshDS()->NbElements();
+        int nb_e03 = pMesh->GetSubMeshContaining(aE(3))->GetSubMeshDS()->NbElements();
+        int nb_e04 = pMesh->GetSubMeshContaining(aE(4))->GetSubMeshDS()->NbElements();
+        int nb_e05 = pMesh->GetSubMeshContaining(aE(5))->GetSubMeshDS()->NbElements();
+        int nb_e06 = pMesh->GetSubMeshContaining(aE(6))->GetSubMeshDS()->NbElements();
+        int nb_e07 = pMesh->GetSubMeshContaining(aE(7))->GetSubMeshDS()->NbElements();
+        int nb_e08 = pMesh->GetSubMeshContaining(aE(8))->GetSubMeshDS()->NbElements();
+        int nb_e09 = pMesh->GetSubMeshContaining(aE(9))->GetSubMeshDS()->NbElements();
+        int nb_e10 = pMesh->GetSubMeshContaining(aE(10))->GetSubMeshDS()->NbElements();
+        int nb_e11 = pMesh->GetSubMeshContaining(aE(11))->GetSubMeshDS()->NbElements();
+        int nb_e12 = pMesh->GetSubMeshContaining(aE(12))->GetSubMeshDS()->NbElements();
+        //
+        int nb_ok = 0 ;
+        //
+        if ( (nb_e01==nb_e03) && (nb_e03==nb_e05) && (nb_e05==nb_e07) ) {
+          if ( has_only_quad_f1 && has_only_quad_f2 && has_only_quad_f3 && has_only_quad_f4 ) {
+            if ( (nb_e09==nb_e10) && (nb_e08==nb_e06) && (nb_e11==nb_e12) && (nb_e04==nb_e02) ) {
+              if (nb_f5==nb_f6) {
+                nb_ok += 1;
+                aFTr = aM(5);
+              }
+            }
+          }
+        }
+        if ( (nb_e02==nb_e04) && (nb_e04==nb_e06) && (nb_e06==nb_e08) ) {
+          if ( has_only_quad_f1 && has_only_quad_f2 && has_only_quad_f5 && has_only_quad_f6 ) {
+            if ( (nb_e01==nb_e03) && (nb_e10==nb_e12) && (nb_e05==nb_e07) && (nb_e09==nb_e11) ) {
+              if (nb_f3==nb_f4) {
+                nb_ok += 1;
+                aFTr = aM(3);
+              }
+            }
+          }
+        }
+        if ( (nb_e09==nb_e10) && (nb_e10==nb_e11) && (nb_e11==nb_e12) ) {
+          if ( has_only_quad_f3 && has_only_quad_f4 && has_only_quad_f5 && has_only_quad_f6 ) {
+            if ( (nb_e01==nb_e05) && (nb_e02==nb_e06) && (nb_e03==nb_e07) && (nb_e04==nb_e08) ) {
+              if (nb_f1==nb_f2) {
+                nb_ok += 1;
+                aFTr = aM(1);
+              }
+            }
+          }
+        }
+        //
+        if ( nb_ok == 1 ) {
+          isOK = 1;
+        }
+        //
+      }
+    }
+    if (!isOK) {
+      myErrorStatus->myName=5; // more than one face has triangulation
+      myErrorStatus->myComment="Incorrect input mesh";
+      return;
+    }
+  }
+  // \end{E.A.}
   // 
   // 1. Vetrices V00, V001;
   //
@@ -1013,10 +1221,12 @@ void StdMeshers_Penta_3D::MakeBlock()
   TopExp::MapShapesAndAncestors(myShape, TopAbs_VERTEX, TopAbs_EDGE, aMVES);
   //
   // 1.1 Base vertex V000
-  iNbE=aME.Extent();
-  if (iNbE!=4){
+  iNbE = aME.Extent();
+  if (iNbE!= NB_WALL_FACES ){
     MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-    myErrorStatus=7; // too few edges are in base face aFTr 
+    myErrorStatus->myName=7; // too few edges are in base face aFTr
+    myErrorStatus->myComment=SMESH_Comment("Not a quadrilateral face #")
+      <<pMesh->GetMeshDS()->ShapeToIndex( aFTr )<<": "<<iNbE<<" edges" ;
     return;
   }
   const TopoDS_Edge& aE1=TopoDS::Edge(aME(1));
@@ -1031,7 +1241,9 @@ void StdMeshers_Penta_3D::MakeBlock()
   iNbEV=aMEV.Extent();
   if (iNbEV!=3){
     MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-    myErrorStatus=7; // too few edges meet in base vertex 
+    myErrorStatus->myName=7; // too few edges meet in base vertex 
+    myErrorStatus->myComment=SMESH_Comment("3 edges must share vertex #")
+      <<pMesh->GetMeshDS()->ShapeToIndex( aV000 )<<" but there are "<<iNbEV<<" edges";
     return;
   }
   //
@@ -1045,18 +1257,20 @@ void StdMeshers_Penta_3D::MakeBlock()
       const TopoDS_Edge& aE=TopoDS::Edge(aEx);
       TopExp::Vertices(aE, aV[0], aV[1]);
       for (i=0; i<2; ++i) {
-       if (!aV[i].IsSame(aV000)) {
-         aV001=aV[i];
-         bFound=!bFound;
-         break;
-       }
+        if (!aV[i].IsSame(aV000)) {
+          aV001=aV[i];
+          bFound=!bFound;
+          break;
+        }
       }
     }
   }
   //
   if (!bFound) {
     MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-    myErrorStatus=8; // can not find reper V001 
+    myErrorStatus->myName=8; // can not find reper V001
+    myErrorStatus->myComment=SMESH_Comment("Can't find opposite vertex for vertex #")
+      <<pMesh->GetMeshDS()->ShapeToIndex( aV000 );
     return;
   }
   //DEB
@@ -1073,17 +1287,18 @@ void StdMeshers_Penta_3D::MakeBlock()
   iNbE=aME.Extent();
   if (iNbE!=1) {
     MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-    myErrorStatus=9; // number of shells in source shape !=1 
+    myErrorStatus->myName=9; // number of shells in source shape !=1
+    myErrorStatus->myComment=SMESH_Comment("Unexpected nb of shells ")<<iNbE;
     return;
   }
   //
   // 2. Load Block
   const TopoDS_Shell& aShell=TopoDS::Shell(aME(1));
   myBlock.Load(aShell, aV000, aV001);
-  iErr=myBlock.ErrorStatus();
+  iErr = myBlock.ErrorStatus();
   if (iErr) {
     MESSAGE("StdMeshers_Penta_3D::MakeBlock() ");
-    myErrorStatus=100; // SMESHBlock: Load operation failed
+    myErrorStatus=myBlock.GetError(); // SMESHBlock: Load operation failed
     return;
   }
 }
@@ -1093,8 +1308,6 @@ void StdMeshers_Penta_3D::MakeBlock()
 //=======================================================================
 void StdMeshers_Penta_3D::CheckData()
 {
-  myErrorStatus=0;
-  //
   int i, iNb;
   int iNbEx[]={8, 12, 6};
   //
@@ -1106,14 +1319,16 @@ void StdMeshers_Penta_3D::CheckData()
   //
   if (myShape.IsNull()){
     MESSAGE("StdMeshers_Penta_3D::CheckData() ");
-    myErrorStatus=2; // null shape
+    myErrorStatus->myName=2; // null shape
+    myErrorStatus->myComment="Null shape";
     return;
   }
   //
   aST=myShape.ShapeType();
   if (!(aST==TopAbs_SOLID || aST==TopAbs_SHELL)) {
     MESSAGE("StdMeshers_Penta_3D::CheckData() ");
-    myErrorStatus=3; // not compatible type of shape
+    myErrorStatus->myName=3; // not compatible type of shape
+    myErrorStatus->myComment=SMESH_Comment("Wrong shape type (TopAbs_ShapeEnum) ")<<aST;
     return;
   }
   //
@@ -1123,7 +1338,8 @@ void StdMeshers_Penta_3D::CheckData()
     iNb=aM.Extent();
     if (iNb!=iNbEx[i]){
       MESSAGE("StdMeshers_Penta_3D::CheckData() ");
-      myErrorStatus=4; // number of subshape is not compatible
+      myErrorStatus->myName=4; // number of sub-shape is not compatible
+      myErrorStatus->myComment="Wrong number of sub-shapes of a block";
       return;
     }
   }
@@ -1154,10 +1370,10 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
   bool rev1, CumOri = false;
   TopExp_Explorer exp( theFace, TopAbs_EDGE );
   int nbEdges = 0;
-  for ( ; exp.More(); exp.Next() )
-  {
-    if ( ++nbEdges > 4 )
+  for ( ; exp.More(); exp.Next() ) {
+    if ( ++nbEdges > NB_WALL_FACES ) {
       return false; // more than 4 edges in theFace
+    }
     TopoDS_Edge e = TopoDS::Edge( exp.Current() );
     if ( theBaseEdge.IsSame( e ))
       continue;
@@ -1174,8 +1390,9 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
     else
       e2 = e;
   }
-  if ( nbEdges < 4 )
-    return false; // lass than 4 edges in theFace
+  if ( nbEdges < NB_WALL_FACES ) {
+    return false; // less than 4 edges in theFace
+  }
 
   // submeshes corresponding to shapes
   SMESHDS_SubMesh* smFace = theMesh->MeshElements( theFace );
@@ -1188,7 +1405,7 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
   SMESHDS_SubMesh* smVft = theMesh->MeshElements( vft );
   if (!smFace || !smb || !smt || !sm1 || !sm2 || !smVfb || !smVlb || !smVft ) {
     MESSAGE( "NULL submesh " <<smFace<<" "<<smb<<" "<<smt<<" "<<
-            sm1<<" "<<sm2<<" "<<smVfb<<" "<<smVlb<<" "<<smVft);
+             sm1<<" "<<sm2<<" "<<smVfb<<" "<<smVlb<<" "<<smVft);
     return false;
   }
   if ( smb->NbNodes() != smt->NbNodes() || sm1->NbNodes() != sm2->NbNodes() ) {
@@ -1200,13 +1417,32 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
     return false;
   }
   if ( sm1->NbNodes() * smb->NbNodes() != smFace->NbNodes() ) {
-    MESSAGE( "Wrong nb face nodes: " <<
-            sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
-    return false;
+    // check quadratic case
+    if ( myCreateQuadratic ) {
+      int n1 = sm1->NbNodes()/2;
+      int n2 = smb->NbNodes()/2;
+      int n3 = sm1->NbNodes() - n1;
+      int n4 = smb->NbNodes() - n2;
+      int nf = sm1->NbNodes()*smb->NbNodes() - n3*n4;
+      if( nf != smFace->NbNodes() ) {
+        MESSAGE( "Wrong nb face nodes: " <<
+                 sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
+        return false;
+      }
+    }
+    else {
+      MESSAGE( "Wrong nb face nodes: " <<
+               sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
+      return false;
+    }
   }
   // IJ size
   int vsize = sm1->NbNodes() + 2;
   int hsize = smb->NbNodes() + 2;
+  if(myCreateQuadratic) {
+    vsize = vsize - sm1->NbNodes()/2 -1;
+    hsize = hsize - smb->NbNodes()/2 -1;
+  }
 
   // load nodes from theBaseEdge
 
@@ -1226,12 +1462,15 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
   double range = l - f;
   SMDS_NodeIteratorPtr nIt = smb->GetNodes();
   const SMDS_MeshNode* node;
-  while ( nIt->more() )
-  {
+  while ( nIt->more() ) {
     node = nIt->next();
+    if(myTool->IsMedium(node))
+      continue;
     const SMDS_EdgePosition* pos =
-      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
-    if ( !pos ) return false;
+      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition() );
+    if ( !pos ) {
+      return false;
+    }
     double u = ( pos->GetUParameter() - f ) / range;
     vector<const SMDS_MeshNode*> & nVec = theIJNodes[ u ];
     nVec.resize( vsize, nullNode );
@@ -1246,19 +1485,21 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
 
   map< double, const SMDS_MeshNode*> sortedNodes; // sort by param on edge
   nIt = sm1->GetNodes();
-  while ( nIt->more() )
-  {
+  while ( nIt->more() ) {
     node = nIt->next();
+    if(myTool->IsMedium(node))
+      continue;
     const SMDS_EdgePosition* pos =
-      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
-    if ( !pos ) return false;
+      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition() );
+    if ( !pos ) {
+      return false;
+    }
     sortedNodes.insert( make_pair( pos->GetUParameter(), node ));
   }
   loadedNodes.insert( nVecf[ vsize - 1 ] = smVft->GetNodes()->next() );
   map< double, const SMDS_MeshNode*>::iterator u_n = sortedNodes.begin();
   int row = rev1 ? vsize - 1 : 0;
-  for ( ; u_n != sortedNodes.end(); u_n++ )
-  {
+  for ( ; u_n != sortedNodes.end(); u_n++ ) {
     if ( rev1 ) row--;
     else        row++;
     loadedNodes.insert( nVecf[ row ] = u_n->second );
@@ -1267,7 +1508,7 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
   // try to load the rest nodes
 
   // get all faces from theFace
-  set<const SMDS_MeshElement*> allFaces, foundFaces;
+  TIDSortedElemSet allFaces, foundFaces;
   SMDS_ElemIteratorPtr eIt = smFace->GetElements();
   while ( eIt->more() ) {
     const SMDS_MeshElement* e = eIt->next();
@@ -1285,8 +1526,7 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
   StdMeshers_IJNodeMap::iterator par_nVec_2 = par_nVec_1;
   // loop on columns
   int col = 0;
-  for ( par_nVec_2++; par_nVec_2 != theIJNodes.end(); par_nVec_1++, par_nVec_2++ )
-  {
+  for ( par_nVec_2++; par_nVec_2 != theIJNodes.end(); par_nVec_1++, par_nVec_2++ ) {
     col++;
     row = 0;
     const SMDS_MeshNode* n1 = par_nVec_1->second[ row ];
@@ -1294,11 +1534,11 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
     const SMDS_MeshElement* face = 0;
     do {
       // look for a face by 2 nodes
-      face = SMESH_MeshEditor::FindFaceInSet( n1, n2, allFaces, foundFaces );
-      if ( face )
-      {
+      face = SMESH_MeshAlgos::FindFaceInSet( n1, n2, allFaces, foundFaces );
+      if ( face ) {
         int nbFaceNodes = face->NbNodes();
-        if ( nbFaceNodes > 4 ) {
+        if ( (!myCreateQuadratic && nbFaceNodes>4) ||
+             (myCreateQuadratic && nbFaceNodes>8) ) {
           MESSAGE(" Too many nodes in a face: " << nbFaceNodes );
           return false;
         }
@@ -1308,6 +1548,8 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
         eIt = face->nodesIterator() ;
         while ( !found && eIt->more() ) {
           node = static_cast<const SMDS_MeshNode*>( eIt->next() );
+          if(myTool->IsMedium(node))
+            continue;
           found = loadedNodes.insert( node ).second;
           if ( !found && node != n1 && node != n2 )
             n3 = node;
@@ -1320,18 +1562,21 @@ bool StdMeshers_Penta_3D::LoadIJNodes(StdMeshers_IJNodeMap & theIJNodes,
           par_nVec_2->second[ row ] = node;
           foundFaces.insert( face );
           n2 = node;
-          if ( nbFaceNodes == 4 )
+          if ( nbFaceNodes==4 || (myCreateQuadratic && nbFaceNodes==8) ) {
             n1 = par_nVec_1->second[ row ];
+          }
         }
-        else if (nbFaceNodes == 3 &&
-                 n3 == par_nVec_1->second[ row ] )
+        else if ( (nbFaceNodes==3 || (myCreateQuadratic && nbFaceNodes==6) )  &&
+                  n3 == par_nVec_1->second[ row ] ) {
           n1 = n3;
+        }
         else {
           MESSAGE( "Not quad mesh, column "<< col );
           return false;
         }
       }
-    } while ( face && n1 && n2 );
+    }
+    while ( face && n1 && n2 );
 
     if ( row < vsize - 1 ) {
       MESSAGE( "Too few nodes in column "<< col <<": "<< row+1);
@@ -1390,24 +1635,47 @@ int StdMeshers_SMESHBlock::ErrorStatus() const
 {
   return myErrorStatus;
 }
+
+//================================================================================
+/*!
+ * \brief Return problem description
+ */
+//================================================================================
+
+SMESH_ComputeErrorPtr StdMeshers_SMESHBlock::GetError() const
+{
+  SMESH_ComputeErrorPtr err = SMESH_ComputeError::New();
+  string & text = err->myComment;
+  switch ( myErrorStatus ) {
+  case 2:
+  case 3: text = "Internal error of StdMeshers_Penta_3D"; break; 
+  case 4: text = "Can't compute normalized parameters of a point inside a block"; break;
+  case 5: text = "Can't compute coordinates by normalized parameters inside a block"; break;
+  case 6: text = "Can't detect block sub-shapes. Not a block?"; break;
+  }
+  if (!text.empty())
+    err->myName = myErrorStatus;
+  return err;
+}
+
 //=======================================================================
 //function : Load
 //purpose  : 
 //=======================================================================
 void StdMeshers_SMESHBlock::Load(const TopoDS_Shell& theShell)
 {
-  
   TopoDS_Vertex aV000, aV001;
   //
   Load(theShell, aV000, aV001);
 }
+
 //=======================================================================
 //function : Load
 //purpose  : 
 //=======================================================================
 void StdMeshers_SMESHBlock::Load(const TopoDS_Shell& theShell,
-                                const TopoDS_Vertex& theV000,
-                                const TopoDS_Vertex& theV001)
+                                 const TopoDS_Vertex& theV000,
+                                 const TopoDS_Vertex& theV001)
 {
   myErrorStatus=0;
   //
@@ -1416,39 +1684,41 @@ void StdMeshers_SMESHBlock::Load(const TopoDS_Shell& theShell,
   bool bOk;
   //
   myShapeIDMap.Clear();  
-  bOk=myTBlock.LoadBlockShapes(myShell, theV000, theV001, myShapeIDMap);
+  bOk = myTBlock.LoadBlockShapes(myShell, theV000, theV001, myShapeIDMap);
   if (!bOk) {
-    myErrorStatus=2;
+    myErrorStatus=6;
     return;
   }
 }
+
 //=======================================================================
 //function : ComputeParameters
 //purpose  : 
 //=======================================================================
 void StdMeshers_SMESHBlock::ComputeParameters(const gp_Pnt& thePnt, 
-                                             gp_XYZ& theXYZ)
+                                              gp_XYZ& theXYZ)
 {
   ComputeParameters(thePnt, myShell, theXYZ);
 }
+
 //=======================================================================
 //function : ComputeParameters
 //purpose  : 
 //=======================================================================
 void StdMeshers_SMESHBlock::ComputeParameters(const gp_Pnt& thePnt,
-                                             const TopoDS_Shape& theShape,
-                                             gp_XYZ& theXYZ)
+                                              const TopoDS_Shape& theShape,
+                                              gp_XYZ& theXYZ)
 {
   myErrorStatus=0;
   //
   int aID;
   bool bOk;
   //
-  aID=ShapeID(theShape);
+  aID = ShapeID(theShape);
   if (myErrorStatus) {
     return;
   }
-  bOk=myTBlock.ComputeParameters(thePnt, theXYZ, aID);
+  bOk = myTBlock.ComputeParameters(thePnt, theXYZ, aID);
   if (!bOk) {
     myErrorStatus=4; // problems with computation Parameters 
     return;
@@ -1469,12 +1739,12 @@ void StdMeshers_SMESHBlock::ComputeParameters(const double& theU,
   int aID;
   bool bOk=false;
   //
-  aID=ShapeID(theShape);
+  aID = ShapeID(theShape);
   if (myErrorStatus) {
     return;
   }
   if ( SMESH_Block::IsEdgeID( aID ))
-      bOk=myTBlock.EdgeParameters( aID, theU, theXYZ );
+    bOk = myTBlock.EdgeParameters( aID, theU, theXYZ );
   if (!bOk) {
     myErrorStatus=4; // problems with computation Parameters 
     return;
@@ -1485,30 +1755,30 @@ void StdMeshers_SMESHBlock::ComputeParameters(const double& theU,
 //function : Point
 //purpose  : 
 //=======================================================================
- void StdMeshers_SMESHBlock::Point(const gp_XYZ& theParams,
-                                  gp_Pnt& aP3D)
+void StdMeshers_SMESHBlock::Point(const gp_XYZ& theParams, gp_Pnt& aP3D)
 {
   TopoDS_Shape aS;
   //
   Point(theParams, aS, aP3D);
 }
+
 //=======================================================================
 //function : Point
 //purpose  : 
 //=======================================================================
- void StdMeshers_SMESHBlock::Point(const gp_XYZ& theParams,
-                                  const TopoDS_Shape& theShape,
-                                  gp_Pnt& aP3D)
+void StdMeshers_SMESHBlock::Point(const gp_XYZ& theParams,
+                                  const TopoDS_Shape& theShape,
+                                  gp_Pnt& aP3D)
 {
-  myErrorStatus=0;
+  myErrorStatus = 0;
   //
   int aID;
-  bool bOk=false;
+  bool bOk = false;
   gp_XYZ aXYZ(99.,99.,99.);
   aP3D.SetXYZ(aXYZ);
   //
   if (theShape.IsNull()) {
-    bOk=myTBlock.ShellPoint(theParams, aXYZ);
+    bOk = myTBlock.ShellPoint(theParams, aXYZ);
   }
   //
   else {
@@ -1518,22 +1788,23 @@ void StdMeshers_SMESHBlock::ComputeParameters(const double& theU,
     }
     //
     if (SMESH_Block::IsVertexID(aID)) {
-      bOk=myTBlock.VertexPoint(aID, aXYZ);
+      bOk = myTBlock.VertexPoint(aID, aXYZ);
     }
     else if (SMESH_Block::IsEdgeID(aID)) {
-      bOk=myTBlock.EdgePoint(aID, theParams, aXYZ);
+      bOk = myTBlock.EdgePoint(aID, theParams, aXYZ);
     }
     //
     else if (SMESH_Block::IsFaceID(aID)) {
-      bOk=myTBlock.FacePoint(aID, theParams, aXYZ);
+      bOk = myTBlock.FacePoint(aID, theParams, aXYZ);
     }
   }
   if (!bOk) {
-    myErrorStatus=4; // problems with point computation 
+    myErrorStatus=5; // problems with point computation 
     return;
   }
   aP3D.SetXYZ(aXYZ);
 }
+
 //=======================================================================
 //function : ShapeID
 //purpose  : 
@@ -1561,6 +1832,7 @@ int StdMeshers_SMESHBlock::ShapeID(const TopoDS_Shape& theShape)
   myErrorStatus=2; // unknown shape;
   return aID;
 }
+
 //=======================================================================
 //function : Shape
 //purpose  : 
@@ -1580,3 +1852,114 @@ const TopoDS_Shape& StdMeshers_SMESHBlock::Shape(const int theID)
   const TopoDS_Shape& aS=myShapeIDMap.FindKey(theID);
   return aS;
 }
+
+
+//=======================================================================
+//function : Evaluate
+//purpose  : 
+//=======================================================================
+bool StdMeshers_Penta_3D::Evaluate(SMESH_Mesh& aMesh, 
+                                   const TopoDS_Shape& aShape,
+                                   MapShapeNbElems& aResMap)
+{
+  MESSAGE("StdMeshers_Penta_3D::Evaluate()");
+
+  // find face contains only triangles
+  vector < SMESH_subMesh * >meshFaces;
+  TopTools_SequenceOfShape aFaces;
+  int NumBase = 0, i = 0;
+  for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
+    i++;
+    aFaces.Append(exp.Current());
+    SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
+    meshFaces.push_back(aSubMesh);
+    MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i]);
+    if( anIt == aResMap.end() ) {
+      NumBase = 0;
+      break;
+    }
+    std::vector<int> aVec = (*anIt).second;
+    int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
+    int nbqua = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
+    if( nbtri>0 && nbqua==0 ) {
+      NumBase = i;
+    }
+  }
+
+  if(NumBase==0) {
+    std::vector<int> aResVec(SMDSEntity_Last);
+    for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
+    SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
+    aResMap.insert(std::make_pair(sm,aResVec));
+    myErrorStatus->myName    = COMPERR_ALGO_FAILED;
+    myErrorStatus->myComment = "Submesh can not be evaluated";
+    return false;
+  }
+
+  // find number of 1d elems for base face
+  int nb1d = 0;
+  TopTools_MapOfShape Edges1;
+  for (TopExp_Explorer exp(aFaces.Value(NumBase), TopAbs_EDGE); exp.More(); exp.Next()) {
+    Edges1.Add(exp.Current());
+    SMESH_subMesh *sm = aMesh.GetSubMesh(exp.Current());
+    if( sm ) {
+      MapShapeNbElemsItr anIt = aResMap.find(sm);
+      if( anIt == aResMap.end() ) continue;
+      std::vector<int> aVec = (*anIt).second;
+      nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
+    }
+  }
+  // find face opposite to base face
+  int OppNum = 0;
+  for(i=1; i<=6; i++) {
+    if(i==NumBase) continue;
+    bool IsOpposite = true;
+    for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
+      if( Edges1.Contains(exp.Current()) ) {
+        IsOpposite = false;
+        break;
+      }
+    }
+    if(IsOpposite) {
+      OppNum = i;
+      break;
+    }
+  }
+  // find number of 2d elems on side faces
+  int nb2d = 0;
+  for(i=1; i<=6; i++) {
+    if( i==OppNum || i==NumBase ) continue;
+    MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
+    if( anIt == aResMap.end() ) continue;
+    std::vector<int> aVec = (*anIt).second;
+    nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
+  }
+
+  MapShapeNbElemsItr anIt = aResMap.find( meshFaces[NumBase-1] );
+  std::vector<int> aVec = (*anIt).second;
+  int nb2d_face0 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
+  int nb0d_face0 = aVec[SMDSEntity_Node];
+
+  anIt = aResMap.find( meshFaces[OppNum-1] );
+  for(i=SMDSEntity_Node; i<SMDSEntity_Last; i++)
+    (*anIt).second[i] = aVec[i];
+
+  SMESH_MesherHelper aTool (aMesh);
+  bool _quadraticMesh = aTool.IsQuadraticSubMesh(aShape);
+
+  std::vector<int> aResVec(SMDSEntity_Last);
+  for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
+  if(_quadraticMesh) {
+    aResVec[SMDSEntity_Quad_Penta] = nb2d_face0 * ( nb2d/nb1d );
+    aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 );
+  }
+  else {
+    aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
+    aResVec[SMDSEntity_Penta] = nb2d_face0 * ( nb2d/nb1d );
+  }
+  SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
+  aResMap.insert(std::make_pair(sm,aResVec));
+
+  return true;
+}
+