Salome HOME
0020128: EDF SMESH 926 : Quadratic conversion of BLSURF mesh
[modules/smesh.git] / src / SMESH / SMESH_MesherHelper.cxx
index ef545254b3798103cd861fadb4819bee32b0c051..b695e1276a5be8fc8310f060e75eaec3f694ba97 100644 (file)
@@ -1,4 +1,4 @@
-//  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
+//  Copyright (C) 2007-2010  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
@@ -19,6 +19,7 @@
 //
 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
+
 // File:      SMESH_MesherHelper.cxx
 // Created:   15.02.06 15:22:41
 // Author:    Sergey KUUL
@@ -57,6 +58,8 @@
 
 #include <limits>
 
+using namespace std;
+
 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
 
 namespace {
@@ -353,7 +356,7 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face&   F,
       static_cast<const SMDS_FacePosition*>(n->GetPosition().get());
     uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
     if ( check )
-      uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
+      uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
   }
   else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
   {
@@ -372,12 +375,12 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face&   F,
     else
       uv.SetCoord(0.,0.);
     if ( check || !validU )
-      uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ),/*force=*/ !validU );
+      uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ),/*force=*/ !validU );
 
     // for a node on a seam edge select one of UVs on 2 pcurves
     if ( n2 && IsSeamShape( edgeID ) )
     {
-      uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
+      uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
     }
     else
     { // adjust uv to period
@@ -488,12 +491,12 @@ bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
       }
       Quantity_Parameter U,V;
       projector.LowerDistanceParameters(U,V);
+      uv.SetCoord( U,V );
       if ( nodePnt.Distance( surface->Value( U, V )) > tol )
       {
         MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
         return false;
       }
-      uv.SetCoord( U,V );
     }
     else if ( uv.Modulus() > numeric_limits<double>::min() )
     {
@@ -594,23 +597,43 @@ gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
 
 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge&   E,
                                     const SMDS_MeshNode* n,
+                                    const SMDS_MeshNode* inEdgeNode,
                                     bool*                check)
 {
   double param = 0;
-  const SMDS_PositionPtr Pos = n->GetPosition();
-  if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) {
-    const SMDS_EdgePosition* epos =
-      static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
+  const SMDS_PositionPtr pos = n->GetPosition();
+  if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
+  {
+    const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos.get() );
     param =  epos->GetUParameter();
   }
-  else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) {
-    SMESHDS_Mesh * meshDS = GetMeshDS();
-    int vertexID = n->GetPosition()->GetShapeId();
-    const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
-    param =  BRep_Tool::Parameter( V, E );
+  else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
+  {
+    if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
+    {
+      Standard_Real f,l;
+      BRep_Tool::Range( E, f,l );
+      double uInEdge = GetNodeU( E, inEdgeNode );
+      param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
+    }
+    else
+    {
+      SMESHDS_Mesh * meshDS = GetMeshDS();
+      int vertexID = pos->GetShapeId();
+      const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
+      param =  BRep_Tool::Parameter( V, E );
+    }
   }
   if ( check )
-    *check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
+  {
+    double tol = BRep_Tool::Tolerance( E );
+    double f,l;  BRep_Tool::Range( E, f,l );
+    bool force = ( param < f-tol || param > l+tol );
+    if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
+      force = ( GetMeshDS()->ShapeToIndex( E ) != pos->GetShapeId() );
+
+    *check = CheckNodeU( E, n, param, 2*tol, force );
+  }
   return param;
 }
 
@@ -624,7 +647,8 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
                                     const SMDS_MeshNode* n,
                                     double&              u,
                                     const double         tol,
-                                    const bool           force) const
+                                    const bool           force,
+                                    double*              distance) const
 {
   if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
   {
@@ -643,7 +667,9 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
     {
       gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
       if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
-      if ( nodePnt.Distance( curve->Value( u )) > tol )
+      double dist = nodePnt.Distance( curve->Value( u ));
+      if ( distance ) *distance = dist;
+      if ( dist > tol )
       {
         // u incorrect, project the node to the curve
         GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
@@ -653,17 +679,34 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
           return false;
         }
         Quantity_Parameter U = projector.LowerDistanceParameter();
-        if ( nodePnt.Distance( curve->Value( U )) > tol )
+        u = double( U );
+        dist = nodePnt.Distance( curve->Value( U ));
+        if ( distance ) *distance = dist;
+        if ( dist > tol )
         {
           MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
           return false;
         }
-        u = double( U );
+        //u = double( U );
       }
       else if ( fabs( u ) > numeric_limits<double>::min() )
       {
         ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
       }
+      if (( u < f-tol || u > l+tol ) && force )
+      {
+        // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
+        try
+        {
+          // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
+          double period = curve->Period();
+          u = ( u < f ) ? u + period : u - period;
+        }
+        catch (Standard_Failure& exc)
+        {
+          return false;
+        }
+      }
     }
   }
   return true;
@@ -730,14 +773,19 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
   }
   else if (edgeID>0 || shapeType == TopAbs_EDGE)
   {
+    if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
+         Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
+         Pos1->GetShapeId() != Pos2->GetShapeId() ) // issue 0021006
+    return getMediumNodeOnComposedWire(n1,n2,force3d);
+
     if( myShape.IsNull() )
       E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
     else {
       E = TopoDS::Edge(myShape);
       edgeID = myShapeID;
     }
-    u[0] = GetNodeU(E,n1, force3d ? 0 : &uvOK[0]);
-    u[1] = GetNodeU(E,n2, force3d ? 0 : &uvOK[1]);
+    u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
+    u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
   }
   if(!force3d)
   {
@@ -797,16 +845,16 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
   if ( !F.IsNull() )
   {
     gp_XY UV = ( uv[0] + uv[1] ) / 2.;
-    CheckNodeUV( F, n12, UV, BRep_Tool::Tolerance( F ), /*force=*/true);
+    CheckNodeUV( F, n12, UV, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
     meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
   }
   else if ( !E.IsNull() )
   {
     double U = ( u[0] + u[1] ) / 2.;
-    CheckNodeU( E, n12, U, BRep_Tool::Tolerance( E ), /*force=*/true);
+    CheckNodeU( E, n12, U, 2*BRep_Tool::Tolerance( E ), /*force=*/true);
     meshDS->SetNodeOnEdge(n12, edgeID, U);
   }
-  else
+  else if ( myShapeID > 0 )
   {
     meshDS->SetNodeInVolume(n12, myShapeID);
   }
@@ -814,6 +862,78 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
   return n12;
 }
 
+//================================================================================
+/*!
+ * \brief Makes a medium node if nodes reside different edges
+ */
+//================================================================================
+
+const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
+                                                                     const SMDS_MeshNode* n2,
+                                                                     bool                 force3d)
+{
+  gp_Pnt middle = 0.5 * XYZ(n1) + 0.5 * XYZ(n2);
+  SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
+
+  // To find position on edge and 3D position for n12,
+  // project <middle> to 2 edges and select projection most close to <middle>
+
+  double u = 0, distMiddleProj = Precision::Infinite();
+  int iOkEdge = 0;
+  TopoDS_Edge edges[2];
+  for ( int is2nd = 0; is2nd < 2; ++is2nd )
+  {
+    // get an edge
+    const SMDS_MeshNode* n = is2nd ? n2 : n1;
+    TopoDS_Shape shape = GetSubShapeByNode( n, GetMeshDS() );
+    if ( shape.IsNull() || shape.ShapeType() != TopAbs_EDGE )
+      continue;
+
+    // project to get U of projection and distance from middle to projection
+    TopoDS_Edge edge = edges[ is2nd ] = TopoDS::Edge( shape );
+    double node2MiddleDist = middle.Distance( XYZ(n) );
+    double foundU = GetNodeU( edge, n ), foundDist = node2MiddleDist;
+    CheckNodeU( edge, n12, foundU, 2*BRep_Tool::Tolerance(edge), /*force=*/true, &foundDist );
+    if ( foundDist < node2MiddleDist )
+    {
+      distMiddleProj = foundDist;
+      u = foundU;
+      iOkEdge = is2nd;
+    }
+  }
+  if ( Precision::IsInfinite( distMiddleProj ))
+  {
+    // both projections failed; set n12 on the edge of n1 with U of a common vertex
+    TopoDS_Vertex vCommon;
+    if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
+      u = BRep_Tool::Parameter( vCommon, edges[0] );
+    else
+    {
+      double f,l, u0 = GetNodeU( edges[0], n1 );
+      BRep_Tool::Range( edges[0],f,l );
+      u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
+    }
+    iOkEdge = 0;
+    distMiddleProj = 0;
+  }
+
+  // move n12 to position of a successfull projection
+  double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
+  if ( !force3d && distMiddleProj > 2*tol )
+  {
+    TopLoc_Location loc; double f,l;
+    Handle(Geom_Curve) curve = BRep_Tool::Curve( edges[iOkEdge],loc,f,l );
+    gp_Pnt p = curve->Value( u );
+    GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
+  }
+
+  GetMeshDS()->SetNodeOnEdge(n12, edges[iOkEdge], u);
+
+  myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
+
+  return n12;
+}
+
 //=======================================================================
 //function : AddNode
 //purpose  : Creates a node
@@ -969,6 +1089,45 @@ SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
   return elem;
 }
 
+//=======================================================================
+//function : AddPolygonalFace
+//purpose  : Creates polygon, with additional nodes in quadratic mesh
+//=======================================================================
+
+SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
+                                                     const int                           id,
+                                                     const bool                          force3d)
+{
+  SMESHDS_Mesh * meshDS = GetMeshDS();
+  SMDS_MeshFace* elem = 0;
+
+  if(!myCreateQuadratic) {
+    if(id)
+      elem = meshDS->AddPolygonalFaceWithID(nodes, id);
+    else
+      elem = meshDS->AddPolygonalFace(nodes);
+  }
+  else {
+    vector<const SMDS_MeshNode*> newNodes;
+    for ( int i = 0; i < nodes.size(); ++i )
+    {
+      const SMDS_MeshNode* n1 = nodes[i];
+      const SMDS_MeshNode* n2 = nodes[(i+1)/nodes.size()];
+      const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
+      newNodes.push_back( n1 );
+      newNodes.push_back( n12 );
+    }
+    if(id)
+      elem = meshDS->AddPolygonalFaceWithID(newNodes, id);
+    else
+      elem = meshDS->AddPolygonalFace(newNodes);
+  }
+  if ( mySetElemOnShape && myShapeID > 0 )
+    meshDS->SetMeshElementOnShape( elem, myShapeID );
+
+  return elem;
+}
+
 //=======================================================================
 //function : AddVolume
 //purpose  : Creates quadratic or linear prism
@@ -1159,6 +1318,62 @@ SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
   return elem;
 }
 
+//=======================================================================
+//function : AddPolyhedralVolume
+//purpose  : Creates polyhedron. In quadratic mesh, adds medium nodes
+//=======================================================================
+
+SMDS_MeshVolume*
+SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
+                                         const std::vector<int>&                  quantities,
+                                         const int                                id,
+                                         const bool                               force3d)
+{
+  SMESHDS_Mesh * meshDS = GetMeshDS();
+  SMDS_MeshVolume* elem = 0;
+  if(!myCreateQuadratic)
+  {
+    if(id)
+      elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
+    else
+      elem = meshDS->AddPolyhedralVolume(nodes, quantities);
+  }
+  else
+  {
+    vector<const SMDS_MeshNode*> newNodes;
+    vector<int> newQuantities;
+    for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
+    {
+      int nbNodesInFace = quantities[iFace];
+      newQuantities.push_back(0);
+      for ( int i = 0; i < nbNodesInFace; ++i )
+      {
+        const SMDS_MeshNode* n1 = nodes[ iN + i ];
+        newNodes.push_back( n1 );
+        newQuantities.back()++;
+        
+        const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
+//         if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
+//              n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
+        {
+          const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
+          newNodes.push_back( n12 );
+          newQuantities.back()++;
+        }
+      }
+      iN += nbNodesInFace;
+    }
+    if(id)
+      elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
+    else
+      elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
+  }
+  if ( mySetElemOnShape && myShapeID > 0 )
+    meshDS->SetMeshElementOnShape( elem, myShapeID );
+
+  return elem;
+}
+
 //=======================================================================
 //function : LoadNodeColumns
 //purpose  : Load nodes bound to face into a map of node columns
@@ -1169,244 +1384,63 @@ bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
                                          const TopoDS_Edge& theBaseEdge,
                                          SMESHDS_Mesh*      theMesh)
 {
-  // get vertices of theBaseEdge
-  TopoDS_Vertex vfb, vlb, vft; // first and last, bottom and top vertices
-  TopoDS_Edge eFrw = TopoDS::Edge( theBaseEdge.Oriented( TopAbs_FORWARD ));
-  TopExp::Vertices( eFrw, vfb, vlb );
-
-  // find the other edges of theFace and orientation of e1
-  TopoDS_Edge e1, e2, eTop;
-  bool rev1, CumOri = false;
-  TopExp_Explorer exp( theFace, TopAbs_EDGE );
-  int nbEdges = 0;
-  for ( ; exp.More(); exp.Next() ) {
-    if ( ++nbEdges > 4 ) {
-      return false; // more than 4 edges in theFace
-    }
-    TopoDS_Edge e = TopoDS::Edge( exp.Current() );
-    if ( theBaseEdge.IsSame( e ))
-      continue;
-    TopoDS_Vertex vCommon;
-    if ( !TopExp::CommonVertex( theBaseEdge, e, vCommon ))
-      eTop = e;
-    else if ( vCommon.IsSame( vfb )) {
-      e1 = e;
-      vft = TopExp::LastVertex( e1, CumOri );
-      rev1 = vfb.IsSame( vft );
-      if ( rev1 )
-        vft = TopExp::FirstVertex( e1, CumOri );
-    }
-    else
-      e2 = e;
-  }
-  if ( nbEdges < 4 ) {
-    return false; // less than 4 edges in theFace
-  }
-  if ( e2.IsNull() && vfb.IsSame( vlb ))
-    e2 = e1;
-
-  // submeshes corresponding to shapes
-  SMESHDS_SubMesh* smFace = theMesh->MeshElements( theFace );
-  SMESHDS_SubMesh* smb = theMesh->MeshElements( theBaseEdge );
-  SMESHDS_SubMesh* smt = theMesh->MeshElements( eTop );
-  SMESHDS_SubMesh* sm1 = theMesh->MeshElements( e1 );
-  SMESHDS_SubMesh* sm2 = theMesh->MeshElements( e2 );
-  SMESHDS_SubMesh* smVfb = theMesh->MeshElements( vfb );
-  SMESHDS_SubMesh* smVlb = theMesh->MeshElements( vlb );
-  SMESHDS_SubMesh* smVft = theMesh->MeshElements( vft );
-  if (!smFace || !smb || !smt || !sm1 || !sm2 || !smVfb || !smVlb || !smVft ) {
-    RETURN_BAD_RESULT( "NULL submesh " <<smFace<<" "<<smb<<" "<<smt<<" "<<
-                       sm1<<" "<<sm2<<" "<<smVfb<<" "<<smVlb<<" "<<smVft);
-  }
-  if ( smb->NbNodes() != smt->NbNodes() || sm1->NbNodes() != sm2->NbNodes() ) {
-    RETURN_BAD_RESULT(" Diff nb of nodes on opposite edges" );
-  }
-  if (smVfb->NbNodes() != 1 || smVlb->NbNodes() != 1 || smVft->NbNodes() != 1) {
-    RETURN_BAD_RESULT("Empty submesh of vertex");
-  }
-  // define whether mesh is quadratic
-  bool isQuadraticMesh = false;
-  SMDS_ElemIteratorPtr eIt = smFace->GetElements();
-  if ( !eIt->more() ) {
-    RETURN_BAD_RESULT("No elements on the face");
-  }
-  const SMDS_MeshElement* e = eIt->next();
-  isQuadraticMesh = e->IsQuadratic();
-  
-  if ( sm1->NbNodes() * smb->NbNodes() != smFace->NbNodes() ) {
-    // check quadratic case
-    if ( isQuadraticMesh ) {
-      // what if there are quadrangles and triangles mixed?
-//       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 {
-      RETURN_BAD_RESULT( "Wrong nb face nodes: " <<
-                         sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
-    }
-  }
-  // IJ size
-  int vsize = sm1->NbNodes() + 2;
-  int hsize = smb->NbNodes() + 2;
-  if(isQuadraticMesh) {
-    vsize = vsize - sm1->NbNodes()/2 -1;
-    hsize = hsize - smb->NbNodes()/2 -1;
-  }
-
-  // load nodes from theBaseEdge
-
-  std::set<const SMDS_MeshNode*> loadedNodes;
-  const SMDS_MeshNode* nullNode = 0;
+  SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
+  if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
+    return false;
 
-  std::vector<const SMDS_MeshNode*> & nVecf = theParam2ColumnMap[ 0.];
-  nVecf.resize( vsize, nullNode );
-  loadedNodes.insert( nVecf[ 0 ] = smVfb->GetNodes()->next() );
+  // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
 
-  std::vector<const SMDS_MeshNode*> & nVecl = theParam2ColumnMap[ 1.];
-  nVecl.resize( vsize, nullNode );
-  loadedNodes.insert( nVecl[ 0 ] = smVlb->GetNodes()->next() );
+  map< double, const SMDS_MeshNode*> sortedBaseNodes;
+  if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
+       || sortedBaseNodes.size() < 2 )
+    return false;
 
-  double f, l;
-  BRep_Tool::Range( eFrw, f, l );
-  double range = l - f;
-  SMDS_NodeIteratorPtr nIt = smb->GetNodes();
-  const SMDS_MeshNode* node;
-  while ( nIt->more() ) {
-    node = nIt->next();
-    if(IsMedium(node, SMDSAbs_Edge))
-      continue;
-    const SMDS_EdgePosition* pos =
-      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
-    if ( !pos ) {
-      return false;
-    }
-    double u = ( pos->GetUParameter() - f ) / range;
-    std::vector<const SMDS_MeshNode*> & nVec = theParam2ColumnMap[ u ];
-    nVec.resize( vsize, nullNode );
-    loadedNodes.insert( nVec[ 0 ] = node );
-  }
-  if ( theParam2ColumnMap.size() != hsize ) {
-    RETURN_BAD_RESULT( "Wrong node positions on theBaseEdge" );
+  int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
+  map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
+  double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
+  for ( ; u_n != sortedBaseNodes.end(); u_n++ )
+  {
+    double par = ( u_n->first - f ) / range;
+    vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
+    nCol.resize( nbRows );
+    nCol[0] = u_n->second;
   }
 
-  // load nodes from e1
+  // fill theParam2ColumnMap column by column by passing from nodes on
+  // theBaseEdge up via mesh faces on theFace
 
-  std::map< double, const SMDS_MeshNode*> sortedNodes; // sort by param on edge
-  nIt = sm1->GetNodes();
-  while ( nIt->more() ) {
-    node = nIt->next();
-    if(IsMedium(node))
-      continue;
-    const SMDS_EdgePosition* pos =
-      dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
-    if ( !pos ) {
-      return false;
-    }
-    sortedNodes.insert( std::make_pair( pos->GetUParameter(), node ));
-  }
-  loadedNodes.insert( nVecf[ vsize - 1 ] = smVft->GetNodes()->next() );
-  std::map< double, const SMDS_MeshNode*>::iterator u_n = sortedNodes.begin();
-  int row  = rev1 ? vsize - 1 : 0;
-  int dRow = rev1 ? -1 : +1;
-  for ( ; u_n != sortedNodes.end(); u_n++ ) {
-    row += dRow;
-    loadedNodes.insert( nVecf[ row ] = u_n->second );
-  }
-
-  // try to load the rest nodes
-
-  // get all faces from theFace
-  TIDSortedElemSet allFaces, foundFaces;
-  eIt = smFace->GetElements();
-  while ( eIt->more() ) {
-    const SMDS_MeshElement* e = eIt->next();
-    if ( e->GetType() == SMDSAbs_Face )
-      allFaces.insert( e );
-  }
-  // Starting from 2 neighbour nodes on theBaseEdge, look for a face
-  // the nodes belong to, and between the nodes of the found face,
-  // look for a not loaded node considering this node to be the next
-  // in a column of the starting second node. Repeat, starting
-  // from nodes next to the previous starting nodes in their columns,
-  // and so on while a face can be found. Then go the the next pair
-  // of nodes on theBaseEdge.
-  TParam2ColumnMap::iterator par_nVec_1 = theParam2ColumnMap.begin();
-  TParam2ColumnMap::iterator par_nVec_2 = par_nVec_1;
-  // loop on columns
-  int col = 0;
-  for ( par_nVec_2++; par_nVec_2 != theParam2ColumnMap.end(); par_nVec_1++, par_nVec_2++ ) {
-    col++;
-    row = 0;
-    const SMDS_MeshNode* n1 = par_nVec_1->second[ row ];
-    const SMDS_MeshNode* n2 = par_nVec_2->second[ row ];
-    const SMDS_MeshElement* face = 0;
-    bool lastColOnClosedFace = ( nVecf[ row ] == n2 );
-    do {
-      // look for a face by 2 nodes
-      face = SMESH_MeshEditor::FindFaceInSet( n1, n2, allFaces, foundFaces );
-      if ( face ) {
-        int nbFaceNodes = face->NbNodes();
-        if ( face->IsQuadratic() )
-          nbFaceNodes /= 2;
-        if ( nbFaceNodes>4 ) {
-          RETURN_BAD_RESULT(" Too many nodes in a face: " << nbFaceNodes );
-        }
-        // look for a not loaded node of the <face>
-        bool found = false;
-        const SMDS_MeshNode* n3 = 0; // a node defferent from n1 and n2
-        for ( int i = 0; i < nbFaceNodes && !found; ++i ) {
-          node = face->GetNode( i );
-          found = loadedNodes.insert( node ).second;
-          if ( !found && node != n1 && node != n2 )
-            n3 = node;
-        }
-        if ( lastColOnClosedFace && row + 1 < vsize ) {
-          node = nVecf[ row + 1 ];
-          found = ( face->GetNodeIndex( node ) >= 0 );
-        }
-        if ( found ) {
-          if ( ++row > vsize - 1 ) {
-            RETURN_BAD_RESULT( "Too many nodes in column "<< col <<": "<< row+1);
-          }
-          par_nVec_2->second[ row ] = node;
-          foundFaces.insert( face );
-          n2 = node;
-          if ( nbFaceNodes==4 ) {
-            n1 = par_nVec_1->second[ row ];
-          }
-        }
-        else if ( nbFaceNodes==3 && n3 == par_nVec_1->second[ row + 1 ] ) {
-          n1 = n3;
-        }
-        else  {
-          RETURN_BAD_RESULT( "Not quad mesh, column "<< col );
-        }
+  TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
+  TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
+  TIDSortedElemSet emptySet, avoidSet;
+  for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
+  {
+    vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
+    vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
+
+    int i1, i2, iRow = 0;
+    const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
+    // find face sharing node n1 and n2 and belonging to faceSubMesh
+    while ( const SMDS_MeshElement* face =
+            SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
+    {
+      if ( faceSubMesh->Contains( face ))
+      {
+        int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
+        if ( nbNodes != 4 )
+          return false;
+        n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
+        n2 = face->GetNode( (i1+2) % 4 );
+        if ( ++iRow >= nbRows )
+          return false;
+        nCol1[ iRow ] = n1;
+        nCol2[ iRow ] = n2;
+        avoidSet.clear();
       }
+      avoidSet.insert( face );
     }
-    while ( face && n1 && n2 );
-
-    if ( row < vsize - 1 ) {
-      MESSAGE( "Too few nodes in column "<< col <<": "<< row+1);
-      MESSAGE( "Base node 1: "<< par_nVec_1->second[0]);
-      MESSAGE( "Base node 2: "<< par_nVec_2->second[0]);
-      if ( n1 ) { MESSAGE( "Current node 1: "<< n1); }
-      else      { MESSAGE( "Current node 1: NULL");  }
-      if ( n2 ) { MESSAGE( "Current node 2: "<< n2); }
-      else      { MESSAGE( "Current node 2: NULL");  }
-      MESSAGE( "first base node: "<< theParam2ColumnMap.begin()->second[0]);
-      MESSAGE( "last base node: "<< theParam2ColumnMap.rbegin()->second[0]);
-      return false;
-    }
-  } // loop on columns
-
+    if ( iRow + 1 < nbRows ) // compact if necessary
+      nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
+  }
   return true;
 }
 
@@ -1487,6 +1521,26 @@ bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMes
     shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape );
 }
 
+//================================================================================
+/*!
+ * \brief Return maximal tolerance of shape
+ */
+//================================================================================
+
+double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
+{
+  double tol = Precision::Confusion();
+  TopExp_Explorer exp;
+  for ( exp.Init( shape, TopAbs_FACE ); exp.More(); exp.Next() )
+    tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Face( exp.Current())));
+  for ( exp.Init( shape, TopAbs_EDGE ); exp.More(); exp.Next() )
+    tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Edge( exp.Current())));
+  for ( exp.Init( shape, TopAbs_VERTEX ); exp.More(); exp.Next() )
+    tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Vertex( exp.Current())));
+
+  return tol;
+}
+
 //=======================================================================
 //function : IsQuadraticMesh
 //purpose  : Check mesh without geometry for: if all elements on this shape are quadratic,
@@ -1532,16 +1586,23 @@ double SMESH_MesherHelper::GetOtherParam(const double param) const
   return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
 }
 
+//#include <Perf_Meter.hxx>
+
 //=======================================================================
 namespace { // Structures used by FixQuadraticElements()
 //=======================================================================
 
 #define __DMP__(txt) \
-//cout << txt
+  //cout << txt
 #define MSG(txt) __DMP__(txt<<endl)
 #define MSGBEG(txt) __DMP__(txt)
 
-  const double straightTol2 = 1e-33; // to detect straing links
+  //const double straightTol2 = 1e-33; // to detect straing links
+  bool isStraightLink(double linkLen2, double middleNodeMove2)
+  {
+    // straight if <node move> < 1/15 * <link length>
+    return middleNodeMove2 < 1/15./15. * linkLen2;
+  }
 
   struct QFace;
   // ---------------------------------------
@@ -1578,8 +1639,10 @@ namespace { // Structures used by FixQuadraticElements()
     { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
     gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
     bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
-    bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
-
+    bool IsStraight() const
+    { return isStraightLink( (XYZ(node1())-XYZ(node2())).SquareModulus(),
+                             _nodeMove.SquareMagnitude());
+    }
     bool operator<(const QLink& other) const {
       return (node1()->GetID() == other.node1()->GetID() ?
               node2()->GetID() < other.node2()->GetID() :
@@ -1601,7 +1664,7 @@ namespace { // Structures used by FixQuadraticElements()
     TChainLink(const QLink* qlink=0):_qlink(qlink) {
       _qfaces[0] = _qfaces[1] = 0;
     }
-    void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
+    void SetFace(const QFace* face) const { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
 
     bool IsBoundary() const { return !_qfaces[1]; }
 
@@ -1760,7 +1823,7 @@ namespace { // Structures used by FixQuadraticElements()
   }
   //================================================================================
   /*!
-   * \brief Make up chain of links
+   * \brief Make up chain of links
    *  \param iSide - link to add first
    *  \param chain - chain to fill in
    *  \param pos   - postion of medium nodes the links should have
@@ -1779,22 +1842,33 @@ namespace { // Structures used by FixQuadraticElements()
 
     if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
       MSGBEG( *this );
+      TLinkSet links;
       list< const QFace* > faces( 1, this );
-      for (list< const QFace* >::iterator fIt = faces.begin(); fIt != faces.end(); ++fIt ) {
-        const QFace* face = *fIt;
+      while ( !faces.empty() ) {
+        const QFace* face = faces.front();
         for ( int i = 0; i < face->_sides.size(); ++i ) {
           if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
             face->_sideIsAdded[i] = true;
-            TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
+            // find a face side in the chain
+            TLinkInSet chLink = links.insert( TChainLink(face->_sides[i])).first;
+//             TChain::iterator chLink = chain.begin();
+//             for ( ; chLink != chain.end(); ++chLink )
+//               if ( chLink->_qlink == face->_sides[i] )
+//                 break;
+//             if ( chLink == chain.end() )
+//               chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
+            // add a face to a chained link and put a continues face in the queue
             chLink->SetFace( face );
             if ( face->_sides[i]->MediumPos() >= pos )
               if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
                 faces.push_back( contFace );
           }
         }
+        faces.pop_front();
       }
       if ( error < ERR_TRI )
         error = ERR_TRI;
+      chain.insert( chain.end(), links.begin(),links.end() );
       return false;
     }
     _sideIsAdded[iSide] = true; // not to add this link to chain again
@@ -1860,7 +1934,7 @@ namespace { // Structures used by FixQuadraticElements()
       TLinkInSet link = links.find( _sides[iL] );
       if ( link == linksEnd ) continue;
       if ( (*link)->MediumPos() > SMDS_TOP_FACE )
-        continue; // We work on faces here, don't go into a volume
+        continue; // We work on faces here, don't go inside a solid
 
       // check link
       if ( link->IsBoundary() ) {
@@ -2003,7 +2077,8 @@ namespace { // Structures used by FixQuadraticElements()
     // propagate to adjacent faces till limit step or boundary
     double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
     double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
-    gp_Vec linkDir1, linkDir2;
+    gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
+    gp_Vec linkDir2(0,0,0);
     try {
       OCC_CATCH_SIGNALS;
       if ( f1 )
@@ -2239,7 +2314,7 @@ namespace { // Structures used by FixQuadraticElements()
 
   enum TSplitTriaResult {
     _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
-    _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK };
+    _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
 
   TSplitTriaResult splitTrianglesIntoChains( TChain &            allLinks,
                                              vector< TChain> &   resultChains,
@@ -2305,6 +2380,8 @@ namespace { // Structures used by FixQuadraticElements()
         const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
         if ( !botTria )
         { // the column ends
+          if ( botLink == startLink )
+            return _TWISTED_CHAIN; // issue 0020951
           linkSet.erase( botLink );
           if ( iRow != rowChains.size() )
             return _FEW_ROWS; // different nb of rows in columns
@@ -2364,8 +2441,10 @@ namespace { // Structures used by FixQuadraticElements()
         // next bottom link ends at the new corner
         linkSet.erase( botLink );
         botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
-        if ( botLink == linksEnd || botLink == (isCase2 ? midQuadLink : sideLink))
+        if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
           return _NO_BOTLINK;
+        if ( midQuadLink == startLink || sideLink == startLink )
+          return _TWISTED_CHAIN; // issue 0020951
         linkSet.erase( midQuadLink );
         linkSet.erase( sideLink );
 
@@ -2409,7 +2488,7 @@ namespace { // Structures used by FixQuadraticElements()
 
     return _OK;
   }
-}
+} //namespace
 
 //=======================================================================
 /*!
@@ -2422,37 +2501,49 @@ namespace { // Structures used by FixQuadraticElements()
 
 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
 {
-  // apply algorithm to solids or geom faces
+  // 0. Apply algorithm to solids or geom faces
   // ----------------------------------------------
   if ( myShape.IsNull() ) {
     if ( !myMesh->HasShapeToMesh() ) return;
     SetSubShape( myMesh->GetShapeToMesh() );
 
+#ifdef _DEBUG_
+    int nbSolids = 0;
+    TopTools_IndexedMapOfShape solids;
+    TopExp::MapShapes(myShape,TopAbs_SOLID,solids);
+    nbSolids = solids.Extent();
+#endif
     TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
     for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
-      faces.Add( f.Current() );
+      faces.Add( f.Current() ); // not in solid
     }
     for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
       if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
         for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
-          faces.Add( f.Current() );
+          faces.Add( f.Current() ); // in not meshed solid
       }
       else { // fix nodes in the solid and its faces
+        MSG("FIX SOLID " << nbSolids-- << " #" << GetMeshDS()->ShapeToIndex(s.Current()));
         SMESH_MesherHelper h(*myMesh);
         h.SetSubShape( s.Current() );
         h.FixQuadraticElements(false);
       }
     }
     // fix nodes on geom faces
+#ifdef _DEBUG_
+    int nbfaces = faces.Extent();
+#endif
     for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
+      MSG("FIX FACE " << nbfaces-- << " #" << GetMeshDS()->ShapeToIndex(fIt.Key()));
       SMESH_MesherHelper h(*myMesh);
       h.SetSubShape( fIt.Key() );
       h.FixQuadraticElements(true);
     }
+    //perf_print_all_meters(1);
     return;
   }
 
-  // Find out type of elements and get iterator on them
+  // 1. Find out type of elements and get iterator on them
   // ---------------------------------------------------
 
   SMDS_ElemIteratorPtr elemIt;
@@ -2471,7 +2562,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
   if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
     return;
 
-  // Fill in auxiliary data structures
+  // 2. Fill in auxiliary data structures
   // ----------------------------------
 
   set< QLink > links;
@@ -2480,7 +2571,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
   set< QFace >::iterator pFace;
 
   bool isCurved = false;
-  bool hasRectFaces = false;
+  //bool hasRectFaces = false;
   set<int> nbElemNodeSet;
 
   if ( elemType == SMDSAbs_Volume )
@@ -2513,9 +2604,9 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
         if ( pFace->NbVolumes() == 0 )
           pFace->AddSelfToLinks();
         pFace->SetVolume( vol );
-        hasRectFaces = hasRectFaces ||
-          ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
-            volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
+//         hasRectFaces = hasRectFaces ||
+//           ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
+//             volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
 #ifdef _DEBUG_
         if ( nbN == 6 )
           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
@@ -2551,13 +2642,13 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
       // store QFace
       pFace = faces.insert( QFace( faceLinks )).first;
       pFace->AddSelfToLinks();
-      hasRectFaces = ( hasRectFaces || nbN == 4 );
+      //hasRectFaces = ( hasRectFaces || nbN == 4 );
     }
   }
   if ( !isCurved )
     return; // no curved edges of faces
 
-  // Compute displacement of medium nodes
+  // 3. Compute displacement of medium nodes
   // -------------------------------------
 
   // two loops on faces: the first is to treat boundary links, the second is for internal ones
@@ -2567,6 +2658,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
   for ( ; isInside < 2; ++isInside ) {
     MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
     SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
+    SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
 
     for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
       if ( bool(isInside) == pFace->IsBoundary() )
@@ -2605,7 +2697,12 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
           TChain& chain = chains[iC];
           if ( chain.empty() ) continue;
           if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
-            MSG("3D straight");
+            MSG("3D straight - ignore");
+            continue;
+          }
+          if ( chain.front()->MediumPos() > bndPos ||
+               chain.back()->MediumPos() > bndPos ) {
+            MSG("Internal chain - ignore");
             continue;
           }
           // mesure chain length and compute link position along the chain
@@ -2642,6 +2739,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
             {
               face = TopoDS::Face( f );
               Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
+              bool isStraight[2];
               for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
               {
                 TChainLink& link = is1 ? chain.back() : chain.front();
@@ -2652,10 +2750,14 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
                 // uvMove = uvm - uv12
                 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
                 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
+                if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
+                  nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
+                isStraight[is1] = isStraightLink( (uv2-uv1).SquareModulus(),uvMove.SquareModulus());
               }
-              if ( move0.SquareMagnitude() < straightTol2 &&
-                   move1.SquareMagnitude() < straightTol2 ) {
-                MSG("2D straight");
+//               if ( move0.SquareMagnitude() < straightTol2 &&
+//                    move1.SquareMagnitude() < straightTol2 ) {
+              if ( isStraight[0] && isStraight[1] ) {
+                MSG("2D straight - ignore");
                 continue; // straight - no need to move nodes of internal links
               }
             }
@@ -2707,7 +2809,8 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
               {
                 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
                 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
-                MSG( "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
+                MSG( "TOO LONG MOVE \t" <<
+                     "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
                      "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
                      "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
                      "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
@@ -2725,7 +2828,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
     } // loop on faces
   }
 
-  // Move nodes
+  // 4. Move nodes
   // -----------
 
   for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {