Salome HOME
Copyrights update 2015.
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
index 4e5e0478af102aa743419d1065ede20a14e6f9e1..e474f5bf7f671fbbd5dcd8087bf8b6bdc2f52bff 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2013  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
@@ -6,7 +6,7 @@
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
 // License as published by the Free Software Foundation; either
-// version 2.1 of the License.
+// version 2.1 of the License, or (at your option) any later version.
 //
 // This library is distributed in the hope that it will be useful,
 // but WITHOUT ANY WARRANTY; without even the implied warranty of
@@ -43,6 +43,7 @@
 #include <BRep_Tool.hxx>
 #include <GCPnts_AbscissaPoint.hxx>
 #include <Geom2dAdaptor_Curve.hxx>
+#include <Geom_Line.hxx>
 #include <TopExp.hxx>
 #include <TopExp_Explorer.hxx>
 #include <TopoDS.hxx>
@@ -151,21 +152,26 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
       myMissingVertexNodes = true;
 
     // check if the edge has a non-uniform parametrization (issue 0020705)
-    if ( !myC2d[i].IsNull() && myEdgeLength[i] > DBL_MIN)
+    if ( !myC2d[i].IsNull() )
     {
-      Geom2dAdaptor_Curve A2dC( myC2d[i],
-                                std::min( myFirst[i], myLast[i] ),
-                                std::max( myFirst[i], myLast[i] ));
-      double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2., p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
-      double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
-      double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
-      //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
-      myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
-      //if ( !myIsUniform[i] ) to implement Value3d(u)
+      if ( myEdgeLength[i] > DBL_MIN)
       {
-        double fp,lp;
-        Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
-        myC3dAdaptor[i].Load( C3d, fp,lp );
+        Geom2dAdaptor_Curve A2dC( myC2d[i],
+                                  std::min( myFirst[i], myLast[i] ),
+                                  std::max( myFirst[i], myLast[i] ));
+        double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2., p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
+        double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
+        double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
+        //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
+        myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
+        Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],d2,d4);
+        myC3dAdaptor[i].Load( C3d, d2,d4 );
+      }
+      else
+      {
+        const TopoDS_Vertex& V = TopoDS::Vertex( vExp.Value() );
+        Handle(Geom_Curve) C3d = new Geom_Line( BRep_Tool::Pnt( V ), gp::DX() );
+        myC3dAdaptor[i].Load( C3d, 0, 0.5 * BRep_Tool::Tolerance( V ));
       }
     }
     // reverse a proxy submesh
@@ -251,7 +257,8 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const StdMeshers_FaceSide*  theSide,
  */
 //================================================================================
 
-StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec& theSideNodes)
+StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
+                                         const TopoDS_Face& theFace)
 {
   myEdge.resize( 1 );
   myEdgeID.resize( 1, -1 );
@@ -272,12 +279,42 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec& theSideNodes)
   if ( !myPoints.empty() )
   {
     myPoints[0].normParam = 0;
-    gp_Pnt pPrev = SMESH_TNodeXYZ( myPoints[0].node );
-    for ( size_t i = 1; i < myPoints.size(); ++i )
+    if ( myPoints[0].node &&
+         myPoints.back().node &&
+         myPoints[ myNbPonits/2 ].node )
     {
-      gp_Pnt p = SMESH_TNodeXYZ( myPoints[i].node );
-      myLength += ( myPoints[i].normParam = p.Distance( pPrev ));
-      pPrev = p;
+      gp_Pnt pPrev = SMESH_TNodeXYZ( myPoints[0].node );
+      for ( size_t i = 1; i < myPoints.size(); ++i )
+      {
+        gp_Pnt p = SMESH_TNodeXYZ( myPoints[i].node );
+        myLength += p.Distance( pPrev );
+        myPoints[i].normParam = myLength;
+        pPrev = p;
+      }
+    }
+    else if ( !theFace.IsNull() )
+    {
+      TopLoc_Location loc;
+      Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace, loc );
+      gp_Pnt pPrev = surf->Value( myPoints[0].u, myPoints[0].v );
+      for ( size_t i = 1; i < myPoints.size(); ++i )
+      {
+        gp_Pnt p = surf->Value( myPoints[i].u, myPoints[i].v );
+        myLength += p.Distance( pPrev );
+        myPoints[i].normParam = myLength;
+        pPrev = p;
+      }
+    }
+    else
+    {
+      gp_Pnt2d pPrev = myPoints[0].UV();
+      for ( size_t i = 1; i < myPoints.size(); ++i )
+      {
+        gp_Pnt2d p = myPoints[i].UV();
+        myLength += p.Distance( pPrev );
+        myPoints[i].normParam = myLength;
+        pPrev = p;
+      }
     }
     if ( myLength > std::numeric_limits<double>::min() )
       for ( size_t i = 1; i < myPoints.size(); ++i )
@@ -694,6 +731,7 @@ void StdMeshers_FaceSide::Reverse()
   }
   for ( size_t i = 0; i < myEdge.size(); ++i )
   {
+    if ( myEdge[i].IsNull() ) continue; // for a side on points only
     double fp,lp;
     Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
     if ( !C3d.IsNull() )
@@ -926,6 +964,18 @@ gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
     return myC2d[ i ]->Value(par);
 
   }
+  else if ( !myPoints.empty() )
+  {
+    int i = U * double( myPoints.size()-1 );
+    while ( i > 0 && myPoints[ i ].normParam > U )
+      --i;
+    while ( i+1 < myPoints.size() && myPoints[ i+1 ].normParam < U )
+      ++i;
+    double r = (( U - myPoints[ i ].normParam ) /
+                ( myPoints[ i+1 ].normParam - myPoints[ i ].normParam ));
+    return ( myPoints[ i   ].UV() * ( 1 - r ) +
+             myPoints[ i+1 ].UV() * r );
+  }
   return myDefaultPnt2d;
 }
 
@@ -944,7 +994,7 @@ gp_Pnt StdMeshers_FaceSide::Value3d(double U) const
   double     r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
 
   double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
-    
+
   // check parametrization of curve
   if( !myIsUniform[i] )
   {
@@ -968,7 +1018,8 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
                                               SMESH_Mesh &         theMesh,
                                               const bool           theIgnoreMediumNodes,
                                               TError &             theError,
-                                              SMESH_ProxyMesh::Ptr theProxyMesh)
+                                              SMESH_ProxyMesh::Ptr theProxyMesh,
+                                              const bool           theCheckVertexNodes)
 {
   list< TopoDS_Edge > edges, internalEdges;
   list< int > nbEdgesInWires;
@@ -993,17 +1044,18 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
     // as StdMeshers_FaceSide::GetUVPtStruct() requires
     if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
     {
-      while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
-                                       theMesh.GetMeshDS()))
-      {
-        wireEdges.splice(wireEdges.end(), wireEdges,
-                         wireEdges.begin(), ++wireEdges.begin());
-        if ( from->IsSame( wireEdges.front() )) {
-          theError = TError
-            ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
-          return TSideVector(0);
+      if ( theCheckVertexNodes )
+        while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
+                                         theMesh.GetMeshDS()))
+        {
+          wireEdges.splice(wireEdges.end(), wireEdges,
+                           wireEdges.begin(), ++wireEdges.begin());
+          if ( from->IsSame( wireEdges.front() )) {
+            theError = TError
+              ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
+            return TSideVector(0);
+          }
         }
-      }
     }
     else if ( *nbE > 1 ) // Issue 0020676 (Face_pb_netgen.brep) - several internal edges in a wire
     {