Salome HOME
0021490: EDF 2114: RadialQuadrangle fails
[modules/smesh.git] / src / StdMeshers / StdMeshers_RadialQuadrangle_1D2D.cxx
index feb4e6d255c797808607f2dc3a5b847d7d7ac525..c1c846646db260d53e136b34b52bef5f7a9f3cfa 100644 (file)
 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
 //
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
-//
 
 //  SMESH SMESH : implementaion of SMESH idl descriptions
 // File      : StdMeshers_RadialQuadrangle_1D2D.cxx
 // Module    : SMESH
 // Created   : Fri Oct 20 11:37:07 2006
 // Author    : Edward AGAPOV (eap)
-//
+
 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
 
 #include "StdMeshers_NumberOfLayers.hxx"
@@ -147,7 +146,8 @@ namespace
    */
   class TEdgeMarker : public SMESH_subMeshEventListener
   {
-    TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false) {}
+    TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false,
+                                              "StdMeshers_RadialQuadrangle_1D2D::TEdgeMarker") {}
   public:
     //!<  Return static listener
     static SMESH_subMeshEventListener* getListener()
@@ -427,7 +427,6 @@ void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMes
 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
                                                const TopoDS_Shape& aShape)
 {
-  TopExp_Explorer exp;
   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
 
   myHelper = new SMESH_MesherHelper( aMesh );
@@ -479,8 +478,8 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
       for(; itn != theNodes.end(); itn++ ) {
         CNodes.push_back( (*itn).second );
         double ang = (*itn).first - fang;
-        if( ang>PI ) ang = ang - 2*PI;
-        if( ang<-PI ) ang = ang + 2*PI;
+        if( ang>M_PI ) ang = ang - 2.*M_PI;
+        if( ang<-M_PI ) ang = ang + 2.*M_PI;
         Angles.Append( ang ); 
       }
     }
@@ -490,8 +489,7 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
     if ( !computeLayerPositions(P0,P1))
       return false;
 
-    exp.Init( CircEdge, TopAbs_VERTEX );
-    TopoDS_Vertex V1 = TopoDS::Vertex( exp.Current() );
+    TopoDS_Vertex V1 = myHelper->IthVertex(0, CircEdge );
     gp_Pnt2d p2dV = BRep_Tool::Parameters( V1, TopoDS::Face(aShape) );
 
     NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
@@ -528,7 +526,7 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
     // a segment of line
     double fp, lp;
     Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
-    if( fabs(fabs(lp-fp)-PI) > Precision::Confusion() ) {
+    if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
       // not half of circle
       return error(COMPERR_BAD_SHAPE);
     }
@@ -552,8 +550,8 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
     for(; itn != theNodes.end(); itn++ ) {
       CNodes.push_back( (*itn).second );
       double ang = (*itn).first - fang;
-      if( ang>PI ) ang = ang - 2*PI;
-      if( ang<-PI ) ang = ang + 2*PI;
+      if( ang>M_PI ) ang = ang - 2.*M_PI;
+      if( ang<-M_PI ) ang = ang + 2.*M_PI;
       Angles.Append( ang );
     }
     const SMDS_MeshNode* NF = theNodes.begin()->second;
@@ -690,14 +688,23 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
     for(; itn != theNodes.end(); itn++ ) {
       CNodes.push_back( (*itn).second );
       double ang = (*itn).first - fang;
-      if( ang>PI ) ang = ang - 2*PI;
-      if( ang<-PI ) ang = ang + 2*PI;
+      if( ang>M_PI ) ang = ang - 2.*M_PI;
+      if( ang<-M_PI ) ang = ang + 2.*M_PI;
       Angles.Append( ang );
     }
     P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
     gp_Pnt P2( NL->X(), NL->Y(), NL->Z() );
     P0 = aCirc->Location();
 
+    // make P1 belong to LinEdge1
+    TopoDS_Vertex V1 = myHelper->IthVertex( 0, LinEdge1 );
+    TopoDS_Vertex V2 = myHelper->IthVertex( 1, LinEdge1 );
+    gp_Pnt PE1 = BRep_Tool::Pnt(V1);
+    gp_Pnt PE2 = BRep_Tool::Pnt(V2);
+    if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
+        ( P1.Distance(PE2) > Precision::Confusion() ) )
+      std::swap( LinEdge1, LinEdge2 );
+
     bool linEdge1Computed, linEdge2Computed;
     if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdge1Computed))
       return false;
@@ -706,23 +713,12 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
     Nodes2.resize( myLayerPositions.size()+1 );
 
     // check that both linear edges have same hypotheses
-    if ( !computeLayerPositions(P0,P1,LinEdge2, &linEdge2Computed))
+    if ( !computeLayerPositions(P0,P2,LinEdge2, &linEdge2Computed))
          return false;
     if ( Nodes1.size() != myLayerPositions.size()+1 )
       return error("Different hypotheses apply to radial edges");
 
-    exp.Init( LinEdge1, TopAbs_VERTEX );
-    TopoDS_Vertex V1 = TopoDS::Vertex( exp.Current() );
-    exp.Next();
-    TopoDS_Vertex V2 = TopoDS::Vertex( exp.Current() );
-    gp_Pnt PE1 = BRep_Tool::Pnt(V1);
-    gp_Pnt PE2 = BRep_Tool::Pnt(V2);
-    if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
-        ( P1.Distance(PE2) > Precision::Confusion() ) )
-    {
-      std::swap( LinEdge1, LinEdge2 );
-      std::swap( linEdge1Computed, linEdge2Computed );
-    }
+    // find the central vertex
     TopoDS_Vertex VC = V2;
     if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
         ( P2.Distance(PE1) > Precision::Confusion() ) )
@@ -1148,7 +1144,7 @@ bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
     // a segment of line
     double fp, lp;
     Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
-    if( fabs(fabs(lp-fp)-PI) > Precision::Confusion() ) {
+    if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
       // not half of circle
       return error(COMPERR_BAD_SHAPE);
     }
@@ -1286,4 +1282,3 @@ bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
   return false;
 
 }
-