mirror of
https://git.dev.opencascade.org/repos/occt.git
synced 2025-08-04 13:13:25 +03:00
0028599: Replacement of old Boolean operations with new ones in BRepProj_Projection algorithm
The usage of *BRepAlgo_Section* has been replaced with the usage of *BRepAlgoAPI_Section* in *BRepProj_Projection* algorithm. The TODO statements have been removed from the failing test case in the "prj" grid as they are working correctly now. The following changes have been made to improve the performance *BRepAlgoAPI_Section*: 1. Revision of the *IntPolyh_Intersection* class to avoid repeated calculation of the deflection of the same triangulation. 2. Small revision of the Edge/Face intersection algorithm to perform Extrema computation on the whole intersection range of the edge instead of discrete ranges. 3. Implementation of the extrema computation for the Circle and Sphere. 4. Correct computation of the parameter of the point on the Circle.
This commit is contained in:
@@ -113,6 +113,9 @@ void Extrema_ExtCS::Perform(const Adaptor3d_Curve& C,
|
||||
NbT = NbU = NbV = 10;
|
||||
GeomAbs_CurveType myCtype = C.GetType();
|
||||
|
||||
myDone = Standard_False;
|
||||
// Try analytic computation of extrema
|
||||
Standard_Boolean isComputeAnalytic = Standard_True;
|
||||
|
||||
switch(myCtype) {
|
||||
|
||||
@@ -233,6 +236,11 @@ void Extrema_ExtCS::Perform(const Adaptor3d_Curve& C,
|
||||
myExtElCS.Perform(C.Circle(), myS->Plane());
|
||||
break;
|
||||
}
|
||||
else if (myStype == GeomAbs_Sphere)
|
||||
{
|
||||
myExtElCS.Perform(C.Circle(), myS->Sphere());
|
||||
break;
|
||||
}
|
||||
}
|
||||
Standard_FALLTHROUGH
|
||||
case GeomAbs_Hyperbola:
|
||||
@@ -246,157 +254,225 @@ void Extrema_ExtCS::Perform(const Adaptor3d_Curve& C,
|
||||
Standard_FALLTHROUGH
|
||||
default:
|
||||
{
|
||||
Extrema_GenExtCS Ext;
|
||||
Ext.Initialize(*myS, NbU, NbV, mytolS);
|
||||
if(myCtype == GeomAbs_Hyperbola) {
|
||||
Standard_Real tmin = Max(-20., C.FirstParameter());
|
||||
Standard_Real tmax = Min(20., C.LastParameter());
|
||||
Ext.Perform(C, NbT, tmin, tmax, mytolC); // to avoid overflow
|
||||
}
|
||||
else {
|
||||
if((myCtype == GeomAbs_Circle && NbT < 13) ||
|
||||
(myCtype == GeomAbs_BSplineCurve && NbT < 13) )
|
||||
{
|
||||
NbT = 13;
|
||||
}
|
||||
Ext.Perform(C, NbT, mytolC);
|
||||
}
|
||||
isComputeAnalytic = Standard_False;
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
myDone = Ext.IsDone();
|
||||
if (myDone) {
|
||||
Standard_Integer NbExt = Ext.NbExt();
|
||||
Standard_Real T,U,V;
|
||||
Extrema_POnCurv PC;
|
||||
Extrema_POnSurf PS;
|
||||
for (i = 1; i <= NbExt; i++) {
|
||||
PC = Ext.PointOnCurve(i);
|
||||
PS = Ext.PointOnSurface(i);
|
||||
T = PC.Parameter();
|
||||
if (isComputeAnalytic)
|
||||
{
|
||||
if (myExtElCS.IsDone())
|
||||
{
|
||||
myDone = Standard_True;
|
||||
myIsPar = myExtElCS.IsParallel();
|
||||
if (myIsPar)
|
||||
{
|
||||
mySqDist.Append(myExtElCS.SquareDistance(1));
|
||||
}
|
||||
else
|
||||
{
|
||||
Standard_Integer NbExt = myExtElCS.NbExt();
|
||||
for (i = 1; i <= NbExt; i++)
|
||||
{
|
||||
Extrema_POnCurv PC;
|
||||
Extrema_POnSurf PS;
|
||||
myExtElCS.Points(i, PC, PS);
|
||||
Standard_Real Ucurve = PC.Parameter();
|
||||
Standard_Real U, V;
|
||||
PS.Parameter(U, V);
|
||||
AddSolution(C, T, U, V, PC.Value(), PS.Value(), Ext.SquareDistance(i));
|
||||
AddSolution(C, Ucurve, U, V, PC.Value(), PS.Value(), myExtElCS.SquareDistance(i));
|
||||
}
|
||||
|
||||
//Add sharp points
|
||||
Standard_Integer SolNumber = mySqDist.Length();
|
||||
Standard_Address CopyC = (Standard_Address)&C;
|
||||
Adaptor3d_Curve& aC = *(Adaptor3d_Curve*)CopyC;
|
||||
Standard_Integer NbIntervals = aC.NbIntervals(GeomAbs_C1);
|
||||
TColStd_Array1OfReal SharpPoints(1, NbIntervals+1);
|
||||
aC.Intervals(SharpPoints, GeomAbs_C1);
|
||||
|
||||
Extrema_ExtPS aProjPS;
|
||||
aProjPS.Initialize (*myS,
|
||||
myS->FirstUParameter(),
|
||||
myS->LastUParameter(),
|
||||
myS->FirstVParameter(),
|
||||
myS->LastVParameter(),
|
||||
mytolS,
|
||||
mytolS);
|
||||
|
||||
for (i = 2; i < SharpPoints.Upper(); ++i)
|
||||
if (mySqDist.Length() == 0 && NbExt > 0)
|
||||
{
|
||||
T = SharpPoints(i);
|
||||
gp_Pnt aPnt = C.Value(T);
|
||||
aProjPS.Perform (aPnt);
|
||||
if (!aProjPS.IsDone())
|
||||
continue;
|
||||
Standard_Integer NbProj = aProjPS.NbExt(), jmin = 0;
|
||||
Standard_Real MinSqDist = RealLast();
|
||||
for (j = 1; j <= NbProj; j++)
|
||||
// Analytical extrema seem to be out of curve/surface boundaries.
|
||||
// Try extremity points of curve.
|
||||
gp_Pnt aPOnC[2], aPOnS[2];
|
||||
Standard_Real aT[2] = { myucinf, myucsup }, U[2], V[2];
|
||||
Standard_Real aDist[2] = { -1, -1 };
|
||||
for (i = 0; i < 2; ++i)
|
||||
{
|
||||
Standard_Real aSqDist = aProjPS.SquareDistance(j);
|
||||
if (aSqDist < MinSqDist)
|
||||
if (Precision::IsInfinite(aT[i]))
|
||||
continue;
|
||||
|
||||
aPOnC[i] = C.Value(aT[i]);
|
||||
switch (myStype)
|
||||
{
|
||||
MinSqDist = aSqDist;
|
||||
jmin = j;
|
||||
case GeomAbs_Plane:
|
||||
{
|
||||
ElSLib::Parameters(myS->Plane(), aPOnC[i], U[i], V[i]);
|
||||
aPOnS[i] = ElSLib::Value(U[i], V[i], myS->Plane());
|
||||
break;
|
||||
}
|
||||
case GeomAbs_Sphere:
|
||||
{
|
||||
ElSLib::Parameters(myS->Sphere(), aPOnC[i], U[i], V[i]);
|
||||
aPOnS[i] = ElSLib::Value(U[i], V[i], myS->Sphere());
|
||||
break;
|
||||
}
|
||||
case GeomAbs_Cylinder:
|
||||
{
|
||||
ElSLib::Parameters(myS->Cylinder(), aPOnC[i], U[i], V[i]);
|
||||
aPOnS[i] = ElSLib::Value(U[i], V[i], myS->Cylinder());
|
||||
break;
|
||||
}
|
||||
case GeomAbs_Torus:
|
||||
{
|
||||
ElSLib::Parameters(myS->Torus(), aPOnC[i], U[i], V[i]);
|
||||
aPOnS[i] = ElSLib::Value(U[i], V[i], myS->Torus());
|
||||
break;
|
||||
}
|
||||
case GeomAbs_Cone:
|
||||
{
|
||||
ElSLib::Parameters(myS->Cone(), aPOnC[i], U[i], V[i]);
|
||||
aPOnS[i] = ElSLib::Value(U[i], V[i], myS->Cone());
|
||||
break;
|
||||
}
|
||||
default:
|
||||
continue;
|
||||
}
|
||||
|
||||
aDist[i] = aPOnC[i].SquareDistance(aPOnS[i]);
|
||||
}
|
||||
if (jmin != 0)
|
||||
|
||||
Standard_Boolean bAdd[2] = {Standard_False, Standard_False};
|
||||
|
||||
// Choose solution to add
|
||||
if (aDist[0] >= 0. && aDist[1] >= 0.)
|
||||
{
|
||||
aProjPS.Point(jmin).Parameter(U,V);
|
||||
AddSolution(C, T, U, V,
|
||||
aPnt, aProjPS.Point(jmin).Value(), MinSqDist);
|
||||
Standard_Real aDiff = aDist[0] - aDist[1];
|
||||
// Both computed -> take only minimal
|
||||
if (Abs(aDiff) < Precision::Confusion())
|
||||
// Add both
|
||||
bAdd[0] = bAdd[1] = Standard_True;
|
||||
else if (aDiff < 0)
|
||||
// Add first
|
||||
bAdd[0] = Standard_True;
|
||||
else
|
||||
// Add second
|
||||
bAdd[1] = Standard_True;
|
||||
}
|
||||
}
|
||||
//Cut sharp solutions to keep only minimum and maximum
|
||||
Standard_Integer imin = SolNumber + 1, imax = mySqDist.Length();
|
||||
for (i = SolNumber + 1; i <= mySqDist.Length(); i++)
|
||||
{
|
||||
if (mySqDist(i) < mySqDist(imin))
|
||||
imin = i;
|
||||
if (mySqDist(i) > mySqDist(imax))
|
||||
imax = i;
|
||||
}
|
||||
if (mySqDist.Length() > SolNumber + 2)
|
||||
{
|
||||
Standard_Real MinSqDist = mySqDist(imin);
|
||||
Standard_Real MaxSqDist = mySqDist(imax);
|
||||
Extrema_POnCurv MinPC = myPOnC(imin);
|
||||
Extrema_POnCurv MaxPC = myPOnC(imax);
|
||||
Extrema_POnSurf MinPS = myPOnS(imin);
|
||||
Extrema_POnSurf MaxPS = myPOnS(imax);
|
||||
else if (aDist[0] >= 0.)
|
||||
// Add first
|
||||
bAdd[0] = Standard_True;
|
||||
else if (aDist[1] >= 0.)
|
||||
// Add second
|
||||
bAdd[1] = Standard_True;
|
||||
|
||||
mySqDist.Remove(SolNumber + 1, mySqDist.Length());
|
||||
myPOnC.Remove(SolNumber + 1, myPOnC.Length());
|
||||
myPOnS.Remove(SolNumber + 1, myPOnS.Length());
|
||||
|
||||
mySqDist.Append(MinSqDist);
|
||||
myPOnC.Append(MinPC);
|
||||
myPOnS.Append(MinPS);
|
||||
mySqDist.Append(MaxSqDist);
|
||||
myPOnC.Append(MaxPC);
|
||||
myPOnS.Append(MaxPS);
|
||||
for (i = 0; i < 2; ++i)
|
||||
{
|
||||
if (bAdd[i])
|
||||
AddSolution(C, aT[i], U[i], V[i], aPOnC[i], aPOnS[i], aDist[i]);
|
||||
}
|
||||
}
|
||||
}
|
||||
return;
|
||||
}
|
||||
break;
|
||||
}
|
||||
|
||||
myDone = myExtElCS.IsDone();
|
||||
if (myDone) {
|
||||
myIsPar = myExtElCS.IsParallel();
|
||||
if (myIsPar) {
|
||||
mySqDist.Append(myExtElCS.SquareDistance(1));
|
||||
// Elementary extrema is not done, try generic solution
|
||||
Extrema_GenExtCS Ext;
|
||||
Ext.Initialize(*myS, NbU, NbV, mytolS);
|
||||
if (myCtype == GeomAbs_Hyperbola) {
|
||||
Standard_Real tmin = Max(-20., C.FirstParameter());
|
||||
Standard_Real tmax = Min(20., C.LastParameter());
|
||||
Ext.Perform(C, NbT, tmin, tmax, mytolC); // to avoid overflow
|
||||
}
|
||||
else {
|
||||
if ((myCtype == GeomAbs_Circle && NbT < 13) ||
|
||||
(myCtype == GeomAbs_BSplineCurve && NbT < 13))
|
||||
{
|
||||
NbT = 13;
|
||||
}
|
||||
else {
|
||||
Standard_Integer NbExt = myExtElCS.NbExt();
|
||||
Standard_Real U, V;
|
||||
for (i = 1; i <= NbExt; i++) {
|
||||
Extrema_POnCurv PC;
|
||||
Extrema_POnSurf PS;
|
||||
myExtElCS.Points(i, PC, PS);
|
||||
Standard_Real Ucurve = PC.Parameter();
|
||||
PS.Parameter(U, V);
|
||||
AddSolution(C, Ucurve, U, V, PC.Value(), PS.Value(), myExtElCS.SquareDistance(i));
|
||||
}
|
||||
if(mySqDist.Length() == 0 && NbExt > 0)
|
||||
Ext.Perform(C, NbT, mytolC);
|
||||
}
|
||||
|
||||
myDone = Ext.IsDone();
|
||||
if (myDone) {
|
||||
Standard_Integer NbExt = Ext.NbExt();
|
||||
Standard_Real T, U, V;
|
||||
Extrema_POnCurv PC;
|
||||
Extrema_POnSurf PS;
|
||||
for (i = 1; i <= NbExt; i++) {
|
||||
PC = Ext.PointOnCurve(i);
|
||||
PS = Ext.PointOnSurface(i);
|
||||
T = PC.Parameter();
|
||||
PS.Parameter(U, V);
|
||||
AddSolution(C, T, U, V, PC.Value(), PS.Value(), Ext.SquareDistance(i));
|
||||
}
|
||||
|
||||
//Add sharp points
|
||||
Standard_Integer SolNumber = mySqDist.Length();
|
||||
Standard_Address CopyC = (Standard_Address)&C;
|
||||
Adaptor3d_Curve& aC = *(Adaptor3d_Curve*)CopyC;
|
||||
Standard_Integer NbIntervals = aC.NbIntervals(GeomAbs_C1);
|
||||
TColStd_Array1OfReal SharpPoints(1, NbIntervals + 1);
|
||||
aC.Intervals(SharpPoints, GeomAbs_C1);
|
||||
|
||||
Extrema_ExtPS aProjPS;
|
||||
aProjPS.Initialize(*myS,
|
||||
myS->FirstUParameter(),
|
||||
myS->LastUParameter(),
|
||||
myS->FirstVParameter(),
|
||||
myS->LastVParameter(),
|
||||
mytolS,
|
||||
mytolS);
|
||||
|
||||
for (i = 2; i < SharpPoints.Upper(); ++i)
|
||||
{
|
||||
T = SharpPoints(i);
|
||||
gp_Pnt aPnt = C.Value(T);
|
||||
aProjPS.Perform(aPnt);
|
||||
if (!aProjPS.IsDone())
|
||||
continue;
|
||||
Standard_Integer NbProj = aProjPS.NbExt(), jmin = 0;
|
||||
Standard_Real MinSqDist = RealLast();
|
||||
for (j = 1; j <= NbProj; j++)
|
||||
{
|
||||
//Analytical extremas seem to be out of curve/surface boundaries.
|
||||
//For plane it is possible to add extremity points of curve
|
||||
if(myStype == GeomAbs_Plane)
|
||||
Standard_Real aSqDist = aProjPS.SquareDistance(j);
|
||||
if (aSqDist < MinSqDist)
|
||||
{
|
||||
gp_Pln aPln = myS->Plane();
|
||||
gp_Pnt PC, PP;
|
||||
if(!Precision::IsInfinite(myucinf))
|
||||
{
|
||||
PC = C.Value(myucinf);
|
||||
ElSLib::PlaneParameters(aPln.Position(), PC, U, V);
|
||||
PP = ElSLib::PlaneValue(U, V, aPln.Position());
|
||||
AddSolution(C, myucinf, U, V, PC, PP, PC.SquareDistance(PP));
|
||||
}
|
||||
if(!Precision::IsInfinite(myucsup))
|
||||
{
|
||||
PC = C.Value(myucsup);
|
||||
ElSLib::PlaneParameters(aPln.Position(), PC, U, V);
|
||||
PP = ElSLib::PlaneValue(U, V, aPln.Position());
|
||||
AddSolution(C, myucsup, U, V, PC, PP, PC.SquareDistance(PP));
|
||||
}
|
||||
MinSqDist = aSqDist;
|
||||
jmin = j;
|
||||
}
|
||||
}
|
||||
if (jmin != 0)
|
||||
{
|
||||
aProjPS.Point(jmin).Parameter(U, V);
|
||||
AddSolution(C, T, U, V,
|
||||
aPnt, aProjPS.Point(jmin).Value(), MinSqDist);
|
||||
}
|
||||
}
|
||||
//Cut sharp solutions to keep only minimum and maximum
|
||||
Standard_Integer imin = SolNumber + 1, imax = mySqDist.Length();
|
||||
for (i = SolNumber + 1; i <= mySqDist.Length(); i++)
|
||||
{
|
||||
if (mySqDist(i) < mySqDist(imin))
|
||||
imin = i;
|
||||
if (mySqDist(i) > mySqDist(imax))
|
||||
imax = i;
|
||||
}
|
||||
if (mySqDist.Length() > SolNumber + 2)
|
||||
{
|
||||
Standard_Real MinSqDist = mySqDist(imin);
|
||||
Standard_Real MaxSqDist = mySqDist(imax);
|
||||
Extrema_POnCurv MinPC = myPOnC(imin);
|
||||
Extrema_POnCurv MaxPC = myPOnC(imax);
|
||||
Extrema_POnSurf MinPS = myPOnS(imin);
|
||||
Extrema_POnSurf MaxPS = myPOnS(imax);
|
||||
|
||||
mySqDist.Remove(SolNumber + 1, mySqDist.Length());
|
||||
myPOnC.Remove(SolNumber + 1, myPOnC.Length());
|
||||
myPOnS.Remove(SolNumber + 1, myPOnS.Length());
|
||||
|
||||
mySqDist.Append(MinSqDist);
|
||||
myPOnC.Append(MinPC);
|
||||
myPOnS.Append(MinPS);
|
||||
mySqDist.Append(MaxSqDist);
|
||||
myPOnC.Append(MaxPC);
|
||||
myPOnS.Append(MaxPS);
|
||||
}
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
|
@@ -33,11 +33,13 @@
|
||||
#include <gp_Vec.hxx>
|
||||
#include <IntAna_IntConicQuad.hxx>
|
||||
#include <IntAna_Quadric.hxx>
|
||||
#include <IntAna_QuadQuadGeo.hxx>
|
||||
#include <Precision.hxx>
|
||||
#include <Standard_NotImplemented.hxx>
|
||||
#include <Standard_OutOfRange.hxx>
|
||||
#include <StdFail_InfiniteSolutions.hxx>
|
||||
#include <StdFail_NotDone.hxx>
|
||||
#include <TColStd_ListOfInteger.hxx>
|
||||
|
||||
Extrema_ExtElCS::Extrema_ExtElCS()
|
||||
{
|
||||
@@ -534,19 +536,137 @@ void Extrema_ExtElCS::Perform(const gp_Circ& ,
|
||||
|
||||
|
||||
|
||||
//=======================================================================
|
||||
//function : Extrema_ExtElCS
|
||||
//purpose : Circle/Sphere
|
||||
//=======================================================================
|
||||
Extrema_ExtElCS::Extrema_ExtElCS(const gp_Circ& C,
|
||||
const gp_Sphere& S)
|
||||
{ Perform(C, S);}
|
||||
|
||||
|
||||
|
||||
//void Extrema_ExtElCS::Perform(const gp_Circ& C,
|
||||
// const gp_Sphere& S)
|
||||
void Extrema_ExtElCS::Perform(const gp_Circ& ,
|
||||
const gp_Sphere& )
|
||||
const gp_Sphere& S)
|
||||
{
|
||||
throw Standard_NotImplemented();
|
||||
Perform(C, S);
|
||||
}
|
||||
|
||||
//=======================================================================
|
||||
//function : Perform
|
||||
//purpose : Circle/Sphere
|
||||
//=======================================================================
|
||||
void Extrema_ExtElCS::Perform(const gp_Circ& C,
|
||||
const gp_Sphere& S)
|
||||
{
|
||||
myDone = Standard_False;
|
||||
|
||||
if (gp_Lin(C.Axis()).SquareDistance(S.Location()) < Precision::SquareConfusion())
|
||||
{
|
||||
// Circle and sphere are parallel
|
||||
myIsPar = Standard_True;
|
||||
myDone = Standard_True;
|
||||
|
||||
// Compute distance from circle to the sphere
|
||||
Standard_Real aSqDistLoc = C.Location().SquareDistance(S.Location());
|
||||
Standard_Real aSqDist = aSqDistLoc + C.Radius() * C.Radius();
|
||||
Standard_Real aDist = sqrt(aSqDist) - S.Radius();
|
||||
mySqDist = new TColStd_HArray1OfReal(1, 1);
|
||||
mySqDist->SetValue(1, aDist * aDist);
|
||||
return;
|
||||
}
|
||||
|
||||
// Intersect sphere with circle's plane
|
||||
gp_Pln CPln(C.Location(), C.Axis().Direction());
|
||||
IntAna_QuadQuadGeo anInter(CPln, S);
|
||||
if (!anInter.IsDone())
|
||||
// not done
|
||||
return;
|
||||
|
||||
if (anInter.TypeInter() != IntAna_Circle)
|
||||
{
|
||||
// Intersection is empty or just a point.
|
||||
// The parallel case has already been considered,
|
||||
// thus, here we have to find only one minimal solution
|
||||
myNbExt = 1;
|
||||
myDone = Standard_True;
|
||||
|
||||
mySqDist = new TColStd_HArray1OfReal(1, 1);
|
||||
myPoint1 = new Extrema_HArray1OfPOnCurv(1, 1);
|
||||
myPoint2 = new Extrema_HArray1OfPOnSurf(1, 1);
|
||||
|
||||
// Compute parameter on circle
|
||||
const Standard_Real aT = ElCLib::Parameter(C, S.Location());
|
||||
// Compute point on circle
|
||||
gp_Pnt aPOnC = ElCLib::Value(aT, C);
|
||||
|
||||
// Compute parameters on sphere
|
||||
Standard_Real aU, aV;
|
||||
ElSLib::Parameters(S, aPOnC, aU, aV);
|
||||
// Compute point on sphere
|
||||
gp_Pnt aPOnS = ElSLib::Value(aU, aV, S);
|
||||
|
||||
// Save solution
|
||||
myPoint1->SetValue(1, Extrema_POnCurv(aT, aPOnC));
|
||||
myPoint2->SetValue(1, Extrema_POnSurf(aU, aV, aPOnS));
|
||||
mySqDist->SetValue(1, aPOnC.SquareDistance(aPOnS));
|
||||
return;
|
||||
}
|
||||
|
||||
// Here, the intersection is a circle
|
||||
|
||||
// Intersection circle
|
||||
gp_Circ aCInt = anInter.Circle(1);
|
||||
|
||||
// Perform intersection of the input circle with the intersection circle
|
||||
Extrema_ExtElC anExtC(C, aCInt);
|
||||
Standard_Boolean isExtremaCircCircValid = anExtC.IsDone() // Check if intersection is done
|
||||
&& !anExtC.IsParallel() // Parallel case has already been considered
|
||||
&& anExtC.NbExt() > 0; // Check that some solutions have been found
|
||||
if (!isExtremaCircCircValid)
|
||||
// not done
|
||||
return;
|
||||
|
||||
myDone = Standard_True;
|
||||
|
||||
// Few solutions
|
||||
Standard_Real aNbExt = anExtC.NbExt();
|
||||
// Find the minimal distance
|
||||
Standard_Real aMinSqDist = ::RealLast();
|
||||
for (Standard_Integer i = 1; i <= aNbExt; ++i)
|
||||
{
|
||||
Standard_Real aSqDist = anExtC.SquareDistance(i);
|
||||
if (aSqDist < aMinSqDist)
|
||||
aMinSqDist = aSqDist;
|
||||
}
|
||||
|
||||
// Collect all solutions close to the minimal one
|
||||
TColStd_ListOfInteger aSols;
|
||||
for (Standard_Integer i = 1; i <= aNbExt; ++i)
|
||||
{
|
||||
Standard_Real aDiff = anExtC.SquareDistance(i) - aMinSqDist;
|
||||
if (aDiff < Precision::SquareConfusion())
|
||||
aSols.Append(i);
|
||||
}
|
||||
|
||||
// Save all minimal solutions
|
||||
myNbExt = aSols.Extent();
|
||||
|
||||
mySqDist = new TColStd_HArray1OfReal(1, myNbExt);
|
||||
myPoint1 = new Extrema_HArray1OfPOnCurv(1, myNbExt);
|
||||
myPoint2 = new Extrema_HArray1OfPOnSurf(1, myNbExt);
|
||||
|
||||
TColStd_ListIteratorOfListOfInteger it(aSols);
|
||||
for (Standard_Integer iSol = 1; it.More(); it.Next(), ++iSol)
|
||||
{
|
||||
Extrema_POnCurv P1, P2;
|
||||
anExtC.Points(it.Value(), P1, P2);
|
||||
|
||||
// Compute parameters on sphere
|
||||
Standard_Real aU, aV;
|
||||
ElSLib::Parameters(S, P1.Value(), aU, aV);
|
||||
// Compute point on sphere
|
||||
gp_Pnt aPOnS = ElSLib::Value(aU, aV, S);
|
||||
|
||||
// Save solution
|
||||
myPoint1->SetValue(iSol, P1);
|
||||
myPoint2->SetValue(iSol, Extrema_POnSurf(aU, aV, aPOnS));
|
||||
mySqDist->SetValue(iSol, P1.Value().SquareDistance(aPOnS));
|
||||
}
|
||||
}
|
||||
|
||||
Extrema_ExtElCS::Extrema_ExtElCS(const gp_Circ& C,
|
||||
|
Reference in New Issue
Block a user