Salome HOME
Re-fix regression of SALOME_TESTS/Grids/smesh/3D_mesh_Extrusion_01/B2
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_2D.cxx
index f8dd7501666ac0edcbb18001285b135c2e1469b3..8413288b8913ee001ac109dd0d310f5bb5a9a4b5 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2016  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
 #include "SMDS_EdgePosition.hxx"
 #include "SMDS_FacePosition.hxx"
 #include "SMESHDS_Hypothesis.hxx"
+#include "SMESHDS_Mesh.hxx"
 #include "SMESHDS_SubMesh.hxx"
 #include "SMESH_Block.hxx"
 #include "SMESH_Comment.hxx"
 #include "SMESH_Gen.hxx"
 #include "SMESH_Mesh.hxx"
+#include "SMESH_MeshAlgos.hxx"
 #include "SMESH_MesherHelper.hxx"
 #include "SMESH_Pattern.hxx"
 #include "SMESH_subMesh.hxx"
 #include "SMESH_subMeshEventListener.hxx"
 
-#include "utilities.h"
+#include <utilities.h>
 
 #include <BRepAdaptor_Surface.hxx>
+#include <BRepMesh_Delaun.hxx>
 #include <BRep_Tool.hxx>
 #include <Bnd_B2d.hxx>
 #include <GeomAPI_ProjectPointOnSurf.hxx>
 #include <TopExp_Explorer.hxx>
 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
 #include <TopTools_ListIteratorOfListOfShape.hxx>
+#include <TopTools_MapOfShape.hxx>
 #include <TopoDS.hxx>
+#include <TopoDS_Solid.hxx>
 #include <gp_Ax2.hxx>
 #include <gp_Ax3.hxx>
+#include <gp_GTrsf.hxx>
 
 
 using namespace std;
@@ -256,7 +262,7 @@ namespace {
   //================================================================================
   /*!
    * \brief find new nodes belonging to one free border of mesh on face
-    * \param sm - submesh on edge or vertex containg nodes to choose from
+    * \param sm - submesh on edge or vertex containing nodes to choose from
     * \param face - the face bound by the submesh
     * \param u2nodes - map to fill with nodes
     * \param seamNodes - set of found nodes
@@ -391,15 +397,16 @@ namespace {
   {
     double f,l;
     Handle(Geom2d_Curve) c1 = BRep_Tool::CurveOnSurface( E1, F, f, l );
-    gp_Pnt2d uvLast1 = c1->Value( E1.Orientation() == TopAbs_REVERSED ? f : l );
+    gp_Pnt2d uvFirst1 = c1->Value( f );
+    gp_Pnt2d uvLast1  = c1->Value( l );
 
     Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( E2, F, f, l );
-    gp_Pnt2d uvFirst2 = c2->Value( f );
-    gp_Pnt2d uvLast2  = c2->Value( l );
-    double tol2 = 1e-5 * uvLast2.SquareDistance( uvFirst2 );
+    gp_Pnt2d uvFirst2 = c2->Value( E2.Orientation() == TopAbs_REVERSED ? l : f );
+    double tol2 = Max( Precision::PConfusion() * Precision::PConfusion(),
+                       1e-5 * uvLast1.SquareDistance( uvFirst1 ));
 
-    return (( uvLast1.SquareDistance( uvFirst2 ) < tol2 ) ||
-            ( uvLast1.SquareDistance( uvLast2 ) < tol2 ));
+    return (( uvFirst2.SquareDistance( uvFirst1 ) < tol2 ) ||
+            ( uvFirst2.SquareDistance( uvLast1  ) < tol2 ));
   }
 
   //================================================================================
@@ -419,23 +426,25 @@ namespace {
                   TAssocTool::TNodeNodeMap&          src2tgtNodes,
                   bool&                              is1DComputed)
   {
-    SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
-    SMESHDS_Mesh* srcMeshDS = srcMesh->GetMeshDS();
-
     src2tgtNodes.clear();
 
     // get ordered src EDGEs
     TError err;
     srcWires = StdMeshers_FaceSide::GetFaceWires( srcFace, *srcMesh,/*skipMediumNodes=*/0, err);
-    if ( err && !err->IsOK() )
+    if (( err && !err->IsOK() ) ||
+        ( srcWires.empty() ))
       return err;
 
+    SMESH_MesherHelper srcHelper( *srcMesh );
+    srcHelper.SetSubShape( srcFace );
+
     // make corresponding sequence of tgt EDGEs
     tgtWires.resize( srcWires.size() );
     for ( size_t iW = 0; iW < srcWires.size(); ++iW )
     {
-      list< TopoDS_Edge > tgtEdges;
       StdMeshers_FaceSidePtr srcWire = srcWires[iW];
+
+      list< TopoDS_Edge > tgtEdges;
       TopTools_IndexedMapOfShape edgeMap; // to detect seam edges
       for ( int iE = 0; iE < srcWire->NbEdges(); ++iE )
       {
@@ -453,7 +462,8 @@ namespace {
           {
             list< TopoDS_Edge >::iterator eIt = tgtEdges.begin();
             std::advance( eIt, index-1 );
-            eIt->Reverse();
+            if ( are2dConnected( tgtEdges.back(), *eIt, tgtFace ))
+              eIt->Reverse();
           }
           else
           {
@@ -469,56 +479,58 @@ namespace {
             tgtE = nE.second;
         }
         tgtEdges.push_back( tgtE );
+      }
 
+      tgtWires[ iW ].reset( new StdMeshers_FaceSide( tgtFace, tgtEdges, tgtMesh,
+                                                     /*theIsForward = */ true,
+                                                     /*theIgnoreMediumNodes = */false));
+      StdMeshers_FaceSidePtr tgtWire = tgtWires[ iW ];
 
-        // Fill map of src to tgt nodes with nodes on edges
+      // Fill map of src to tgt nodes with nodes on edges
 
-        if ( srcMesh->GetSubMesh( srcE )->IsEmpty() ||
-             tgtMesh->GetSubMesh( tgtE )->IsEmpty() )
+      for ( int iE = 0; iE < srcWire->NbEdges(); ++iE )
+      {
+        if ( srcMesh->GetSubMesh( srcWire->Edge(iE) )->IsEmpty() ||
+             tgtMesh->GetSubMesh( tgtWire->Edge(iE) )->IsEmpty() )
         {
           // add nodes on VERTEXes for a case of not meshes EDGEs
-          const TopoDS_Shape&  srcV = SMESH_MesherHelper::IthVertex( 0, srcE );
-          const TopoDS_Shape&  tgtV = shape2ShapeMap( srcV, /*isSrc=*/true );
-          const SMDS_MeshNode* srcN = SMESH_Algo::VertexNode( TopoDS::Vertex( srcV ), srcMeshDS );
-          const SMDS_MeshNode* tgtN = SMESH_Algo::VertexNode( TopoDS::Vertex( tgtV ), tgtMeshDS );
+          const SMDS_MeshNode* srcN = srcWire->VertexNode( iE );
+          const SMDS_MeshNode* tgtN = tgtWire->VertexNode( iE );
           if ( srcN && tgtN )
             src2tgtNodes.insert( make_pair( srcN, tgtN ));
         }
         else
         {
-          const bool skipMediumNodes = true;
-          map< double, const SMDS_MeshNode* > srcNodes, tgtNodes;
-          if ( !SMESH_Algo::GetSortedNodesOnEdge( srcMeshDS, srcE, skipMediumNodes, srcNodes) ||
-               !SMESH_Algo::GetSortedNodesOnEdge( tgtMeshDS, tgtE, skipMediumNodes, tgtNodes ))
-            return SMESH_ComputeError::New( COMPERR_BAD_INPUT_MESH,
-                                            "Invalid node parameters on edges");
+          const bool skipMedium = true, isFwd = true;
+          StdMeshers_FaceSide srcEdge( srcFace, srcWire->Edge(iE), srcMesh, isFwd, skipMedium);
+          StdMeshers_FaceSide tgtEdge( tgtFace, tgtWire->Edge(iE), tgtMesh, isFwd, skipMedium);
+          
+          vector< const SMDS_MeshNode* > srcNodes = srcEdge.GetOrderedNodes();
+          vector< const SMDS_MeshNode* > tgtNodes = tgtEdge.GetOrderedNodes();
 
           if (( srcNodes.size() != tgtNodes.size() ) && tgtNodes.size() > 0 )
             return SMESH_ComputeError::New( COMPERR_BAD_INPUT_MESH,
                                             "Different number of nodes on edges");
           if ( !tgtNodes.empty() )
           {
-            map< double, const SMDS_MeshNode* >::iterator u_tn = tgtNodes.begin();
-            if ( srcE.Orientation() == tgtE.Orientation() )
-            {
-              map< double, const SMDS_MeshNode* >::iterator u_sn = srcNodes.begin();
-              for ( ; u_tn != tgtNodes.end(); ++u_tn, ++u_sn)
-                src2tgtNodes.insert( make_pair( u_sn->second, u_tn->second ));
-            }
-            else
+            vector< const SMDS_MeshNode* >::iterator tn = tgtNodes.begin();
+            //if ( srcWire->Edge(iE).Orientation() == tgtWire->Edge(iE).Orientation() )
             {
-              map< double, const SMDS_MeshNode* >::reverse_iterator u_sn = srcNodes.rbegin();
-              for ( ; u_tn != tgtNodes.end(); ++u_tn, ++u_sn)
-                src2tgtNodes.insert( make_pair( u_sn->second, u_tn->second ));
+              vector< const SMDS_MeshNode* >::iterator sn = srcNodes.begin();
+              for ( ; tn != tgtNodes.end(); ++tn, ++sn)
+                src2tgtNodes.insert( make_pair( *sn, *tn ));
             }
+            // else
+            // {
+            //   vector< const SMDS_MeshNode* >::reverse_iterator sn = srcNodes.rbegin();
+            //   for ( ; tn != tgtNodes.end(); ++tn, ++sn)
+            //     src2tgtNodes.insert( make_pair( *sn, *tn ));
+            // }
             is1DComputed = true;
           }
         }
       } // loop on EDGEs of a WIRE
 
-      tgtWires[ iW ].reset( new StdMeshers_FaceSide( tgtFace, tgtEdges, tgtMesh,
-                                                     /*theIsForward = */ true,
-                                                     /*theIgnoreMediumNodes = */false));
     } // loop on WIREs
 
     return TError();
@@ -549,18 +561,24 @@ namespace {
 
     // transformation to get location of target nodes from source ones
     StdMeshers_ProjectionUtils::TrsfFinder3D trsf;
+    bool trsfIsOK = false;
     if ( tgtFace.IsPartner( srcFace ))
     {
-      gp_Trsf srcTrsf = srcFace.Location();
-      gp_Trsf tgtTrsf = tgtFace.Location();
-      trsf.Set( srcTrsf.Inverted() * tgtTrsf );
+      gp_GTrsf srcTrsf = srcFace.Location().Transformation();
+      gp_GTrsf tgtTrsf = tgtFace.Location().Transformation();
+      gp_GTrsf t = srcTrsf.Inverted().Multiplied( tgtTrsf );
+      trsf.Set( t );
       // check
       gp_Pnt srcP = BRep_Tool::Pnt( srcWires[0]->FirstVertex() );
       gp_Pnt tgtP = BRep_Tool::Pnt( tgtWires[0]->FirstVertex() );
-      if ( tgtP.Distance( trsf.Transform( srcP )) > tol )
-        trsf.Set( tgtTrsf.Inverted() * srcTrsf );
+      trsfIsOK = ( tgtP.Distance( trsf.Transform( srcP )) < tol );
+      if ( !trsfIsOK )
+      {
+        trsf.Set( tgtTrsf.Inverted().Multiplied( srcTrsf ));
+        trsfIsOK = ( tgtP.Distance( trsf.Transform( srcP )) < tol );
+      }
     }
-    else
+    if ( !trsfIsOK )
     {
       // Try to find the 3D transformation
 
@@ -568,12 +586,13 @@ namespace {
       vector< gp_XYZ > srcPnts, tgtPnts;
       srcPnts.reserve( totNbSeg );
       tgtPnts.reserve( totNbSeg );
+      gp_XYZ srcBC( 0,0,0 ), tgtBC( 0,0,0 );
       for ( size_t iW = 0; iW < srcWires.size(); ++iW )
       {
         const double minSegLen = srcWires[iW]->Length() / totNbSeg;
         for ( int iE = 0; iE < srcWires[iW]->NbEdges(); ++iE )
         {
-          int nbSeg    = Max( 1, int( srcWires[iW]->EdgeLength( iE ) / minSegLen ));
+          size_t nbSeg = Max( 1, int( srcWires[iW]->EdgeLength( iE ) / minSegLen ));
           double srcU  = srcWires[iW]->FirstParameter( iE );
           double tgtU  = tgtWires[iW]->FirstParameter( iE );
           double srcDu = ( srcWires[iW]->LastParameter( iE )- srcU ) / nbSeg;
@@ -584,6 +603,8 @@ namespace {
             tgtPnts.push_back( tgtWires[iW]->Value3d( tgtU ).XYZ() );
             srcU += srcDu;
             tgtU += tgtDu;
+            srcBC += srcPnts.back();
+            tgtBC += tgtPnts.back();
           }
         }
       }
@@ -592,10 +613,11 @@ namespace {
 
       // check trsf
 
-      bool trsfIsOK = true;
       const int nbTestPnt = 20;
       const size_t  iStep = Max( 1, int( srcPnts.size() / nbTestPnt ));
       // check boundary
+      gp_Pnt trsfTgt = trsf.Transform( srcBC / srcPnts.size() );
+      trsfIsOK = ( trsfTgt.SquareDistance( tgtBC / tgtPnts.size() ) < tol*tol );
       for ( size_t i = 0; ( i < srcPnts.size() && trsfIsOK ); i += iStep )
       {
         gp_Pnt trsfTgt = trsf.Transform( srcPnts[i] );
@@ -606,8 +628,8 @@ namespace {
       {
         BRepAdaptor_Surface srcSurf( srcFace );
         gp_Pnt srcP =
-          srcSurf.Value( 0.5 * ( srcSurf.FirstUParameter() + srcSurf.LastUParameter() ),
-                         0.5 * ( srcSurf.FirstVParameter() + srcSurf.LastVParameter() ));
+          srcSurf.Value( 0.321 * ( srcSurf.FirstUParameter() + srcSurf.LastUParameter() ),
+                         0.123 * ( srcSurf.FirstVParameter() + srcSurf.LastVParameter() ));
         gp_Pnt tgtTrsfP = trsf.Transform( srcP );
         TopLoc_Location loc;
         GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( tgtFace, loc, 0.1*tol );
@@ -625,7 +647,7 @@ namespace {
     // Make new faces
 
     // prepare the helper to adding quadratic elements if necessary
-    helper.SetSubShape( tgtFace );
+    //helper.SetSubShape( tgtFace );
     helper.IsQuadraticSubMesh( tgtFace );
 
     SMESHDS_SubMesh* srcSubDS = srcMeshDS->MeshElements( srcFace );
@@ -721,13 +743,13 @@ namespace {
         {
         case SMDS_TOP_FACE:
         {
-          if ( nbOkPos < 10 ) break;
+          if ( nbOkPos > 10 ) break;
           gp_XY uv = helper.GetNodeUV( tgtFace, n ), uvBis = uv;
           if (( helper.CheckNodeUV( tgtFace, n, uv, tol )) &&
               (( uv - uvBis ).SquareModulus() < tol2d ))
             ++nbOkPos;
           else
-            nbOkPos = 0;
+            nbOkPos = -((int) src2tgtNodes.size() );
           break;
         }
         case SMDS_TOP_EDGE:
@@ -802,14 +824,14 @@ namespace {
         // find trsf
         const int totNbSeg = 50;
         vector< gp_XY > srcPnts, tgtPnts;
-        srcPnts.resize( totNbSeg );
-        tgtPnts.resize( totNbSeg );
+        srcPnts.reserve( totNbSeg );
+        tgtPnts.reserve( totNbSeg );
         for ( size_t iW = 0; iW < srcWires.size(); ++iW )
         {
           const double minSegLen = srcWires[iW]->Length() / totNbSeg;
           for ( int iE = 0; iE < srcWires[iW]->NbEdges(); ++iE )
           {
-            int nbSeg    = Max( 1, int( srcWires[iW]->EdgeLength( iE ) / minSegLen ));
+            size_t nbSeg = Max( 1, int( srcWires[iW]->EdgeLength( iE ) / minSegLen ));
             double srcU  = srcWires[iW]->FirstParameter( iE );
             double tgtU  = tgtWires[iW]->FirstParameter( iE );
             double srcDu = ( srcWires[iW]->LastParameter( iE )- srcU ) / nbSeg;
@@ -888,10 +910,9 @@ namespace {
           case SMDS_TOP_EDGE: {
             TopoDS_Shape srcEdge = srcHelper.GetSubShapeByNode( srcNode, srcHelper.GetMeshDS() );
             TopoDS_Edge  tgtEdge = TopoDS::Edge( shape2ShapeMap( srcEdge, /*isSrc=*/true ));
-            tgtMeshDS->SetNodeOnEdge( n, TopoDS::Edge( tgtEdge ));
-            double U = srcHelper.GetNodeU( TopoDS::Edge( srcEdge ), srcNode );
+            double U = Precision::Infinite();
             helper.CheckNodeU( tgtEdge, n, U, Precision::PConfusion());
-            n->SetPosition(SMDS_PositionPtr(new SMDS_EdgePosition( U )));
+            tgtMeshDS->SetNodeOnEdge( n, TopoDS::Edge( tgtEdge ), U );
             break;
           }
           case SMDS_TOP_VERTEX: {
@@ -900,6 +921,7 @@ namespace {
             tgtMeshDS->SetNodeOnVertex( n, TopoDS::Vertex( tgtV ));
             break;
           }
+          default:;
           }
           srcN_tgtN->second = n;
         }
@@ -911,10 +933,165 @@ namespace {
       case 3: helper.AddFace(tgtNodes[0], tgtNodes[2], tgtNodes[1]); break;
       case 4: helper.AddFace(tgtNodes[0], tgtNodes[3], tgtNodes[2], tgtNodes[1]); break;
       }
-    }
+    }  // loop on all mesh faces on srcFace
+
     return true;
+  }
+
+  //================================================================================
+  /*!
+   * \brief Preform projection in case of quadrilateral faces
+   */
+  //================================================================================
 
-  } // bool projectBy2DSimilarity(...)
+  bool projectQuads(const TopoDS_Face&                 tgtFace,
+                    const TopoDS_Face&                 srcFace,
+                    const TSideVector&                 tgtWires,
+                    const TSideVector&                 srcWires,
+                    const TAssocTool::TShapeShapeMap&  shape2ShapeMap,
+                    TAssocTool::TNodeNodeMap&          src2tgtNodes,
+                    const bool                         is1DComputed)
+  {
+    // SMESH_Mesh * tgtMesh = tgtWires[0]->GetMesh();
+    // SMESH_Mesh * srcMesh = srcWires[0]->GetMesh();
+    // //SMESHDS_Mesh * tgtMeshDS = tgtMesh->GetMeshDS();
+    // SMESHDS_Mesh * srcMeshDS = srcMesh->GetMeshDS();
+
+    // if ( srcWires[0]->NbEdges() != 4 )
+    //   return false;
+    // if ( !is1DComputed )
+    //   return false;
+    // for ( int iE = 0; iE < 4; ++iE )
+    // {
+    //   SMESHDS_SubMesh* sm = srcMeshDS->MeshElements( srcWires[0]->Edge( iE ));
+    //   if ( !sm ) return false;
+    //   if ( sm->NbNodes() + sm->NbElements() == 0 ) return false;
+    // }
+    // if ( BRepAdaptor_Surface( tgtFace ).GetType() != GeomAbs_Plane )
+    //   return false;
+    // // if ( BRepAdaptor_Surface( tgtFace ).GetType() == GeomAbs_Plane &&
+    // //      BRepAdaptor_Surface( srcFace ).GetType() == GeomAbs_Plane )
+    // //   return false; // too easy
+
+    // // load EDGEs to SMESH_Block
+
+    // SMESH_Block block;
+    // TopTools_IndexedMapOfOrientedShape blockSubShapes;
+    // {
+    //   const TopoDS_Solid& box = srcMesh->PseudoShape();
+    //   TopoDS_Shell shell = TopoDS::Shell( TopExp_Explorer( box, TopAbs_SHELL ).Current() );
+    //   TopoDS_Vertex v;
+    //   block.LoadBlockShapes( shell, v, v, blockSubShapes ); // fill all since operator[] is missing
+    // }
+    // const SMESH_Block::TShapeID srcFaceBID = SMESH_Block::ID_Fxy0;
+    // const SMESH_Block::TShapeID tgtFaceBID = SMESH_Block::ID_Fxy1;
+    // vector< int > edgeBID;
+    // block.GetFaceEdgesIDs( srcFaceBID, edgeBID ); // u0, u1, 0v, 1v
+    // blockSubShapes.Substitute( edgeBID[0], srcWires[0]->Edge(0) );
+    // blockSubShapes.Substitute( edgeBID[1], srcWires[0]->Edge(2) );
+    // blockSubShapes.Substitute( edgeBID[2], srcWires[0]->Edge(3) );
+    // blockSubShapes.Substitute( edgeBID[3], srcWires[0]->Edge(1) );
+    // block.GetFaceEdgesIDs( tgtFaceBID, edgeBID ); // u0, u1, 0v, 1v
+    // blockSubShapes.Substitute( edgeBID[0], tgtWires[0]->Edge(0) );
+    // blockSubShapes.Substitute( edgeBID[1], tgtWires[0]->Edge(2) );
+    // blockSubShapes.Substitute( edgeBID[2], tgtWires[0]->Edge(3) );
+    // blockSubShapes.Substitute( edgeBID[3], tgtWires[0]->Edge(1) );
+    // block.LoadFace( srcFace, srcFaceBID, blockSubShapes );
+    // block.LoadFace( tgtFace, tgtFaceBID, blockSubShapes );
+
+    // // remember connectivity of new faces in terms of ( node-or-XY )
+
+    // typedef std::pair< const SMDS_MeshNode*, gp_XYZ > TNodeOrXY; // node-or-XY
+    // typedef std::vector< TNodeOrXY* >                 TFaceConn; // face connectivity
+    // std::vector< TFaceConn >                    newFacesVec;     // connectivity of all faces
+    // std::map< const SMDS_MeshNode*, TNodeOrXY > srcNode2tgtNXY;  // src node -> node-or-XY
+
+    // TAssocTool::TNodeNodeMap::iterator                                       srcN_tgtN;
+    // std::map< const SMDS_MeshNode*, TNodeOrXY >::iterator                    srcN_tgtNXY;
+    // std::pair< std::map< const SMDS_MeshNode*, TNodeOrXY >::iterator, bool > n2n_isNew;
+    // TNodeOrXY nullNXY( (SMDS_MeshNode*)NULL, gp_XYZ(0,0,0) );
+
+    // SMESHDS_SubMesh* srcSubDS = srcMeshDS->MeshElements( srcFace );
+    // newFacesVec.resize( srcSubDS->NbElements() );
+    // int iFaceSrc = 0;
+
+    // SMDS_ElemIteratorPtr elemIt = srcSubDS->GetElements();
+    // while ( elemIt->more() ) // loop on all mesh faces on srcFace
+    // {
+    //   const SMDS_MeshElement* elem = elemIt->next();
+    //   TFaceConn& tgtNodes = newFacesVec[ iFaceSrc++ ];
+
+    //   const int nbN = elem->NbCornerNodes(); 
+    //   tgtNodes.resize( nbN );
+    //   for ( int i = 0; i < nbN; ++i ) // loop on nodes of the source element
+    //   {
+    //     const SMDS_MeshNode* srcNode = elem->GetNode(i);
+    //     n2n_isNew = srcNode2tgtNXY.insert( make_pair( srcNode, nullNXY ));
+    //     TNodeOrXY & tgtNodeOrXY = n2n_isNew.first->second;
+    //     if ( n2n_isNew.second ) // new src node encounters
+    //     {
+    //       srcN_tgtN = src2tgtNodes.find( srcNode );
+    //       if ( srcN_tgtN != src2tgtNodes.end() )
+    //       {
+    //         tgtNodeOrXY.first = srcN_tgtN->second; // tgt node exists
+    //       }
+    //       else 
+    //       {
+    //         // find XY of src node withing the quadrilateral srcFace
+    //         if ( !block.ComputeParameters( SMESH_TNodeXYZ( srcNode ),
+    //                                        tgtNodeOrXY.second, srcFaceBID ))
+    //           return false;
+    //       }
+    //     }
+    //     tgtNodes[ i ] = & tgtNodeOrXY;
+    //   }
+    // }
+
+    // // as all XY are computed, create tgt nodes and faces
+
+    // SMESH_MesherHelper helper( *tgtMesh );
+    // helper.SetSubShape( tgtFace );
+    // if ( is1DComputed )
+    //   helper.IsQuadraticSubMesh( tgtFace );
+    // else
+    //   helper.SetIsQuadratic( srcSubDS->GetElements()->next()->IsQuadratic() );
+    // helper.SetElementsOnShape( true );
+    // Handle(Geom_Surface) tgtSurface = BRep_Tool::Surface( tgtFace );
+
+    // SMESH_MesherHelper srcHelper( *srcMesh );
+    // srcHelper.SetSubShape( srcFace );
+
+    // vector< const SMDS_MeshNode* > tgtNodes;
+    // gp_XY uv;
+
+    // for ( size_t iFaceTgt = 0; iFaceTgt < newFacesVec.size(); ++iFaceTgt )
+    // {
+    //   TFaceConn& tgtConn = newFacesVec[ iFaceTgt ];
+    //   tgtNodes.resize( tgtConn.size() );
+    //   for ( size_t iN = 0; iN < tgtConn.size(); ++iN )
+    //   {
+    //     const SMDS_MeshNode* & tgtN = tgtConn[ iN ]->first;
+    //     if ( !tgtN ) // create a node
+    //     {
+    //       if ( !block.FaceUV( tgtFaceBID, tgtConn[iN]->second, uv ))
+    //         return false;
+    //       gp_Pnt p = tgtSurface->Value( uv.X(), uv.Y() );
+    //       tgtN = helper.AddNode( p.X(), p.Y(), p.Z(), uv.X(), uv.Y() );
+    //     }
+    //     tgtNodes[ tgtNodes.size() - iN - 1] = tgtN; // reversed orientation
+    //   }
+    //   switch ( tgtNodes.size() )
+    //   {
+    //   case 3: helper.AddFace(tgtNodes[0], tgtNodes[1], tgtNodes[2]); break;
+    //   case 4: helper.AddFace(tgtNodes[0], tgtNodes[1], tgtNodes[2], tgtNodes[3]); break;
+    //   default:
+    //     if ( tgtNodes.size() > 4 )
+    //       helper.AddPolygonalFace( tgtNodes );
+    //   }
+    // }
+    return false; //true;
+
+  } // bool projectQuads(...)
 
   //================================================================================
   /*!
@@ -922,12 +1099,12 @@ namespace {
    */
   //================================================================================
 
-  void fixDistortedFaces( SMESH_MesherHelper& helper,
+  bool fixDistortedFaces( SMESH_MesherHelper& helper,
                           TSideVector&        tgtWires )
   {
     SMESH_subMesh* faceSM = helper.GetMesh()->GetSubMesh( helper.GetSubShape() );
 
-    if ( helper.IsDistorted2D( faceSM ))
+    if ( helper.IsDistorted2D( faceSM, /*checkUV=*/true ))
     {
       SMESH_MeshEditor editor( helper.GetMesh() );
       SMESHDS_SubMesh* smDS = faceSM->GetSubMeshDS();
@@ -965,7 +1142,290 @@ namespace {
       set<const SMDS_MeshNode*> fixedNodes;
       editor.Smooth( faces, fixedNodes, algo, /*nbIterations=*/ 10,
                      /*theTgtAspectRatio=*/1.0, /*the2D=*/!isPlanar);
+
+      helper.ToFixNodeParameters( true );
+
+      return !helper.IsDistorted2D( faceSM, /*checkUV=*/true );
+    }
+    return true;
+  }
+
+  typedef list< pair< const SMDS_MeshNode*, const BRepMesh_Triangle* > > TNodeTriaList;
+
+  //================================================================================
+  /*!
+   * \brief Add in-FACE nodes surrounding a given node to a queue
+   */
+  //================================================================================
+
+  void addCloseNodes( const SMDS_MeshNode*     srcNode,
+                      const BRepMesh_Triangle* bmTria,
+                      const int                srcFaceID,
+                      TNodeTriaList &          noTriQueue )
+  {
+    // find in-FACE nodes
+    SMDS_ElemIteratorPtr elems = srcNode->GetInverseElementIterator(SMDSAbs_Face);
+    while ( elems->more() )
+    {
+      const SMDS_MeshElement* elem = elems->next();
+      if ( elem->getshapeId() == srcFaceID )
+      {
+        for ( int i = 0, nb = elem->NbNodes(); i < nb; ++i )
+        {
+          const SMDS_MeshNode* n = elem->GetNode( i );
+          if ( !n->isMarked() )
+            noTriQueue.push_back( make_pair( n, bmTria ));
+        }
+      }
+    }
+  }
+
+  //================================================================================
+  /*!
+   * \brief Find a delauney triangle containing a given 2D point and return
+   *        barycentric coordinates within the found triangle
+   */
+  //================================================================================
+
+  const BRepMesh_Triangle* findTriangle( const gp_XY&                            uv,
+                                         const BRepMesh_Triangle*                bmTria,
+                                         Handle(BRepMesh_DataStructureOfDelaun)& triaDS,
+                                         double                                  bc[3] )
+  {
+    int   nodeIDs[3];
+    gp_XY nodeUVs[3];
+    int   linkIDs[3];
+    Standard_Boolean ori[3];
+
+    while ( bmTria )
+    {
+      // check bmTria
+
+      triaDS->ElementNodes( *bmTria, nodeIDs );
+      nodeUVs[0] = triaDS->GetNode( nodeIDs[0] ).Coord();
+      nodeUVs[1] = triaDS->GetNode( nodeIDs[1] ).Coord();
+      nodeUVs[2] = triaDS->GetNode( nodeIDs[2] ).Coord();
+
+      SMESH_MeshAlgos::GetBarycentricCoords( uv,
+                                             nodeUVs[0], nodeUVs[1], nodeUVs[2],
+                                             bc[0], bc[1] );
+      if ( bc[0] >= 0 && bc[1] >= 0 && bc[0] + bc[1] <= 1 )
+      {
+        bc[2] = 1 - bc[0] - bc[1];
+        return bmTria;
+      }
+
+      // look for a neighbor triangle, which is adjacent to a link intersected
+      // by a segment( triangle center -> uv )
+
+      gp_XY gc = ( nodeUVs[0] + nodeUVs[1] + nodeUVs[2] ) / 3.;
+      gp_XY seg = uv - gc;
+
+      bmTria->Edges( linkIDs, ori );
+      int triaID = triaDS->IndexOf( *bmTria );
+      bmTria = 0;
+
+      for ( int i = 0; i < 3; ++i )
+      {
+        const BRepMesh_PairOfIndex & triIDs = triaDS->ElementsConnectedTo( linkIDs[i] );
+        if ( triIDs.Extent() < 2 )
+          continue; // no neighbor triangle
+
+        // check if a link intersects gc2uv
+        const BRepMesh_Edge & link = triaDS->GetLink( linkIDs[i] );
+        const BRepMesh_Vertex & n1 = triaDS->GetNode( link.FirstNode() );
+        const BRepMesh_Vertex & n2 = triaDS->GetNode( link.LastNode() );
+        gp_XY uv1 = n1.Coord();
+        gp_XY lin = n2.Coord() - uv1; // link direction
+
+        double crossSegLin = seg ^ lin;
+        if ( Abs( crossSegLin ) < std::numeric_limits<double>::min() )
+          continue; // parallel
+
+        double uSeg = ( uv1 - gc ) ^ lin / crossSegLin;
+        if ( 0. <= uSeg && uSeg <= 1. )
+        {
+          bmTria = & triaDS->GetElement( triIDs.Index( 1 + ( triIDs.Index(1) == triaID )));
+          break;
+        }
+      }
     }
+    return bmTria;
+  }
+
+  //================================================================================
+  /*!
+   * \brief Morph mesh on the target face to lie within FACE boundary w/o distortion
+   *
+   * algo:
+   * - make a CDT on the src FACE
+   * - find a triangle containing a src node and get its barycentric coordinates
+   * - move the node to a point with the same barycentric coordinates in a corresponding
+   *   tgt triangle
+   */
+  //================================================================================
+
+  bool morph( SMESH_MesherHelper&             tgtHelper,
+              const TopoDS_Face&              tgtFace,
+              const TopoDS_Face&              srcFace,
+              const TSideVector&              tgtWires,
+              const TSideVector&              srcWires,
+              const TAssocTool::TNodeNodeMap& src2tgtNodes )
+  {
+    if ( srcWires.size() != tgtWires.size() ) return false;
+
+    // count boundary points
+    int iP = 1, nbP = 0;
+    for ( size_t iW = 0; iW < srcWires.size(); ++iW )
+      nbP += srcWires[iW]->NbPoints() - 1; // 1st and last points coincide
+
+    // fill boundary points
+    BRepMesh::Array1OfVertexOfDelaun srcVert( 1, 1 + nbP ), tgtVert( 1, 1 + nbP );
+    vector< const SMDS_MeshNode* > bndSrcNodes( nbP + 1 ); bndSrcNodes[0] = 0;
+    BRepMesh_Vertex v( 0, 0, BRepMesh_Frontier );
+    for ( size_t iW = 0; iW < srcWires.size(); ++iW )
+    {
+      const UVPtStructVec& srcPnt = srcWires[iW]->GetUVPtStruct();
+      const UVPtStructVec& tgtPnt = tgtWires[iW]->GetUVPtStruct();
+      if ( srcPnt.size() != tgtPnt.size() ) return false;
+
+      for ( int i = 0, nb = srcPnt.size() - 1;  i < nb;  ++i, ++iP )
+      {
+        bndSrcNodes[ iP ]  = srcPnt[i].node;
+        srcPnt[i].node->setIsMarked( true );
+
+        v.ChangeCoord() = srcPnt[i].UV();
+        srcVert( iP )   = v;
+        v.ChangeCoord() = tgtPnt[i].UV();
+        tgtVert( iP )   = v;
+      }
+    }
+    // triangulate the srcFace in 2D
+    BRepMesh_Delaun delauney( srcVert );
+    Handle(BRepMesh_DataStructureOfDelaun) triaDS = delauney.Result();
+
+    Handle(ShapeAnalysis_Surface) tgtSurface = tgtHelper.GetSurface( tgtFace );
+    SMESHDS_Mesh* srcMesh = srcWires[0]->GetMesh()->GetMeshDS();
+    SMESHDS_Mesh* tgtMesh = tgtHelper.GetMeshDS();
+    const SMDS_MeshNode *srcNode, *tgtNode;
+    const BRepMesh_Triangle *bmTria;
+
+    // un-mark internal src nodes; later we will mark moved nodes
+    SMDS_NodeIteratorPtr nIt = srcMesh->MeshElements( srcFace )->GetNodes();
+    if ( !nIt || !nIt->more() ) return true;
+    while ( nIt->more() )
+      ( srcNode = nIt->next() )->setIsMarked( false );
+
+    // initialize a queue of nodes with starting triangles
+    const int srcFaceID = srcNode->getshapeId();
+    TNodeTriaList noTriQueue;
+    size_t iBndSrcN = 1;
+    for ( ; iBndSrcN < bndSrcNodes.size() &&  noTriQueue.empty();  ++iBndSrcN )
+    {
+      // get a triangle
+      const BRepMesh::ListOfInteger & linkIds = triaDS->LinksConnectedTo( iBndSrcN );
+      const BRepMesh_PairOfIndex &    triaIds = triaDS->ElementsConnectedTo( linkIds.First() );
+      const BRepMesh_Triangle&           tria = triaDS->GetElement( triaIds.Index(1) );
+
+      addCloseNodes( bndSrcNodes[ iBndSrcN ], &tria, srcFaceID, noTriQueue );
+    }
+
+    // Move tgt nodes
+
+    double bc[3]; // barycentric coordinates
+    int    nodeIDs[3];
+    bool   checkUV = true;
+    const SMDS_FacePosition* pos;
+
+    while ( !noTriQueue.empty() )
+    {
+      srcNode = noTriQueue.front().first;
+      bmTria  = noTriQueue.front().second;
+      noTriQueue.pop_front();
+      if ( srcNode->isMarked() )
+        continue;
+      srcNode->setIsMarked( true );
+
+      // find a delauney triangle containing the src node
+      gp_XY uv = tgtHelper.GetNodeUV( srcFace, srcNode, NULL, &checkUV );
+      bmTria = findTriangle( uv, bmTria, triaDS, bc );
+      if ( !bmTria )
+        continue;
+
+      // compute new coordinates for a corresponding tgt node
+      gp_XY uvNew( 0., 0. ), nodeUV;
+      triaDS->ElementNodes( *bmTria, nodeIDs );
+      for ( int i = 0; i < 3; ++i )
+        uvNew += bc[i] * tgtVert( nodeIDs[i]).Coord();
+      gp_Pnt xyz = tgtSurface->Value( uvNew );
+
+      // find and move tgt node
+      TAssocTool::TNodeNodeMap::const_iterator n2n = src2tgtNodes.find( srcNode );
+      if ( n2n == src2tgtNodes.end() ) continue;
+      tgtNode = n2n->second;
+      tgtMesh->MoveNode( tgtNode, xyz.X(), xyz.Y(), xyz.Z() );
+
+      if (( pos = dynamic_cast< const SMDS_FacePosition* >( tgtNode->GetPosition() )))
+        const_cast<SMDS_FacePosition*>( pos )->SetParameters( uvNew.X(), uvNew.Y() );
+
+      addCloseNodes( srcNode, bmTria, srcFaceID, noTriQueue );
+
+      // assure that all src nodes are visited
+      for ( ; iBndSrcN < bndSrcNodes.size() &&  noTriQueue.empty();  ++iBndSrcN )
+      {
+        const BRepMesh::ListOfInteger & linkIds = triaDS->LinksConnectedTo( iBndSrcN );
+        const BRepMesh_PairOfIndex &    triaIds = triaDS->ElementsConnectedTo( linkIds.First() );
+        const BRepMesh_Triangle&           tria = triaDS->GetElement( triaIds.Index(1) );
+        addCloseNodes( bndSrcNodes[ iBndSrcN ], &tria, srcFaceID, noTriQueue );
+      }
+    }
+
+    return true;
+  }
+
+  //=======================================================================
+  /*
+   * Set initial association of VERTEXes for the case of projection
+   * from a quadrangle FACE to a closed FACE, where opposite src EDGEs
+   * have different nb of segments
+   */
+  //=======================================================================
+
+  void initAssoc4Quad2Closed(const TopoDS_Shape&          tgtFace,
+                             SMESH_MesherHelper&          tgtHelper,
+                             const TopoDS_Shape&          srcFace,
+                             SMESH_Mesh*                  srcMesh,
+                             TAssocTool::TShapeShapeMap & assocMap)
+  {
+    if ( !tgtHelper.HasRealSeam() || srcFace.ShapeType() != TopAbs_FACE )
+      return; // no seam edge
+    list< TopoDS_Edge > tgtEdges, srcEdges;
+    list< int > tgtNbEW, srcNbEW;
+    int tgtNbW = SMESH_Block::GetOrderedEdges( TopoDS::Face( tgtFace ), tgtEdges, tgtNbEW );
+    int srcNbW = SMESH_Block::GetOrderedEdges( TopoDS::Face( srcFace ), srcEdges, srcNbEW );
+    if ( tgtNbW != 1 || srcNbW != 1 ||
+         tgtNbEW.front() != 4 || srcNbEW.front() != 4 )
+      return; // not quads
+
+    int srcNbSeg[4];
+    list< TopoDS_Edge >::iterator edgeS = srcEdges.begin(), edgeT = tgtEdges.begin();
+    for ( int i = 0; edgeS != srcEdges.end(); ++i, ++edgeS )
+      if ( SMESHDS_SubMesh* sm = srcMesh->GetMeshDS()->MeshElements( *edgeS ))
+        srcNbSeg[ i ] = sm->NbNodes();
+      else
+        return; // not meshed
+    if ( srcNbSeg[0] == srcNbSeg[2] && srcNbSeg[1] == srcNbSeg[3] )
+      return; // same nb segments
+    if ( srcNbSeg[0] != srcNbSeg[2] && srcNbSeg[1] != srcNbSeg[3] )
+      return; // all different nb segments
+
+    edgeS = srcEdges.begin();
+    if ( srcNbSeg[0] != srcNbSeg[2] )
+      ++edgeS;
+    TAssocTool::InsertAssociation( tgtHelper.IthVertex( 0,*edgeT ),
+                                   tgtHelper.IthVertex( 0,*edgeS ), assocMap );
+    TAssocTool::InsertAssociation( tgtHelper.IthVertex( 1,*edgeT ),
+                                   tgtHelper.IthVertex( 1,*edgeS ), assocMap );
   }
 
 } // namespace
@@ -980,7 +1440,6 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
 {
   _src2tgtNodes.clear();
 
-  MESSAGE("Projection_2D Compute");
   if ( !_sourceHypo )
     return false;
 
@@ -990,6 +1449,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     srcMesh = tgtMesh;
 
   SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
+  SMESH_MesherHelper helper( theMesh );
 
   // ---------------------------
   // Make sub-shapes association
@@ -998,16 +1458,20 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
   TopoDS_Face   tgtFace = TopoDS::Face( theShape.Oriented(TopAbs_FORWARD));
   TopoDS_Shape srcShape = _sourceHypo->GetSourceFace().Oriented(TopAbs_FORWARD);
 
+  helper.SetSubShape( tgtFace );
+
   TAssocTool::TShapeShapeMap shape2ShapeMap;
   TAssocTool::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
+  if ( shape2ShapeMap.IsEmpty() )
+    initAssoc4Quad2Closed( tgtFace, helper, srcShape, srcMesh, shape2ShapeMap );
   if ( !TAssocTool::FindSubShapeAssociation( tgtFace, tgtMesh, srcShape, srcMesh,
                                              shape2ShapeMap)  ||
        !shape2ShapeMap.IsBound( tgtFace ))
   {
     if ( srcShape.ShapeType() == TopAbs_FACE )
     {
-      int nbE1 = SMESH_MesherHelper::Count( tgtFace, TopAbs_EDGE, /*ignoreSame=*/true );
-      int nbE2 = SMESH_MesherHelper::Count( srcShape, TopAbs_EDGE, /*ignoreSame=*/true );
+      int nbE1 = helper.Count( tgtFace, TopAbs_EDGE, /*ignoreSame=*/true );
+      int nbE2 = helper.Count( srcShape, TopAbs_EDGE, /*ignoreSame=*/true );
       if ( nbE1 != nbE2 )
         return error(COMPERR_BAD_SHAPE,
                      SMESH_Comment("Different number of edges in source and target faces: ")
@@ -1048,25 +1512,32 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
   if ( err && !err->IsOK() )
     return error( err );
 
-  bool done = false;
+  bool projDone = false;
 
-  if ( !done )
+  if ( !projDone )
   {
     // try to project from the same face with different location
-    done = projectPartner( tgtFace, srcFace, tgtWires, srcWires,
-                           shape2ShapeMap, _src2tgtNodes, is1DComputed );
+    projDone = projectPartner( tgtFace, srcFace, tgtWires, srcWires,
+                               shape2ShapeMap, _src2tgtNodes, is1DComputed );
   }
-  if ( !done )
+  if ( !projDone )
   {
     // projection in case if the faces are similar in 2D space
-    done = projectBy2DSimilarity( tgtFace, srcFace, tgtWires, srcWires,
-                                  shape2ShapeMap, _src2tgtNodes, is1DComputed);
+    projDone = projectBy2DSimilarity( tgtFace, srcFace, tgtWires, srcWires,
+                                      shape2ShapeMap, _src2tgtNodes, is1DComputed );
+  }
+  if ( !projDone )
+  {
+    // projection in case of quadrilateral faces
+    // NOT IMPLEMENTED, returns false
+    projDone = projectQuads( tgtFace, srcFace, tgtWires, srcWires,
+                             shape2ShapeMap, _src2tgtNodes, is1DComputed);
   }
 
-  SMESH_MesherHelper helper( theMesh );
-  helper.SetSubShape( tgtFace );
+  // it will remove mesh built on edges and vertices in failure case
+  MeshCleaner cleaner( tgtSubMesh );
 
-  if ( !done )
+  if ( !projDone )
   {
     _src2tgtNodes.clear();
     // --------------------
@@ -1076,20 +1547,22 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     // Check if node projection to a face is needed
     Bnd_B2d uvBox;
     SMDS_ElemIteratorPtr faceIt = srcSubMesh->GetSubMeshDS()->GetElements();
-    int nbFaceNodes = 0;
-    for ( ; nbFaceNodes < 3 && faceIt->more();  ) {
+    set< const SMDS_MeshNode* > faceNodes;
+    for ( ; faceNodes.size() < 3 && faceIt->more();  ) {
       const SMDS_MeshElement* face = faceIt->next();
       SMDS_ElemIteratorPtr nodeIt = face->nodesIterator();
       while ( nodeIt->more() ) {
         const SMDS_MeshNode* node = static_cast<const SMDS_MeshNode*>( nodeIt->next() );
-        if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
-          nbFaceNodes++;
+        if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE &&
+             faceNodes.insert( node ).second )
           uvBox.Add( helper.GetNodeUV( srcFace, node ));
-        }
       }
     }
-    const bool toProjectNodes =
-      ( nbFaceNodes > 0 && ( uvBox.IsVoid() || uvBox.SquareExtent() < DBL_MIN ));
+    bool toProjectNodes = false;
+    if ( faceNodes.size() == 1 )
+      toProjectNodes = ( uvBox.IsVoid() || uvBox.CornerMin().IsEqual( gp_XY(0,0), 1e-12 ));
+    else if ( faceNodes.size() > 1 )
+      toProjectNodes = ( uvBox.IsVoid() || uvBox.SquareExtent() < DBL_MIN );
 
     // Find the corresponding source and target vertex
     // and <theReverse> flag needed to call mapper.Apply()
@@ -1097,13 +1570,10 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     TopoDS_Vertex srcV1, tgtV1;
     bool reverse = false;
 
-    if ( _sourceHypo->HasVertexAssociation() ) {
-      srcV1 = _sourceHypo->GetSourceVertex(1);
-      tgtV1 = _sourceHypo->GetTargetVertex(1);
-    } else {
-      srcV1 = TopoDS::Vertex( TopExp_Explorer( srcFace, TopAbs_VERTEX ).Current() );
-      tgtV1 = TopoDS::Vertex( shape2ShapeMap( srcV1, /*isSrc=*/true ));
-    }
+    TopExp_Explorer vSrcExp( srcFace, TopAbs_VERTEX );
+    srcV1 = TopoDS::Vertex( vSrcExp.Current() );
+    tgtV1 = TopoDS::Vertex( shape2ShapeMap( srcV1, /*isSrc=*/true ));
+
     list< TopoDS_Edge > tgtEdges, srcEdges;
     list< int > nbEdgesInWires;
     SMESH_Block::GetOrderedEdges( tgtFace, tgtEdges, nbEdgesInWires, tgtV1 );
@@ -1115,7 +1585,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
       TopoDS_Shape srcE1bis = shape2ShapeMap( tgtE1 );
       reverse = ( ! srcE1.IsSame( srcE1bis ));
       if ( reverse &&
-           _sourceHypo->HasVertexAssociation() &&
+           //_sourceHypo->HasVertexAssociation() &&
            nbEdgesInWires.front() > 2 &&
            helper.IsRealSeam( tgtEdges.front() ))
       {
@@ -1124,19 +1594,83 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
         // we can't use only theReverse flag to correctly associate source
         // and target faces in the mapper. Thus we select srcV1 so that
         // GetOrderedEdges() to return EDGEs in a needed order
-        list< TopoDS_Edge >::iterator edge = srcEdges.begin();
-        for ( ; edge != srcEdges.end(); ++edge ) {
-          if ( srcE1bis.IsSame( *edge )) {
-            srcV1 = helper.IthVertex( 0, *edge );
-            break;
+        TopoDS_Face tgtFaceBis = tgtFace;
+        TopTools_MapOfShape checkedVMap( tgtEdges.size() );
+        checkedVMap.Add ( srcV1 );
+        for ( vSrcExp.Next(); vSrcExp.More(); )
+        {
+          tgtFaceBis.Reverse();
+          tgtEdges.clear();
+          SMESH_Block::GetOrderedEdges( tgtFaceBis, tgtEdges, nbEdgesInWires, tgtV1 );
+          bool ok = true;
+          list< TopoDS_Edge >::iterator edgeS = srcEdges.begin(), edgeT = tgtEdges.begin();
+          for ( ; edgeS != srcEdges.end() && ok ; ++edgeS, ++edgeT )
+            ok = edgeT->IsSame( shape2ShapeMap( *edgeS, /*isSrc=*/true ));
+          if ( ok )
+            break; // FOUND!
+
+          reverse = !reverse;
+          if ( reverse )
+          {
+            vSrcExp.Next();
+            while ( vSrcExp.More() && !checkedVMap.Add( vSrcExp.Current() ))
+              vSrcExp.Next();
+          }
+          else
+          {
+            srcV1 = TopoDS::Vertex( vSrcExp.Current() );
+            tgtV1 = TopoDS::Vertex( shape2ShapeMap( srcV1, /*isSrc=*/true ));
+            srcEdges.clear();
+            SMESH_Block::GetOrderedEdges( srcFace, srcEdges, nbEdgesInWires, srcV1 );
+          }
+        }
+      }
+      // for the case: project to a closed face from a non-closed face w/o vertex assoc;
+      // avoid projecting to a seam from two EDGEs with different nb nodes on them
+      // ( test mesh_Projection_2D_01/B1 )
+      if ( !_sourceHypo->HasVertexAssociation() &&
+           nbEdgesInWires.front() > 2 &&
+           helper.IsRealSeam( tgtEdges.front() ))
+      {
+        TopoDS_Shape srcEdge1 = shape2ShapeMap( tgtEdges.front() );
+        list< TopoDS_Edge >::iterator srcEdge2 =
+          std::find( srcEdges.begin(), srcEdges.end(), srcEdge1);
+        list< TopoDS_Edge >::iterator srcEdge3 =
+          std::find( srcEdges.begin(), srcEdges.end(), srcEdge1.Reversed());
+        if ( srcEdge2 == srcEdges.end() || srcEdge3 == srcEdges.end() ) // srcEdge1 is not a seam
+        {
+          // find srcEdge2 which also will be projected to tgtEdges.front()
+          for ( srcEdge2 = srcEdges.begin(); srcEdge2 != srcEdges.end(); ++srcEdge2 )
+            if ( !srcEdge1.IsSame( *srcEdge2 ) &&
+                 tgtEdges.front().IsSame( shape2ShapeMap( *srcEdge2, /*isSrc=*/true )))
+              break;
+          // compare nb nodes on srcEdge1 and srcEdge2
+          if ( srcEdge2 != srcEdges.end() )
+          {
+            int nbN1 = 0, nbN2 = 0;
+            if ( SMESHDS_SubMesh* sm = srcMesh->GetMeshDS()->MeshElements( srcEdge1 ))
+              nbN1 = sm->NbNodes();
+            if ( SMESHDS_SubMesh* sm = srcMesh->GetMeshDS()->MeshElements( *srcEdge2 ))
+              nbN2 = sm->NbNodes();
+            if ( nbN1 != nbN2 )
+              srcV1 = helper.IthVertex( 1, srcEdges.front() );
           }
         }
       }
     }
-    else if ( nbEdgesInWires.front() == 1 )
+    else if ( nbEdgesInWires.front() == 1 ) // a sole edge in a wire
     {
-      // TODO::Compare orientation of curves in a sole edge
-      //RETURN_BAD_RESULT("Not implemented case");
+      TopoDS_Edge srcE1 = srcEdges.front(), tgtE1 = tgtEdges.front();
+      for ( size_t iW = 0; iW < srcWires.size(); ++iW )
+      {
+        StdMeshers_FaceSidePtr srcWire = srcWires[iW];
+        for ( int iE = 0; iE < srcWire->NbEdges(); ++iE )
+          if ( srcE1.IsSame( srcWire->Edge( iE )))
+          {
+            reverse = ( tgtE1.Orientation() != tgtWires[iW]->Edge( iE ).Orientation() );
+            break;
+          }
+      }
     }
     else
     {
@@ -1145,7 +1679,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
 
     // Load pattern from the source face
     SMESH_Pattern mapper;
-    mapper.Load( srcMesh, srcFace, toProjectNodes, srcV1 );
+    mapper.Load( srcMesh, srcFace, toProjectNodes, srcV1, /*keepNodes=*/true );
     if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK )
       return error(COMPERR_BAD_INPUT_MESH,"Can't load mesh pattern from the source face");
 
@@ -1156,8 +1690,11 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     // Compute mesh on a target face
 
     mapper.Apply( tgtFace, tgtV1, reverse );
-    if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK )
+    if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK ) {
+      // std::ofstream file("/tmp/Pattern.smp" );
+      // mapper.Save( file );
       return error("Can't apply source mesh pattern to the face");
+    }
 
     // Create the mesh
 
@@ -1166,10 +1703,18 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK )
       return error("Can't make mesh by source mesh pattern");
 
-    // it will remove mesh built by pattern mapper on edges and vertices
-    // in failure case
-    MeshCleaner cleaner( tgtSubMesh );
+    // fill _src2tgtNodes
+    std::vector< const SMDS_MeshNode* > *srcNodes, *tgtNodes;
+    mapper.GetInOutNodes( srcNodes, tgtNodes );
+    size_t nbN = std::min( srcNodes->size(), tgtNodes->size() );
+    for ( size_t i = 0; i < nbN; ++i )
+      if ( (*srcNodes)[i] && (*tgtNodes)[i] )
+        _src2tgtNodes.insert( make_pair( (*srcNodes)[i], (*tgtNodes)[i] ));
 
+
+  } // end of projection using Pattern mapping
+
+  {
     // -------------------------------------------------------------------------
     // mapper doesn't take care of nodes already existing on edges and vertices,
     // so we must merge nodes created by it with existing ones
@@ -1209,7 +1754,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
         continue; // do not treat sm of degen VERTEX
       }
 
-      // Sort new and old nodes of a submesh separately
+      // Sort new and old nodes of a sub-mesh separately
 
       bool isSeam = helper.IsRealSeam( sm->GetId() );
 
@@ -1279,8 +1824,8 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
              u2nodesOnSeam.size()            > 0 &&
              seam.ShapeType() == TopAbs_EDGE )
         {
-          int nbE1 = SMESH_MesherHelper::Count( tgtFace, TopAbs_EDGE, /*ignoreSame=*/true );
-          int nbE2 = SMESH_MesherHelper::Count( srcFace, TopAbs_EDGE, /*ignoreSame=*/true );
+          int nbE1 = helper.Count( tgtFace, TopAbs_EDGE, /*ignoreSame=*/true );
+          int nbE2 = helper.Count( srcFace, TopAbs_EDGE, /*ignoreSame=*/true );
           if ( nbE1 != nbE2 ) // 2 EDGEs are mapped to a seam EDGE
           {
             // find the 2 EDGEs of srcFace
@@ -1329,18 +1874,11 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     if ( nbFaceBeforeMerge != nbFaceAtferMerge && !helper.HasDegeneratedEdges() )
       return error(COMPERR_BAD_INPUT_MESH, "Probably invalid node parameters on geom faces");
 
-
-    // ----------------------------------------------------------------
-    // The mapper can create distorted faces by placing nodes out of the FACE
-    // boundary -- fix bad faces by smoothing
-    // ----------------------------------------------------------------
-
-    fixDistortedFaces( helper, tgtWires );
-
     // ----------------------------------------------------------------
     // The mapper can't create quadratic elements, so convert if needed
     // ----------------------------------------------------------------
 
+    SMDS_ElemIteratorPtr faceIt;
     faceIt         = srcSubMesh->GetSubMeshDS()->GetElements();
     bool srcIsQuad = faceIt->next()->IsQuadratic();
     faceIt         = tgtSubMesh->GetSubMeshDS()->GetElements();
@@ -1354,12 +1892,22 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
 
       editor.ConvertToQuadratic(/*theForce3d=*/false, tgtFaces, false);
     }
+  } // end of coincident nodes and quadratic elements treatment
 
-    cleaner.Release(); // not to remove mesh
-
-  } // end of projection using Pattern mapping
 
+  if ( !projDone || is1DComputed )
+    // ----------------------------------------------------------------
+    // The mapper can create distorted faces by placing nodes out of the FACE
+    // boundary, also bad face can be created if EDGEs already discretized
+    // --> fix bad faces by smoothing
+    // ----------------------------------------------------------------
+    if ( helper.IsDistorted2D( tgtSubMesh, /*checkUV=*/false ))
+    {
+      morph( helper, tgtFace, srcFace, tgtWires, srcWires, _src2tgtNodes );
 
+      if ( !fixDistortedFaces( helper, tgtWires ))
+        return error("Invalid mesh generated");
+    }
   // ---------------------------
   // Check elements orientation
   // ---------------------------
@@ -1405,6 +1953,8 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
     }
   }
 
+  cleaner.Release(); // not to remove mesh
+
   return true;
 }