Salome HOME
IMP 22635: EDF 8345 - Creation of group based on groups
[modules/smesh.git] / src / SMESH / SMESH_MeshEditor.cxx
index 931caa43414d2a0aea00e1bffb4e2a288bf62963..b3401c601fc5c2fd06f21761a7852e98b3a43223 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2015  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
@@ -50,6 +50,7 @@
 #include <Basics_OCCTVersion.hxx>
 
 #include "utilities.h"
+#include "chrono.hxx"
 
 #include <BRepAdaptor_Surface.hxx>
 #include <BRepBuilderAPI_MakeEdge.hxx>
@@ -130,7 +131,7 @@ SMESH_MeshEditor::SMESH_MeshEditor( SMESH_Mesh* theMesh )
  */
 //================================================================================
 
-void SMESH_MeshEditor::CrearLastCreated()
+void SMESH_MeshEditor::ClearLastCreated()
 {
   myLastCreatedNodes.Clear();
   myLastCreatedElems.Clear();
@@ -518,14 +519,12 @@ int SMESH_MeshEditor::FindShape (const SMDS_MeshElement * theElem)
   }
   else
   {
-    const map<int,SMESHDS_SubMesh*>& id2sm = GetMeshDS()->SubMeshes();
-    map<int,SMESHDS_SubMesh*>::const_iterator id_sm = id2sm.begin();
-    for ( ; id_sm != id2sm.end(); ++id_sm )
-      if ( id_sm->second->Contains( theElem ))
-        return id_sm->first;
+    SMESHDS_SubMeshIteratorPtr smIt = GetMeshDS()->SubMeshes();
+    while ( const SMESHDS_SubMesh* sm = smIt->next() )
+      if ( sm->Contains( theElem ))
+        return sm->GetID();
   }
 
-  //MESSAGE ("::FindShape() - SHAPE NOT FOUND")
   return 0;
 }
 
@@ -1174,7 +1173,7 @@ int SMESH_MeshEditor::Reorient2D (TIDSortedElemSet &       theFaces,
     avoidSet.clear();
     avoidSet.insert(theFace);
 
-    NLink link( theFace->GetNode( 0 ), 0 );
+    NLink link( theFace->GetNode( 0 ), (SMDS_MeshNode *) 0 );
 
     const int nbNodes = theFace->NbCornerNodes();
     for ( int i = 0; i < nbNodes; ++i ) // loop on links of theFace
@@ -1245,6 +1244,92 @@ int SMESH_MeshEditor::Reorient2D (TIDSortedElemSet &       theFaces,
   return nbReori;
 }
 
+//================================================================================
+/*!
+ * \brief Reorient faces basing on orientation of adjacent volumes.
+ * \param theFaces - faces to reorient. If empty, all mesh faces are treated.
+ * \param theVolumes - reference volumes.
+ * \param theOutsideNormal - to orient faces to have their normal
+ *        pointing either \a outside or \a inside the adjacent volumes.
+ * \return number of reoriented faces.
+ */
+//================================================================================
+
+int SMESH_MeshEditor::Reorient2DBy3D (TIDSortedElemSet & theFaces,
+                                      TIDSortedElemSet & theVolumes,
+                                      const bool         theOutsideNormal)
+{
+  int nbReori = 0;
+
+  SMDS_ElemIteratorPtr faceIt;
+  if ( theFaces.empty() )
+    faceIt = GetMeshDS()->elementsIterator( SMDSAbs_Face );
+  else
+    faceIt = elemSetIterator( theFaces );
+
+  vector< const SMDS_MeshNode* > faceNodes;
+  TIDSortedElemSet checkedVolumes;
+  set< const SMDS_MeshNode* > faceNodesSet;
+  SMDS_VolumeTool volumeTool;
+
+  while ( faceIt->more() ) // loop on given faces
+  {
+    const SMDS_MeshElement* face = faceIt->next();
+    if ( face->GetType() != SMDSAbs_Face )
+      continue;
+
+    const int nbCornersNodes = face->NbCornerNodes();
+    faceNodes.assign( face->begin_nodes(), face->end_nodes() );
+
+    checkedVolumes.clear();
+    SMDS_ElemIteratorPtr vIt = faceNodes[ 0 ]->GetInverseElementIterator( SMDSAbs_Volume );
+    while ( vIt->more() )
+    {
+      const SMDS_MeshElement* volume = vIt->next();
+
+      if ( !checkedVolumes.insert( volume ).second )
+        continue;
+      if ( !theVolumes.empty() && !theVolumes.count( volume ))
+        continue;
+
+      // is volume adjacent?
+      bool allNodesCommon = true;
+      for ( int iN = 1; iN < nbCornersNodes && allNodesCommon; ++iN )
+        allNodesCommon = ( volume->GetNodeIndex( faceNodes[ iN ]) > -1 );
+      if ( !allNodesCommon )
+        continue;
+
+      // get nodes of a corresponding volume facet
+      faceNodesSet.clear();
+      faceNodesSet.insert( faceNodes.begin(), faceNodes.end() );
+      volumeTool.Set( volume );
+      int facetID = volumeTool.GetFaceIndex( faceNodesSet );
+      if ( facetID < 0 ) continue;
+      volumeTool.SetExternalNormal();
+      const SMDS_MeshNode** facetNodes = volumeTool.GetFaceNodes( facetID );
+
+      // compare order of faceNodes and facetNodes
+      const int iQ = 1 + ( nbCornersNodes < faceNodes.size() );
+      int iNN[2];
+      for ( int i = 0; i < 2; ++i )
+      {
+        const SMDS_MeshNode* n = facetNodes[ i*iQ ];
+        for ( int iN = 0; iN < nbCornersNodes; ++iN )
+          if ( faceNodes[ iN ] == n )
+          {
+            iNN[ i ] = iN;
+            break;
+          }
+      }
+      bool isOutside = Abs( iNN[0]-iNN[1] ) == 1 ? iNN[0] < iNN[1] : iNN[0] > iNN[1];
+      if ( isOutside != theOutsideNormal )
+        nbReori += Reorient( face );
+    }
+  }  // loop on given faces
+
+  return nbReori;
+}
+
 //=======================================================================
 //function : getBadRate
 //purpose  :
@@ -1615,7 +1700,7 @@ namespace
     };
   const int theHexTo4Prisms_FB[6*4+1] = // front-back
     {
-      0, 3, 8, 1, 2, 9,    3, 7, 8, 2, 6, 9,    7, 4, 8, 6, 5, 9,   4, 0, 8, 5, 1, 9,    -1
+      0, 3, 9, 1, 2, 8,    3, 7, 9, 2, 6, 8,    7, 4, 9, 6, 5, 8,   4, 0, 9, 5, 1, 8,    -1
     };
 
   const int theHexTo2Prisms_BT_1[6*2+1] =
@@ -2095,10 +2180,6 @@ namespace
 void SMESH_MeshEditor::SplitVolumes (const TFacetOfElem & theElems,
                                      const int            theMethodFlags)
 {
-  // std-like iterator on coordinates of nodes of mesh element
-  typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > NXyzIterator;
-  NXyzIterator xyzEnd;
-
   SMDS_VolumeTool    volTool;
   SMESH_MesherHelper helper( *GetMesh()), fHelper(*GetMesh());
   fHelper.ToFixNodeParameters( true );
@@ -3617,13 +3698,8 @@ static bool getClosestUV (Extrema_GenExtPS& projector,
   if ( projector.IsDone() ) {
     double u, v, minVal = DBL_MAX;
     for ( int i = projector.NbExt(); i > 0; i-- )
-#if OCC_VERSION_LARGE > 0x06040000 // Porting to OCCT6.5.1
       if ( projector.SquareDistance( i ) < minVal ) {
         minVal = projector.SquareDistance( i );
-#else
-      if ( projector.Value( i ) < minVal ) {
-        minVal = projector.Value( i );
-#endif
         projector.Point( i ).Parameter( u, v );
       }
     result.SetCoord( u, v );
@@ -4259,6 +4335,30 @@ void SMESH_MeshEditor::sweepElement(const SMDS_MeshElement*               elem,
       }
     }
   }
+  else if ( elem->GetType() == SMDSAbs_Edge )
+  {
+    // orient a new face same as adjacent one
+    int i1, i2;
+    const SMDS_MeshElement* e;
+    TIDSortedElemSet dummy;
+    if (( e = SMESH_MeshAlgos::FindFaceInSet( nextNod[0], prevNod[0], dummy,dummy, &i1, &i2 )) ||
+        ( e = SMESH_MeshAlgos::FindFaceInSet( prevNod[1], nextNod[1], dummy,dummy, &i1, &i2 )) ||
+        ( e = SMESH_MeshAlgos::FindFaceInSet( prevNod[0], prevNod[1], dummy,dummy, &i1, &i2 )))
+    {
+      // there is an adjacent face, check order of nodes in it
+      bool sameOrder = ( Abs( i2 - i1 ) == 1 ) ? ( i2 > i1 ) : ( i2 < i1 );
+      if ( sameOrder )
+      {
+        std::swap( itNN[0],    itNN[1] );
+        std::swap( prevNod[0], prevNod[1] );
+        std::swap( nextNod[0], nextNod[1] );
+        isSingleNode.swap( isSingleNode[0], isSingleNode[1] );
+        if ( nbSame > 0 )
+          sames[0] = 1 - sames[0];
+        iNotSameNode = 1 - iNotSameNode;
+      }
+    }
+  }
 
   int iSameNode = 0, iBeforeSame = 0, iAfterSame = 0, iOpposSame = 0;
   if ( nbSame > 0 ) {
@@ -4364,11 +4464,11 @@ void SMESH_MeshEditor::sweepElement(const SMDS_MeshElement*               elem,
                 return; // medium node on axis
               }
               else if(sames[0]==0)
-                aNewElem = aMesh->AddFace(prevNod[0], nextNod[1], prevNod[1],
-                                          nextNod[2], midlNod[1], prevNod[2]);
+                aNewElem = aMesh->AddFace(prevNod[0], prevNod[1], nextNod[1],
+                                          prevNod[2], midlNod[1], nextNod[2] );
               else // sames[0]==1
-                aNewElem = aMesh->AddFace(prevNod[0], nextNod[0], prevNod[1],
-                                          midlNod[0], nextNod[2], prevNod[2]);
+                aNewElem = aMesh->AddFace(prevNod[0], prevNod[1], nextNod[0],
+                                          prevNod[2], nextNod[2], midlNod[0]);
             }
           }
           else if ( nbDouble == 3 )
@@ -4534,8 +4634,8 @@ void SMESH_MeshEditor::sweepElement(const SMDS_MeshElement*               elem,
 
       default:
         break;
-      }
-    }
+      } // switch ( baseType )
+    } // scope
 
     if ( !aNewElem && elem->GetType() == SMDSAbs_Face ) // try to create a polyherdal prism
     {
@@ -4589,7 +4689,8 @@ void SMESH_MeshEditor::sweepElement(const SMDS_MeshElement*               elem,
           polyedre_nodes.resize( prevNbNodes );
       }
       aNewElem = aMesh->AddPolyhedralVolume (polyedre_nodes, quantities);
-    }
+
+    } //  // try to create a polyherdal prism
 
     if ( aNewElem ) {
       newElems.push_back( aNewElem );
@@ -4601,7 +4702,7 @@ void SMESH_MeshEditor::sweepElement(const SMDS_MeshElement*               elem,
     for ( iNode = 0; iNode < nbNodes; iNode++ )
       prevNod[ iNode ] = nextNod[ iNode ];
 
-  } // for steps
+  } // loop on steps
 }
 
 //=======================================================================
@@ -5138,50 +5239,328 @@ SMESH_MeshEditor::RotationSweep(TIDSortedElemSet & theElems,
   return newGroupIDs;
 }
 
+//=======================================================================
+//function : ExtrusParam
+//purpose  : standard construction
+//=======================================================================
+
+SMESH_MeshEditor::ExtrusParam::ExtrusParam( const gp_Vec&  theStep,
+                                            const int      theNbSteps,
+                                            const int      theFlags,
+                                            const double   theTolerance):
+  myDir( theStep ),
+  myFlags( theFlags ),
+  myTolerance( theTolerance ),
+  myElemsToUse( NULL )
+{
+  mySteps = new TColStd_HSequenceOfReal;
+  const double stepSize = theStep.Magnitude();
+  for (int i=1; i<=theNbSteps; i++ )
+    mySteps->Append( stepSize );
+
+  if (( theFlags & EXTRUSION_FLAG_SEW ) &&
+      ( theTolerance > 0 ))
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByDirAndSew;
+  }
+  else
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByDir;
+  }
+}
 
 //=======================================================================
-//function : CreateNode
-//purpose  :
+//function : ExtrusParam
+//purpose  : steps are given explicitly
+//=======================================================================
+
+SMESH_MeshEditor::ExtrusParam::ExtrusParam( const gp_Dir&                   theDir,
+                                            Handle(TColStd_HSequenceOfReal) theSteps,
+                                            const int                       theFlags,
+                                            const double                    theTolerance):
+  myDir( theDir ),
+  mySteps( theSteps ),
+  myFlags( theFlags ),
+  myTolerance( theTolerance ),
+  myElemsToUse( NULL )
+{
+  if (( theFlags & EXTRUSION_FLAG_SEW ) &&
+      ( theTolerance > 0 ))
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByDirAndSew;
+  }
+  else
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByDir;
+  }
+}
+
+//=======================================================================
+//function : ExtrusParam
+//purpose  : for extrusion by normal
+//=======================================================================
+
+SMESH_MeshEditor::ExtrusParam::ExtrusParam( const double theStepSize,
+                                            const int    theNbSteps,
+                                            const int    theFlags,
+                                            const int    theDim ):
+  myDir( 1,0,0 ),
+  mySteps( new TColStd_HSequenceOfReal ),
+  myFlags( theFlags ),
+  myTolerance( 0 ),
+  myElemsToUse( NULL )
+{
+  for (int i = 0; i < theNbSteps; i++ )
+    mySteps->Append( theStepSize );
+
+  if ( theDim == 1 )
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByNormal1D;
+  }
+  else
+  {
+    myMakeNodesFun = & SMESH_MeshEditor::ExtrusParam::makeNodesByNormal2D;
+  }
+}
+
 //=======================================================================
-const SMDS_MeshNode* SMESH_MeshEditor::CreateNode(const double x,
-                                                  const double y,
-                                                  const double z,
-                                                  const double tolnode,
-                                                  SMESH_SequenceOfNode& aNodes)
+//function : ExtrusParam::SetElementsToUse
+//purpose  : stores elements to use for extrusion by normal, depending on
+//           state of EXTRUSION_FLAG_USE_INPUT_ELEMS_ONLY flag
+//=======================================================================
+
+void SMESH_MeshEditor::ExtrusParam::SetElementsToUse( const TIDSortedElemSet& elems )
 {
-  // myLastCreatedElems.Clear();
-  // myLastCreatedNodes.Clear();
+  myElemsToUse = ToUseInpElemsOnly() ? & elems : 0;
+}
 
-  gp_Pnt P1(x,y,z);
-  SMESHDS_Mesh * aMesh = myMesh->GetMeshDS();
+//=======================================================================
+//function : ExtrusParam::beginStepIter
+//purpose  : prepare iteration on steps
+//=======================================================================
 
-  // try to search in sequence of existing nodes
-  // if aNodes.Length()>0 we 'nave to use given sequence
-  // else - use all nodes of mesh
-  if(aNodes.Length()>0) {
-    int i;
-    for(i=1; i<=aNodes.Length(); i++) {
-      gp_Pnt P2(aNodes.Value(i)->X(),aNodes.Value(i)->Y(),aNodes.Value(i)->Z());
-      if(P1.Distance(P2)<tolnode)
-        return aNodes.Value(i);
+void SMESH_MeshEditor::ExtrusParam::beginStepIter( bool withMediumNodes )
+{
+  myWithMediumNodes = withMediumNodes;
+  myNextStep = 1;
+  myCurSteps.clear();
+}
+//=======================================================================
+//function : ExtrusParam::moreSteps
+//purpose  : are there more steps?
+//=======================================================================
+
+bool SMESH_MeshEditor::ExtrusParam::moreSteps()
+{
+  return myNextStep <= mySteps->Length() || !myCurSteps.empty();
+}
+//=======================================================================
+//function : ExtrusParam::nextStep
+//purpose  : returns the next step
+//=======================================================================
+
+double SMESH_MeshEditor::ExtrusParam::nextStep()
+{
+  double res = 0;
+  if ( !myCurSteps.empty() )
+  {
+    res = myCurSteps.back();
+    myCurSteps.pop_back();
+  }
+  else if ( myNextStep <= mySteps->Length() )
+  {
+    myCurSteps.push_back( mySteps->Value( myNextStep ));
+    ++myNextStep;
+    if ( myWithMediumNodes )
+    {
+      myCurSteps.back() /= 2.;
+      myCurSteps.push_back( myCurSteps.back() );
     }
+    res = nextStep();
   }
-  else {
-    SMDS_NodeIteratorPtr itn = aMesh->nodesIterator();
-    while(itn->more()) {
-      const SMDS_MeshNode* aN = static_cast<const SMDS_MeshNode*> (itn->next());
-      gp_Pnt P2(aN->X(),aN->Y(),aN->Z());
-      if(P1.Distance(P2)<tolnode)
-        return aN;
+  return res;
+}
+
+//=======================================================================
+//function : ExtrusParam::makeNodesByDir
+//purpose  : create nodes for standard extrusion
+//=======================================================================
+
+int SMESH_MeshEditor::ExtrusParam::
+makeNodesByDir( SMESHDS_Mesh*                     mesh,
+                const SMDS_MeshNode*              srcNode,
+                std::list<const SMDS_MeshNode*> & newNodes,
+                const bool                        makeMediumNodes)
+{
+  gp_XYZ p = SMESH_TNodeXYZ( srcNode );
+
+  int nbNodes = 0;
+  for ( beginStepIter( makeMediumNodes ); moreSteps(); ++nbNodes ) // loop on steps
+  {
+    p += myDir.XYZ() * nextStep();
+    const SMDS_MeshNode * newNode = mesh->AddNode( p.X(), p.Y(), p.Z() );
+    newNodes.push_back( newNode );
+  }
+  return nbNodes;
+}
+
+//=======================================================================
+//function : ExtrusParam::makeNodesByDirAndSew
+//purpose  : create nodes for standard extrusion with sewing
+//=======================================================================
+
+int SMESH_MeshEditor::ExtrusParam::
+makeNodesByDirAndSew( SMESHDS_Mesh*                     mesh,
+                      const SMDS_MeshNode*              srcNode,
+                      std::list<const SMDS_MeshNode*> & newNodes,
+                      const bool                        makeMediumNodes)
+{
+  gp_XYZ P1 = SMESH_TNodeXYZ( srcNode );
+
+  int nbNodes = 0;
+  for ( beginStepIter( makeMediumNodes ); moreSteps(); ++nbNodes ) // loop on steps
+  {
+    P1 += myDir.XYZ() * nextStep();
+
+    // try to search in sequence of existing nodes
+    // if myNodes.Length()>0 we 'nave to use given sequence
+    // else - use all nodes of mesh
+    const SMDS_MeshNode * node = 0;
+    if ( myNodes.Length() > 0 ) {
+      int i;
+      for(i=1; i<=myNodes.Length(); i++) {
+        gp_XYZ P2 = SMESH_TNodeXYZ( myNodes.Value(i) );
+        if (( P1 - P2 ).SquareModulus() < myTolerance * myTolerance )
+        {
+          node = myNodes.Value(i);
+          break;
+        }
+      }
+    }
+    else {
+      SMDS_NodeIteratorPtr itn = mesh->nodesIterator();
+      while(itn->more()) {
+        SMESH_TNodeXYZ P2( itn->next() );
+        if (( P1 - P2 ).SquareModulus() < myTolerance * myTolerance )
+        {
+          node = P2._node;
+          break;
+        }
+      }
+    }
+
+    if ( !node )
+      node = mesh->AddNode( P1.X(), P1.Y(), P1.Z() );
+
+    newNodes.push_back( node );
+
+  } // loop on steps
+
+  return nbNodes;
+}
+
+//=======================================================================
+//function : ExtrusParam::makeNodesByNormal2D
+//purpose  : create nodes for extrusion using normals of faces
+//=======================================================================
+
+int SMESH_MeshEditor::ExtrusParam::
+makeNodesByNormal2D( SMESHDS_Mesh*                     mesh,
+                     const SMDS_MeshNode*              srcNode,
+                     std::list<const SMDS_MeshNode*> & newNodes,
+                     const bool                        makeMediumNodes)
+{
+  const bool alongAvgNorm = ( myFlags & EXTRUSION_FLAG_BY_AVG_NORMAL );
+
+  gp_XYZ p = SMESH_TNodeXYZ( srcNode );
+
+  // get normals to faces sharing srcNode
+  vector< gp_XYZ > norms, baryCenters;
+  gp_XYZ norm, avgNorm( 0,0,0 );
+  SMDS_ElemIteratorPtr faceIt = srcNode->GetInverseElementIterator( SMDSAbs_Face );
+  while ( faceIt->more() )
+  {
+    const SMDS_MeshElement* face = faceIt->next();
+    if ( myElemsToUse && !myElemsToUse->count( face ))
+      continue;
+    if ( SMESH_MeshAlgos::FaceNormal( face, norm, /*normalized=*/true ))
+    {
+      norms.push_back( norm );
+      avgNorm += norm;
+      if ( !alongAvgNorm )
+      {
+        gp_XYZ bc(0,0,0);
+        int nbN = 0;
+        for ( SMDS_ElemIteratorPtr nIt = face->nodesIterator(); nIt->more(); ++nbN )
+          bc += SMESH_TNodeXYZ( nIt->next() );
+        baryCenters.push_back( bc / nbN );
+      }
     }
   }
 
-  // create new node and return it
-  const SMDS_MeshNode* NewNode = aMesh->AddNode(x,y,z);
-  //myLastCreatedNodes.Append(NewNode);
-  return NewNode;
+  if ( norms.empty() ) return 0;
+
+  double normSize = avgNorm.Modulus();
+  if ( normSize < std::numeric_limits<double>::min() )
+    return 0;
+
+  if ( myFlags & EXTRUSION_FLAG_BY_AVG_NORMAL ) // extrude along avgNorm
+  {
+    myDir = avgNorm;
+    return makeNodesByDir( mesh, srcNode, newNodes, makeMediumNodes );
+  }
+
+  avgNorm /= normSize;
+
+  int nbNodes = 0;
+  for ( beginStepIter( makeMediumNodes ); moreSteps(); ++nbNodes ) // loop on steps
+  {
+    gp_XYZ pNew = p;
+    double stepSize = nextStep();
+
+    if ( norms.size() > 1 )
+    {
+      for ( size_t iF = 0; iF < norms.size(); ++iF ) // loop on faces
+      {
+        // translate plane of a face
+        baryCenters[ iF ] += norms[ iF ] * stepSize;
+
+        // find point of intersection of the face plane located at baryCenters[ iF ]
+        // and avgNorm located at pNew
+        double d    = -( norms[ iF ] * baryCenters[ iF ]); // d of plane equation ax+by+cz+d=0
+        double dot  = ( norms[ iF ] * avgNorm );
+        if ( dot < std::numeric_limits<double>::min() )
+          dot = stepSize * 1e-3;
+        double step = -( norms[ iF ] * pNew + d ) / dot;
+        pNew += step * avgNorm;
+      }
+    }
+    else
+    {
+      pNew += stepSize * avgNorm;
+    }
+    p = pNew;
+
+    const SMDS_MeshNode * newNode = mesh->AddNode( p.X(), p.Y(), p.Z() );
+    newNodes.push_back( newNode );
+  }
+  return nbNodes;
 }
 
+//=======================================================================
+//function : ExtrusParam::makeNodesByNormal1D
+//purpose  : create nodes for extrusion using normals of edges
+//=======================================================================
+
+int SMESH_MeshEditor::ExtrusParam::
+makeNodesByNormal1D( SMESHDS_Mesh*                     mesh,
+                     const SMDS_MeshNode*              srcNode,
+                     std::list<const SMDS_MeshNode*> & newNodes,
+                     const bool                        makeMediumNodes)
+{
+  throw SALOME_Exception("Extrusion 1D by Normal not implemented");
+  return 0;
+}
 
 //=======================================================================
 //function : ExtrusionSweep
@@ -5193,20 +5572,11 @@ SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
                                   const gp_Vec&        theStep,
                                   const int            theNbSteps,
                                   TTElemOfElemListMap& newElemsMap,
-                                  const bool           theMakeGroups,
                                   const int            theFlags,
                                   const double         theTolerance)
 {
-  ExtrusParam aParams;
-  aParams.myDir = gp_Dir(theStep);
-  aParams.myNodes.Clear();
-  aParams.mySteps = new TColStd_HSequenceOfReal;
-  int i;
-  for(i=1; i<=theNbSteps; i++)
-    aParams.mySteps->Append(theStep.Magnitude());
-
-  return
-    ExtrusionSweep(theElems,aParams,newElemsMap,theMakeGroups,theFlags,theTolerance);
+  ExtrusParam aParams( theStep, theNbSteps, theFlags, theTolerance );
+  return ExtrusionSweep( theElems, aParams, newElemsMap );
 }
 
 
@@ -5218,10 +5588,7 @@ SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
 SMESH_MeshEditor::PGroupIDs
 SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
                                   ExtrusParam&         theParams,
-                                  TTElemOfElemListMap& newElemsMap,
-                                  const bool           theMakeGroups,
-                                  const int            theFlags,
-                                  const double         theTolerance)
+                                  TTElemOfElemListMap& newElemsMap)
 {
   myLastCreatedElems.Clear();
   myLastCreatedNodes.Clear();
@@ -5231,7 +5598,8 @@ SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
 
   SMESHDS_Mesh* aMesh = GetMeshDS();
 
-  int nbsteps = theParams.mySteps->Length();
+  const int nbSteps = theParams.NbSteps();
+  theParams.SetElementsToUse( theElems );
 
   TNodeOfNodeListMap mapNewNodes;
   //TNodeOfNodeVecMap mapNewNodes;
@@ -5243,14 +5611,16 @@ SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
                                      myMesh->NbVolumes(ORDER_QUADRATIC) );
   // loop on theElems
   TIDSortedElemSet::iterator itElem;
-  for ( itElem = theElems.begin(); itElem != theElems.end(); itElem++ ) {
+  for ( itElem = theElems.begin(); itElem != theElems.end(); itElem++ )
+  {
     // check element type
     const SMDS_MeshElement* elem = *itElem;
     if ( !elem  || elem->GetType() == SMDSAbs_Volume )
       continue;
 
+    const size_t nbNodes = elem->NbNodes();
     vector<TNodeOfNodeListMapItr> & newNodesItVec = mapElemNewNodes[ elem ];
-    newNodesItVec.reserve( elem->NbNodes() );
+    newNodesItVec.reserve( nbNodes );
 
     // loop on elem nodes
     SMDS_ElemIteratorPtr itN = elem->nodesIterator();
@@ -5279,55 +5649,33 @@ SMESH_MeshEditor::ExtrusionSweep (TIDSortedElemSet &   theElems,
               needMediumNodes = true;
           }
         }
-
-        double coord[] = { node->X(), node->Y(), node->Z() };
-        for ( int i = 0; i < nbsteps; i++ )
+        // create nodes for all steps
+        if ( theParams.MakeNodes( GetMeshDS(), node, listNewNodes, needMediumNodes ))
         {
-          if ( needMediumNodes ) // create a medium node
+          list<const SMDS_MeshNode*>::iterator newNodesIt = listNewNodes.begin();
+          for ( ; newNodesIt != listNewNodes.end(); ++newNodesIt )
           {
-            double x = coord[0] + theParams.myDir.X()*theParams.mySteps->Value(i+1)/2.;
-            double y = coord[1] + theParams.myDir.Y()*theParams.mySteps->Value(i+1)/2.;
-            double z = coord[2] + theParams.myDir.Z()*theParams.mySteps->Value(i+1)/2.;
-            if( theFlags & EXTRUSION_FLAG_SEW ) {
-              const SMDS_MeshNode * newNode = CreateNode(x, y, z,
-                                                         theTolerance, theParams.myNodes);
-              listNewNodes.push_back( newNode );
-            }
-            else {
-              const SMDS_MeshNode * newNode = aMesh->AddNode(x, y, z);
-              myLastCreatedNodes.Append(newNode);
-              srcNodes.Append( node );
-              listNewNodes.push_back( newNode );
-            }
-          }
-          // create a corner node
-          coord[0] = coord[0] + theParams.myDir.X()*theParams.mySteps->Value(i+1);
-          coord[1] = coord[1] + theParams.myDir.Y()*theParams.mySteps->Value(i+1);
-          coord[2] = coord[2] + theParams.myDir.Z()*theParams.mySteps->Value(i+1);
-          if( theFlags & EXTRUSION_FLAG_SEW ) {
-            const SMDS_MeshNode * newNode = CreateNode(coord[0], coord[1], coord[2],
-                                                       theTolerance, theParams.myNodes);
-            listNewNodes.push_back( newNode );
-          }
-          else {
-            const SMDS_MeshNode * newNode = aMesh->AddNode( coord[0], coord[1], coord[2] );
-            myLastCreatedNodes.Append(newNode);
+            myLastCreatedNodes.Append( *newNodesIt );
             srcNodes.Append( node );
-            listNewNodes.push_back( newNode );
           }
         }
+        else
+        {
+          break; // newNodesItVec will be shorter than nbNodes
+        }
       }
       newNodesItVec.push_back( nIt );
     }
     // make new elements
-    sweepElement( elem, newNodesItVec, newElemsMap[elem], nbsteps, srcElems );
+    if ( newNodesItVec.size() == nbNodes )
+      sweepElement( elem, newNodesItVec, newElemsMap[elem], nbSteps, srcElems );
   }
 
-  if( theFlags & EXTRUSION_FLAG_BOUNDARY ) {
-    makeWalls( mapNewNodes, newElemsMap, mapElemNewNodes, theElems, nbsteps, srcElems );
+  if ( theParams.ToMakeBoundary() ) {
+    makeWalls( mapNewNodes, newElemsMap, mapElemNewNodes, theElems, nbSteps, srcElems );
   }
   PGroupIDs newGroupIDs;
-  if ( theMakeGroups )
+  if ( theParams.ToMakeGroups() )
     newGroupIDs = generateGroups( srcNodes, srcElems, "extruded");
 
   return newGroupIDs;
@@ -6428,7 +6776,7 @@ SMESH_MeshEditor::Transform (TIDSortedElemSet & theElems,
 
   if ( ( theMakeGroups && theCopy ) ||
        ( theMakeGroups && theTargetMesh ) )
-    newGroupIDs = generateGroups( srcNodes, srcElems, groupPostfix, theTargetMesh );
+    newGroupIDs = generateGroups( srcNodes, srcElems, groupPostfix, theTargetMesh, false );
 
   return newGroupIDs;
 }
@@ -6436,9 +6784,11 @@ SMESH_MeshEditor::Transform (TIDSortedElemSet & theElems,
 //=======================================================================
 /*!
  * \brief Create groups of elements made during transformation
- * \param nodeGens - nodes making corresponding myLastCreatedNodes
- * \param elemGens - elements making corresponding myLastCreatedElems
- * \param postfix - to append to names of new groups
+ *  \param nodeGens - nodes making corresponding myLastCreatedNodes
+ *  \param elemGens - elements making corresponding myLastCreatedElems
+ *  \param postfix - to append to names of new groups
+ *  \param targetMesh - mesh to create groups in
+ *  \param topPresent - is there "top" elements that are created by sweeping
  */
 //=======================================================================
 
@@ -6446,14 +6796,17 @@ SMESH_MeshEditor::PGroupIDs
 SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
                                  const SMESH_SequenceOfElemPtr& elemGens,
                                  const std::string&             postfix,
-                                 SMESH_Mesh*                    targetMesh)
+                                 SMESH_Mesh*                    targetMesh,
+                                 const bool                     topPresent)
 {
   PGroupIDs newGroupIDs( new list<int> );
   SMESH_Mesh* mesh = targetMesh ? targetMesh : GetMesh();
 
   // Sort existing groups by types and collect their names
 
-  // to store an old group and a generated new ones
+  // containers to store an old group and generated new ones;
+  // 1st new group is for result elems of different type than a source one;
+  // 2nd new group is for same type result elems ("top" group at extrusion)
   using boost::tuple;
   using boost::make_tuple;
   typedef tuple< SMESHDS_GroupBase*, SMESHDS_Group*, SMESHDS_Group* > TOldNewGroup;
@@ -6483,6 +6836,7 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
 
   // Loop on nodes and elements to add them in new groups
 
+  vector< const SMDS_MeshElement* > resultElems;
   for ( int isNodes = 0; isNodes < 2; ++isNodes )
   {
     const SMESH_SequenceOfElemPtr& gens  = isNodes ? nodeGens : elemGens;
@@ -6505,7 +6859,7 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
         continue;
       }
       // collect all elements made by the iElem-th sourceElem
-      list< const SMDS_MeshElement* > resultElems;
+      resultElems.clear();
       if ( const SMDS_MeshElement* resElem = elems( iElem ))
         if ( resElem != sourceElem )
           resultElems.push_back( resElem );
@@ -6514,25 +6868,23 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
           if ( resElem != sourceElem )
             resultElems.push_back( resElem );
 
-      // there must be a top element
       const SMDS_MeshElement* topElem = 0;
-      if ( isNodes )
+      if ( isNodes ) // there must be a top element
       {
         topElem = resultElems.back();
         resultElems.pop_back();
       }
       else
       {
-        list< const SMDS_MeshElement* >::reverse_iterator resElemIt = resultElems.rbegin();
+        vector< const SMDS_MeshElement* >::reverse_iterator resElemIt = resultElems.rbegin();
         for ( ; resElemIt != resultElems.rend() ; ++resElemIt )
           if ( (*resElemIt)->GetType() == sourceElem->GetType() )
           {
             topElem = *resElemIt;
-            resultElems.erase( --(resElemIt.base()) ); // erase *resElemIt
+            *resElemIt = 0; // erase *resElemIt
             break;
           }
       }
-
       // add resultElems to groups originted from ones the sourceElem belongs to
       list< TOldNewGroup >::iterator gOldNew, gLast = groupsOldNew.end();
       for ( gOldNew = groupsOldNew.begin(); gOldNew != gLast; ++gOldNew )
@@ -6542,16 +6894,17 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
         {
           // fill in a new group
           SMDS_MeshGroup & newGroup = gOldNew->get<1>()->SMDSGroup();
-          list< const SMDS_MeshElement* >::iterator resLast = resultElems.end(), resElemIt;
+          vector< const SMDS_MeshElement* >::iterator resLast = resultElems.end(), resElemIt;
           for ( resElemIt = resultElems.begin(); resElemIt != resLast; ++resElemIt )
-            newGroup.Add( *resElemIt );
+            if ( *resElemIt )
+              newGroup.Add( *resElemIt );
 
           // fill a "top" group
           if ( topElem )
           {
             SMDS_MeshGroup & newTopGroup = gOldNew->get<2>()->SMDSGroup();
             newTopGroup.Add( topElem );
-          }
+         }
         }
       }
     } // loop on created elements
@@ -6565,7 +6918,6 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
     SMESHDS_GroupBase* oldGroupDS =   orderedOldNewGroups[i]->get<0>();
     SMESHDS_Group*   newGroups[2] = { orderedOldNewGroups[i]->get<1>(),
                                       orderedOldNewGroups[i]->get<2>() };
-    const int nbNewGroups = !newGroups[0]->IsEmpty() + !newGroups[1]->IsEmpty();
     for ( int is2nd = 0; is2nd < 2; ++is2nd )
     {
       SMESHDS_Group* newGroupDS = newGroups[ is2nd ];
@@ -6579,11 +6931,21 @@ SMESH_MeshEditor::generateGroups(const SMESH_SequenceOfElemPtr& nodeGens,
         newGroupDS->SetType( newGroupDS->GetElements()->next()->GetType() );
 
         // make a name
-        const bool isTop = ( nbNewGroups == 2 &&
+        const bool isTop = ( topPresent &&
                              newGroupDS->GetType() == oldGroupDS->GetType() &&
                              is2nd );
 
         string name = oldGroupDS->GetStoreName();
+        { // remove trailing whitespaces (issue 22599)
+          size_t size = name.size();
+          while ( size > 1 && isspace( name[ size-1 ]))
+            --size;
+          if ( size != name.size() )
+          {
+            name.resize( size );
+            oldGroupDS->SetStoreName( name.c_str() );
+          }
+        }
         if ( !targetMesh ) {
           string suffix = ( isTop ? "top": postfix.c_str() );
           name += "_";
@@ -8751,6 +9113,7 @@ void SMESH_MeshEditor::ConvertToQuadratic(const bool theForce3d, const bool theT
   aHelper.SetIsQuadratic( true );
   aHelper.SetIsBiQuadratic( theToBiQuad );
   aHelper.SetElementsOnShape(true);
+  aHelper.ToFixNodeParameters( true );
 
   // convert elements assigned to sub-meshes
   int nbCheckedElems = 0;
@@ -9975,7 +10338,7 @@ SMESH_MeshEditor::FindMatchingNodes(set<const SMDS_MeshElement*>& theSide1,
 
 void SMESH_MeshEditor::DoubleElements( const TIDSortedElemSet& theElements )
 {
-  CrearLastCreated();
+  ClearLastCreated();
   SMESHDS_Mesh* mesh = GetMeshDS();
 
   // get an element type and an iterator over elements
@@ -10134,6 +10497,7 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh*     theMeshDS,
       {
         // duplicate node
         aNewNode = theMeshDS->AddNode( aCurrNode->X(), aCurrNode->Y(), aCurrNode->Z() );
+        copyPosition( aCurrNode, aNewNode );
         theNodeNodeMap[ aCurrNode ] = aNewNode;
         myLastCreatedNodes.Append( aNewNode );
       }
@@ -10146,10 +10510,8 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh*     theMeshDS,
     if ( theIsDoubleElem )
       AddElement(newNodes, anElem->GetType(), anElem->IsPoly());
     else
-      {
-      MESSAGE("ChangeElementNodes");
       theMeshDS->ChangeElementNodes( anElem, &newNodes[ 0 ], anElem->NbNodes() );
-      }
+
     res = true;
   }
   return res;
@@ -10160,8 +10522,8 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh*     theMeshDS,
   \brief Creates a hole in a mesh by doubling the nodes of some particular elements
   \param theNodes - identifiers of nodes to be doubled
   \param theModifiedElems - identifiers of elements to be updated by the new (doubled)
-         nodes. If list of element identifiers is empty then nodes are doubled but
-         they not assigned to elements
+  nodes. If list of element identifiers is empty then nodes are doubled but
+  they not assigned to elements
   \return TRUE if operation has been completed successfully, FALSE otherwise
 */
 //================================================================================
@@ -10197,6 +10559,7 @@ bool SMESH_MeshEditor::DoubleNodes( const std::list< int >& theListOfNodes,
     const SMDS_MeshNode* aNewNode = aMeshDS->AddNode( aNode->X(), aNode->Y(), aNode->Z() );
     if ( aNewNode )
     {
+      copyPosition( aNode, aNewNode );
       anOldNodeToNewNode[ aNode ] = aNewNode;
       myLastCreatedNodes.Append( aNewNode );
     }
@@ -10298,15 +10661,12 @@ namespace {
     }
     void Perform(const gp_Pnt& aPnt, double theTol)
     {
+      theTol *= theTol;
       _state = TopAbs_OUT;
       _extremum.Perform(aPnt);
       if ( _extremum.IsDone() )
         for ( int iSol = 1; iSol <= _extremum.NbExt() && _state == TopAbs_OUT; ++iSol)
-#if OCC_VERSION_LARGE > 0x06040000 // Porting to OCCT6.5.1
           _state = ( _extremum.SquareDistance(iSol) <= theTol ? TopAbs_IN : TopAbs_OUT );
-#else
-          _state = ( _extremum.Value(iSol) <= theTol ? TopAbs_IN : TopAbs_OUT );
-#endif
     }
     TopAbs_State State() const
     {
@@ -10818,6 +11178,7 @@ bool SMESH_MeshEditor::DoubleNodesOnGroupBoundaries( const std::vector<TIDSorted
                         }
                       double *coords = grid->GetPoint(oldId);
                       SMDS_MeshNode *newNode = meshDS->AddNode(coords[0], coords[1], coords[2]);
+                      copyPosition( meshDS->FindNodeVtk( oldId ), newNode );
                       int newId = newNode->getVtkId();
                       nodeDomains[oldId][idom] = newId; // cloned node for other domains
                       //MESSAGE("-+-+-c     oldNode " << oldId << " domain " << idomain << " newNode " << newId << " domain " << idom << " size=" <<nodeDomains[oldId].size());
@@ -11156,6 +11517,14 @@ bool SMESH_MeshEditor::DoubleNodesOnGroupBoundaries( const std::vector<TIDSorted
         }
     }
 
+  // Remove empty groups (issue 0022812)
+  std::map<std::string, SMESH_Group*>::iterator name_group = mapOfJunctionGroups.begin();
+  for ( ; name_group != mapOfJunctionGroups.end(); ++name_group )
+  {
+    if ( name_group->second && name_group->second->GetGroupDS()->IsEmpty() )
+      myMesh->RemoveGroup( name_group->second->GetGroupDS()->GetID() );
+  }
+
   meshDS->CleanDownWardConnectivity(); // Mesh has been modified, downward connectivity is no more usable, free memory
   grid->BuildLinks();
 
@@ -11223,6 +11592,7 @@ bool SMESH_MeshEditor::CreateFlatElementsOnFacesGroups(const std::vector<TIDSort
               if (!clonedNodes.count(node))
                 {
                   clone = meshDS->AddNode(node->X(), node->Y(), node->Z());
+                  copyPosition( node, clone );
                   clonedNodes[node] = clone;
                 }
               else
@@ -11239,6 +11609,7 @@ bool SMESH_MeshEditor::CreateFlatElementsOnFacesGroups(const std::vector<TIDSort
                   if (!intermediateNodes.count(node))
                     {
                       inter = meshDS->AddNode(node->X(), node->Y(), node->Z());
+                      copyPosition( node, inter );
                       intermediateNodes[node] = inter;
                     }
                   else
@@ -12183,3 +12554,45 @@ int SMESH_MeshEditor::MakeBoundaryMesh(const TIDSortedElemSet& elements,
   }
   return nbAddedBnd;
 }
+
+//================================================================================
+/*!
+ * \brief Copy node position and set \a to node on the same geometry
+ */
+//================================================================================
+
+void SMESH_MeshEditor::copyPosition( const SMDS_MeshNode* from,
+                                     const SMDS_MeshNode* to )
+{
+  if ( !from || !to ) return;
+
+  SMDS_PositionPtr pos = from->GetPosition();
+  if ( !pos || from->getshapeId() < 1 ) return;
+
+  switch ( pos->GetTypeOfPosition() )
+  {
+  case SMDS_TOP_3DSPACE: break;
+
+  case SMDS_TOP_FACE:
+  {
+    const SMDS_FacePosition* fPos = static_cast< const SMDS_FacePosition* >( pos );
+    GetMeshDS()->SetNodeOnFace( to, from->getshapeId(),
+                                fPos->GetUParameter(), fPos->GetVParameter() );
+    break;
+  }
+  case SMDS_TOP_EDGE:
+  {
+    // WARNING: it is dangerous to set equal nodes on one EDGE!!!!!!!!
+    const SMDS_EdgePosition* ePos = static_cast< const SMDS_EdgePosition* >( pos );
+    GetMeshDS()->SetNodeOnEdge( to, from->getshapeId(), ePos->GetUParameter() );
+    break;
+  }
+  case SMDS_TOP_VERTEX:
+  {
+    GetMeshDS()->SetNodeOnVertex( to, from->getshapeId() );
+    break;
+  }
+  case SMDS_TOP_UNSPEC:
+  default:;
+  }
+}