Salome HOME
SALOME Forum bug: http://www.salome-platform.org/forum/forum_10/967838025
authoreap <eap@opencascade.com>
Fri, 23 Aug 2013 07:44:07 +0000 (07:44 +0000)
committereap <eap@opencascade.com>
Fri, 23 Aug 2013 07:44:07 +0000 (07:44 +0000)
Don't use BRep_Tool::Degenerated() which sometimes gives a wrong answer,
use SMESH_Algo::isDegenerated() instead

src/SMESH/SMESH_Algo.cxx
src/SMESH/SMESH_Algo.hxx
src/SMESH/SMESH_MeshEditor.cxx
src/SMESH/SMESH_MesherHelper.cxx
src/SMESH/SMESH_subMesh.cxx
src/StdMeshers/StdMeshers_Prism_3D.cxx
src/StdMeshers/StdMeshers_Quadrangle_2D.cxx
src/StdMeshers/StdMeshers_ViscousLayers.cxx

index f08542dc0f73c1ce3e9c48387d055c853a36fdbf..518c0d6688f29623cd9cc014ddf4a346e8fe7358 100644 (file)
@@ -310,10 +310,10 @@ SMESH_Algo::GetAppliedHypothesis(SMESH_Mesh &         aMesh,
 double SMESH_Algo::EdgeLength(const TopoDS_Edge & E)
 {
   double UMin = 0, UMax = 0;
-  if (BRep_Tool::Degenerated(E))
-    return 0;
   TopLoc_Location L;
   Handle(Geom_Curve) C = BRep_Tool::Curve(E, L, UMin, UMax);
+  if ( C.IsNull() )
+    return 0.;
   GeomAdaptor_Curve AdaptCurve(C, UMin, UMax); //range is important for periodic curves
   double length = GCPnts_AbscissaPoint::Length(AdaptCurve, UMin, UMax);
   return length;
@@ -580,6 +580,20 @@ bool SMESH_Algo::isStraight( const TopoDS_Edge & E,
   return true;
 }
 
+//================================================================================
+/*!
+ * \brief Return true if an edge has no 3D curve
+ */
+//================================================================================
+
+bool SMESH_Algo::isDegenerated( const TopoDS_Edge & E )
+{
+  double f,l;
+  TopLoc_Location loc;
+  Handle(Geom_Curve) C = BRep_Tool::Curve( E, loc, f,l );
+  return C.IsNull();
+}
+
 //================================================================================
 /*!
  * \brief Return the node built on a vertex
index 4d22e655dc432cfa12bc636e1f309734cd2ccd7b..0208934b862ade5b5c9d3c50e2a8d94435cea35c 100644 (file)
@@ -361,6 +361,10 @@ public:
    * \brief Return true if an edge can be considered straight
    */
   static bool isStraight( const TopoDS_Edge & E, const bool degenResult=false );
+  /*!
+   * \brief Return true if an edge has no 3D curve
+   */
+  static bool isDegenerated( const TopoDS_Edge & E );
 
   /*!
    * \brief Return the node built on a vertex
index 9078ce9f0a56f9cafd8d224763a6f79b4f6d85b6..00521f94cbc208e971718ac26932319bfe764c0a 100644 (file)
@@ -4976,7 +4976,7 @@ SMESH_MeshEditor::ExtrusionAlongTrack (TIDSortedElemSet &   theElements,
   if( aS.ShapeType() == TopAbs_EDGE ) {
     aTrackEdge = TopoDS::Edge( aS );
     // the Edge must not be degenerated
-    if ( BRep_Tool::Degenerated( aTrackEdge ) )
+    if ( SMESH_Algo::isDegenerated( aTrackEdge ) )
       return EXTR_BAD_PATH_SHAPE;
     TopExp::Vertices( aTrackEdge, aV1, aV2 );
     aItN = theTrack->GetFather()->GetSubMesh( aV1 )->GetSubMeshDS()->GetNodes();
@@ -5264,7 +5264,7 @@ SMESH_MeshEditor::ExtrusionAlongTrack (TIDSortedElemSet &   theElements,
   else if( aS.ShapeType() == TopAbs_EDGE ) {
     aTrackEdge = TopoDS::Edge( aS );
     // the Edge must not be degenerated
-    if ( BRep_Tool::Degenerated( aTrackEdge ) )
+    if ( SMESH_Algo::isDegenerated( aTrackEdge ) )
       return EXTR_BAD_PATH_SHAPE;
     TopExp::Vertices( aTrackEdge, aV1, aV2 );
     const SMDS_MeshNode* aN1 = SMESH_Algo::VertexNode( aV1, pMeshDS );
@@ -5290,7 +5290,7 @@ SMESH_MeshEditor::ExtrusionAlongTrack (TIDSortedElemSet &   theElements,
     TopExp_Explorer eExp(aS, TopAbs_EDGE);
     for(; eExp.More(); eExp.Next()) {
       TopoDS_Edge E = TopoDS::Edge( eExp.Current() );
-      if( BRep_Tool::Degenerated(E) ) continue;
+      if( SMESH_Algo::isDegenerated(E) ) continue;
       SMESH_subMesh* SM = theTrack->GetSubMesh(E);
       if(SM) {
         LSM.push_back(SM);
index ea08fda34fe3c0a2efab714cd04d9ead8e1e5972..ad330f22ae71d6c4414b5a11dda53e9c6b886345 100644 (file)
@@ -281,7 +281,7 @@ void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
         }
 
         // look for a degenerated edge
-        if ( BRep_Tool::Degenerated( edge )) {
+        if ( SMESH_Algo::isDegenerated( edge )) {
           myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
           for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
             myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
@@ -3825,7 +3825,7 @@ namespace { // Structures used by FixQuadraticElements()
       {
         // check if the EDGE needs checking
         const TopoDS_Edge& edge = TopoDS::Edge( edgeIt.Current() );
-        if ( BRep_Tool::Degenerated( edge ) )
+        if ( SMESH_Algo::isDegenerated( edge ) )
           continue;
         if ( theHelper.IsRealSeam( edge ) &&
              edge.Orientation() == TopAbs_REVERSED )
index cf05b50bc4defa46c123074ca2ac5310520ad652..7f96a76fdb0505717158b8b207d4f0c7e13f8058 100644 (file)
@@ -1965,7 +1965,7 @@ bool SMESH_subMesh::checkComputeError(SMESH_Algo*         theAlgo,
       if ( _computeState != COMPUTE_OK )
       {
         if ( _subShape.ShapeType() == TopAbs_EDGE &&
-             BRep_Tool::Degenerated( TopoDS::Edge( _subShape )) )
+             SMESH_Algo::isDegenerated( TopoDS::Edge( _subShape )) )
           _computeState = COMPUTE_OK;
         else if ( theComputeOK )
           _computeError = SMESH_ComputeError::New(COMPERR_NO_MESH_ON_SHAPE,"",theAlgo);
index 6ba56c39c926bf4fd20095540b8f854342e29363..22398fd87a94005a82806cd9d450d68a15f696eb 100644 (file)
@@ -2088,7 +2088,7 @@ bool StdMeshers_Prism_3D::initPrism(Prism_3D::TPrismTopo& thePrism,
     if ( botSM ) {
       if ( ! botSM->GetSubShape().IsSame( thePrism.myBottom )) {
         std::swap( botSM, topSM );
-        if ( ! botSM->GetSubShape().IsSame( thePrism.myBottom ))
+        if ( !botSM || ! botSM->GetSubShape().IsSame( thePrism.myBottom ))
           return toSM( error( COMPERR_BAD_INPUT_MESH,
                               "Incompatible non-structured sub-meshes"));
       }
@@ -2571,7 +2571,7 @@ bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> &
     list< TopoDS_Edge >::const_iterator edgeIt = prism.myBottomEdges.begin();
     for ( int iE = 0; iE < prism.myNbEdgesInWires.front(); ++iE, ++edgeIt )
     {
-      if ( BRep_Tool::Degenerated( *edgeIt )) continue;
+      if ( SMESH_Algo::isDegenerated( *edgeIt )) continue;
       const TParam2ColumnMap* u2colMap =
         GetParam2ColumnMap( MeshDS()->ShapeToIndex( *edgeIt ), isReverse );
       if ( !u2colMap ) return false;
index ef7e6100b4c06a87a2e0f30062a111b934d4425b..de0523143bd3ebe30604009bf7d70f669d16bf2f 100644 (file)
@@ -901,7 +901,7 @@ FaceQuadStruct::Ptr StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh &
             sideEdges.splice(sideEdges.begin(), edges, --edges.end());
         }
       }
-      if ( sideEdges.size() == 1 && BRep_Tool::Degenerated( sideEdges.front() ))
+      if ( sideEdges.size() == 1 && SMESH_Algo::isDegenerated( sideEdges.front() ))
         degenSides.push_back( nbSides );
 
       quad->side.push_back(new StdMeshers_FaceSide(F, sideEdges, &aMesh, nbSides < QUAD_TOP_SIDE,
index 5ec3c8fa19de78041c42cc09802e152f5a0cc28e..7c9588671f657f7ccd4e9a03b2d3e25387e720eb 100644 (file)
@@ -786,7 +786,7 @@ namespace
     for ( ; eExp.More(); eExp.Next() )
     {
       const TopoDS_Edge& E = TopoDS::Edge( eExp.Current() );
-      if ( BRep_Tool::Degenerated( E )) continue;
+      if ( SMESH_Algo::isDegenerated( E )) continue;
       // check if 2D curve is concave
       BRepAdaptor_Curve2d curve( E, F );
       const int nbIntervals = curve.NbIntervals( GeomAbs_C2 );