1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-08-04 13:13:25 +03:00

0025416: Wrong section curve

1. Restriction line is processed in IntTools_FaceFace with using methods of GeomInt_IntSS class.
2. Check, if Restriction- and Walking-lines (or Restriction-Restriction lines) are coincided, has been added in IntPatch_ImpPrmIntersection.cxx (at that RLine is considered to be isoline only).
3. Check, if RLine and GLine are coincided, has been added in IntPatch_ImpImpIntersection.cxx.
4. Create new class IntPatch_PointLine, which is inherited from IntPatch_Line.
5. The reason of exception (in DEBUG MODE) has been eliminated.

New test cases for issue #25416 were added.

tests/bugs/modalg_5/bug24650 was modified.
This commit is contained in:
nbv
2015-03-25 09:52:07 +03:00
committed by apn
parent 2a8523acca
commit d4b867e617
29 changed files with 1526 additions and 418 deletions

View File

@@ -1 +1,2 @@
GeomInt_IntSS_1.cxx
GeomInt_VectorOfReal.hxx

View File

@@ -35,11 +35,13 @@ uses StdFail,
IntSurf,
IntPatch,
ApproxInt
ApproxInt,
Bnd
is
class IntSS;
imported VectorOfReal;
class LineConstructor;

View File

@@ -30,8 +30,11 @@ uses Intersection from IntPatch,
Surface from Geom,
HSurface from GeomAdaptor,
TopolTool from Adaptor3d,
Line from IntPatch
Line from IntPatch,
RLine from IntPatch,
Box2d from Bnd,
VectorOfReal from GeomInt
raises NotDone from StdFail,
OutOfRange from Standard
@@ -191,6 +194,27 @@ is
aTolCheck:out Real from Standard;
aTolAngCheck:out Real from Standard);
TreatRLine(myclass; theRL: RLine from IntPatch;
theHS1, theHS2: HSurface from GeomAdaptor;
theC3d: out Curve from Geom;
theC2d1, theC2d2: out Curve from Geom2d;
theTolReached: out Real from Standard);
---Purpose: converts RLine to Geom(2d)_Curve.
BuildPCurves(myclass; f, l : Real from Standard;
Tol: out Real from Standard;
S: Surface from Geom;
C: Curve from Geom;
C2d: out Curve from Geom2d);
---Purpose: creates 2D-curve on given surface from given 3D-curve
TrimILineOnSurfBoundaries(myclass; theC2d1, theC2d2: Curve from Geom2d;
theBound1, theBound2: Box2d from Bnd;
theArrayOfParameters: out VectorOfReal from GeomInt);
---Purpose: puts into theArrayOfParameters the parameters of intersection
-- points of given theC2d1 and theC2d2 curves with the boundaries
-- of the source surface.
--- Private methods
--
--

View File

@@ -14,8 +14,12 @@
// Alternatively, this file may be used under the terms of Open CASCADE
// commercial license or contractual agreement.
#include <algorithm>
#include <GeomInt_IntSS.ixx>
#include <GeomInt.hxx>
#include <Precision.hxx>
#include <gp_Pnt.hxx>
#include <gp_Pnt2d.hxx>
@@ -75,13 +79,22 @@
#include <GeomInt_WLApprox.hxx>
#include <Extrema_ExtPS.hxx>
static
void BuildPCurves (Standard_Real f,
Standard_Real l,
Standard_Real& Tol,
const Handle (Geom_Surface)& S,
const Handle (Geom_Curve)& C,
Handle (Geom2d_Curve)& C2d);
#include <GeomAdaptor_HSurface.hxx>
#include <gp_Lin2d.hxx>
#include <Geom2d_Line.hxx>
#include <IntPatch_RLine.hxx>
#include <Geom2dAdaptor.hxx>
#include <Adaptor2d_HCurve2d.hxx>
#include <Approx_CurveOnSurface.hxx>
#include <GeomAdaptor.hxx>
#include <GeomProjLib.hxx>
#include <GeomAPI_ProjectPointOnSurf.hxx>
#include <TColStd_Array1OfReal.hxx>
#include <TColgp_Array1OfPnt2d.hxx>
#include <TColStd_Array1OfInteger.hxx>
#include <Geom2d_BSplineCurve.hxx>
#include <Geom2dAPI_InterCurveCurve.hxx>
#include <IntRes2d_IntersectionPoint.hxx>
static
void Parameters(const Handle(GeomAdaptor_HSurface)& HS1,
@@ -145,20 +158,6 @@ static
void GetQuadric(const Handle(GeomAdaptor_HSurface)& HS1,
IntSurf_Quadric& quad1);
static
void TreatRLine(const Handle(IntPatch_RLine)& theRL,
const Handle(GeomAdaptor_HSurface)& theHS1,
const Handle(GeomAdaptor_HSurface)& theHS2,
const Standard_Boolean theApprox1,
const Standard_Boolean theApprox2,
Handle(Geom_Curve)& theC3d,
Handle(Geom2d_Curve)& theC2d1,
Handle(Geom2d_Curve)& theC2d2,
Standard_Real& theTolReached);
//
class ProjectPointOnSurf
{
public:
@@ -251,13 +250,13 @@ Standard_Real ProjectPointOnSurf::LowerDistance() const
//function : MakeCurve
//purpose :
//=======================================================================
void GeomInt_IntSS::MakeCurve(const Standard_Integer Index,
const Handle(Adaptor3d_TopolTool) & dom1,
const Handle(Adaptor3d_TopolTool) & dom2,
const Standard_Real Tol,
const Standard_Boolean Approx,
const Standard_Boolean ApproxS1,
const Standard_Boolean ApproxS2)
void GeomInt_IntSS::MakeCurve(const Standard_Integer Index,
const Handle(Adaptor3d_TopolTool) & dom1,
const Handle(Adaptor3d_TopolTool) & dom2,
const Standard_Real Tol,
const Standard_Boolean Approx,
const Standard_Boolean ApproxS1,
const Standard_Boolean ApproxS2)
{
Standard_Boolean myApprox1, myApprox2, myApprox;
@@ -1036,96 +1035,92 @@ Standard_Real ProjectPointOnSurf::LowerDistance() const
default: break;
}
if(!isAnalS1 || !isAnalS2)
Handle(IntPatch_RLine) RL =
Handle(IntPatch_RLine)::DownCast(L);
Handle(Geom_Curve) aC3d;
Handle(Geom2d_Curve) aC2d1, aC2d2;
Standard_Real aTolReached;
TreatRLine(RL, myHS1, myHS2, aC3d,
aC2d1, aC2d2, aTolReached);
if(aC3d.IsNull())
break;
Bnd_Box2d aBox1, aBox2;
const Standard_Real aU1f = myHS1->FirstUParameter(),
aV1f = myHS1->FirstVParameter(),
aU1l = myHS1->LastUParameter(),
aV1l = myHS1->LastVParameter();
const Standard_Real aU2f = myHS2->FirstUParameter(),
aV2f = myHS2->FirstVParameter(),
aU2l = myHS2->LastUParameter(),
aV2l = myHS2->LastVParameter();
aBox1.Add(gp_Pnt2d(aU1f, aV1f));
aBox1.Add(gp_Pnt2d(aU1l, aV1l));
aBox2.Add(gp_Pnt2d(aU2f, aV2f));
aBox2.Add(gp_Pnt2d(aU2l, aV2l));
GeomInt_VectorOfReal anArrayOfParameters;
//We consider here that the intersection line is same-parameter-line
anArrayOfParameters.Append(aC3d->FirstParameter());
anArrayOfParameters.Append(aC3d->LastParameter());
TrimILineOnSurfBoundaries(aC2d1, aC2d2, aBox1, aBox2, anArrayOfParameters);
const Standard_Integer aNbIntersSolutionsm1 = anArrayOfParameters.Length() - 1;
//Trim RLine found.
for(Standard_Integer anInd = 0; anInd < aNbIntersSolutionsm1; anInd++)
{
Handle(IntPatch_RLine) RL =
Handle(IntPatch_RLine)::DownCast(L);
Handle(Geom_Curve) aC3d;
Handle(Geom2d_Curve) aC2d1, aC2d2;
Standard_Real aTolReached;
TreatRLine(RL, myHS1, myHS2, myApprox1, myApprox2, aC3d,
aC2d1, aC2d2, aTolReached);
sline.Append(aC3d);
slineS1.Append(aC2d1);
slineS2.Append(aC2d2);
const Standard_Real aParF = anArrayOfParameters(anInd),
aParL = anArrayOfParameters(anInd+1);
if((aParL - aParF) <= Precision::PConfusion())
continue;
const Standard_Real aPar = 0.5*(aParF + aParL);
gp_Pnt2d aPt;
Handle(Geom2d_Curve) aCurv2d1, aCurv2d2;
if(!aC2d1.IsNull())
{
aC2d1->D0(aPar, aPt);
if(aBox1.IsOut(aPt))
continue;
if(myApprox1)
aCurv2d1 = new Geom2d_TrimmedCurve(aC2d1, aParF, aParL);
}
if(!aC2d2.IsNull())
{
aC2d2->D0(aPar, aPt);
if(aBox2.IsOut(aPt))
continue;
if(myApprox2)
aCurv2d2 = new Geom2d_TrimmedCurve(aC2d2, aParF, aParL);
}
Handle(Geom_Curve) aCurv3d = new Geom_TrimmedCurve(aC3d, aParF, aParL);
sline.Append(aCurv3d);
slineS1.Append(aCurv2d1);
slineS2.Append(aCurv2d2);
}
}
break;
}
}
//
//=======================================================================
//function : AdjustUPeriodic
//purpose :
//=======================================================================
void TreatRLine(const Handle(IntPatch_RLine)& theRL,
const Handle(GeomAdaptor_HSurface)& theHS1,
const Handle(GeomAdaptor_HSurface)& theHS2,
const Standard_Boolean theApprox1,
const Standard_Boolean theApprox2,
Handle(Geom_Curve)& theC3d,
Handle(Geom2d_Curve)& theC2d1,
Handle(Geom2d_Curve)& theC2d2,
Standard_Real& theTolReached)
{
Handle(GeomAdaptor_HSurface) aGAHS;
Handle(Adaptor2d_HCurve2d) anAHC2d;
Standard_Real tf, tl;
gp_Lin2d aL;
// It is assumed that 2d curve is 2d line (rectangular surface domain)
if(theRL->IsArcOnS1())
{
aGAHS = theHS1;
anAHC2d = theRL->ArcOnS1();
theRL->ParamOnS1(tf, tl);
theC2d1 = Geom2dAdaptor::MakeCurve(theRL->ArcOnS1()->Curve2d());
}
else if (theRL->IsArcOnS2())
{
aGAHS = theHS2;
anAHC2d = theRL->ArcOnS2();
theRL->ParamOnS2(tf, tl);
theC2d2 = Geom2dAdaptor::MakeCurve(theRL->ArcOnS2()->Curve2d());
}
else
{
return;
}
//
//To provide sameparameter it is necessary to get 3d curve as
//approximation of curve on surface.
Standard_Integer aMaxDeg = 8;
Standard_Integer aMaxSeg = 1000;
Approx_CurveOnSurface anApp(anAHC2d, aGAHS, tf, tl, Precision::Confusion(),
GeomAbs_C1, aMaxDeg, aMaxSeg,
Standard_True, Standard_False);
if(anApp.HasResult())
{
theC3d = anApp.Curve3d();
theTolReached = anApp.MaxError3d();
Standard_Real aTol = Precision::Confusion();
if(theRL->IsArcOnS1() && theApprox2)
{
Handle(Geom_Surface) aS = GeomAdaptor::MakeSurface(theHS1->Surface());
BuildPCurves (tf, tl, aTol,
aS, theC3d, theC2d2);
}
if(theRL->IsArcOnS2() && theApprox1)
{
Handle(Geom_Surface) aS = GeomAdaptor::MakeSurface(theHS2->Surface());
BuildPCurves (tf, tl, aTol,
aS, theC3d, theC2d1);
}
theTolReached = Max(theTolReached, aTol);
}
//
}
//
//=======================================================================
//function : AdjustUPeriodic
//purpose :
//=======================================================================
void AdjustUPeriodic (const Handle(Geom_Surface)& aS, Handle(Geom2d_Curve)& aC2D)
{
@@ -1168,53 +1163,6 @@ void TreatRLine(const Handle(IntPatch_RLine)& theRL,
//
//
//=======================================================================
//function : BuildPCurves
//purpose :
//=======================================================================
void BuildPCurves (Standard_Real f, Standard_Real l,
Standard_Real& Tol,
const Handle(Geom_Surface)& S,
const Handle(Geom_Curve)& C,
Handle(Geom2d_Curve)& C2d)
{
if (C2d.IsNull())
{
C2d = GeomProjLib::Curve2d(C,f,l,S,Tol);
if (S->IsUPeriodic())
{
Standard_Real umin,umax,vmin,vmax;
S->Bounds(umin,umax,vmin,vmax);
// Recadre dans le domaine UV de la face
const Standard_Real period = S->UPeriod();
gp_Pnt2d Pf=C2d->Value(f);
//
Standard_Real U0=Pf.X();
//
const Standard_Real aEpsilon=Epsilon(10.);//1.77e-15
if (fabs(U0)<aEpsilon)
U0=0.;
if (fabs(U0-period)<aEpsilon)
U0=period;
//
const Standard_Real aEps=Precision::PConfusion();//1.e-9
Standard_Real du=0.;
while(U0 <(umin-aEps)) {
U0+=period;
du+=period;
}
while(U0>(umax+aEps)) {
U0-=period;
du-=period;
}
if (du!=0.) {
gp_Vec2d T1(du, 0.);
C2d->Translate(T1);
}
}
}
}
//
//=======================================================================
//function : GetQuadric
//purpose :
//=======================================================================
@@ -2106,26 +2054,312 @@ Standard_Boolean FindPoint(const gp_Pnt2d& theFirstPoint,
return Standard_False;
}
/*
static
void DumpWLine (const Handle(IntPatch_WLine)& theWLine);
//=======================================================================
//function : DumpWLine
//function : TreatRLine
//purpose : Approx of Restriction line
//=======================================================================
void GeomInt_IntSS::TreatRLine(const Handle(IntPatch_RLine)& theRL,
const Handle(GeomAdaptor_HSurface)& theHS1,
const Handle(GeomAdaptor_HSurface)& theHS2,
Handle(Geom_Curve)& theC3d,
Handle(Geom2d_Curve)& theC2d1,
Handle(Geom2d_Curve)& theC2d2,
Standard_Real& theTolReached)
{
Handle(GeomAdaptor_HSurface) aGAHS;
Handle(Adaptor2d_HCurve2d) anAHC2d;
Standard_Real tf, tl;
gp_Lin2d aL;
// It is assumed that 2d curve is 2d line (rectangular surface domain)
if(theRL->IsArcOnS1())
{
aGAHS = theHS1;
anAHC2d = theRL->ArcOnS1();
theRL->ParamOnS1(tf, tl);
theC2d1 = Geom2dAdaptor::MakeCurve(anAHC2d->Curve2d());
}
else if (theRL->IsArcOnS2())
{
aGAHS = theHS2;
anAHC2d = theRL->ArcOnS2();
theRL->ParamOnS2(tf, tl);
theC2d2 = Geom2dAdaptor::MakeCurve(anAHC2d->Curve2d());
}
else
{
return;
}
//
//To provide sameparameter it is necessary to get 3d curve as
//approximation of curve on surface.
Standard_Integer aMaxDeg = 8;
Standard_Integer aMaxSeg = 1000;
Approx_CurveOnSurface anApp(anAHC2d, aGAHS, tf, tl, Precision::Confusion(),
GeomAbs_C1, aMaxDeg, aMaxSeg,
Standard_True, Standard_False);
if(!anApp.HasResult())
return;
theC3d = anApp.Curve3d();
theTolReached = anApp.MaxError3d();
Standard_Real aTol = Precision::Confusion();
if(theRL->IsArcOnS1())
{
Handle(Geom_Surface) aS = GeomAdaptor::MakeSurface(theHS2->Surface());
BuildPCurves (tf, tl, aTol,
aS, theC3d, theC2d2);
}
if(theRL->IsArcOnS2())
{
Handle(Geom_Surface) aS = GeomAdaptor::MakeSurface(theHS1->Surface());
BuildPCurves (tf, tl, aTol,
aS, theC3d, theC2d1);
}
theTolReached = Max(theTolReached, aTol);
}
//=======================================================================
//function : BuildPCurves
//purpose :
//=======================================================================
void DumpWLine (const Handle(IntPatch_WLine)& theWLine)
void GeomInt_IntSS::BuildPCurves (Standard_Real f,
Standard_Real l,
Standard_Real& Tol,
const Handle (Geom_Surface)& S,
const Handle (Geom_Curve)& C,
Handle (Geom2d_Curve)& C2d)
{
Standard_Integer i, nbp;
Standard_Real u1,v1,u2,v2;
if (!C2d.IsNull()) {
return;
}
//
nbp=theWLine->NbPnts();
printf(" nbp=%d\n", nbp);
for(i=1;i<=nbp;i++) {
const IntSurf_PntOn2S& aPoint = theWLine->Point(i);
const gp_Pnt& aP=aPoint.Value();
aPoint.Parameters(u1,v1,u2,v2);
printf("point i%d %lf %lf %lf [ %lf %lf ] [ %lf %lf ]\n",
i, aP.X(), aP.Y(), aP.Z(), u1, v1, u2, v2);
Standard_Real umin,umax,vmin,vmax;
//
S->Bounds(umin, umax, vmin, vmax);
// in class ProjLib_Function the range of parameters is shrank by 1.e-09
if((l - f) > 2.e-09) {
C2d = GeomProjLib::Curve2d(C,f,l,S,umin,umax,vmin,vmax,Tol);
//
if (C2d.IsNull()) {
// proj. a circle that goes through the pole on a sphere to the sphere
Tol += Precision::Confusion();
C2d = GeomProjLib::Curve2d(C,f,l,S,Tol);
}
}
else {
if((l - f) > Epsilon(Abs(f))) {
GeomAPI_ProjectPointOnSurf aProjector1, aProjector2;
gp_Pnt P1 = C->Value(f);
gp_Pnt P2 = C->Value(l);
aProjector1.Init(P1, S);
aProjector2.Init(P2, S);
if(aProjector1.IsDone() && aProjector2.IsDone()) {
Standard_Real U=0., V=0.;
aProjector1.LowerDistanceParameters(U, V);
gp_Pnt2d p1(U, V);
aProjector2.LowerDistanceParameters(U, V);
gp_Pnt2d p2(U, V);
if(p1.Distance(p2) > gp::Resolution()) {
TColgp_Array1OfPnt2d poles(1,2);
TColStd_Array1OfReal knots(1,2);
TColStd_Array1OfInteger mults(1,2);
poles(1) = p1;
poles(2) = p2;
knots(1) = f;
knots(2) = l;
mults(1) = mults(2) = 2;
C2d = new Geom2d_BSplineCurve(poles,knots,mults,1);
// compute reached tolerance.begin
gp_Pnt PMid = C->Value((f + l) * 0.5);
aProjector1.Perform(PMid);
if(aProjector1.IsDone()) {
aProjector1.LowerDistanceParameters(U, V);
gp_Pnt2d pmidproj(U, V);
gp_Pnt2d pmidcurve2d = C2d->Value((f + l) * 0.5);
Standard_Real adist = pmidcurve2d.Distance(pmidproj);
Tol = (adist > Tol) ? adist : Tol;
}
// compute reached tolerance.end
}
}
}
}
//
if (S->IsUPeriodic() && !C2d.IsNull()) {
// Recadre dans le domaine UV de la face
Standard_Real aTm, U0, aEps, period, du, U0x;
Standard_Boolean bAdjust;
//
aEps = Precision::PConfusion();
period = S->UPeriod();
//
aTm = .5*(f + l);
gp_Pnt2d pm = C2d->Value(aTm);
U0 = pm.X();
//
bAdjust =
GeomInt::AdjustPeriodic(U0, umin, umax, period, U0x, du, aEps);
if (bAdjust) {
gp_Vec2d T1(du, 0.);
C2d->Translate(T1);
}
}
}
*/
//=======================================================================
//function : TrimILineOnSurfBoundaries
//purpose : This function finds intersection points of given curves with
// surface boundaries and fills theArrayOfParameters by parameters
// along the given curves corresponding of these points.
//=======================================================================
void GeomInt_IntSS::TrimILineOnSurfBoundaries(const Handle(Geom2d_Curve)& theC2d1,
const Handle(Geom2d_Curve)& theC2d2,
const Bnd_Box2d& theBound1,
const Bnd_Box2d& theBound2,
GeomInt_VectorOfReal& theArrayOfParameters)
{
//Rectangular boundaries of two surfaces: [0]:U=Ufirst, [1]:U=Ulast,
// [2]:V=Vfirst, [3]:V=Vlast
const Standard_Integer aNumberOfCurves = 4;
Handle(Geom2d_Curve) aCurS1Bounds[aNumberOfCurves];
Handle(Geom2d_Curve) aCurS2Bounds[aNumberOfCurves];
Standard_Real aU1f=0.0, aV1f=0.0, aU1l=0.0, aV1l=0.0;
Standard_Real aU2f=0.0, aV2f=0.0, aU2l=0.0, aV2l=0.0;
theBound1.Get(aU1f, aV1f, aU1l, aV1l);
theBound2.Get(aU2f, aV2f, aU2l, aV2l);
Standard_Real aDelta = aV1l-aV1f;
if(Abs(aDelta) > RealSmall())
{
if(!Precision::IsInfinite(aU1f))
{
aCurS1Bounds[0] = new Geom2d_Line(gp_Pnt2d(aU1f, aV1f), gp_Dir2d(0.0, 1.0));
if(!Precision::IsInfinite(aDelta))
aCurS1Bounds[0] = new Geom2d_TrimmedCurve(aCurS1Bounds[0], 0, aDelta);
}
if(!Precision::IsInfinite(aU1l))
{
aCurS1Bounds[1] = new Geom2d_Line(gp_Pnt2d(aU1l, aV1f), gp_Dir2d(0.0, 1.0));
if(!Precision::IsInfinite(aDelta))
aCurS1Bounds[1] = new Geom2d_TrimmedCurve(aCurS1Bounds[1], 0, aDelta);
}
}
aDelta = aU1l-aU1f;
if(Abs(aDelta) > RealSmall())
{
if(!Precision::IsInfinite(aV1f))
{
aCurS1Bounds[2] = new Geom2d_Line(gp_Pnt2d(aU1f, aV1f), gp_Dir2d(1.0, 0.0));
if(!Precision::IsInfinite(aDelta))
aCurS1Bounds[2] = new Geom2d_TrimmedCurve(aCurS1Bounds[2], 0, aDelta);
}
if(!Precision::IsInfinite(aV1l))
{
aCurS1Bounds[3] = new Geom2d_Line(gp_Pnt2d(aU1l, aV1l), gp_Dir2d(1.0, 0.0));
if(!Precision::IsInfinite(aDelta))
aCurS1Bounds[3] = new Geom2d_TrimmedCurve(aCurS1Bounds[3], 0, aDelta);
}
}
aDelta = aV2l-aV2f;
if(Abs(aDelta) > RealSmall())
{
if(!Precision::IsInfinite(aU2f))
{
aCurS2Bounds[0] = new Geom2d_Line(gp_Pnt2d(aU2f, aV2f), gp_Dir2d(0.0, 1.0));
if(!Precision::IsInfinite(aDelta))
aCurS2Bounds[0] = new Geom2d_TrimmedCurve(aCurS2Bounds[0], 0, aDelta);
}
if(!Precision::IsInfinite(aU2l))
{
aCurS2Bounds[1] = new Geom2d_Line(gp_Pnt2d(aU2l, aV2f), gp_Dir2d(0.0, 1.0));
if(!Precision::IsInfinite(aDelta))
aCurS2Bounds[1] = new Geom2d_TrimmedCurve(aCurS2Bounds[1], 0, aDelta);
}
}
aDelta = aU2l-aU2f;
if(Abs(aDelta) > RealSmall())
{
if(!Precision::IsInfinite(aV2f))
{
aCurS2Bounds[2] = new Geom2d_Line(gp_Pnt2d(aU2f, aV2f), gp_Dir2d(1.0, 0.0));
if(!Precision::IsInfinite(aDelta))
aCurS2Bounds[2] = new Geom2d_TrimmedCurve(aCurS2Bounds[2], 0, aDelta);
}
if(!Precision::IsInfinite(aV2l))
{
aCurS2Bounds[3] = new Geom2d_Line(gp_Pnt2d(aU2l, aV2l), gp_Dir2d(1.0, 0.0));
if(!Precision::IsInfinite(aDelta))
aCurS2Bounds[3] = new Geom2d_TrimmedCurve(aCurS2Bounds[3], 0, aDelta);
}
}
Geom2dAPI_InterCurveCurve anIntCC2d;
if(!theC2d1.IsNull())
{
for(Standard_Integer aCurID = 0; aCurID < aNumberOfCurves; aCurID++)
{
if(aCurS1Bounds[aCurID].IsNull())
continue;
anIntCC2d.Init(theC2d1, aCurS1Bounds[aCurID]);
for (Standard_Integer aPntID = 1; aPntID <= anIntCC2d.NbPoints(); aPntID++)
{
const Standard_Real aParam = anIntCC2d.Intersector().Point(aPntID).ParamOnFirst();
theArrayOfParameters.Append(aParam);
}
for (Standard_Integer aSegmID = 1; aSegmID <= anIntCC2d.NbSegments(); aSegmID++)
{
Handle(Geom2d_Curve) aCS, aCSTemp;
anIntCC2d.Segment(aSegmID, aCS, aCSTemp);
theArrayOfParameters.Append(aCS->FirstParameter());
theArrayOfParameters.Append(aCS->LastParameter());
}
}
}
if(!theC2d2.IsNull())
{
for(Standard_Integer aCurID = 0; aCurID < aNumberOfCurves; aCurID++)
{
if(aCurS2Bounds[aCurID].IsNull())
continue;
anIntCC2d.Init(theC2d2, aCurS2Bounds[aCurID]);
for (Standard_Integer aPntID = 1; aPntID <= anIntCC2d.NbPoints(); aPntID++)
{
const Standard_Real aParam = anIntCC2d.Intersector().Point(aPntID).ParamOnFirst();
theArrayOfParameters.Append(aParam);
}
for (Standard_Integer aSegmID = 1; aSegmID <= anIntCC2d.NbSegments(); aSegmID++)
{
Handle(Geom2d_Curve) aCS, aCSTemp;
anIntCC2d.Segment(aSegmID, aCS, aCSTemp);
theArrayOfParameters.Append(aCS->FirstParameter());
theArrayOfParameters.Append(aCS->LastParameter());
}
}
}
std::sort(theArrayOfParameters.begin(), theArrayOfParameters.end());
}

View File

@@ -0,0 +1,23 @@
// Created on: 2015-02-18
// Created by: Nikolai BUKHALOV
// Copyright (c) 2003-2015 OPEN CASCADE SAS
//
// This file is part of Open CASCADE Technology software library.
//
// This library is free software; you can redistribute it and/or modify it under
// the terms of the GNU Lesser General Public License version 2.1 as published
// by the Free Software Foundation, with special exception defined in the file
// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
// distribution for complete text of the license and disclaimer of any warranty.
//
// Alternatively, this file may be used under the terms of Open CASCADE
// commercial license or contractual agreement.
#ifndef GeomInt_VectorOfReal_HeaderFile
#define GeomInt_VectorOfReal_HeaderFile
#include <NCollection_Vector.hxx>
typedef NCollection_Vector<Standard_Real> GeomInt_VectorOfReal;
#endif