Improve projection of ellipse and circle on a plane in case of the same parametrization of the original curve and the projected one is not necessary. Now the projection is a canonical curve instead of B-spline.
#include <Geom_Parabola.hxx>
#include <Geom_Hyperbola.hxx>
#include <Geom_Ellipse.hxx>
+#include <GeomLib_Tool.hxx>
+#include <math_Jacobi.hxx>
+#include <math_Matrix.hxx>
gp_Ax2 Axis;
Standard_Real R1 =0., R2 =0.;
- if ( Type != GeomAbs_Line) // on garde le parametrage
- myKeepParam = Standard_True;
- else // on prend le choix utilisateur.
- myKeepParam = KeepParametrization;
+ myKeepParam = KeepParametrization;
switch ( Type) {
case GeomAbs_Line:
Standard_Real Tol2 = myTolerance*myTolerance;
if (VDx.SquareMagnitude() < Tol2 ||
- VDy.SquareMagnitude() < Tol2 ) {
- myIsApprox = Standard_True;
+ VDy.SquareMagnitude() < Tol2 ||
+ VDx.CrossSquareMagnitude(VDy) < Tol2)
+ {
+ myIsApprox = Standard_True;
}
- if (!myIsApprox &&
- gp_Dir(VDx).IsNormal(gp_Dir(VDy),Precision::Angular())) {
+ if (!myIsApprox)
+ {
Dx = gp_Dir(VDx);
Dy = gp_Dir(VDy);
gp_Pnt O = Axis.Location();
gp_Pnt Py = ProjectPnt(myPlane,myDirection,O.Translated(R2*gp_Vec(Y)));
Standard_Real Major = P.Distance(Px);
Standard_Real Minor = P.Distance(Py);
- gp_Ax2 Axe( P, Dx^Dy,Dx);
- if ( Abs( Major - Minor) < Precision::Confusion()) {
- myType = GeomAbs_Circle;
- gp_Circ Circ(Axe, Major);
- GeomCirclePtr = new Geom_Circle(Circ);
-// Modified by Sergey KHROMOV - Tue Jan 29 16:57:29 2002 Begin
- GeomAdaptor_Curve aGACurve(GeomCirclePtr);
- myResult = new GeomAdaptor_HCurve(aGACurve);
-// Modified by Sergey KHROMOV - Tue Jan 29 16:57:30 2002 End
- }
- else if ( Major > Minor) {
- myType = GeomAbs_Ellipse;
- Elips = gp_Elips( Axe, Major, Minor);
-
- GeomEllipsePtr = new Geom_Ellipse(Elips) ;
+ if (myKeepParam)
+ {
+ myIsApprox = !gp_Dir(VDx).IsNormal(gp_Dir(VDy), Precision::Angular());
+ }
+ else
+ {
+ // Since it is not necessary to keep the same parameter for the point on the original and on the projected curves,
+ // we will use the following approach to find axes of the projected ellipse and provide the canonical curve:
+ // https://www.geometrictools.com/Documentation/ParallelProjectionEllipse.pdf
+ math_Matrix aMatrA(1, 2, 1, 2);
+ // A = Jp^T * Pr(Je), where
+ // Pr(Je) - projection of axes of original ellipse to the target plane
+ // Jp - X and Y axes of the target plane
+ aMatrA(1, 1) = myPlane.XDirection().XYZ().Dot(VDx.XYZ());
+ aMatrA(1, 2) = myPlane.XDirection().XYZ().Dot(VDy.XYZ());
+ aMatrA(2, 1) = myPlane.YDirection().XYZ().Dot(VDx.XYZ());
+ aMatrA(2, 2) = myPlane.YDirection().XYZ().Dot(VDy.XYZ());
+
+ math_Matrix aMatrDelta2(1, 2, 1, 2, 0.0);
+ // | 1/MajorRad^2 0 |
+ // Delta^2 = | |
+ // | 0 1/MajorRad^2 |
+ aMatrDelta2(1, 1) = 1.0 / (R1 * R1);
+ aMatrDelta2(2, 2) = 1.0 / (R2 * R2);
+
+ math_Matrix aMatrAInv = aMatrA.Inverse();
+ math_Matrix aMatrM = aMatrAInv.Transposed() * aMatrDelta2 * aMatrAInv;
+
+ // perform eigenvalues calculation
+ math_Jacobi anEigenCalc(aMatrM);
+ if (anEigenCalc.IsDone())
+ {
+ // radii of the projected ellipse
+ Minor = 1.0 / Sqrt(anEigenCalc.Value(1));
+ Major = 1.0 / Sqrt(anEigenCalc.Value(2));
+
+ // calculate the rotation angle for the plane axes to meet the correct axes of the projected ellipse
+ // (swap eigenvectors in respect to major and minor axes)
+ const math_Matrix& anEigenVec = anEigenCalc.Vectors();
+ gp_Trsf2d aTrsfInPlane;
+ aTrsfInPlane.SetValues(anEigenVec(1, 2), anEigenVec(1, 1), 0.0,
+ anEigenVec(2, 2), anEigenVec(2, 1), 0.0);
+ gp_Trsf aRot;
+ aRot.SetRotation(gp_Ax1(P, myPlane.Direction()), aTrsfInPlane.RotationPart());
+
+ Dx = myPlane.XDirection().Transformed(aRot);
+ Dy = myPlane.YDirection().Transformed(aRot);
+ }
+ else
+ {
+ myIsApprox = Standard_True;
+ }
+ }
+
+ if (!myIsApprox)
+ {
+ gp_Ax2 Axe(P, Dx^Dy, Dx);
+
+ if (Abs(Major - Minor) < Precision::Confusion()) {
+ myType = GeomAbs_Circle;
+ gp_Circ Circ(Axe, Major);
+ GeomCirclePtr = new Geom_Circle(Circ);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:29 2002 Begin
- GeomAdaptor_Curve aGACurve(GeomEllipsePtr);
- myResult = new GeomAdaptor_HCurve(aGACurve);
+ GeomAdaptor_Curve aGACurve(GeomCirclePtr);
+ myResult = new GeomAdaptor_HCurve(aGACurve);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:30 2002 End
- }
- else {
- myIsApprox = Standard_True;
- myType = GeomAbs_BSplineCurve;
- PerformApprox(myCurve,myPlane,myDirection,ApproxCurve);
+ }
+ else if ( Major > Minor) {
+ myType = GeomAbs_Ellipse;
+ Elips = gp_Elips( Axe, Major, Minor);
+
+ GeomEllipsePtr = new Geom_Ellipse(Elips);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:29 2002 Begin
- GeomAdaptor_Curve aGACurve(ApproxCurve);
- myResult = new GeomAdaptor_HCurve(aGACurve);
+ GeomAdaptor_Curve aGACurve(GeomEllipsePtr);
+ myResult = new GeomAdaptor_HCurve(aGACurve);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:30 2002 End
- }
+ }
+ else {
+ myIsApprox = Standard_True;
+ }
+ }
}
- else {
- myIsApprox = Standard_True;
+
+ // No way to build the canonical curve, approximate as B-spline
+ if (myIsApprox)
+ {
myType = GeomAbs_BSplineCurve;
PerformApprox(myCurve,myPlane,myDirection,ApproxCurve);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:29 2002 Begin
myResult = new GeomAdaptor_HCurve(aGACurve);
// Modified by Sergey KHROMOV - Tue Jan 29 16:57:30 2002 End
}
+ else if (GeomCirclePtr || GeomEllipsePtr)
+ {
+ Handle(Geom_Curve) aResultCurve = GeomCirclePtr;
+ if (aResultCurve.IsNull())
+ aResultCurve = GeomEllipsePtr;
+ // start and end parameters of the projected curve
+ Standard_Real aParFirst = myCurve->FirstParameter();
+ Standard_Real aParLast = myCurve->LastParameter();
+ gp_Pnt aPntFirst = ProjectPnt(myPlane, myDirection, myCurve->Value(aParFirst));
+ gp_Pnt aPntLast = ProjectPnt(myPlane, myDirection, myCurve->Value(aParLast));
+ GeomLib_Tool::Parameter(aResultCurve, aPntFirst, Precision::Confusion(), myFirstPar);
+ GeomLib_Tool::Parameter(aResultCurve, aPntLast, Precision::Confusion(), myLastPar);
+ while (myLastPar <= myFirstPar)
+ myLastPar += myResult->Period();
+ }
}
break;
case GeomAbs_Parabola:
{
+ myKeepParam = Standard_True;
// P(u) = O + (u*u)/(4*f) * Xc + u * Yc
// ==> Q(u) = f(P(u))
// = f(O) + (u*u)/(4*f) * f(Xc) + u * f(Yc)
break;
case GeomAbs_Hyperbola:
{
+ myKeepParam = Standard_True;
// P(u) = O + R1 * Cosh(u) * Xc + R2 * Sinh(u) * Yc
// ==> Q(u) = f(P(u))
// = f(O) + R1 * Cosh(u) * f(Xc) + R2 * Sinh(u) * f(Yc)
Handle(Geom_BezierCurve) ProjCu =
Handle(Geom_BezierCurve)::DownCast(BezierCurvePtr->Copy());
+ myKeepParam = Standard_True;
myIsApprox = Standard_False;
myType = Type;
for ( Standard_Integer i = 1; i <= NbPoles; i++) {
Handle(Geom_BSplineCurve) ProjectedBSplinePtr =
Handle(Geom_BSplineCurve)::DownCast(BSplineCurvePtr->Copy()) ;
+ myKeepParam = Standard_True;
myIsApprox = Standard_False;
myType = Type;
for ( Standard_Integer i = 1; i <= BSplineCurvePtr->NbPoles(); i++) {
break;
default:
{
+ myKeepParam = Standard_True;
myIsApprox = Standard_True;
myType = GeomAbs_BSplineCurve;
PerformApprox(myCurve,myPlane,myDirection,ApproxCurve);
uplevel #0 offsetperform result
}
+proc ProjectCurvePointToPlaneAlongDir {curve param pln {dir {}}} {
+ upvar $pln p
+ upvar $curve c
+ cvalue c $param x y z
+ if {[llength $dir] == 0 } {
+ # project to plane along the normal
+ regexp {Axis :([-0-9.+eE]+), ([-0-9.+eE]+), ([-0-9.+eE]+)} [dump p] full dx dy dz
+ lappend dir $dx $dy $dz
+ }
+ line ln x y z [lindex $dir 0] [lindex $dir 1] [lindex $dir 2]
+ intersect pt ln p
+ regexp {Point : ([-0-9.+eE]+), ([-0-9.+eE]+), ([-0-9.+eE]+)} [dump pt] full x y z
+
+ set pntOnPlane {}
+ lappend pntOnPlane $x $y $z
+ return $pntOnPlane
+}
-
-
-
-
-
-
-
+proc CheckProjectionToPlane {nbSamples origCurve origParam0 origParam1 projCurve projParam0 projParam1 pln {dir {}} {tolerance 1.e-7}} {
+ upvar $pln p
+ upvar $origCurve origC
+ upvar $projCurve projC
+
+ set isOk 1
+ for {set i 0} {$i <= $nbSamples} {incr i} {
+ set parOrig [expr $origParam0 + ($origParam1 - $origParam0) * $i / $nbSamples]
+ set parProj [expr $projParam0 + ($projParam1 - $projParam0) * $i / $nbSamples]
+
+ set pnt [ProjectCurvePointToPlaneAlongDir origC $parOrig p $dir]
+ cvalue projC $parProj X Y Z
+
+ set dx [expr [lindex $pnt 0]-[dval X]]
+ set dy [expr [lindex $pnt 1]-[dval Y]]
+ set dz [expr [lindex $pnt 2]-[dval Z]]
+
+ if {[expr $dx*$dx + $dy*$dy + $dz*$dz] < [expr $tolerance*$tolerance]} {
+ puts "OK: Projection correct"
+ } else {
+ puts "ERROR: Projection incorrect"
+ set isOk 0
+ }
+ }
+ return $isOk
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 1 10
+projonplane r c p 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 1 10
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p $projDir]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+circle c 0 0 0 10
+plane p 0 0 0 1 0 1
+projonplane r c p 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+circle c 0 0 0 10
+plane p 0 0 0 1 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p $projDir]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+ellipse c 0 0 0 0.866025403784439 0 0.5 0.5 0 -0.866025403784439 20 10
+plane p 0 0 0 0 0 1
+projonplane r c p 0
+
+if {![regexp {Circle} [dump r]]} {
+ puts "ERROR: Projected curve is not a circle"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c 0 p $projDir]
+parameters r [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r [dval shift] [expr [dval shift]+2*[dval pi]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 1 10
+projonplane r c p 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 1 10
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+circle c 0 0 0 10
+plane p 0 0 0 1 0 1
+projonplane r c p 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+circle c 0 0 0 10
+plane p 0 0 0 1 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+ellipse c 0 0 0 0.866025403784439 0 0.5 0.5 0 -0.866025403784439 20 10
+plane p 0 0 0 0 0 1
+projonplane r c p 1
+
+if {![regexp {Circle} [dump r]]} {
+ puts "ERROR: Projected curve is not a circle"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+plane p 0 0 0 0 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+if {[CheckProjectionToPlane 100 c 0 [expr 2*[dval pi]] r 0 [expr 2*[dval pi]] p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 1 10
+projonplane r c p 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 1 10
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p $projDir]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+circle c 0 0 0 10
+trim c c $startPar $endPar
+plane p 0 0 0 1 0 1
+projonplane r c p 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+circle c 0 0 0 10
+trim c c $startPar $endPar
+plane p 0 0 0 1 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p $projDir]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0.866025403784439 0 0.5 0.5 0 -0.866025403784439 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 0 1
+projonplane r c p 0
+
+if {![regexp {Circle} [dump r]]} {
+ puts "ERROR: Projected curve is not a circle"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 0
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+# calculate a parametric shift on the projected curve
+set pnt [ProjectCurvePointToPlaneAlongDir c $startPar p $projDir]
+trim ru r
+parameters ru [lindex $pnt 0] [lindex $pnt 1] [lindex $pnt 2] 0.1 shift
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r [dval shift] [expr $endPar-$startPar+[dval shift]] p $projDir]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 1 10
+projonplane r c p 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 1 10
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+circle c 0 0 0 10
+trim c c $startPar $endPar
+plane p 0 0 0 1 0 1
+projonplane r c p 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+circle c 0 0 0 10
+trim c c $startPar $endPar
+plane p 0 0 0 1 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {BSplineCurve} [dump r]]} {
+ puts "ERROR: Projected curve is not a B-spline curve"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0.866025403784439 0 0.5 0.5 0 -0.866025403784439 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 0 1
+projonplane r c p 1
+
+if {![regexp {Circle} [dump r]]} {
+ puts "ERROR: Projected curve is not a circle"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p {} 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}
--- /dev/null
+puts ""
+puts "=========================================================================="
+puts "OCC31016: Projection of an ellipse or a circle is a B-spline in some cases"
+puts "=========================================================================="
+puts ""
+
+set projDir { 1 1 1 }
+set startPar 1
+set endPar 4.5
+
+ellipse c 0 0 0 0 0 1 2 1 0 20 10
+trim c c $startPar $endPar
+plane p 0 0 0 0 0 1
+projonplane r c p [lindex $projDir 0] [lindex $projDir 1] [lindex $projDir 2] 1
+
+if {![regexp {Ellipse} [dump r]]} {
+ puts "ERROR: Projected curve is not an ellipse"
+}
+
+if {[CheckProjectionToPlane 100 c $startPar $endPar r $startPar $endPar p $projDir 1.e-6]} {
+ puts ""
+ puts "OK: All sample points are projected correctly"
+ puts ""
+} else {
+ puts ""
+ puts "ERROR: Projection is incorrect for some points"
+ puts ""
+}