]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Regression of smesh/mesh_Projection_2D_01/B7
authoreap <eap@opencascade.com>
Wed, 13 May 2015 13:03:36 +0000 (16:03 +0300)
committereap <eap@opencascade.com>
Wed, 13 May 2015 13:03:36 +0000 (16:03 +0300)
   Fix SMESH_MeshEditor::Smooth()

Optimize a bit SMDS_MeshNode::GetInverseElementIterator()

Minimize usage of BRepClass_FaceClassifier in StdMeshers_Import_1D2D::Compute()

src/SMDS/SMDS_MeshNode.cxx
src/SMESH/SMESH_MeshEditor.cxx
src/StdMeshers/StdMeshers_Import_1D2D.cxx

index dd7bf6a423d1fa22038dfb91031b927992543276..3bd528bb766e2a5087596d7a1a5808e8605da430 100644 (file)
@@ -126,7 +126,7 @@ void SMDS_MeshNode::SetPosition(const SMDS_PositionPtr& aPos)
 
 const SMDS_PositionPtr& SMDS_MeshNode::GetPosition() const
 {
-        return myPosition;
+  return myPosition;
 }
 
 //=======================================================================
@@ -149,30 +149,26 @@ public:
   SMDS_MeshNode_MyInvIterator(SMDS_Mesh *mesh, vtkIdType* cells, int ncells, SMDSAbs_ElementType type) :
     myMesh(mesh), myCells(cells), myNcells(ncells), myType(type), iter(0)
   {
-    //MESSAGE("SMDS_MeshNode_MyInvIterator : ncells " << myNcells);
-    cellList.clear();
+    cellList.reserve( ncells );
     if (type == SMDSAbs_All)
+      cellList.assign( cells, cells + ncells );
+    else
       for (int i = 0; i < ncells; i++)
-        cellList.push_back(cells[i]);
-    else for (int i = 0; i < ncells; i++)
       {
         int vtkId = cells[i];
         int smdsId = myMesh->fromVtkToSmds(vtkId);
         const SMDS_MeshElement* elem = myMesh->FindElement(smdsId);
         if (elem->GetType() == type)
-          {
-            //MESSAGE("Add element vtkId " << vtkId << " " << elem->GetType())
-            cellList.push_back(vtkId);
-          }
+        {
+          cellList.push_back(vtkId);
+        }
       }
     myCells = cellList.empty() ? 0 : &cellList[0];
     myNcells = cellList.size();
-    //MESSAGE("myNcells="<<myNcells);
   }
 
   bool more()
   {
-    //MESSAGE("iter " << iter << " ncells " << myNcells);
     return (iter < myNcells);
   }
 
index 602ade3a3d2959cb633167ee3044b8f4e896adc9..b45e0117bb235707b4499f60e482f6619c487794 100644 (file)
@@ -3768,8 +3768,11 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
   // smooth elements on each TopoDS_Face separately
   // ===============================================
 
-  set< int >::reverse_iterator fId = faceIdSet.rbegin(); // treate 0 fId at the end
-  for ( ; fId != faceIdSet.rend(); ++fId ) {
+  SMESH_MesherHelper helper( *GetMesh() );
+
+  set< int >::reverse_iterator fId = faceIdSet.rbegin(); // treat 0 fId at the end
+  for ( ; fId != faceIdSet.rend(); ++fId )
+  {
     // get face surface and submesh
     Handle(Geom_Surface) surface;
     SMESHDS_SubMesh* faceSubMesh = 0;
@@ -3777,7 +3780,8 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
     double fToler2 = 0, f,l;
     double u1 = 0, u2 = 0, v1 = 0, v2 = 0;
     bool isUPeriodic = false, isVPeriodic = false;
-    if ( *fId ) {
+    if ( *fId )
+    {
       face = TopoDS::Face( aMesh->IndexToShape( *fId ));
       surface = BRep_Tool::Surface( face );
       faceSubMesh = aMesh->MeshElements( *fId );
@@ -3790,6 +3794,7 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
       if ( isVPeriodic )
         surface->VPeriod();
       surface->Bounds( u1, u2, v1, v2 );
+      helper.SetSubShape( face );
     }
     // ---------------------------------------------------------
     // for elements on a face, find movable and fixed nodes and
@@ -3811,7 +3816,8 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
     int nbElemOnFace = 0;
     itElem = theElems.begin();
     // loop on not yet smoothed elements: look for elems on a face
-    while ( itElem != theElems.end() ) {
+    while ( itElem != theElems.end() )
+    {
       if ( faceSubMesh && nbElemOnFace == faceSubMesh->NbElements() )
         break; // all elements found
 
@@ -3867,12 +3873,15 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
 
       // get nodes to check UV
       list< const SMDS_MeshNode* > uvCheckNodes;
+      const SMDS_MeshNode* nodeInFace = 0;
       itN = elem->nodesIterator();
       nn = 0; nbn =  elem->NbNodes();
       if(elem->IsQuadratic())
         nbn = nbn/2;
       while ( nn++ < nbn ) {
         node = static_cast<const SMDS_MeshNode*>( itN->next() );
+        if ( node->GetPosition()->GetDim() == 2 )
+          nodeInFace = node;
         if ( uvMap.find( node ) == uvMap.end() )
           uvCheckNodes.push_back( node );
         // add nodes of elems sharing node
@@ -3898,41 +3907,21 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
         const SMDS_PositionPtr& pos = node->GetPosition();
         posType = pos ? pos->GetTypeOfPosition() : SMDS_TOP_3DSPACE;
         // get existing UV
-        switch ( posType ) {
-        case SMDS_TOP_FACE: {
-          SMDS_FacePosition* fPos = ( SMDS_FacePosition* ) pos;
-          uv.SetCoord( fPos->GetUParameter(), fPos->GetVParameter() );
-          break;
-        }
-        case SMDS_TOP_EDGE: {
-          TopoDS_Shape S = aMesh->IndexToShape( node->getshapeId() );
-          Handle(Geom2d_Curve) pcurve;
-          if ( !S.IsNull() && S.ShapeType() == TopAbs_EDGE )
-            pcurve = BRep_Tool::CurveOnSurface( TopoDS::Edge( S ), face, f,l );
-          if ( !pcurve.IsNull() ) {
-            double u = (( SMDS_EdgePosition* ) pos )->GetUParameter();
-            uv = pcurve->Value( u ).XY();
-          }
-          break;
-        }
-        case SMDS_TOP_VERTEX: {
-          TopoDS_Shape S = aMesh->IndexToShape( node->getshapeId() );
-          if ( !S.IsNull() && S.ShapeType() == TopAbs_VERTEX )
-            uv = BRep_Tool::Parameters( TopoDS::Vertex( S ), face ).XY();
-          break;
-        }
-        default:;
-        }
-        // check existing UV
-        bool project = true;
-        gp_Pnt pNode ( node->X(), node->Y(), node->Z() );
-        double dist1 = DBL_MAX, dist2 = 0;
-        if ( posType != SMDS_TOP_3DSPACE ) {
-          dist1 = pNode.SquareDistance( surface->Value( uv.X(), uv.Y() ));
-          project = dist1 > fToler2;
-        }
+        if ( pos )
+        {
+          bool toCheck = true;
+          uv = helper.GetNodeUV( face, node, nodeInFace, &toCheck );
+        }
+        // compute not existing UV
+        bool project = ( posType == SMDS_TOP_3DSPACE );
+        // double dist1 = DBL_MAX, dist2 = 0;
+        // if ( posType != SMDS_TOP_3DSPACE ) {
+        //   dist1 = pNode.SquareDistance( surface->Value( uv.X(), uv.Y() ));
+        //   project = dist1 > fToler2;
+        // }
         if ( project ) { // compute new UV
           gp_XY newUV;
+          gp_Pnt pNode = SMESH_TNodeXYZ( node );
           if ( !getClosestUV( projector, pNode, newUV )) {
             MESSAGE("Node Projection Failed " << node);
           }
@@ -3942,9 +3931,9 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
             if ( isVPeriodic )
               newUV.SetY( ElCLib::InPeriod( newUV.Y(), v1, v2 ));
             // check new UV
-            if ( posType != SMDS_TOP_3DSPACE )
-              dist2 = pNode.SquareDistance( surface->Value( newUV.X(), newUV.Y() ));
-            if ( dist2 < dist1 )
+            // if ( posType != SMDS_TOP_3DSPACE )
+            //   dist2 = pNode.SquareDistance( surface->Value( newUV.X(), newUV.Y() ));
+            // if ( dist2 < dist1 )
               uv = newUV;
           }
         }
@@ -4012,9 +4001,8 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
         uv2 = pcurve->Value( f );
         int iPar = Abs( uv1.X() - uv2.X() ) > Abs( uv1.Y() - uv2.Y() ) ? 1 : 2;
         // assure uv1 < uv2
-        if ( uv1.Coord( iPar ) > uv2.Coord( iPar )) {
-          gp_Pnt2d tmp = uv1; uv1 = uv2; uv2 = tmp;
-        }
+        if ( uv1.Coord( iPar ) > uv2.Coord( iPar ))
+          std::swap( uv1, uv2 );
         // get nodes on seam and its vertices
         list< const SMDS_MeshNode* > seamNodes;
         SMDS_NodeIteratorPtr nSeamIt = sm->GetNodes();
@@ -4064,12 +4052,14 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
                   setMovableNodes.find( n ) == setMovableNodes.end() )
                 continue;
               // add only nodes being closer to uv2 than to uv1
-              gp_Pnt pMid (0.5 * ( n->X() + nSeam->X() ),
-                           0.5 * ( n->Y() + nSeam->Y() ),
-                           0.5 * ( n->Z() + nSeam->Z() ));
-              gp_XY uv;
-              getClosestUV( projector, pMid, uv );
-              if ( uv.Coord( iPar ) > uvMap[ n ]->Coord( iPar ) ) {
+              // gp_Pnt pMid (0.5 * ( n->X() + nSeam->X() ),
+              //              0.5 * ( n->Y() + nSeam->Y() ),
+              //              0.5 * ( n->Z() + nSeam->Z() ));
+              // gp_XY uv;
+              // getClosestUV( projector, pMid, uv );
+              double x = uvMap[ n ]->Coord( iPar );
+              if ( Abs( uv1.Coord( iPar ) - x ) >
+                   Abs( uv2.Coord( iPar ) - x )) {
                 nodesNearSeam.insert( n );
                 nbUseMap2++;
               }
@@ -4178,8 +4168,6 @@ void SMESH_MeshEditor::Smooth (TIDSortedElemSet &          theElems,
     // move medium nodes of quadratic elements
     if ( isQuadratic )
     {
-      SMESH_MesherHelper helper( *GetMesh() );
-      helper.SetSubShape( face );
       vector<const SMDS_MeshNode*> nodes;
       bool checkUV;
       list< const SMDS_MeshElement* >::iterator elemIt = elemsOnFace.begin();
index 7d3097997dd2a04d230a566dbf3af73fe41f677e..196814b789a73d927d0a3d711c8ca531abd95af2 100644 (file)
@@ -192,11 +192,13 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
   const int shapeID = tgtMesh->ShapeToIndex( geomFace );
   const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
 
+
   Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
   const bool reverse =
     ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace ) == TopAbs_REVERSED );
   gp_Pnt p; gp_Vec du, dv;
 
+  // BRepClass_FaceClassifier is most time consuming, so minimize its usage
   BRepClass_FaceClassifier classifier;
   Bnd_B2d bndBox2d;
   Bnd_Box bndBox3d;
@@ -217,7 +219,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
     bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
 
     BRepBndLib::Add( geomFace, bndBox3d );
-    bndBox3d.Enlarge( 1e-4 * sqrt( bndBox3d.SquareExtent() ));
+    bndBox3d.Enlarge( 1e-5 * sqrt( bndBox3d.SquareExtent() ));
   }
 
   set<int> subShapeIDs;
@@ -277,6 +279,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
   vector<TopAbs_State>         nodeState;
   vector<const SMDS_MeshNode*> newNodes; // of a face
   set   <const SMDS_MeshNode*> bndNodes; // nodes classified ON
+  vector<bool>                 isNodeIn; // nodes classified IN, by node ID
 
   for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
   {
@@ -330,6 +333,9 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
         {
           if ( !subShapeIDs.count( newNode->getshapeId() ))
             break; // node is Imported onto other FACE
+          if ( newNode->GetID() < (int) isNodeIn.size() &&
+               isNodeIn[ newNode->GetID() ])
+            isIn = true;
           if ( !isIn && bndNodes.count( *node ))
             nodeState[ i ] = TopAbs_ON;
         }
@@ -362,10 +368,15 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
             newNode = tgtMesh->AddNode( nXYZ.X(), nXYZ.Y(), nXYZ.Z());
             tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
             nbCreatedNodes++;
+            if ( newNode->GetID() >= (int) isNodeIn.size() )
+            {
+              isNodeIn.push_back( false ); // allow allocate more than newNode->GetID()
+              isNodeIn.resize( newNode->GetID() + 1, false );
+            }
             if ( nodeState[i] == TopAbs_ON )
               bndNodes.insert( *node );
             else
-              isIn = true;
+              isNodeIn[ newNode->GetID() ] = isIn = true;
           }
         }
         if ( !(newNodes[i] = newNode ) || isOut )