Salome HOME
merge from branch V6_main tag mergeto_BR_SMDS_MEMIMP_29nov10 BR_SMDS_MEMIMP
authorprascle <prascle>
Mon, 29 Nov 2010 13:13:50 +0000 (13:13 +0000)
committerprascle <prascle>
Mon, 29 Nov 2010 13:13:50 +0000 (13:13 +0000)
doc/docutils/conf.py
src/SMESH/SMESH_MesherHelper.cxx
src/SMESH/SMESH_MesherHelper.hxx

index 84a74a5704b8517f6cef804b703a7f0ef27b2d13..aebd4b16314bb44ce830687b6d6a02d8d72cf36b 100644 (file)
@@ -54,9 +54,9 @@ copyright = '2010 EDF R&D'
 # built documents.
 #
 # The short X.Y version.
-version = '5.1.4'
+version = '6.2.0'
 # The full version, including alpha/beta/rc tags.
-release = '5.1.4'
+release = '6.2.0'
 
 # The language for content autogenerated by Sphinx. Refer to documentation
 # for a list of supported languages.
index 78c47ae1277a4f043fe6552f22d82676ebe8d9a4..bcdd2096a0e5f22528ef9bec48357005c07124e7 100644 (file)
@@ -315,6 +315,31 @@ void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
   myTLinkNodeMap.insert( make_pair(link,n12));
 }
 
+//================================================================================
+/*!
+ * \brief Return true if position of nodes on the shape hasn't yet been checked or
+ * the positions proved to be invalid
+ */
+//================================================================================
+
+bool SMESH_MesherHelper::toCheckPosOnShape(int shapeID ) const
+{
+  map< int,bool >::const_iterator id_ok = myNodePosShapesValidity.find( shapeID );
+  return ( id_ok == myNodePosShapesValidity.end() || !id_ok->second );
+}
+
+//================================================================================
+/*!
+ * \brief Set validity of positions of nodes on the shape.
+ * Once set, validity is not changed
+ */
+//================================================================================
+
+void SMESH_MesherHelper::setPosOnShapeValidity(int shapeID, bool ok ) const
+{
+  ((SMESH_MesherHelper*)this)->myNodePosShapesValidity.insert( make_pair( shapeID, ok));
+}
+
 //=======================================================================
 //function : GetUVOnSeam
 //purpose  : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
@@ -477,7 +502,8 @@ bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
                                      const double         tol,
                                      const bool           force) const
 {
-  if ( force || !myOkNodePosShapes.count( n->getshapeId() ))
+  int shapeID = n->getshapeId();
+  if ( force || toCheckPosOnShape( shapeID ))
   {
     double toldis = tol;
     double tolmin = 1.e-7*myMesh->GetMeshDS()->getMaxDim(); // nodes coordinates are stored in float format
@@ -491,6 +517,7 @@ bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
          Precision::IsInfinite( uv.Y() ) ||
          nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > toldis )
     {
+      setPosOnShapeValidity( shapeID, false );
       // uv incorrect, project the node to surface
       GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
       projector.Perform( nodePnt );
@@ -507,10 +534,14 @@ bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
         MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
         return false;
       }
+      // store the fixed UV on the face
+      if ( myShape.IsSame(F) && shapeID == myShapeID )
+        const_cast<SMDS_MeshNode*>(n)->SetPosition
+          ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
     }
     else if ( uv.Modulus() > numeric_limits<double>::min() )
     {
-      ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->getshapeId() );
+      setPosOnShapeValidity( shapeID, true );
     }
   }
   return true;
@@ -660,7 +691,8 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
                                     const bool           force,
                                     double*              distance) const
 {
-  if ( force || !myOkNodePosShapes.count( n->getshapeId() ))
+  int shapeID = n->getshapeId();
+  if ( force || toCheckPosOnShape( shapeID ))
   {
     double toldis = tol;
     double tolmin = 1.e-7*myMesh->GetMeshDS()->getMaxDim(); // nodes coordinates are stored in float format
@@ -684,6 +716,7 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
       if ( distance ) *distance = dist;
       if ( dist > toldis )
       {
+        setPosOnShapeValidity( shapeID, false );
         // u incorrect, project the node to the curve
         int edgeID = GetMeshDS()->ShapeToIndex( E );
         TID2ProjectorOnCurve& i2proj = const_cast< TID2ProjectorOnCurve&>( myEdge2Projector );
@@ -711,11 +744,14 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
           MESSAGE("distance " << nodePnt.Distance(curve->Value( U )) << " " << toldis);
           return false;
         }
-        //u = double( U );
+        // store the fixed U on the edge
+        if ( myShape.IsSame(E) && shapeID == myShapeID )
+          const_cast<SMDS_MeshNode*>(n)->SetPosition
+            ( SMDS_PositionPtr( new SMDS_EdgePosition( U )));
       }
       else if ( fabs( u ) > numeric_limits<double>::min() )
       {
-        ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->getshapeId() );
+        setPosOnShapeValidity( shapeID, true );
       }
       if (( u < f-tol || u > l+tol ) && force )
       {
@@ -758,6 +794,10 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
   SMDS_MeshNode* n12;
   SMESHDS_Mesh* meshDS = GetMeshDS();
 
+  if ( IsSeamShape( n1->getshapeId() ))
+    // to get a correct UV of a node on seam, the second node must have checked UV
+    std::swap( n1, n2 );
+
   // get type of shape for the new medium node
   int faceID = -1, edgeID = -1;
   const SMDS_PositionPtr Pos1 = n1->GetPosition();
@@ -1676,9 +1716,9 @@ namespace { // Structures used by FixQuadraticElements()
               node2()->GetID() < other.node2()->GetID() :
               node1()->GetID() < other.node1()->GetID());
     }
-    struct PtrComparator {
-      bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
-    };
+//     struct PtrComparator {
+//       bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
+//     };
   };
   // ---------------------------------------------------------
   /*!
@@ -1725,7 +1765,7 @@ namespace { // Structures used by FixQuadraticElements()
   /*!
    * \brief Face shared by two volumes and bound by QLinks
    */
-  struct QFace: public TIDSortedElemSet
+  struct QFace: public TIDSortedNodeSet
   {
     mutable const SMDS_MeshElement* _volumes[2];
     mutable vector< const QLink* >  _sides;
@@ -1764,6 +1804,8 @@ namespace { // Structures used by FixQuadraticElements()
 
     bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
 
+    bool IsSpoiled(const QLink* bentLink ) const;
+
     TLinkInSet GetBoundaryLink( const TLinkSet&      links,
                                 const TChainLink&    avoidLink,
                                 TLinkInSet *         notBoundaryLink = 0,
@@ -1808,7 +1850,7 @@ namespace { // Structures used by FixQuadraticElements()
   ostream& operator << (ostream& out, const QFace& f)
   {
     out <<"QFace nodes: "/*<< &f << "  "*/;
-    for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
+    for ( TIDSortedNodeSet::const_iterator n = f.begin(); n != f.end(); ++n )
       out << (*n)->GetID() << " ";
     out << " \tvolumes: "
         << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
@@ -2150,6 +2192,37 @@ namespace { // Structures used by FixQuadraticElements()
     return fullLen;
   }
 
+  //================================================================================
+  /*!
+   * \brief Checks if the face is distorted due to bentLink
+   */
+  //================================================================================
+
+  bool QFace::IsSpoiled(const QLink* bentLink ) const
+  {
+    // code is valid for convex faces only
+    gp_XYZ gc(0,0,0);
+    for ( TIDSortedNodeSet::const_iterator n = begin(); n!=end(); ++n)
+      gc += XYZ( *n ) / size();
+    for (unsigned i = 0; i < _sides.size(); ++i )
+    {
+      if ( _sides[i] == bentLink ) continue;
+      gp_Vec linkNorm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
+      gp_Vec vecOut( gc, _sides[i]->MiddlePnt() );
+      if ( linkNorm * vecOut < 0 )
+        linkNorm.Reverse();
+      double mag2 = linkNorm.SquareMagnitude();
+      if ( mag2 > numeric_limits<double>::min() )
+        linkNorm /= sqrt( mag2 );
+      gp_Vec vecBent    ( _sides[i]->MiddlePnt(), bentLink->MediumPnt());
+      gp_Vec vecStraight( _sides[i]->MiddlePnt(), bentLink->MiddlePnt());
+      if ( vecBent * linkNorm > -0.1*vecStraight.Magnitude() )
+        return true;
+    }
+    return false;
+    
+  }
+
   //================================================================================
   /*!
    * \brief Find pairs of continues faces 
@@ -2295,41 +2368,20 @@ namespace { // Structures used by FixQuadraticElements()
     TLinkSet linkSet( allLinks.begin(), allLinks.end());
     TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
 
-    // move in 2d if we are on geom face
-//     TopoDS_Face face;
-//     TopLoc_Location loc;
-//     SMESH_MesherHelper faceHelper( *helper.GetMesh());
-//     while ( linkIt->IsBoundary()) ++linkIt;
-//     if ( linkIt == linksEnd ) return;
-//     if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
-//       bool checkPos = true;
-//       TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
-//       if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
-//         face = TopoDS::Face( f );
-//         helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
-//         if (checkPos)
-//           face.Nullify();
-//         else
-//           faceHelper.SetSubShape( face );
-//       }
-//     }
     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
     {
       if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
       {
-//         if ( !face.IsNull() ) {
-//           const SMDS_MeshNode* inFaceNode =
-//             faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
-//           gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
-//           gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
-//           gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
-//           gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
-//           gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
-//           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
-//         }
-//         else {
-          linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
-          //}
+        // move iff a boundary link is bent towards inside of a face (issue 0021084)
+        const QFace* face = linkIt->_qfaces[0];
+        gp_XYZ pIn = ( face->_sides[0]->MiddlePnt() +
+                       face->_sides[1]->MiddlePnt() +
+                       face->_sides[2]->MiddlePnt() ) / 3.;
+        gp_XYZ insideDir( pIn - (*linkIt)->MiddlePnt());
+        bool linkBentInside = ((*linkIt)->_nodeMove.Dot( insideDir ) > 0 );
+        //if ( face->IsSpoiled( linkIt->_qlink ))
+        if ( linkBentInside )
+          face->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
       }
     }
   }
@@ -2600,7 +2652,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
 
   bool isCurved = false;
   //bool hasRectFaces = false;
-  set<int> nbElemNodeSet;
+  //set<int> nbElemNodeSet;
 
   if ( elemType == SMDSAbs_Volume )
   {
@@ -2613,7 +2665,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
       for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
       {
         int nbN = volTool.NbFaceNodes( iF );
-        nbElemNodeSet.insert( nbN );
+        //nbElemNodeSet.insert( nbN );
         const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
         vector< const QLink* > faceLinks( nbN/2 );
         for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
@@ -2655,7 +2707,7 @@ void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
       const SMDS_MeshElement* face = elemIt->next();
       if ( !face->IsQuadratic() )
         continue;
-      nbElemNodeSet.insert( face->NbNodes() );
+      //nbElemNodeSet.insert( face->NbNodes() );
       int nbN = face->NbNodes()/2;
       vector< const QLink* > faceLinks( nbN );
       for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
index 2f89ce7899239f2abc9617f95c7e8c864f8032b5..d904f95ded3491c0ec6ef67954b8cadf812ed355 100644 (file)
@@ -501,7 +501,10 @@ protected:
   // to create quadratic elements
   bool            myCreateQuadratic;
   bool            mySetElemOnShape;
-  std::set< int > myOkNodePosShapes;
+
+  std::map< int,bool > myNodePosShapesValidity;
+  bool toCheckPosOnShape(int shapeID ) const;
+  void setPosOnShapeValidity(int shapeID, bool ok ) const;
 
 };