mirror of
https://git.dev.opencascade.org/repos/occt.git
synced 2025-04-10 18:51:21 +03:00
0029162: Geom2dInt_GInter algorithm does not find intersection of ellipse and line
Analytical intersection algorithm is implemented for ellipse-line intersection
This commit is contained in:
parent
a89a630e2a
commit
69f87d091e
@ -6,7 +6,6 @@ IntCurve_IntConicConic.cxx
|
|||||||
IntCurve_IntConicConic.hxx
|
IntCurve_IntConicConic.hxx
|
||||||
IntCurve_IntConicConic.lxx
|
IntCurve_IntConicConic.lxx
|
||||||
IntCurve_IntConicConic_1.cxx
|
IntCurve_IntConicConic_1.cxx
|
||||||
IntCurve_IntConicConic_1.hxx
|
|
||||||
IntCurve_IntConicConic_Tool.cxx
|
IntCurve_IntConicConic_Tool.cxx
|
||||||
IntCurve_IntConicConic_Tool.hxx
|
IntCurve_IntConicConic_Tool.hxx
|
||||||
IntCurve_IntConicCurveGen.gxx
|
IntCurve_IntConicCurveGen.gxx
|
||||||
|
@ -28,7 +28,6 @@
|
|||||||
#include <IntAna2d_IntPoint.hxx>
|
#include <IntAna2d_IntPoint.hxx>
|
||||||
#include <IntCurve_IConicTool.hxx>
|
#include <IntCurve_IConicTool.hxx>
|
||||||
#include <IntCurve_IntConicConic.hxx>
|
#include <IntCurve_IntConicConic.hxx>
|
||||||
#include <IntCurve_IntConicConic_1.hxx>
|
|
||||||
#include <IntCurve_PConic.hxx>
|
#include <IntCurve_PConic.hxx>
|
||||||
#include <IntRes2d_Domain.hxx>
|
#include <IntRes2d_Domain.hxx>
|
||||||
#include <Precision.hxx>
|
#include <Precision.hxx>
|
||||||
@ -37,7 +36,6 @@
|
|||||||
//=======================================================================
|
//=======================================================================
|
||||||
// Perform() for
|
// Perform() for
|
||||||
// Line - Parabola
|
// Line - Parabola
|
||||||
// Line - Elipse
|
|
||||||
// Line - Hyperbola
|
// Line - Hyperbola
|
||||||
// Circle - Parabola
|
// Circle - Parabola
|
||||||
// Circle - Elipse
|
// Circle - Elipse
|
||||||
@ -190,34 +188,6 @@ void IntCurve_IntConicConic::Perform(const gp_Lin2d& L,
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
//=======================================================================
|
|
||||||
//function : Perform
|
|
||||||
//purpose : Line - Elipse
|
|
||||||
//=======================================================================
|
|
||||||
void IntCurve_IntConicConic::Perform(const gp_Lin2d& L,
|
|
||||||
const IntRes2d_Domain& DL,
|
|
||||||
const gp_Elips2d& E,
|
|
||||||
const IntRes2d_Domain& DE,
|
|
||||||
const Standard_Real TolConf,
|
|
||||||
const Standard_Real Tol)
|
|
||||||
{
|
|
||||||
|
|
||||||
this->ResetFields();
|
|
||||||
IntCurve_IConicTool ITool(L);
|
|
||||||
IntCurve_PConic PCurve(E);
|
|
||||||
PCurve.SetAccuracy(20);
|
|
||||||
|
|
||||||
Inter.SetReversedParameters(ReversedParameters());
|
|
||||||
if(! DE.IsClosed()) {
|
|
||||||
IntRes2d_Domain D(DE);
|
|
||||||
D.SetEquivalentParameters(DE.FirstParameter(),DE.FirstParameter()+M_PI+M_PI);
|
|
||||||
Inter.Perform(ITool,DL,PCurve,D,TolConf,Tol);
|
|
||||||
}
|
|
||||||
else {
|
|
||||||
Inter.Perform(ITool,DL,PCurve,DE,TolConf,Tol);
|
|
||||||
}
|
|
||||||
this->SetValues(Inter);
|
|
||||||
}
|
|
||||||
|
|
||||||
//=======================================================================
|
//=======================================================================
|
||||||
//function : Perform
|
//function : Perform
|
||||||
|
@ -27,7 +27,6 @@
|
|||||||
#include <gp_Vec2d.hxx>
|
#include <gp_Vec2d.hxx>
|
||||||
#include <IntCurve_IConicTool.hxx>
|
#include <IntCurve_IConicTool.hxx>
|
||||||
#include <IntCurve_IntConicConic.hxx>
|
#include <IntCurve_IntConicConic.hxx>
|
||||||
#include <IntCurve_IntConicConic_1.hxx>
|
|
||||||
#include <IntCurve_IntConicConic_Tool.hxx>
|
#include <IntCurve_IntConicConic_Tool.hxx>
|
||||||
#include <IntCurve_PConic.hxx>
|
#include <IntCurve_PConic.hxx>
|
||||||
#include <IntImpParGen.hxx>
|
#include <IntImpParGen.hxx>
|
||||||
@ -37,6 +36,7 @@
|
|||||||
#include <IntRes2d_TypeTrans.hxx>
|
#include <IntRes2d_TypeTrans.hxx>
|
||||||
#include <Precision.hxx>
|
#include <Precision.hxx>
|
||||||
#include <Standard_ConstructionError.hxx>
|
#include <Standard_ConstructionError.hxx>
|
||||||
|
#include <Extrema_ExtElC2d.hxx>
|
||||||
|
|
||||||
Standard_Boolean Affichage=Standard_False;
|
Standard_Boolean Affichage=Standard_False;
|
||||||
Standard_Boolean AffichageGraph=Standard_True;
|
Standard_Boolean AffichageGraph=Standard_True;
|
||||||
@ -2245,3 +2245,476 @@ const IntRes2d_IntersectionPoint SegmentToPoint( const IntRes2d_IntersectionPoin
|
|||||||
}
|
}
|
||||||
return(IntRes2d_IntersectionPoint(Pa.Value(),u1,u2,t1,t2,Standard_False));
|
return(IntRes2d_IntersectionPoint(Pa.Value(),u1,u2,t1,t2,Standard_False));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
//=======================================================================
|
||||||
|
//function : LineEllipseGeometricIntersection
|
||||||
|
//purpose :
|
||||||
|
//=======================================================================
|
||||||
|
void LineEllipseGeometricIntersection(const gp_Lin2d& Line,
|
||||||
|
const gp_Elips2d& Ellipse,
|
||||||
|
const Standard_Real ,
|
||||||
|
const Standard_Real TolTang,
|
||||||
|
PeriodicInterval& EInt1,
|
||||||
|
PeriodicInterval& EInt2,
|
||||||
|
Standard_Integer& nbsol)
|
||||||
|
{
|
||||||
|
|
||||||
|
const gp_Ax22d& anElAxis = Ellipse.Axis();
|
||||||
|
gp_Trsf2d aTr;
|
||||||
|
aTr.SetTransformation(anElAxis.XAxis());
|
||||||
|
gp_Elips2d aTEllipse = Ellipse.Transformed(aTr);
|
||||||
|
gp_Lin2d aTLine = Line.Transformed(aTr);
|
||||||
|
Standard_Real aDY = aTLine.Position().Direction().Y();
|
||||||
|
Standard_Boolean IsVert = Abs(aDY) > 1. - 2. * Epsilon(1.);
|
||||||
|
//
|
||||||
|
Standard_Real a = aTEllipse.MajorRadius();
|
||||||
|
Standard_Real b = aTEllipse.MinorRadius();
|
||||||
|
Standard_Real a2 = a * a;
|
||||||
|
Standard_Real b2 = b * b;
|
||||||
|
//
|
||||||
|
Standard_Real eps0 = 1.e-12;
|
||||||
|
if (b / a < 1.e-5)
|
||||||
|
{
|
||||||
|
eps0 = 1.e-6;
|
||||||
|
}
|
||||||
|
//
|
||||||
|
Standard_Real anA, aB, aC;
|
||||||
|
aTLine.Coefficients(anA, aB, aC);
|
||||||
|
if (IsVert)
|
||||||
|
{
|
||||||
|
aC += aB * aTLine.Position().Location().Y();
|
||||||
|
aB = 0.;
|
||||||
|
}
|
||||||
|
//
|
||||||
|
Standard_Real x1 = 0., y1 = 0., x2 = 0., y2 = 0.;
|
||||||
|
if (Abs(aB) > eps0 )
|
||||||
|
{
|
||||||
|
Standard_Real m = -anA / aB;
|
||||||
|
Standard_Real m2 = m * m;
|
||||||
|
Standard_Real c = -aC / aB;
|
||||||
|
Standard_Real c2 = c * c;
|
||||||
|
Standard_Real D = a2 * m2 + b2 - c2;
|
||||||
|
if (D < 0.)
|
||||||
|
{
|
||||||
|
Extrema_ExtElC2d anExt(aTLine, aTEllipse);
|
||||||
|
Standard_Integer i, imin = 0;
|
||||||
|
Standard_Real dmin = RealLast();
|
||||||
|
for (i = 1; i <= anExt.NbExt(); ++i)
|
||||||
|
{
|
||||||
|
if (anExt.SquareDistance(i) < dmin)
|
||||||
|
{
|
||||||
|
dmin = anExt.SquareDistance(i);
|
||||||
|
imin = i;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
if (imin > 0 && dmin <= TolTang * TolTang)
|
||||||
|
{
|
||||||
|
nbsol = 1;
|
||||||
|
Extrema_POnCurv2d aP1, aP2;
|
||||||
|
anExt.Points(imin, aP1, aP2);
|
||||||
|
Standard_Real pe1 = aP2.Parameter();
|
||||||
|
EInt1.SetValues(pe1, pe1);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
nbsol = 0;
|
||||||
|
}
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
D = Sqrt(D);
|
||||||
|
Standard_Real n = a2 * m2 + b2;
|
||||||
|
Standard_Real k = a * b * D / n;
|
||||||
|
Standard_Real l = -a2 * m * c / n;
|
||||||
|
x1 = l + k;
|
||||||
|
y1 = m * x1 + c;
|
||||||
|
x2 = l - k;
|
||||||
|
y2 = m * x2 + c;
|
||||||
|
nbsol = 2;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
x1 = -aC / anA;
|
||||||
|
if (Abs(x1) > a + TolTang)
|
||||||
|
{
|
||||||
|
nbsol = 0;
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
else if (Abs(x1) >= a - Epsilon(1. + a))
|
||||||
|
{
|
||||||
|
nbsol = 1;
|
||||||
|
y1 = 0.;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
y1 = b * Sqrt(1. - x1 * x1 / a2);
|
||||||
|
x2 = x1;
|
||||||
|
y2 = -y1;
|
||||||
|
nbsol = 2;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
gp_Pnt2d aP1(x1, y1);
|
||||||
|
gp_Pnt2d aP2(x2, y2);
|
||||||
|
Standard_Real pe1 = 0., pe2 = 0.;
|
||||||
|
pe1 = ElCLib::Parameter(aTEllipse, aP1);
|
||||||
|
if (nbsol > 1)
|
||||||
|
{
|
||||||
|
pe2 = ElCLib::Parameter(aTEllipse, aP2);
|
||||||
|
if (pe2 < pe1)
|
||||||
|
{
|
||||||
|
Standard_Real t = pe1;
|
||||||
|
pe1 = pe2;
|
||||||
|
pe2 = t;
|
||||||
|
}
|
||||||
|
EInt2.SetValues(pe2, pe2);
|
||||||
|
}
|
||||||
|
EInt1.SetValues(pe1, pe1);
|
||||||
|
|
||||||
|
|
||||||
|
}
|
||||||
|
//=======================================================================
|
||||||
|
//function : ProjectOnLAndIntersectWithLDomain
|
||||||
|
//purpose :
|
||||||
|
//=======================================================================
|
||||||
|
void ProjectOnLAndIntersectWithLDomain(const gp_Elips2d& Ellipse
|
||||||
|
, const gp_Lin2d& Line
|
||||||
|
, PeriodicInterval& EDomainAndRes
|
||||||
|
, Interval& LDomain
|
||||||
|
, PeriodicInterval* EllipseSolution
|
||||||
|
, Interval* LineSolution
|
||||||
|
, Standard_Integer &NbSolTotal
|
||||||
|
, const IntRes2d_Domain& RefLineDomain
|
||||||
|
, const IntRes2d_Domain&)
|
||||||
|
{
|
||||||
|
|
||||||
|
if (EDomainAndRes.IsNull()) return;
|
||||||
|
//-------------------------------------------------------------------------
|
||||||
|
//-- On cherche l intervalle correspondant sur C2
|
||||||
|
//-- Puis on intersecte l intervalle avec le domaine de C2
|
||||||
|
//-- Enfin, on cherche l intervalle correspondant sur C1
|
||||||
|
//--
|
||||||
|
|
||||||
|
Standard_Real Linf = ElCLib::Parameter(Line
|
||||||
|
, ElCLib::Value(EDomainAndRes.Binf, Ellipse));
|
||||||
|
Standard_Real Lsup = ElCLib::Parameter(Line
|
||||||
|
, ElCLib::Value(EDomainAndRes.Bsup, Ellipse));
|
||||||
|
|
||||||
|
Interval LInter(Linf, Lsup); //-- Necessairement Borne
|
||||||
|
|
||||||
|
Interval LInterAndDomain = LDomain.IntersectionWithBounded(LInter);
|
||||||
|
|
||||||
|
if (!LInterAndDomain.IsNull) {
|
||||||
|
|
||||||
|
Standard_Real DomLinf = (RefLineDomain.HasFirstPoint()) ? RefLineDomain.FirstParameter() : -Precision::Infinite();
|
||||||
|
Standard_Real DomLsup = (RefLineDomain.HasLastPoint()) ? RefLineDomain.LastParameter() : Precision::Infinite();
|
||||||
|
|
||||||
|
Linf = LInterAndDomain.Binf;
|
||||||
|
Lsup = LInterAndDomain.Bsup;
|
||||||
|
|
||||||
|
if (Linf<DomLinf) {
|
||||||
|
Linf = DomLinf;
|
||||||
|
}
|
||||||
|
if (Lsup<DomLinf) {
|
||||||
|
Lsup = DomLinf;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (Linf>DomLsup) {
|
||||||
|
Linf = DomLsup;
|
||||||
|
}
|
||||||
|
if (Lsup>DomLsup) {
|
||||||
|
Lsup = DomLsup;
|
||||||
|
}
|
||||||
|
|
||||||
|
LInterAndDomain.Binf = Linf;
|
||||||
|
LInterAndDomain.Bsup = Lsup;
|
||||||
|
|
||||||
|
|
||||||
|
Standard_Real Einf = EDomainAndRes.Binf;
|
||||||
|
Standard_Real Esup = EDomainAndRes.Bsup;
|
||||||
|
|
||||||
|
if (Einf >= Esup) { Einf = EDomainAndRes.Binf; Esup = EDomainAndRes.Bsup; }
|
||||||
|
EllipseSolution[NbSolTotal] = PeriodicInterval(Einf, Esup);
|
||||||
|
if (EllipseSolution[NbSolTotal].Length() > M_PI)
|
||||||
|
EllipseSolution[NbSolTotal].Complement();
|
||||||
|
|
||||||
|
LineSolution[NbSolTotal] = LInterAndDomain;
|
||||||
|
NbSolTotal++;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
//=======================================================================
|
||||||
|
//function : Perform
|
||||||
|
//purpose : Line - Elipse
|
||||||
|
//=======================================================================
|
||||||
|
void IntCurve_IntConicConic::Perform(const gp_Lin2d& L, const
|
||||||
|
IntRes2d_Domain& DL, const gp_Elips2d& E,
|
||||||
|
const IntRes2d_Domain& DE, const Standard_Real TolConf,
|
||||||
|
const Standard_Real Tol)
|
||||||
|
{
|
||||||
|
Standard_Boolean TheReversedParameters = ReversedParameters();
|
||||||
|
this->ResetFields();
|
||||||
|
this->SetReversedParameters(TheReversedParameters);
|
||||||
|
|
||||||
|
Standard_Integer nbsol = 0;
|
||||||
|
PeriodicInterval EInt1, EInt2;
|
||||||
|
|
||||||
|
LineEllipseGeometricIntersection(L, E, TolConf, Tol, EInt1, EInt2, nbsol);
|
||||||
|
done = Standard_True;
|
||||||
|
if (nbsol == 0)
|
||||||
|
{
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
//
|
||||||
|
if (nbsol == 2 && EInt2.Bsup == EInt1.Binf + PIpPI) {
|
||||||
|
Standard_Real FirstBound = DE.FirstParameter();
|
||||||
|
Standard_Real LastBound = DE.LastParameter();
|
||||||
|
Standard_Real FirstTol = DE.FirstTolerance();
|
||||||
|
Standard_Real LastTol = DE.LastTolerance();
|
||||||
|
if (EInt1.Binf == 0 && FirstBound - FirstTol > EInt1.Bsup)
|
||||||
|
{
|
||||||
|
nbsol = 1;
|
||||||
|
EInt1.SetValues(EInt2.Binf, EInt2.Bsup);
|
||||||
|
}
|
||||||
|
else if (EInt2.Bsup == PIpPI && LastBound + LastTol < EInt2.Binf)
|
||||||
|
{
|
||||||
|
nbsol = 1;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
//
|
||||||
|
PeriodicInterval EDomain(DE);
|
||||||
|
Standard_Real deltat = EDomain.Bsup - EDomain.Binf;
|
||||||
|
while (EDomain.Binf >= PIpPI) EDomain.Binf -= PIpPI;
|
||||||
|
while (EDomain.Binf < 0.0) EDomain.Binf += PIpPI;
|
||||||
|
EDomain.Bsup = EDomain.Binf + deltat;
|
||||||
|
//
|
||||||
|
Standard_Real BinfModif = EDomain.Binf;
|
||||||
|
Standard_Real BsupModif = EDomain.Bsup;
|
||||||
|
BinfModif -= DE.FirstTolerance() / E.MinorRadius();
|
||||||
|
BsupModif += DE.LastTolerance() / E.MinorRadius();
|
||||||
|
deltat = BsupModif - BinfModif;
|
||||||
|
if (deltat <= PIpPI) {
|
||||||
|
EDomain.Binf = BinfModif;
|
||||||
|
EDomain.Bsup = BsupModif;
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
Standard_Real t = PIpPI - deltat;
|
||||||
|
t *= 0.5;
|
||||||
|
EDomain.Binf = BinfModif + t;
|
||||||
|
EDomain.Bsup = BsupModif - t;
|
||||||
|
}
|
||||||
|
deltat = EDomain.Bsup - EDomain.Binf;
|
||||||
|
while (EDomain.Binf >= PIpPI) EDomain.Binf -= PIpPI;
|
||||||
|
while (EDomain.Binf < 0.0) EDomain.Binf += PIpPI;
|
||||||
|
EDomain.Bsup = EDomain.Binf + deltat;
|
||||||
|
//
|
||||||
|
Interval LDomain(DL);
|
||||||
|
|
||||||
|
Standard_Integer NbSolTotal = 0;
|
||||||
|
|
||||||
|
PeriodicInterval SolutionEllipse[4];
|
||||||
|
Interval SolutionLine[4];
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
//----------- Treatment of first geometric interval EInt1 ----
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
PeriodicInterval EDomainAndRes = EDomain.FirstIntersection(EInt1);
|
||||||
|
|
||||||
|
ProjectOnLAndIntersectWithLDomain(E, L, EDomainAndRes, LDomain, SolutionEllipse
|
||||||
|
, SolutionLine, NbSolTotal, DL, DE);
|
||||||
|
|
||||||
|
EDomainAndRes = EDomain.SecondIntersection(EInt1);
|
||||||
|
|
||||||
|
ProjectOnLAndIntersectWithLDomain(E, L, EDomainAndRes, LDomain, SolutionEllipse
|
||||||
|
, SolutionLine, NbSolTotal, DL, DE);
|
||||||
|
|
||||||
|
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
//----------- Treatment of second geometric interval EInt2 ----
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
if (nbsol == 2)
|
||||||
|
{
|
||||||
|
EDomainAndRes = EDomain.FirstIntersection(EInt2);
|
||||||
|
|
||||||
|
ProjectOnLAndIntersectWithLDomain(E, L, EDomainAndRes, LDomain, SolutionEllipse
|
||||||
|
, SolutionLine, NbSolTotal, DL, DE);
|
||||||
|
|
||||||
|
EDomainAndRes = EDomain.SecondIntersection(EInt2);
|
||||||
|
|
||||||
|
ProjectOnLAndIntersectWithLDomain(E, L, EDomainAndRes, LDomain, SolutionEllipse
|
||||||
|
, SolutionLine, NbSolTotal, DL, DE);
|
||||||
|
}
|
||||||
|
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
//-- Calculation of Transitions at Positions.
|
||||||
|
//----------------------------------------------------------------------
|
||||||
|
Standard_Real R = E.MinorRadius();
|
||||||
|
Standard_Integer i;
|
||||||
|
Standard_Real MaxTol = TolConf;
|
||||||
|
if (MaxTol<Tol) MaxTol = Tol;
|
||||||
|
if (MaxTol<1.0e-10) MaxTol = 1.0e-10;
|
||||||
|
|
||||||
|
for (i = 0; i<NbSolTotal; i++) {
|
||||||
|
if ((R * SolutionEllipse[i].Length())<MaxTol
|
||||||
|
&& (SolutionLine[i].Length())<MaxTol) {
|
||||||
|
|
||||||
|
Standard_Real t = (SolutionEllipse[i].Binf + SolutionEllipse[i].Bsup)*0.5;
|
||||||
|
SolutionEllipse[i].Binf = SolutionEllipse[i].Bsup = t;
|
||||||
|
|
||||||
|
t = (SolutionLine[i].Binf + SolutionLine[i].Bsup)*0.5;
|
||||||
|
SolutionLine[i].Binf = SolutionLine[i].Bsup = t;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
//
|
||||||
|
if (NbSolTotal) {
|
||||||
|
gp_Ax22d EllipseAxis = E.Axis();
|
||||||
|
gp_Ax2d LineAxis = L.Position();
|
||||||
|
gp_Pnt2d P1a, P2a, P1b, P2b;
|
||||||
|
gp_Vec2d Tan1, Tan2, Norm1;
|
||||||
|
gp_Vec2d Norm2(0.0, 0.0);
|
||||||
|
IntRes2d_Transition T1a, T2a, T1b, T2b;
|
||||||
|
IntRes2d_Position Pos1a, Pos1b, Pos2a, Pos2b;
|
||||||
|
|
||||||
|
ElCLib::EllipseD1(SolutionEllipse[0].Binf, EllipseAxis, E.MajorRadius(), E.MinorRadius(), P1a, Tan1);
|
||||||
|
ElCLib::LineD1(SolutionLine[0].Binf, LineAxis, P2a, Tan2);
|
||||||
|
|
||||||
|
Standard_Boolean isOpposite = (Tan1.Dot(Tan2) < 0.0);
|
||||||
|
for (i = 0; i<NbSolTotal; i++)
|
||||||
|
{
|
||||||
|
Standard_Real p1 = SolutionEllipse[i].Binf;
|
||||||
|
Standard_Real p2 = SolutionEllipse[i].Bsup;
|
||||||
|
Standard_Real q1 = DE.FirstParameter();
|
||||||
|
Standard_Real q2 = DE.LastParameter();
|
||||||
|
|
||||||
|
if (p1>q2) {
|
||||||
|
do {
|
||||||
|
p1 -= PIpPI;
|
||||||
|
p2 -= PIpPI;
|
||||||
|
} while ((p1>q2));
|
||||||
|
}
|
||||||
|
else if (p2<q1) {
|
||||||
|
do {
|
||||||
|
p1 += PIpPI;
|
||||||
|
p2 += PIpPI;
|
||||||
|
} while ((p2<q1));
|
||||||
|
}
|
||||||
|
if (p1<q1 && p2>q1) {
|
||||||
|
p1 = q1;
|
||||||
|
}
|
||||||
|
if (p1<q2 && p2>q2) {
|
||||||
|
p2 = q2;
|
||||||
|
}
|
||||||
|
|
||||||
|
SolutionEllipse[i].Binf = p1;
|
||||||
|
SolutionEllipse[i].Bsup = p2;
|
||||||
|
|
||||||
|
Standard_Real Linf = isOpposite ? SolutionLine[i].Bsup : SolutionLine[i].Binf;
|
||||||
|
Standard_Real Lsup = isOpposite ? SolutionLine[i].Binf : SolutionLine[i].Bsup;
|
||||||
|
|
||||||
|
if (Linf > Lsup) {
|
||||||
|
Standard_Real T = SolutionEllipse[i].Binf;
|
||||||
|
SolutionEllipse[i].Binf = SolutionEllipse[i].Bsup;
|
||||||
|
SolutionEllipse[i].Bsup = T;
|
||||||
|
T = Linf; Linf = Lsup; Lsup = T;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
ElCLib::EllipseD2(SolutionEllipse[i].Binf, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1a, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Linf, LineAxis, P2a, Tan2);
|
||||||
|
|
||||||
|
IntImpParGen::DeterminePosition(Pos1a, DE, P1a, SolutionEllipse[i].Binf);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2a, DL, P2a, Linf);
|
||||||
|
Determine_Transition_LC(Pos1a, Tan1, Norm1, T1a, Pos2a, Tan2, Norm2, T2a, Tol);
|
||||||
|
Standard_Real Einf;
|
||||||
|
if (Pos1a == IntRes2d_End) {
|
||||||
|
Einf = DE.LastParameter();
|
||||||
|
P1a = DE.LastPoint();
|
||||||
|
Linf = ElCLib::Parameter(L, P1a);
|
||||||
|
|
||||||
|
ElCLib::EllipseD2(Einf, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1a, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Linf, LineAxis, P2a, Tan2);
|
||||||
|
IntImpParGen::DeterminePosition(Pos1a, DE, P1a, Einf);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2a, DL, P2a, Linf);
|
||||||
|
Determine_Transition_LC(Pos1a, Tan1, Norm1, T1a, Pos2a, Tan2, Norm2, T2a, Tol);
|
||||||
|
}
|
||||||
|
else if (Pos1a == IntRes2d_Head) {
|
||||||
|
Einf = DE.FirstParameter();
|
||||||
|
P1a = DE.FirstPoint();
|
||||||
|
Linf = ElCLib::Parameter(L, P1a);
|
||||||
|
|
||||||
|
ElCLib::EllipseD2(Einf, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1a, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Linf, LineAxis, P2a, Tan2);
|
||||||
|
IntImpParGen::DeterminePosition(Pos1a, DE, P1a, Einf);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2a, DL, P2a, Linf);
|
||||||
|
Determine_Transition_LC(Pos1a, Tan1, Norm1, T1a, Pos2a, Tan2, Norm2, T2a, Tol);
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
Einf = NormalizeOnCircleDomain(SolutionEllipse[i].Binf, DE);
|
||||||
|
}
|
||||||
|
|
||||||
|
IntRes2d_IntersectionPoint NewPoint1(P1a, Linf, Einf, T2a, T1a, ReversedParameters());
|
||||||
|
|
||||||
|
if ((SolutionLine[i].Length() + SolutionEllipse[i].Length()) >0.0) {
|
||||||
|
|
||||||
|
ElCLib::EllipseD2(SolutionEllipse[i].Binf, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1b, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Lsup, LineAxis, P2b, Tan2);
|
||||||
|
|
||||||
|
IntImpParGen::DeterminePosition(Pos1b, DE, P1b, SolutionEllipse[i].Bsup);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2b, DL, P2b, Lsup);
|
||||||
|
Determine_Transition_LC(Pos1b, Tan1, Norm1, T1b, Pos2b, Tan2, Norm2, T2b, Tol);
|
||||||
|
Standard_Real Esup;
|
||||||
|
if (Pos1b == IntRes2d_End) {
|
||||||
|
Esup = DL.LastParameter();
|
||||||
|
P1b = DE.LastPoint();
|
||||||
|
Lsup = ElCLib::Parameter(L, P1b);
|
||||||
|
ElCLib::EllipseD2(Esup, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1b, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Lsup, LineAxis, P2b, Tan2);
|
||||||
|
|
||||||
|
IntImpParGen::DeterminePosition(Pos1b, DE, P1b, Esup);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2b, DL, P2b, Lsup);
|
||||||
|
Determine_Transition_LC(Pos1b, Tan1, Norm1, T1b, Pos2b, Tan2, Norm2, T2b, Tol);
|
||||||
|
}
|
||||||
|
else if (Pos1b == IntRes2d_Head) {
|
||||||
|
Esup = DE.FirstParameter();
|
||||||
|
P1b = DE.FirstPoint();
|
||||||
|
Lsup = ElCLib::Parameter(L, P1b);
|
||||||
|
ElCLib::EllipseD2(Esup, EllipseAxis, E.MajorRadius(),
|
||||||
|
E.MinorRadius(), P1b, Tan1, Norm1);
|
||||||
|
ElCLib::LineD1(Lsup, LineAxis, P2b, Tan2);
|
||||||
|
|
||||||
|
IntImpParGen::DeterminePosition(Pos1b, DE, P1b, Esup);
|
||||||
|
IntImpParGen::DeterminePosition(Pos2b, DL, P2b, Lsup);
|
||||||
|
Determine_Transition_LC(Pos1b, Tan1, Norm1, T1b, Pos2b, Tan2, Norm2, T2b, Tol);
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
Esup = NormalizeOnCircleDomain(SolutionEllipse[i].Bsup, DE);
|
||||||
|
}
|
||||||
|
|
||||||
|
IntRes2d_IntersectionPoint NewPoint2(P1b, Lsup, Esup, T2b, T1b, ReversedParameters());
|
||||||
|
|
||||||
|
if (((Abs(Esup - Einf)*R > MaxTol) && (Abs(Lsup - Linf) > MaxTol))
|
||||||
|
|| (T1a.TransitionType() != T2a.TransitionType())) {
|
||||||
|
IntRes2d_IntersectionSegment NewSeg(NewPoint1, NewPoint2, isOpposite, ReversedParameters());
|
||||||
|
Append(NewSeg);
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
if (Pos1a != IntRes2d_Middle || Pos2a != IntRes2d_Middle) {
|
||||||
|
Insert(NewPoint1);
|
||||||
|
}
|
||||||
|
if (Pos1b != IntRes2d_Middle || Pos2b != IntRes2d_Middle) {
|
||||||
|
Insert(NewPoint2);
|
||||||
|
}
|
||||||
|
|
||||||
|
}
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
Insert(NewPoint1);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
@ -1,92 +0,0 @@
|
|||||||
// Created on: 1992-05-06
|
|
||||||
// Created by: Laurent BUCHARD
|
|
||||||
// Copyright (c) 1992-1999 Matra Datavision
|
|
||||||
// Copyright (c) 1999-2014 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 IntCurve_IntConicConic_1_HeaderFile
|
|
||||||
#define IntCurve_IntConicConic_1_HeaderFile
|
|
||||||
|
|
||||||
#include <IntCurve_IntConicConic_Tool.hxx>
|
|
||||||
|
|
||||||
|
|
||||||
//======================================================================
|
|
||||||
//=== I n t e r s e c t i o n C e r c l e C e r c l e =====
|
|
||||||
//======================================================================
|
|
||||||
|
|
||||||
//----------------------------------------------------------------------
|
|
||||||
void CircleCircleGeometricIntersection(const gp_Circ2d& C1
|
|
||||||
,const gp_Circ2d& C2
|
|
||||||
,const Standard_Real Tol
|
|
||||||
,PeriodicInterval& C1_Res1
|
|
||||||
,PeriodicInterval& C1_Res2
|
|
||||||
,Standard_Integer& nbsol);
|
|
||||||
//----------------------------------------------------------------------
|
|
||||||
void CircleCircleDomainIntersection(const gp_Circ2d& C1
|
|
||||||
,const gp_Circ2d& C2
|
|
||||||
,const Standard_Real Tol
|
|
||||||
,PeriodicInterval& Res1
|
|
||||||
,PeriodicInterval& C1_Res2
|
|
||||||
,Standard_Integer& nbsol);
|
|
||||||
//----------------------------------------------------------------------
|
|
||||||
void ProjectOnC2AndIntersectWithC2Domain(const gp_Circ2d& Circle1
|
|
||||||
,const gp_Circ2d& Circle2
|
|
||||||
,PeriodicInterval& C1DomainAndRes
|
|
||||||
,PeriodicInterval& C2Domain
|
|
||||||
,PeriodicInterval* SolutionC1
|
|
||||||
,PeriodicInterval* SolutionC2
|
|
||||||
,Standard_Integer &NbSolTotal
|
|
||||||
,const Standard_Boolean IdentCircles);
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
//======================================================================
|
|
||||||
//=== I n t e r s e c t i o n L i g n e C e r c l e =====
|
|
||||||
//======================================================================
|
|
||||||
void LineCircleGeometricIntersection(const gp_Lin2d& Line
|
|
||||||
,const gp_Circ2d& Circle
|
|
||||||
,const Standard_Real Tol
|
|
||||||
,PeriodicInterval& C1Int
|
|
||||||
,PeriodicInterval& C2Int
|
|
||||||
,Standard_Integer& nbsol);
|
|
||||||
|
|
||||||
|
|
||||||
void ProjectOnLAndIntersectWithLDomain(const gp_Circ2d& Circle
|
|
||||||
,const gp_Lin2d& Line
|
|
||||||
,PeriodicInterval& CDomainAndRes
|
|
||||||
,Interval& LDomain
|
|
||||||
,PeriodicInterval* CircleSolution
|
|
||||||
,Interval* LineSolution
|
|
||||||
,Standard_Integer &NbSolTotal);
|
|
||||||
|
|
||||||
|
|
||||||
//======================================================================
|
|
||||||
//=== I n t e r s e c t i o n L i g n e L i g n e =====
|
|
||||||
//======================================================================
|
|
||||||
|
|
||||||
void DomainIntersection(const IntRes2d_Domain& Domain
|
|
||||||
,const Standard_Real U1inf
|
|
||||||
,const Standard_Real U1sup
|
|
||||||
,Standard_Real& Res1inf
|
|
||||||
,Standard_Real& Res1sup);
|
|
||||||
|
|
||||||
void LineLineGeometricIntersection(const gp_Lin2d& L1
|
|
||||||
,const gp_Lin2d& L2
|
|
||||||
,const Standard_Real Tol
|
|
||||||
,Standard_Real& U1
|
|
||||||
,Standard_Real& U2
|
|
||||||
,Standard_Real& SinDemiAngle
|
|
||||||
,Standard_Integer& nbsol);
|
|
||||||
|
|
||||||
|
|
||||||
#endif
|
|
22
tests/bugs/modalg_7/bug29162
Normal file
22
tests/bugs/modalg_7/bug29162
Normal file
@ -0,0 +1,22 @@
|
|||||||
|
puts "========"
|
||||||
|
puts "OCC29162"
|
||||||
|
puts "========"
|
||||||
|
puts ""
|
||||||
|
##################################################
|
||||||
|
# Geom2dInt_GInter algorithm does not find
|
||||||
|
# intersection of ellipse and line
|
||||||
|
##################################################
|
||||||
|
|
||||||
|
ellipse e -610.348096534595 -710.720096056787 0.999999902285153 0.000442074298181498 15.0066332711999 0.291884102212871
|
||||||
|
trim e e 3.09462291909258 9.37780822627216
|
||||||
|
|
||||||
|
line l -625.34430362036 -680.713463264921 -0.000252749178714602 -0.999999968058926
|
||||||
|
trim l l 0 30.0132665523925
|
||||||
|
|
||||||
|
set cx 0
|
||||||
|
set cy 0
|
||||||
|
|
||||||
|
regexp {Intersection point 1 : ([-0-9.+eE]+) ([-0-9.+eE]+)} [2dintersect e l ] full cx cy
|
||||||
|
checkreal n1 $cx -625.35188801291508 1.e-7 1.e-7
|
||||||
|
checkreal n2 $cy -710.72104765746838 1.e-7 1.e-7
|
||||||
|
|
@ -1,4 +1,6 @@
|
|||||||
# !!!! This file is generated automatically, do not edit manually! See end script
|
# !!!! This file is generated automatically, do not edit manually! See end script
|
||||||
|
|
||||||
|
set LinuxDiff 3
|
||||||
set filename Inventor_scissor.stp
|
set filename Inventor_scissor.stp
|
||||||
|
|
||||||
set ref_data {
|
set ref_data {
|
||||||
@ -7,7 +9,7 @@ TPSTAT : Faulties = 0 ( 0 ) Warnings = 2 ( 2 ) Summary = 2 ( 2 )
|
|||||||
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
|
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
|
||||||
NBSHAPES : Solid = 3 ( 3 ) Shell = 3 ( 3 ) Face = 127 ( 127 ) Summary = 838 ( 838 )
|
NBSHAPES : Solid = 3 ( 3 ) Shell = 3 ( 3 ) Face = 127 ( 127 ) Summary = 838 ( 838 )
|
||||||
STATSHAPE : Solid = 3 ( 3 ) Shell = 3 ( 3 ) Face = 127 ( 127 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 346 ( 346 )
|
STATSHAPE : Solid = 3 ( 3 ) Shell = 3 ( 3 ) Face = 127 ( 127 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 346 ( 346 )
|
||||||
TOLERANCE : MaxTol = 0.0005363486823 ( 0.0005363486821 ) AvgTol = 1.61499994e-005 ( 1.636295799e-005 )
|
TOLERANCE : MaxTol = 0.0005363486823 ( 0.0005363486821 ) AvgTol = 1.614962503e-005 ( 1.636256008e-005 )
|
||||||
LABELS : N0Labels = 4 ( 4 ) N1Labels = 3 ( 3 ) N2Labels = 0 ( 0 ) TotalLabels = 7 ( 7 ) NameLabels = 7 ( 7 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
LABELS : N0Labels = 4 ( 4 ) N1Labels = 3 ( 3 ) N2Labels = 0 ( 0 ) TotalLabels = 7 ( 7 ) NameLabels = 7 ( 7 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
||||||
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
||||||
NCOLORS : NColors = 0 ( 0 )
|
NCOLORS : NColors = 0 ( 0 )
|
||||||
|
@ -1,14 +1,15 @@
|
|||||||
# !!!! This file is generated automatically, do not edit manually! See end script
|
# !!!! This file is generated automatically, do not edit manually! See end script
|
||||||
|
|
||||||
|
set LinuxDiff 3
|
||||||
set filename PRO8843.stp
|
set filename PRO8843.stp
|
||||||
|
|
||||||
set ref_data {
|
set ref_data {
|
||||||
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
||||||
TPSTAT : Faulties = 3 ( 0 ) Warnings = 10 ( 32 ) Summary = 13 ( 32 )
|
TPSTAT : Faulties = 0 ( 0 ) Warnings = 11 ( 32 ) Summary = 11 ( 32 )
|
||||||
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
|
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
|
||||||
NBSHAPES : Solid = 25 ( 25 ) Shell = 45 ( 45 ) Face = 198 ( 198 ) Summary = 1343 ( 1343 )
|
NBSHAPES : Solid = 25 ( 25 ) Shell = 45 ( 45 ) Face = 198 ( 198 ) Summary = 1343 ( 1343 )
|
||||||
STATSHAPE : Solid = 25 ( 25 ) Shell = 45 ( 45 ) Face = 198 ( 198 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 496 ( 496 )
|
STATSHAPE : Solid = 25 ( 25 ) Shell = 45 ( 45 ) Face = 198 ( 198 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 496 ( 496 )
|
||||||
TOLERANCE : MaxTol = 0.02054327471 ( 0.0293421074 ) AvgTol = 0.0005819945737 ( 0.001390253412 )
|
TOLERANCE : MaxTol = 0.02054327471 ( 0.02934181398 ) AvgTol = 0.000584176674 ( 0.00139221264 )
|
||||||
LABELS : N0Labels = 3 ( 3 ) N1Labels = 67 ( 67 ) N2Labels = 0 ( 0 ) TotalLabels = 70 ( 70 ) NameLabels = 5 ( 5 ) ColorLabels = 45 ( 45 ) LayerLabels = 45 ( 45 )
|
LABELS : N0Labels = 3 ( 3 ) N1Labels = 67 ( 67 ) N2Labels = 0 ( 0 ) TotalLabels = 70 ( 70 ) NameLabels = 5 ( 5 ) ColorLabels = 45 ( 45 ) LayerLabels = 45 ( 45 )
|
||||||
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
||||||
NCOLORS : NColors = 2 ( 2 )
|
NCOLORS : NColors = 2 ( 2 )
|
||||||
|
@ -1,15 +1,13 @@
|
|||||||
# !!!! This file is generated automatically, do not edit manually! See end script
|
# !!!! This file is generated automatically, do not edit manually! See end script
|
||||||
|
|
||||||
|
|
||||||
set filename 612319029MB-HEAD-CYLINDER.stp
|
set filename 612319029MB-HEAD-CYLINDER.stp
|
||||||
|
|
||||||
set ref_data {
|
set ref_data {
|
||||||
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
||||||
TPSTAT : Faulties = 0 ( 0 ) Warnings = 1869 ( 2152 ) Summary = 1869 ( 2152 )
|
TPSTAT : Faulties = 0 ( 0 ) Warnings = 1660 ( 2188 ) Summary = 1660 ( 2188 )
|
||||||
CHECKSHAPE : Wires = 36 ( 36 ) Faces = 36 ( 36 ) Shells = 0 ( 0 ) Solids = 1 ( 1 )
|
CHECKSHAPE : Wires = 36 ( 36 ) Faces = 36 ( 36 ) Shells = 0 ( 0 ) Solids = 1 ( 1 )
|
||||||
NBSHAPES : Solid = 1 ( 1 ) Shell = 17 ( 17 ) Face = 6173 ( 6173 ) Summary = 41176 ( 41161 )
|
NBSHAPES : Solid = 1 ( 1 ) Shell = 17 ( 17 ) Face = 6173 ( 6173 ) Summary = 41176 ( 41161 )
|
||||||
STATSHAPE : Solid = 1 ( 1 ) Shell = 17 ( 17 ) Face = 6173 ( 6173 ) FreeWire = 0 ( 0 ) FreeEdge = 4 ( 4 ) SharedEdge = 17483 ( 17472 )
|
STATSHAPE : Solid = 1 ( 1 ) Shell = 17 ( 17 ) Face = 6173 ( 6173 ) FreeWire = 0 ( 0 ) FreeEdge = 4 ( 4 ) SharedEdge = 17483 ( 17472 )
|
||||||
TOLERANCE : MaxTol = 0.07753055119 ( 0.08114875624 ) AvgTol = 0.002252828596 ( 0.008808842178 )
|
TOLERANCE : MaxTol = 0.07752745013 ( 0.08114469901 ) AvgTol = 0.002255222432 ( 0.008779543621 )
|
||||||
LABELS : N0Labels = 3 ( 3 ) N1Labels = 2 ( 2 ) N2Labels = 0 ( 0 ) TotalLabels = 5 ( 5 ) NameLabels = 5 ( 5 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
LABELS : N0Labels = 3 ( 3 ) N1Labels = 2 ( 2 ) N2Labels = 0 ( 0 ) TotalLabels = 5 ( 5 ) NameLabels = 5 ( 5 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
||||||
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
||||||
NCOLORS : NColors = 0 ( 0 )
|
NCOLORS : NColors = 0 ( 0 )
|
||||||
|
@ -5,11 +5,11 @@ set filename Z8INV5.stp
|
|||||||
|
|
||||||
set ref_data {
|
set ref_data {
|
||||||
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
|
||||||
TPSTAT : Faulties = 0 ( 0 ) Warnings = 118 ( 622 ) Summary = 118 ( 622 )
|
TPSTAT : Faulties = 0 ( 0 ) Warnings = 119 ( 620 ) Summary = 119 ( 620 )
|
||||||
CHECKSHAPE : Wires = 16 ( 16 ) Faces = 17 ( 19 ) Shells = 1 ( 1 ) Solids = 0 ( 0 )
|
CHECKSHAPE : Wires = 15 ( 16 ) Faces = 17 ( 17 ) Shells = 1 ( 1 ) Solids = 0 ( 0 )
|
||||||
NBSHAPES : Solid = 22 ( 22 ) Shell = 24 ( 24 ) Face = 1519 ( 1519 ) Summary = 11220 ( 11209 )
|
NBSHAPES : Solid = 22 ( 22 ) Shell = 24 ( 24 ) Face = 1519 ( 1519 ) Summary = 11224 ( 11212 )
|
||||||
STATSHAPE : Solid = 22 ( 22 ) Shell = 24 ( 24 ) Face = 1519 ( 1519 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 4793 ( 4785 )
|
STATSHAPE : Solid = 22 ( 22 ) Shell = 24 ( 24 ) Face = 1519 ( 1519 ) FreeWire = 0 ( 0 ) FreeEdge = 0 ( 0 ) SharedEdge = 4794 ( 4787 )
|
||||||
TOLERANCE : MaxTol = 7.159520237 ( 7.159520237 ) AvgTol = 0.03414100054 ( 0.03415856171 )
|
TOLERANCE : MaxTol = 7.499684301 ( 7.499684301 ) AvgTol = 0.03452487556 ( 0.03544461059 )
|
||||||
LABELS : N0Labels = 25 ( 25 ) N1Labels = 23 ( 23 ) N2Labels = 0 ( 0 ) TotalLabels = 48 ( 48 ) NameLabels = 48 ( 48 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
LABELS : N0Labels = 25 ( 25 ) N1Labels = 23 ( 23 ) N2Labels = 0 ( 0 ) TotalLabels = 48 ( 48 ) NameLabels = 48 ( 48 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
|
||||||
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
|
||||||
NCOLORS : NColors = 0 ( 0 )
|
NCOLORS : NColors = 0 ( 0 )
|
||||||
|
Loading…
x
Reference in New Issue
Block a user