X-Git-Url: http://git.salome-platform.org/gitweb/?p=modules%2Fsmesh.git;a=blobdiff_plain;f=src%2FStdMeshers%2FStdMeshers_FaceSide.cxx;h=b543928dcabd653c2e7adfda7f4f26fca61b11c5;hp=4e5e0478af102aa743419d1065ede20a14e6f9e1;hb=d9faba6c847c1c1a4d4f501ca5ac5725a25a8236;hpb=06236fdb5a448ef7546f484b59fffcf22878d7df diff --git a/src/StdMeshers/StdMeshers_FaceSide.cxx b/src/StdMeshers/StdMeshers_FaceSide.cxx index 4e5e0478a..b543928dc 100644 --- a/src/StdMeshers/StdMeshers_FaceSide.cxx +++ b/src/StdMeshers/StdMeshers_FaceSide.cxx @@ -1,4 +1,4 @@ -// Copyright (C) 2007-2013 CEA/DEN, EDF R&D, OPEN CASCADE +// Copyright (C) 2007-2014 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 @@ -251,7 +251,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 +273,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::min() ) for ( size_t i = 1; i < myPoints.size(); ++i ) @@ -694,6 +725,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 +958,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; }