1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-04-05 18:16:23 +03:00

0023620: Follow up of 0022939 - make Bezier curve/surface evaluation thread-safe

1. Remove cache from Geom_BezierCurve, Geom2d_BezierCurve and Geom_BezierSurface
2. Add cache for Bezier curves into GeomAdaptor_Curve, Geom2dAdaptor_Curve and GeomAdaptor_Surface
3. Update comments in corresponding cache classes
4. Avoid frequent down-casting to B-splines in adaptors
This commit is contained in:
azv 2015-11-02 09:33:04 +03:00 committed by bugmaster
parent 525ec87c53
commit c8b5b3d89e
11 changed files with 555 additions and 883 deletions

View File

@ -31,9 +31,9 @@
#include <TColgp_Array1OfPnt.hxx>
#include <TColgp_Array1OfPnt2d.hxx>
//! \brief A cache class for B-spline curves.
//! \brief A cache class for Bezier and B-spline curves.
//!
//! Defines all data, that can be cached on a span of B-spline curve.
//! Defines all data, that can be cached on a span of a curve.
//! The data should be recalculated in going from span to span.
class BSplCLib_Cache : public Standard_Transient
{
@ -41,10 +41,10 @@ public:
//! Default constructor
Standard_EXPORT BSplCLib_Cache();
//! Constructor for caching of 2D curves
//! \param theDegree degree of the B-spline
//! \param thePeriodic identify the B-spline is periodic
//! \param theFlatKnots knots of B-spline curve (with repetitions)
//! \param thePoles2d array of poles of 2D B-spline
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param theFlatKnots knots of Bezier/B-spline curve (with repetitions)
//! \param thePoles2d array of poles of 2D curve
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT BSplCLib_Cache(const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
@ -52,10 +52,10 @@ public:
const TColgp_Array1OfPnt2d& thePoles2d,
const TColStd_Array1OfReal* theWeights = NULL);
//! Constructor for caching of 3D curves
//! \param theDegree degree of the B-spline
//! \param thePeriodic identify the B-spline is periodic
//! \param theFlatKnots knots of B-spline curve (with repetitions)
//! \param thePoles array of poles of 3D B-spline
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param theFlatKnots knots of Bezier/B-spline curve (with repetitions)
//! \param thePoles array of poles of 3D curve
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT BSplCLib_Cache(const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
@ -69,10 +69,10 @@ public:
//! Recomputes the cache data for 2D curves. Does not verify validity of the cache
//! \param theParameter the value on the knot's axis to identify the span
//! \param theDegree degree of the B-spline
//! \param thePeriodic identify the B-spline is periodic
//! \param theFlatKnots knots of B-spline curve (with repetitions)
//! \param thePoles2d array of poles of 2D B-spline
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param theFlatKnots knots of Bezier/B-spline curve (with repetitions)
//! \param thePoles2d array of poles of 2D curve
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT void BuildCache(const Standard_Real& theParameter,
const Standard_Integer& theDegree,
@ -82,10 +82,10 @@ public:
const TColStd_Array1OfReal* theWeights = NULL);
//! Recomputes the cache data for 3D curves. Does not verify validity of the cache
//! \param theParameter the value on the knot's axis to identify the span
//! \param theDegree degree of the B-spline
//! \param thePeriodic identify the B-spline is periodic
//! \param theFlatKnots knots of B-spline curve (with repetitions)
//! \param thePoles array of poles of 3D B-spline
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param theFlatKnots knots of Bezier/B-spline curve (with repetitions)
//! \param thePoles array of poles of 3D curve
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT void BuildCache(const Standard_Real& theParameter,
const Standard_Integer& theDegree,
@ -94,24 +94,24 @@ public:
const TColgp_Array1OfPnt& thePoles,
const TColStd_Array1OfReal* theWeights = NULL);
//! Calculates the point on B-spline in the selected point
//! Calculates the point on the curve in the specified parameter
//! \param[in] theParameter parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on B-spline)
//! \param[out] thePoint the result of calculation (the point on the curve)
Standard_EXPORT void D0(const Standard_Real& theParameter, gp_Pnt2d& thePoint) const;
Standard_EXPORT void D0(const Standard_Real& theParameter, gp_Pnt& thePoint) const;
//! Calculates the point on B-spline and its first derivative in the selected point
//! Calculates the point on the curve and its first derivative in the specified parameter
//! \param[in] theParameter parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on B-spline)
//! \param[out] theTangent tangent vector (first derivatives) for B-spline in the calculated point
//! \param[out] thePoint the result of calculation (the point on the curve)
//! \param[out] theTangent tangent vector (first derivatives) for the curve in the calculated point
Standard_EXPORT void D1(const Standard_Real& theParameter, gp_Pnt2d& thePoint, gp_Vec2d& theTangent) const;
Standard_EXPORT void D1(const Standard_Real& theParameter, gp_Pnt& thePoint, gp_Vec& theTangent) const;
//! Calculates the point on B-spline and two derivatives in the selected point
//! Calculates the point on the curve and two derivatives in the specified parameter
//! \param[in] theParameter parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on B-spline)
//! \param[out] theTangent tangent vector (1st derivatives) for B-spline in the calculated point
//! \param[out] theCurvature curvature vector (2nd derivatives) for B-spline in the calculated point
//! \param[out] thePoint the result of calculation (the point on the curve)
//! \param[out] theTangent tangent vector (1st derivatives) for the curve in the calculated point
//! \param[out] theCurvature curvature vector (2nd derivatives) for the curve in the calculated point
Standard_EXPORT void D2(const Standard_Real& theParameter,
gp_Pnt2d& thePoint,
gp_Vec2d& theTangent,
@ -121,12 +121,12 @@ public:
gp_Vec& theTangent,
gp_Vec& theCurvature) const;
//! Calculates the point on B-spline and three derivatives in the selected point
//! Calculates the point on the curve and three derivatives in the specified parameter
//! \param[in] theParameter parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on B-spline)
//! \param[out] theTangent tangent vector (1st derivatives) for B-spline in the calculated point
//! \param[out] theCurvature curvature vector (2nd derivatives) for B-spline in the calculated point
//! \param[out] theTorsion second curvature vector (3rd derivatives) for B-spline in the calculated point
//! \param[out] thePoint the result of calculation (the point on the curve)
//! \param[out] theTangent tangent vector (1st derivatives) for the curve in the calculated point
//! \param[out] theCurvature curvature vector (2nd derivatives) for the curve in the calculated point
//! \param[out] theTorsion second curvature vector (3rd derivatives) for the curve in the calculated point
Standard_EXPORT void D3(const Standard_Real& theParameter,
gp_Pnt2d& thePoint,
gp_Vec2d& theTangent,
@ -142,16 +142,16 @@ public:
DEFINE_STANDARD_RTTI(BSplCLib_Cache, Standard_Transient)
protected:
//! Normalizes the parameter for periodical B-splines
//! Normalizes the parameter for periodical curves
//! \param theFlatKnots knots with repetitions
//! \param theParameter the value to be normalized into the knots array
void PeriodicNormalization(const TColStd_Array1OfReal& theFlatKnots, Standard_Real& theParameter) const;
//! Fills array of derivatives in the selected point of the B-spline
//! Fills array of derivatives in the selected point of the curve
//! \param[in] theParameter parameter of the calculation
//! \param[in] theDerivative maximal derivative to be calculated (computes all derivatives lesser than specified)
//! \param[out] theDerivArray result array of derivatives (with size (theDerivative+1)*(PntDim+1),
//! where PntDim = 2 or 3 is a dimension of B-spline curve)
//! where PntDim = 2 or 3 is a dimension of the curve)
void CalculateDerivative(const Standard_Real& theParameter,
const Standard_Integer& theDerivative,
Standard_Real& theDerivArray) const;
@ -163,13 +163,13 @@ private:
// x2 y2 [z2] [w2] etc
// for 2D-curves there is no z conponent, for non-rational curves there is no weight
Standard_Boolean myIsRational; ///< identifies the rationality of B-spline
Standard_Boolean myIsRational; ///< identifies the rationality of Bezier/B-spline curve
Standard_Real mySpanStart; ///< parameter for the first point of the span
Standard_Real mySpanLength; ///< length of the span
Standard_Integer mySpanIndex; ///< index of the span on B-spline curve
Standard_Integer mySpanIndexMax; ///< maximal number of spans on B-spline curve
Standard_Integer myDegree; ///< degree of B-spline
Handle(TColStd_HArray1OfReal) myFlatKnots; ///< knots of B-spline (used for periodic normalization of parameters, exists only for periodical splines)
Standard_Integer mySpanIndex; ///< index of the span on Bezier/B-spline curve
Standard_Integer mySpanIndexMax; ///< maximal number of spans on Bezier/B-spline curve
Standard_Integer myDegree; ///< degree of Bezier/B-spline
Handle(TColStd_HArray1OfReal) myFlatKnots; ///< knots of Bezier/B-spline (used for periodic normalization of parameters, exists only for periodical splines)
};
DEFINE_STANDARD_HANDLE(BSplCLib_Cache, Standard_Transient)

View File

@ -29,23 +29,23 @@
#include <TColStd_Array1OfReal.hxx>
#include <TColStd_Array2OfReal.hxx>
//! \brief A cache class for B-spline surfaces.
//! \brief A cache class for Bezier and B-spline surfaces.
//!
//! Defines all data, that can be cached on a span of B-spline surface.
//! Defines all data, that can be cached on a span of the surface.
//! The data should be recalculated in going from span to span.
class BSplSLib_Cache : public Standard_Transient
{
public:
//! Default constructor
Standard_EXPORT BSplSLib_Cache();
//! Constructor for caching of the span for B-spline surface
//! \param theDegreeU degree along the first parameter (U) of the B-spline
//! \param thePeriodicU identify the B-spline is periodical along U axis
//! \param theFlatKnotsU knots of B-spline curve (with repetition) along U axis
//! \param theDegreeV degree alogn the second parameter (V) of the B-spline
//! \param thePeriodicV identify the B-spline is periodical along V axis
//! \param theFlatKnotsV knots of B-spline curve (with repetition) along V axis
//! \param thePoles array of poles of the B-spline surface
//! Constructor for caching of the span for the surface
//! \param theDegreeU degree along the first parameter (U) of the surface
//! \param thePeriodicU identify the surface is periodical along U axis
//! \param theFlatKnotsU knots of the surface (with repetition) along U axis
//! \param theDegreeV degree alogn the second parameter (V) of the surface
//! \param thePeriodicV identify the surface is periodical along V axis
//! \param theFlatKnotsV knots of the surface (with repetition) along V axis
//! \param thePoles array of poles of the surface
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT BSplSLib_Cache(const Standard_Integer& theDegreeU,
const Standard_Boolean& thePeriodicU,
@ -65,13 +65,13 @@ public:
//! Recomputes the cache data. Does not verify validity of the cache
//! \param theParameterU the parametric value on the U axis to identify the span
//! \param theParameterV the parametric value on the V axis to identify the span
//! \param theDegreeU degree of the B-spline along U axis
//! \param thePeriodicU identify the B-spline is periodic along U axis
//! \param theFlatKnotsU flat knots of B-spline surface along U axis
//! \param theDegreeV degree of the B-spline along V axis
//! \param thePeriodicV identify the B-spline is periodic along V axis
//! \param theFlatKnotsV flat knots of B-spline surface along V axis
//! \param thePoles array of poles of B-spline
//! \param theDegreeU degree along U axis
//! \param thePeriodicU identify whether the surface is periodic along U axis
//! \param theFlatKnotsU flat knots of the surface along U axis
//! \param theDegreeV degree along V axis
//! \param thePeriodicV identify whether the surface is periodic along V axis
//! \param theFlatKnotsV flat knots of the surface along V axis
//! \param thePoles array of poles of the surface
//! \param theWeights array of weights of corresponding poles
Standard_EXPORT void BuildCache(const Standard_Real& theParameterU,
const Standard_Real& theParameterV,
@ -84,16 +84,16 @@ public:
const TColgp_Array2OfPnt& thePoles,
const TColStd_Array2OfReal* theWeights = NULL);
//! Calculates the point on B-spline for specified parameters
//! Calculates the point on the surface for specified parameters
//! \param[in] theU first parameter for calculation of the value
//! \param[in] theV second parameter for calculation of the value
//! \param[out] thePoint the result of calculation (the point on the B-spline)
//! \param[out] thePoint the result of calculation (the point on the surface)
Standard_EXPORT void D0(const Standard_Real& theU, const Standard_Real& theV, gp_Pnt& thePoint) const;
//! Calculates the point on B-spline and its first derivative
//! Calculates the point on the surface and its first derivative
//! \param[in] theU first parameter of calculation of the value
//! \param[in] theV second parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on the B-spline)
//! \param[out] thePoint the result of calculation (the point on the surface)
//! \param[out] theTangentU tangent vector along U axis in the calculated point
//! \param[out] theTangentV tangent vector along V axis in the calculated point
Standard_EXPORT void D1(const Standard_Real& theU,
@ -102,10 +102,10 @@ public:
gp_Vec& theTangentU,
gp_Vec& theTangentV) const;
//! Calculates the point on B-spline and derivatives till second order
//! Calculates the point on the surface and derivatives till second order
//! \param[in] theU first parameter of calculation of the value
//! \param[in] theV second parameter of calculation of the value
//! \param[out] thePoint the result of calculation (the point on B-spline)
//! \param[out] thePoint the result of calculation (the point on the surface)
//! \param[out] theTangentU tangent vector along U axis in the calculated point
//! \param[out] theTangentV tangent vector along V axis in the calculated point
//! \param[out] theCurvatureU curvature vector (2nd derivative on U) along U axis
@ -124,8 +124,8 @@ public:
DEFINE_STANDARD_RTTI(BSplSLib_Cache, Standard_Transient)
protected:
//! Normalizes the parameter for periodical B-splines
//! \param[in] theDegree degree of B-spline along selected direction
//! Normalizes the parameter for periodical surfaces
//! \param[in] theDegree degree along selected direction
//! \param[in] theFlatKnots knots with repetitions along selected direction
//! \param[in,out] theParameter the value to be normalized into the knots array
void PeriodicNormalization(const Standard_Integer& theDegree,
@ -140,13 +140,13 @@ private:
// for non-rational surfaces there is no weight;
// size of array: (max(myDegree)+1) * A*(min(myDegree)+1), where A = 4 or 3
Standard_Boolean myIsRational; ///< identifies the rationality of B-spline
Standard_Boolean myIsRational; ///< identifies the rationality of Bezier/B-spline surface
Standard_Real mySpanStart[2]; ///< parameters (u, v) for the frst point of the span
Standard_Real mySpanLength[2]; ///< lengths of the span along corresponding parameter
Standard_Integer mySpanIndex[2]; ///< indexes of the span on B-spline surface
Standard_Integer mySpanIndex[2]; ///< indexes of the span on Bezier/B-spline surface
Standard_Integer mySpanIndexMax[2]; ///< maximal indexes of span
Standard_Integer myDegree[2]; ///< degrees of B-spline for each parameter
Handle(TColStd_HArray1OfReal) myFlatKnots[2]; ///< arrays of knots of B-spline
Standard_Integer myDegree[2]; ///< degrees of Bezier/B-spline for each parameter
Handle(TColStd_HArray1OfReal) myFlatKnots[2]; ///< arrays of knots of Bezier/B-spline
// (used for periodic normalization of parameters, Null for non-periodical splines)
};

View File

@ -27,7 +27,6 @@
#define No_Standard_DimensionError
#include <BSplCLib.hxx>
#include <Geom_BezierCurve.hxx>
#include <Geom_Geometry.hxx>
#include <gp.hxx>
@ -64,8 +63,7 @@ static Standard_Boolean Rational(const TColStd_Array1OfReal& W)
//purpose :
//=======================================================================
Geom_BezierCurve::Geom_BezierCurve(const TColgp_Array1OfPnt& Poles):
validcache(0), parametercache(0.), spanlenghtcache(1.)
Geom_BezierCurve::Geom_BezierCurve(const TColgp_Array1OfPnt& Poles)
{
Standard_Integer nbpoles = Poles.Length();
if(nbpoles < 2 || nbpoles > (Geom_BezierCurve::MaxDegree() + 1))
@ -87,8 +85,7 @@ validcache(0), parametercache(0.), spanlenghtcache(1.)
//=======================================================================
Geom_BezierCurve::Geom_BezierCurve(const TColgp_Array1OfPnt& Poles,
const TColStd_Array1OfReal& Weights):
validcache(0), parametercache(0.), spanlenghtcache(1.)
const TColStd_Array1OfReal& Weights)
{
// copy the poles
Standard_Integer nbpoles = Poles.Length();
@ -360,8 +357,6 @@ void Geom_BezierCurve::Reverse ()
cweights(nbpoles-i+1) = w;
}
}
UpdateCoefficients();
}
//=======================================================================
@ -383,21 +378,21 @@ void Geom_BezierCurve::Segment(const Standard_Real U1, const Standard_Real U2)
{
closed = (Abs(Value(U1).Distance (Value(U2))) <= Precision::Confusion());
if(!CoefficientsOK(0.)) UpdateCoefficients(0.);
TColStd_Array1OfReal bidflatknots(BSplCLib::FlatBezierKnots(Degree()), 1, 2 * (Degree() + 1));
TColgp_HArray1OfPnt coeffs(1, poles->Size());
if (IsRational()) {
PLib::Trimming(U1,U2,coeffs->ChangeArray1(),
&wcoeffs->ChangeArray1());
PLib::CoefficientsPoles(coeffs->Array1(),
&wcoeffs->Array1(),
poles->ChangeArray1(),
&weights->ChangeArray1());
TColStd_Array1OfReal wcoeffs(1, poles->Size());
BSplCLib::BuildCache(0.0, 1.0, 0, Degree(), bidflatknots,
poles->Array1(), &weights->Array1(), coeffs, &wcoeffs);
PLib::Trimming(U1, U2, coeffs, &wcoeffs);
PLib::CoefficientsPoles(coeffs, &wcoeffs, poles->ChangeArray1(), &weights->ChangeArray1());
}
else {
PLib::Trimming(U1,U2,coeffs->ChangeArray1(), PLib::NoWeights());
PLib::CoefficientsPoles(coeffs->Array1(), PLib::NoWeights(),
poles->ChangeArray1(), PLib::NoWeights());
BSplCLib::BuildCache(0.0, 1.0, 0, Degree(), bidflatknots,
poles->Array1(), BSplCLib::NoWeights(), coeffs, BSplCLib::NoWeights());
PLib::Trimming(U1, U2, coeffs, PLib::NoWeights());
PLib::CoefficientsPoles(coeffs, PLib::NoWeights(), poles->ChangeArray1(), PLib::NoWeights());
}
UpdateCoefficients();
}
//=======================================================================
@ -417,7 +412,6 @@ void Geom_BezierCurve::SetPole (const Standard_Integer Index,
if (Index == 1 || Index == cpoles.Length()) {
closed = (cpoles(1).Distance(cpoles(NbPoles())) <= Precision::Confusion());
}
UpdateCoefficients();
}
//=======================================================================
@ -456,7 +450,6 @@ void Geom_BezierCurve::SetWeight(const Standard_Integer Index,
// set weights of 1.
weights = new TColStd_HArray1OfReal(1,nbpoles);
wcoeffs = new TColStd_HArray1OfReal(1,nbpoles);
weights->Init(1.);
}
@ -464,13 +457,8 @@ void Geom_BezierCurve::SetWeight(const Standard_Integer Index,
cweights(Index) = Weight;
// is it turning into non rational
if (wasrat) {
if (!Rational(cweights)) {
weights.Nullify();
wcoeffs.Nullify();
}
}
UpdateCoefficients();
if (wasrat && !Rational(cweights))
weights.Nullify();
}
//=======================================================================
@ -540,18 +528,7 @@ Standard_Integer Geom_BezierCurve::Degree () const
void Geom_BezierCurve::D0 (const Standard_Real U, gp_Pnt& P ) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD0(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
&wcoeffs->Array1(),
P);
else
BSplCLib::CacheD0(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(), BSplCLib::NoWeights(), P);
BSplCLib::D0(U, Poles(), Weights(), P);
}
//=======================================================================
@ -561,20 +538,7 @@ void Geom_BezierCurve::D0 (const Standard_Real U, gp_Pnt& P ) const
void Geom_BezierCurve::D1(const Standard_Real U, gp_Pnt& P, gp_Vec& V1) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD1(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
&wcoeffs->Array1(),
P,V1);
else
BSplCLib::CacheD1(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
BSplCLib::NoWeights(),
P,V1);
BSplCLib::D1(U, Poles(), Weights(), P, V1);
}
//=======================================================================
@ -587,20 +551,7 @@ void Geom_BezierCurve::D2 (const Standard_Real U,
gp_Vec& V1,
gp_Vec& V2) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD2(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
&wcoeffs->Array1(),
P,V1,V2);
else
BSplCLib::CacheD2(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
BSplCLib::NoWeights(),
P,V1,V2);
BSplCLib::D2(U, Poles(), Weights(), P, V1, V2);
}
//=======================================================================
@ -614,18 +565,7 @@ void Geom_BezierCurve::D3 (const Standard_Real U,
gp_Vec& V2,
gp_Vec& V3) const
{
if(!CoefficientsOK(U))
((Geom_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD3(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
&wcoeffs->Array1(),
P,V1,V2,V3);
else
BSplCLib::CacheD3(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),
BSplCLib::NoWeights(),
P,V1,V2,V3);
BSplCLib::D3(U, Poles(), Weights(), P, V1, V2, V3);
}
//=======================================================================
@ -792,8 +732,6 @@ void Geom_BezierCurve::Transform (const gp_Trsf& T)
for (Standard_Integer i = 1; i <= nbpoles; i++)
cpoles (i).Transform(T);
UpdateCoefficients();
}
//=======================================================================
@ -865,58 +803,10 @@ void Geom_BezierCurve::Init
// set fields
poles = Poles;
coeffs = new TColgp_HArray1OfPnt (1,nbpoles);
if (rational) {
if (rational)
weights = Weights;
wcoeffs = new TColStd_HArray1OfReal (1, nbpoles, 0.0);
}
else {
weights.Nullify();
wcoeffs.Nullify();
}
UpdateCoefficients();
}
//=======================================================================
//function : CoefficientsOK
//purpose :
//=======================================================================
Standard_Boolean Geom_BezierCurve::CoefficientsOK(const Standard_Real U)const
{
return (validcache && ((parametercache == 0. && U < 1.) ||
(parametercache == 1. && U >= 1.)));
}
//=======================================================================
//function : UpdateCoefficients
//purpose :
//=======================================================================
//void Geom_BezierCurve::UpdateCoefficients(const Standard_Real U)
void Geom_BezierCurve::UpdateCoefficients(const Standard_Real )
{
maxderivinvok = 0;
parametercache = 0.;
//
// Idee lumineuse sacrifiee sur l autel des performances.
// if (U >= 1.) parametercache = 1.;
TColStd_Array1OfReal bidflatknots(BSplCLib::FlatBezierKnots(Degree()),
1, 2*(Degree()+1));
if (IsRational())
BSplCLib::BuildCache(parametercache,spanlenghtcache,0,Degree(),
bidflatknots,poles->Array1(),
&weights->Array1(),
coeffs->ChangeArray1(),
&wcoeffs->ChangeArray1());
else
BSplCLib::BuildCache(parametercache,spanlenghtcache,0,Degree(),
bidflatknots,poles->Array1(),
BSplCLib::NoWeights(),
coeffs->ChangeArray1(),
BSplCLib::NoWeights());
validcache = 1;
weights.Nullify();
}

View File

@ -29,6 +29,8 @@
#include <TColgp_Array1OfPnt.hxx>
#include <TColStd_Array1OfReal.hxx>
#include <GeomAbs_Shape.hxx>
#include <BSplCLib.hxx>
class Standard_ConstructionError;
class Standard_DimensionError;
class Standard_RangeError;
@ -304,7 +306,15 @@ public:
//!
//! Raised if the length of W is not equal to the number of poles.
Standard_EXPORT void Weights (TColStd_Array1OfReal& W) const;
//! Returns all the weights of the curve.
const TColStd_Array1OfReal* Weights() const
{
if (!weights.IsNull())
return &weights->Array1();
return BSplCLib::NoWeights();
}
//! Applies the transformation T to this Bezier curve.
Standard_EXPORT void Transform (const gp_Trsf& T) Standard_OVERRIDE;
@ -344,25 +354,12 @@ private:
//! Update rational and closed.
//!
//! if nbpoles < 2 or nbboles > MaDegree + 1
Standard_EXPORT void Init (const Handle(TColgp_HArray1OfPnt)& Poles, const Handle(TColStd_HArray1OfReal)& Weights);
//! returns true if the coefficients have been
//! computed with the right value of cacheparameter
//! for the given U value.
Standard_EXPORT Standard_Boolean CoefficientsOK (const Standard_Real U) const;
//! Recompute the coeficients.
Standard_EXPORT void UpdateCoefficients (const Standard_Real U = 0.0);
void Init (const Handle(TColgp_HArray1OfPnt)& Poles, const Handle(TColStd_HArray1OfReal)& Weights);
Standard_Boolean rational;
Standard_Boolean closed;
Handle(TColgp_HArray1OfPnt) poles;
Handle(TColStd_HArray1OfReal) weights;
Handle(TColgp_HArray1OfPnt) coeffs;
Handle(TColStd_HArray1OfReal) wcoeffs;
Standard_Integer validcache;
Standard_Real parametercache;
Standard_Real spanlenghtcache;
Standard_Real maxderivinv;
Standard_Boolean maxderivinvok;

View File

@ -28,7 +28,6 @@
#include <BSplCLib.hxx>
#include <BSplSLib.hxx>
#include <Geom_BezierCurve.hxx>
#include <Geom_BezierSurface.hxx>
#include <Geom_Curve.hxx>
@ -362,11 +361,6 @@ static void DeleteRatPoleRow
Geom_BezierSurface::Geom_BezierSurface
(const TColgp_Array2OfPnt& SurfacePoles):
ucacheparameter(0.),
vcacheparameter(0.),
ucachespanlenght(1.),
vcachespanlenght(1.),
validcache(0),
maxderivinvok(Standard_False)
{
Standard_Integer NbUPoles = SurfacePoles.ColLength();
@ -397,11 +391,6 @@ Geom_BezierSurface::Geom_BezierSurface
Geom_BezierSurface::Geom_BezierSurface
(const TColgp_Array2OfPnt& SurfacePoles,
const TColStd_Array2OfReal& PoleWeights ):
ucacheparameter(0.),
vcacheparameter(0.),
ucachespanlenght(1.),
vcachespanlenght(1.),
validcache(0),
maxderivinvok(Standard_False)
{
Standard_Integer NbUPoles = SurfacePoles.ColLength();
@ -472,20 +461,13 @@ Geom_BezierSurface::Geom_BezierSurface
Geom_BezierSurface::Geom_BezierSurface
(const Handle(TColgp_HArray2OfPnt)& SurfacePoles,
const Handle(TColgp_HArray2OfPnt)& SurfaceCoefs,
const Handle(TColStd_HArray2OfReal)& PoleWeights,
const Handle(TColStd_HArray2OfReal)& CoefWeights,
const Standard_Boolean IsURational,
const Standard_Boolean IsVRational)
:maxderivinvok(Standard_False)
{
urational = IsURational;
vrational = IsVRational;
ucachespanlenght = 1.;
vcachespanlenght = 1.;
validcache = 1;
ucacheparameter =0.;
vcacheparameter = 0.;
Standard_Integer NbUPoles = SurfacePoles->ColLength();
Standard_Integer NbVPoles = SurfacePoles->RowLength();
@ -493,17 +475,9 @@ Geom_BezierSurface::Geom_BezierSurface
1,NbVPoles) ;
poles->ChangeArray2() = SurfacePoles->Array2();
coeffs = new TColgp_HArray2OfPnt (1,SurfaceCoefs->ColLength(),
1,SurfaceCoefs->RowLength()) ;
coeffs->ChangeArray2() = SurfaceCoefs->Array2();
if ( urational || vrational) {
weights = new TColStd_HArray2OfReal (1,NbUPoles,1,NbVPoles);
weights->ChangeArray2() = PoleWeights->Array2();
wcoeffs = new TColStd_HArray2OfReal (1,SurfaceCoefs->ColLength(),
1,SurfaceCoefs->RowLength()) ;
wcoeffs->ChangeArray2() = CoefWeights->Array2();
}
}
@ -554,10 +528,6 @@ void Geom_BezierSurface::ExchangeUV ()
Standard_Boolean temp = urational;
urational = vrational;
vrational = temp;
coeffs = new TColgp_HArray2OfPnt (LC, UC, LR, UR);
wcoeffs = new TColStd_HArray2OfReal (LC, UC, LR, UR);
UpdateCoefficients();
}
//=======================================================================
@ -680,11 +650,6 @@ void Geom_BezierSurface::InsertPoleColAfter
}
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
UpdateCoefficients();
}
//=======================================================================
@ -723,14 +688,8 @@ void Geom_BezierSurface::InsertPoleColAfter
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
Rational(weights->Array2(), urational, vrational);
UpdateCoefficients();
}
//=======================================================================
@ -796,12 +755,6 @@ void Geom_BezierSurface::InsertPoleRowAfter (const Standard_Integer UIndex,
}
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
UpdateCoefficients();
}
//=======================================================================
@ -840,14 +793,8 @@ void Geom_BezierSurface::InsertPoleRowAfter
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
Rational(weights->Array2(), urational, vrational);
UpdateCoefficients();
}
//=======================================================================
@ -907,11 +854,6 @@ void Geom_BezierSurface::RemovePoleCol (const Standard_Integer VIndex)
}
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
UpdateCoefficients();
}
//=======================================================================
@ -948,11 +890,6 @@ void Geom_BezierSurface::RemovePoleRow (const Standard_Integer UIndex)
}
poles = npoles;
weights = nweights;
coeffs = new TColgp_HArray2OfPnt(1,poles->ColLength(),
1,poles->RowLength());
wcoeffs = new TColStd_HArray2OfReal(1,poles->ColLength(),
1,poles->RowLength());
UpdateCoefficients();
}
//=======================================================================
@ -970,11 +907,46 @@ void Geom_BezierSurface::Segment
Handle(TColgp_HArray2OfPnt) Coefs;
Handle(TColStd_HArray2OfReal) WCoefs;
if (validcache == 0) UpdateCoefficients(0., 0.);
Standard_Integer aMinDegree = UDegree() <= VDegree() ? UDegree() : VDegree();
Standard_Integer aMaxDegree = UDegree() > VDegree() ? UDegree() : VDegree();
Coefs = new TColgp_HArray2OfPnt(1, aMaxDegree + 1, 1, aMinDegree + 1);
if (rat)
WCoefs = new TColStd_HArray2OfReal(1, aMaxDegree + 1, 1, aMinDegree + 1);
TColStd_Array1OfReal biduflatknots(BSplCLib::FlatBezierKnots(UDegree()), 1, 2 * (UDegree() + 1));
TColStd_Array1OfReal bidvflatknots(BSplCLib::FlatBezierKnots(VDegree()), 1, 2 * (VDegree() + 1));
Standard_Real uparameter_11 = 0.5;
Standard_Real uspanlenght_11 = 0.5;
Standard_Real vparameter_11 = 0.5;
Standard_Real vspanlenght_11 = 0.5;
if (urational || vrational) {
BSplSLib::BuildCache(uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11, 0, 0,
UDegree(), VDegree(), 0, 0,
biduflatknots, bidvflatknots,
poles->Array2(),
&weights->Array2(),
Coefs->ChangeArray2(),
&WCoefs->ChangeArray2());
}
else {
BSplSLib::BuildCache(uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11, 0, 0,
UDegree(), VDegree(), 0, 0,
biduflatknots, bidvflatknots,
poles->Array2(),
BSplSLib::NoWeights(),
Coefs->ChangeArray2(),
BSplSLib::NoWeights());
}
// Attention si udeg <= vdeg u et v sont intervertis
// dans les coeffs, il faut donc tout transposer.
if(UDegree() <= VDegree()) {
Handle(TColgp_HArray2OfPnt) coeffs = Coefs;
Handle(TColStd_HArray2OfReal) wcoeffs = WCoefs;
Standard_Integer ii, jj;
Coefs = new (TColgp_HArray2OfPnt)(1,UDegree()+1,1,VDegree()+1);
if (rat) {
@ -986,10 +958,6 @@ void Geom_BezierSurface::Segment
if (rat) WCoefs->SetValue(ii, jj, wcoeffs->Value(jj,ii));
}
}
else {
Coefs = coeffs;
if (rat) {WCoefs = wcoeffs;}
}
// Trim dans la base cannonique et Update des Poles et Coeffs
@ -1014,7 +982,6 @@ void Geom_BezierSurface::Segment
PLib::CoefficientsPoles (Coefs->Array2(), PLib::NoWeights2(),
poles->ChangeArray2(), PLib::NoWeights2());
}
UpdateCoefficients();
}
//=======================================================================
@ -1034,7 +1001,6 @@ void Geom_BezierSurface::SetPole
VIndex > Poles.RowLength() ) Standard_OutOfRange::Raise();
Poles (UIndex, VIndex) = P;
UpdateCoefficients();
}
//=======================================================================
@ -1111,8 +1077,6 @@ void Geom_BezierSurface::SetPoleCol (const Standard_Integer VIndex,
for (Standard_Integer I = CPoles.Lower(); I <= CPoles.Upper(); I++) {
Poles (I, VIndex) = CPoles (I);
}
UpdateCoefficients();
}
//=======================================================================
@ -1134,7 +1098,6 @@ void Geom_BezierSurface::SetPoleRow (const Standard_Integer UIndex,
for (Standard_Integer I = CPoles.Lower(); I <= CPoles.Upper(); I++) {
Poles (UIndex, I) = CPoles (I);
}
UpdateCoefficients();
}
//=======================================================================
@ -1181,16 +1144,12 @@ void Geom_BezierSurface::SetWeight (const Standard_Integer UIndex,
Standard_Boolean wasrat = (urational||vrational);
if (!wasrat) {
// a weight of 1. does not turn to rational
if (Abs(Weight - 1.) <= gp::Resolution()) {
UpdateCoefficients(); //Pour l'appel via SetPole
if (Abs(Weight - 1.) <= gp::Resolution())
return;
}
// set weights of 1.
weights = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength(), 1.);
wcoeffs = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength());
}
TColStd_Array2OfReal & Weights = weights->ChangeArray2();
@ -1208,14 +1167,8 @@ void Geom_BezierSurface::SetWeight (const Standard_Integer UIndex,
}
// is it turning into non rational
if (wasrat) {
if (!(urational || vrational)) {
weights.Nullify();
wcoeffs.Nullify();
}
}
UpdateCoefficients(); //Dans tous cas :Attention a SetPoleCol !
if (wasrat && !(urational || vrational))
weights.Nullify();
}
//=======================================================================
@ -1234,8 +1187,6 @@ void Geom_BezierSurface::SetWeightCol
// set weights of 1.
weights = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength(), 1.);
wcoeffs = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength());
}
TColStd_Array2OfReal & Weights = weights->ChangeArray2();
@ -1257,14 +1208,8 @@ void Geom_BezierSurface::SetWeightCol
Rational(Weights, urational, vrational);
// is it turning into non rational
if (wasrat) {
if (!(urational || vrational)) {
weights.Nullify();
wcoeffs.Nullify();
}
}
UpdateCoefficients();
if (wasrat && !(urational || vrational))
weights.Nullify();
}
//=======================================================================
@ -1283,8 +1228,6 @@ void Geom_BezierSurface::SetWeightRow
// set weights of 1.
weights = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength(), 1.);
wcoeffs = new TColStd_HArray2OfReal (1, poles->ColLength(),
1, poles->RowLength());
}
TColStd_Array2OfReal & Weights = weights->ChangeArray2();
@ -1309,14 +1252,8 @@ void Geom_BezierSurface::SetWeightRow
Rational(Weights, urational, vrational);
// is it turning into non rational
if (wasrat) {
if (!(urational || vrational)) {
weights.Nullify();
wcoeffs.Nullify();
}
}
UpdateCoefficients();
if (wasrat && !(urational || vrational))
weights.Nullify();
}
//=======================================================================
@ -1352,7 +1289,6 @@ void Geom_BezierSurface::UReverse ()
}
}
}
UpdateCoefficients();
}
//=======================================================================
@ -1399,7 +1335,6 @@ void Geom_BezierSurface::VReverse ()
}
}
}
UpdateCoefficients();
}
//=======================================================================
@ -1448,60 +1383,30 @@ void Geom_BezierSurface::D0 (const Standard_Real U,
const Standard_Real V,
gp_Pnt& P ) const
{
if (validcache == 1) {
//
// XAB : cet algorithme devient instable pour les hauts degres
// RBD : Beaucoup moins maintenant avec le calcul d'un nouveau cache
// sur [-1,1].
//
Standard_Real uparameter_11 = (2*ucacheparameter + ucachespanlenght)/2,
uspanlenght_11 = ucachespanlenght/2,
vparameter_11 = (2*vcacheparameter + vcachespanlenght)/2,
vspanlenght_11 = vcachespanlenght/2 ;
if (urational || vrational) {
BSplSLib::CacheD0(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
&wcoeffs->Array2(),
P);
}
else {
BSplSLib::CacheD0(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
BSplSLib::NoWeights(),
P);
}
Standard_Real array_u[2];
Standard_Real array_v[2];
Standard_Integer mult_u[2];
Standard_Integer mult_v[2];
TColStd_Array1OfReal biduknots(array_u[0], 1, 2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0], 1, 2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0], 1, 2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0], 1, 2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
BSplSLib::D0(U, V, 1, 1, poles->Array2(),
&weights->Array2(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P);
}
else {
Standard_Real array_u[2] ;
Standard_Real array_v[2] ;
Standard_Integer mult_u[2] ;
Standard_Integer mult_v[2] ;
TColStd_Array1OfReal biduknots(array_u[0],1,2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0],1,2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0],1,2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0],1,2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
BSplSLib::D0(U, V, 1,1,poles->Array2(),
&weights->Array2(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P) ;
}
else {
BSplSLib::D0(U, V, 1,1,poles->Array2(),
BSplSLib::NoWeights(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P) ;
}
BSplSLib::D0(U, V, 1, 1, poles->Array2(),
BSplSLib::NoWeights(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P);
}
}
@ -1517,59 +1422,29 @@ void Geom_BezierSurface::D1
gp_Vec& D1U,
gp_Vec& D1V ) const
{
if (validcache == 1) {
//
// XAB : cet algorithme devient instable pour les hauts degres
// RBD : Beaucoup moins maintenant avec le calcul d'un nouveau cache
// sur [-1,1].
//
Standard_Real uparameter_11 = (2*ucacheparameter + ucachespanlenght)/2,
uspanlenght_11 = ucachespanlenght/2,
vparameter_11 = (2*vcacheparameter + vcachespanlenght)/2,
vspanlenght_11 = vcachespanlenght/2 ;
if (urational || vrational) {
BSplSLib::CacheD1(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
&wcoeffs->Array2(),
P, D1U, D1V);
}
else {
BSplSLib::CacheD1(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
BSplSLib::NoWeights(),
P, D1U, D1V);
}
Standard_Real array_u[2];
Standard_Real array_v[2];
Standard_Integer mult_u[2];
Standard_Integer mult_v[2];
TColStd_Array1OfReal biduknots(array_u[0], 1, 2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0], 1, 2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0], 1, 2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0], 1, 2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
BSplSLib::D1(U, V, 1, 1, poles->Array2(),
&weights->Array2(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P, D1U, D1V);
}
else {
Standard_Real array_u[2] ;
Standard_Real array_v[2] ;
Standard_Integer mult_u[2] ;
Standard_Integer mult_v[2] ;
TColStd_Array1OfReal biduknots(array_u[0],1,2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0],1,2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0],1,2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0],1,2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
BSplSLib::D1(U, V, 1,1,poles->Array2(),
&weights->Array2(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P,D1U, D1V) ;
}
else {
BSplSLib::D1(U, V, 1,1,poles->Array2(),
BSplSLib::NoWeights(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P,D1U, D1V) ;
}
BSplSLib::D1(U, V, 1, 1, poles->Array2(),
BSplSLib::NoWeights(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P, D1U, D1V);
}
}
@ -1585,63 +1460,31 @@ void Geom_BezierSurface::D2
gp_Vec& D1U, gp_Vec& D1V,
gp_Vec& D2U, gp_Vec& D2V, gp_Vec& D2UV ) const
{
if (validcache == 1) {
//
// XAB : cet algorithme devient instable pour les hauts degres
// RBD : Beaucoup moins maintenant avec le calcul d'un nouveau cache
// sur [-1,1].
//
Standard_Real uparameter_11 = (2*ucacheparameter + ucachespanlenght)/2,
uspanlenght_11 = ucachespanlenght/2,
vparameter_11 = (2*vcacheparameter + vcachespanlenght)/2,
vspanlenght_11 = vcachespanlenght/2 ;
if (urational || vrational) {
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::CacheD2(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
&wcoeffs->Array2(),
P, D1U, D1V, D2U, D2UV , D2V);
}
else {
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::CacheD2(U, V, UDegree(), VDegree(),
uparameter_11, vparameter_11,
uspanlenght_11, vspanlenght_11,
coeffs->Array2(),
BSplSLib::NoWeights(),
P, D1U, D1V, D2U, D2UV , D2V);
}
Standard_Real array_u[2];
Standard_Real array_v[2];
Standard_Integer mult_u[2];
Standard_Integer mult_v[2];
TColStd_Array1OfReal biduknots(array_u[0], 1, 2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0], 1, 2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0], 1, 2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0], 1, 2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::D2(U, V, 1, 1, poles->Array2(),
&weights->Array2(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P, D1U, D1V, D2U, D2V, D2UV);
}
else {
Standard_Real array_u[2] ;
Standard_Real array_v[2] ;
Standard_Integer mult_u[2] ;
Standard_Integer mult_v[2] ;
TColStd_Array1OfReal biduknots(array_u[0],1,2); biduknots(1) = 0.; biduknots(2) = 1.;
TColStd_Array1OfInteger bidumults(mult_u[0],1,2); bidumults.Init(UDegree() + 1);
TColStd_Array1OfReal bidvknots(array_v[0],1,2); bidvknots(1) = 0.; bidvknots(2) = 1.;
TColStd_Array1OfInteger bidvmults(mult_v[0],1,2); bidvmults.Init(VDegree() + 1);
if (urational || vrational) {
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::D2(U, V, 1,1,poles->Array2(),
&weights->Array2(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P,D1U, D1V, D2U, D2V , D2UV) ;
}
else {
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::D2(U, V, 1,1,poles->Array2(),
BSplSLib::NoWeights(),
biduknots,bidvknots,&bidumults,&bidvmults,
UDegree(),VDegree(),
urational,vrational,Standard_False,Standard_False,
P,D1U, D1V, D2U, D2V, D2UV ) ;
}
//-- ATTENTION a l'ORDRE d'appel ds BSPLSLIB
BSplSLib::D2(U, V, 1, 1, poles->Array2(),
BSplSLib::NoWeights(),
biduknots, bidvknots, &bidumults, &bidvmults,
UDegree(), VDegree(),
urational, vrational, Standard_False, Standard_False,
P, D1U, D1V, D2U, D2V, D2UV);
}
}
@ -1946,7 +1789,6 @@ void Geom_BezierSurface::Transform (const gp_Trsf& T)
Poles (I, J).Transform (T);
}
}
UpdateCoefficients();
}
//=======================================================================
@ -2072,7 +1914,7 @@ void Geom_BezierSurface::Resolution(const Standard_Real Tolerance3D,
Handle(Geom_Geometry) Geom_BezierSurface::Copy() const
{
Handle(Geom_BezierSurface) S = new Geom_BezierSurface
(poles, coeffs, weights, wcoeffs, urational, vrational);
(poles, weights, urational, vrational);
return S;
}
@ -2085,70 +1927,11 @@ void Geom_BezierSurface::Init
(const Handle(TColgp_HArray2OfPnt)& Poles,
const Handle(TColStd_HArray2OfReal)& Weights)
{
Standard_Integer NbUPoles = Poles->ColLength();
Standard_Integer NbVPoles = Poles->RowLength();
Standard_Integer maxcls = Max(NbUPoles, NbVPoles);
Standard_Integer mincls = Min(NbUPoles, NbVPoles);
// set fields
poles = Poles;
coeffs = new TColgp_HArray2OfPnt (1,maxcls,1,mincls);
if (urational || vrational) {
poles = Poles;
if (urational || vrational)
weights = Weights;
wcoeffs = new TColStd_HArray2OfReal (1,maxcls,1,mincls);
}
else {
else
weights.Nullify();
wcoeffs.Nullify();
}
UpdateCoefficients();
}
//=======================================================================
//function : UpdateCoefficients
//purpose :
//=======================================================================
void Geom_BezierSurface::UpdateCoefficients(const Standard_Real ,
const Standard_Real )
{
maxderivinvok = Standard_False;
ucacheparameter = 0.;
TColStd_Array1OfReal biduflatknots(BSplCLib::FlatBezierKnots(UDegree()),
1, 2*(UDegree()+1));
vcacheparameter = 0.;
TColStd_Array1OfReal bidvflatknots(BSplCLib::FlatBezierKnots(VDegree()),
1, 2*(VDegree()+1));
Standard_Real uparameter_11 = (2*ucacheparameter + ucachespanlenght)/2,
uspanlenght_11 = ucachespanlenght/2,
vparameter_11 = (2*vcacheparameter + vcachespanlenght)/2,
vspanlenght_11 = vcachespanlenght/2 ;
if ( urational || vrational ) {
BSplSLib::BuildCache(uparameter_11,vparameter_11,
uspanlenght_11,vspanlenght_11,0,0,
UDegree(),VDegree(),0,0,
biduflatknots,bidvflatknots,
poles->Array2(),
&weights->Array2(),
coeffs->ChangeArray2(),
&wcoeffs->ChangeArray2());
}
else {
BSplSLib::BuildCache(uparameter_11,vparameter_11,
uspanlenght_11,vspanlenght_11,0,0,
UDegree(),VDegree(),0,0,
biduflatknots,bidvflatknots,
poles->Array2(),
BSplSLib::NoWeights(),
coeffs->ChangeArray2(),
BSplSLib::NoWeights());
}
validcache = 1;
}

View File

@ -31,6 +31,8 @@
#include <TColgp_Array1OfPnt.hxx>
#include <TColStd_Array1OfReal.hxx>
#include <GeomAbs_Shape.hxx>
#include <BSplSLib.hxx>
class Standard_ConstructionError;
class Standard_DimensionError;
class Standard_RangeError;
@ -471,7 +473,12 @@ public:
//! Raised if the length of P in the U an V direction is not equal to
//! NbUPoles and NbVPoles.
Standard_EXPORT void Poles (TColgp_Array2OfPnt& P) const;
//! Returns the poles of the Bezier surface.
const TColgp_Array2OfPnt& Poles() const
{
return poles->Array2();
}
//! Returns the degree of the surface in the U direction it is
//! NbUPoles - 1
@ -504,6 +511,13 @@ public:
//! equal to NbUPoles and NbVPoles.
Standard_EXPORT void Weights (TColStd_Array2OfReal& W) const;
//! Returns the weights of the Bezier surface.
const TColStd_Array2OfReal* Weights() const
{
if (!weights.IsNull())
return &weights->Array2();
return BSplSLib::NoWeights();
}
//! Returns True if the first control points row and the
//! last control points row are identical. The tolerance
@ -582,7 +596,7 @@ protected:
private:
Standard_EXPORT Geom_BezierSurface(const Handle(TColgp_HArray2OfPnt)& SurfacePoles, const Handle(TColgp_HArray2OfPnt)& SurfaceCoefficients, const Handle(TColStd_HArray2OfReal)& PoleWeights, const Handle(TColStd_HArray2OfReal)& CoefficientWeights, const Standard_Boolean IsURational, const Standard_Boolean IsVRational);
Geom_BezierSurface(const Handle(TColgp_HArray2OfPnt)& SurfacePoles, const Handle(TColStd_HArray2OfReal)& PoleWeights, const Standard_Boolean IsURational, const Standard_Boolean IsVRational);
//! Set poles to Poles, weights to Weights (not
//! copied).
@ -591,22 +605,13 @@ private:
//! coefficient 1.
//!
//! if nbpoles < 2 or nbpoles > MaDegree
Standard_EXPORT void Init (const Handle(TColgp_HArray2OfPnt)& Poles, const Handle(TColStd_HArray2OfReal)& Weights);
void Init (const Handle(TColgp_HArray2OfPnt)& Poles, const Handle(TColStd_HArray2OfReal)& Weights);
//! Recompute the coeficients.
Standard_EXPORT void UpdateCoefficients (const Standard_Real U = 0.0, const Standard_Real V = 0.0);
Standard_Boolean urational;
Standard_Boolean vrational;
Handle(TColgp_HArray2OfPnt) poles;
Handle(TColStd_HArray2OfReal) weights;
Handle(TColgp_HArray2OfPnt) coeffs;
Handle(TColStd_HArray2OfReal) wcoeffs;
Standard_Real ucacheparameter;
Standard_Real vcacheparameter;
Standard_Real ucachespanlenght;
Standard_Real vcachespanlenght;
Standard_Integer validcache;
Standard_Real umaxderivinv;
Standard_Real vmaxderivinv;
Standard_Boolean maxderivinvok;

View File

@ -28,7 +28,6 @@
#define No_Standard_DimensionError
#include <BSplCLib.hxx>
#include <Geom2d_BezierCurve.hxx>
#include <Geom2d_Geometry.hxx>
#include <gp.hxx>
@ -67,8 +66,7 @@ static Standard_Boolean Rational(const TColStd_Array1OfReal& W)
//=======================================================================
Geom2d_BezierCurve::Geom2d_BezierCurve
(const TColgp_Array1OfPnt2d& Poles):
validcache(0), parametercache(0.), spanlenghtcache(1.)
(const TColgp_Array1OfPnt2d& Poles)
{
// copy the poles
@ -89,9 +87,8 @@ Geom2d_BezierCurve::Geom2d_BezierCurve
//=======================================================================
Geom2d_BezierCurve::Geom2d_BezierCurve
(const TColgp_Array1OfPnt2d& Poles,
const TColStd_Array1OfReal& Weights):
validcache(0), parametercache(0.), spanlenghtcache(1.)
(const TColgp_Array1OfPnt2d& Poles,
const TColStd_Array1OfReal& Weights)
{
// copy the poles
@ -346,8 +343,6 @@ void Geom2d_BezierCurve::Reverse ()
cweights(nbpoles-i+1) = w;
}
}
UpdateCoefficients();
}
@ -376,18 +371,21 @@ void Geom2d_BezierCurve::Segment
// WARNING : when calling trimming be carefull that the cache
// is computed regarding 0.0e0 and not 1.0e0
//
TColStd_Array1OfReal bidflatknots(BSplCLib::FlatBezierKnots(Degree()), 1, 2 * (Degree() + 1));
TColgp_Array1OfPnt2d coeffs(1, poles->Size());
if (IsRational()) {
PLib::Trimming(U1,U2,coeffs->ChangeArray1(),&wcoeffs->ChangeArray1());
PLib::CoefficientsPoles(coeffs->Array1(),&wcoeffs->Array1(),
poles->ChangeArray1(),&weights->ChangeArray1());
TColStd_Array1OfReal wcoeffs(1, poles->Size());
BSplCLib::BuildCache(0.0, 1.0, 0, Degree(), bidflatknots,
poles->Array1(), &weights->Array1(), coeffs, &wcoeffs);
PLib::Trimming(U1, U2, coeffs, &wcoeffs);
PLib::CoefficientsPoles(coeffs, &wcoeffs, poles->ChangeArray1(), &weights->ChangeArray1());
}
else {
PLib::Trimming(U1,U2,coeffs->ChangeArray1(), PLib::NoWeights());
PLib::CoefficientsPoles(coeffs->Array1(), PLib::NoWeights(),
poles->ChangeArray1(), PLib::NoWeights());
BSplCLib::BuildCache(0.0, 1.0, 0, Degree(), bidflatknots,
poles->Array1(), BSplCLib::NoWeights(), coeffs, BSplCLib::NoWeights());
PLib::Trimming(U1, U2, coeffs, PLib::NoWeights());
PLib::CoefficientsPoles(coeffs, PLib::NoWeights(), poles->ChangeArray1(), PLib::NoWeights());
}
UpdateCoefficients();
}
@ -409,8 +407,6 @@ void Geom2d_BezierCurve::SetPole
if (Index == 1 || Index == cpoles.Length()) {
closed = (cpoles(1).Distance(cpoles(NbPoles())) <= gp::Resolution());
}
UpdateCoefficients();
}
@ -456,7 +452,6 @@ void Geom2d_BezierCurve::SetWeight
// set weights of 1.
weights = new TColStd_HArray1OfReal(1,nbpoles);
wcoeffs = new TColStd_HArray1OfReal(1,nbpoles);
weights->Init(1.);
}
@ -464,14 +459,8 @@ void Geom2d_BezierCurve::SetWeight
cweights(Index) = Weight;
// is it turning into non rational
if (wasrat) {
if (!Rational(cweights)) {
weights.Nullify();
wcoeffs.Nullify();
}
}
UpdateCoefficients();
if (wasrat && !Rational(cweights))
weights.Nullify();
}
@ -548,16 +537,7 @@ Standard_Integer Geom2d_BezierCurve::Degree () const
void Geom2d_BezierCurve::D0 (const Standard_Real U, gp_Pnt2d& P ) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom2d_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD0(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),&wcoeffs->Array1(),P);
else
BSplCLib::CacheD0(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(), BSplCLib::NoWeights(), P);
BSplCLib::D0(U, Poles(), Weights(), P);
}
//=======================================================================
@ -569,16 +549,7 @@ void Geom2d_BezierCurve::D1(const Standard_Real U,
gp_Pnt2d& P,
gp_Vec2d& V1) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom2d_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD1(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),&wcoeffs->Array1(),P,V1);
else
BSplCLib::CacheD1(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(), BSplCLib::NoWeights(), P,V1);
BSplCLib::D1(U, Poles(), Weights(), P, V1);
}
//=======================================================================
@ -591,16 +562,7 @@ void Geom2d_BezierCurve::D2 (const Standard_Real U,
gp_Vec2d& V1,
gp_Vec2d& V2) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom2d_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD2(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),&wcoeffs->Array1(),P,V1,V2);
else
BSplCLib::CacheD2(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(), BSplCLib::NoWeights(), P,V1,V2);
BSplCLib::D2(U, Poles(), Weights(), P, V1, V2);
}
//=======================================================================
@ -614,16 +576,7 @@ void Geom2d_BezierCurve::D3 (const Standard_Real U,
gp_Vec2d& V2,
gp_Vec2d& V3) const
{
// Idee lumineuse sacrifiee sur l autel des performances.
//
// if(!CoefficientsOK(U))
// ((Geom2d_BezierCurve*)(void*)this)->UpdateCoefficients(U);
if (IsRational())
BSplCLib::CacheD3(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(),&wcoeffs->Array1(),P,V1,V2,V3);
else
BSplCLib::CacheD3(U,Degree(),parametercache,spanlenghtcache,
coeffs->Array1(), BSplCLib::NoWeights(), P,V1,V2,V3);
BSplCLib::D3(U, Poles(), Weights(), P, V1, V2, V3);
}
//=======================================================================
@ -784,8 +737,6 @@ void Geom2d_BezierCurve::Transform (const gp_Trsf2d& T)
for (Standard_Integer i = 1; i <= nbpoles; i++)
cpoles (i).Transform(T);
UpdateCoefficients();
}
@ -862,60 +813,10 @@ void Geom2d_BezierCurve::Init
rational = !Weights.IsNull();
// set fields
poles = Poles;
coeffs = new TColgp_HArray1OfPnt2d (1,nbpoles);
if (rational) {
poles = Poles;
if (rational)
weights = Weights;
wcoeffs = new TColStd_HArray1OfReal (1, nbpoles, 0.0);
}
else {
weights.Nullify();
wcoeffs.Nullify();
}
UpdateCoefficients();
}
//=======================================================================
//function : CoefficientsOK
//purpose :
//=======================================================================
//Standard_Boolean Geom2d_BezierCurve::CoefficientsOK(const Standard_Real U)const
Standard_Boolean Geom2d_BezierCurve::CoefficientsOK(const Standard_Real )const
{
return (validcache) ;
}
//=======================================================================
//function : UpdateCoefficients
//purpose :
//=======================================================================
//void Geom2d_BezierCurve::UpdateCoefficients(const Standard_Real U)
void Geom2d_BezierCurve::UpdateCoefficients(const Standard_Real )
{
maxderivinvok = 0;
parametercache = 0.;
//
// Idee lumineuse sacrifiee sur l autel des performances.
// if (U >= 1.) parametercache = 1.;
TColStd_Array1OfReal bidflatknots(BSplCLib::FlatBezierKnots(Degree()),
1, 2*(Degree()+1));
if (IsRational())
BSplCLib::BuildCache(parametercache,spanlenghtcache,0,Degree(),
bidflatknots,poles->Array1(),&weights->Array1(),
coeffs->ChangeArray1(),&wcoeffs->ChangeArray1());
else
BSplCLib::BuildCache(parametercache,spanlenghtcache,0,Degree(),
bidflatknots,poles->Array1(),
BSplCLib::NoWeights(),
coeffs->ChangeArray1(),
BSplCLib::NoWeights());
validcache = 1;
weights.Nullify();
}

View File

@ -29,6 +29,8 @@
#include <TColgp_Array1OfPnt2d.hxx>
#include <TColStd_Array1OfReal.hxx>
#include <GeomAbs_Shape.hxx>
#include <BSplCLib.hxx>
class Standard_ConstructionError;
class Standard_DimensionError;
class Standard_RangeError;
@ -271,6 +273,11 @@ public:
//! Raised if the length of P is not equal to the number of poles.
Standard_EXPORT void Poles (TColgp_Array1OfPnt2d& P) const;
//! Returns all the poles of the curve.
const TColgp_Array1OfPnt2d& Poles() const
{
return poles->Array1();
}
//! Returns Value (U=1), it is the first control point
//! of the curve.
@ -284,7 +291,15 @@ public:
//!
//! Raised if the length of W is not equal to the number of poles.
Standard_EXPORT void Weights (TColStd_Array1OfReal& W) const;
//! Returns all the weights of the curve.
const TColStd_Array1OfReal* Weights() const
{
if (!weights.IsNull())
return &weights->Array1();
return BSplCLib::NoWeights();
}
//! Applies the transformation T to this Bezier curve.
Standard_EXPORT void Transform (const gp_Trsf2d& T) Standard_OVERRIDE;
@ -327,25 +342,13 @@ private:
//! Update rational and closed.
//!
//! if nbpoles < 2 or nbboles > MaDegree + 1
Standard_EXPORT void Init (const Handle(TColgp_HArray1OfPnt2d)& Poles, const Handle(TColStd_HArray1OfReal)& Weights);
//! returns true if the coefficients have been
//! computed with the right value of cacheparameter
//! for the given U value.
Standard_EXPORT Standard_Boolean CoefficientsOK (const Standard_Real U) const;
//! Recompute the coeficients.
Standard_EXPORT void UpdateCoefficients (const Standard_Real U = 0.0);
void Init (const Handle(TColgp_HArray1OfPnt2d)& Poles, const Handle(TColStd_HArray1OfReal)& Weights);
Standard_Boolean rational;
Standard_Boolean closed;
Handle(TColgp_HArray1OfPnt2d) poles;
Handle(TColStd_HArray1OfReal) weights;
Handle(TColgp_HArray1OfPnt2d) coeffs;
Handle(TColStd_HArray1OfReal) wcoeffs;
Standard_Integer validcache;
Standard_Real parametercache;
Standard_Real spanlenghtcache;
Standard_Real maxderivinv;
Standard_Boolean maxderivinvok;

View File

@ -62,10 +62,9 @@
#include <TColStd_HArray1OfInteger.hxx>
//#include <Geom2dConvert_BSplineCurveKnotSplitting.hxx>
#define myBspl Handle(Geom2d_BSplineCurve)::DownCast (myCurve)
#define PosTol Precision::PConfusion()/2
static const Standard_Real PosTol = Precision::PConfusion() / 2;
static const int maxDerivOrder = 3;
static const Standard_Integer maxDerivOrder = 3;
static const Standard_Real MinStep = 1e-7;
static gp_Vec2d dummyDerivative; // used as empty value for unused derivatives in AdjustDerivative
@ -90,17 +89,18 @@ GeomAbs_Shape Geom2dAdaptor_Curve::LocalContinuity(const Standard_Real U1,
const {
Standard_NoSuchObject_Raise_if(myTypeCurve!=GeomAbs_BSplineCurve," ");
Standard_Integer Nb = myBspl->NbKnots();
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,U1,myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,U1,aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,U2,myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,U2,aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
if ( Abs(newFirst-TK(Index1+1))<Precision::PConfusion()) {
if (Index1 < Nb)Index1++;
@ -109,7 +109,7 @@ GeomAbs_Shape Geom2dAdaptor_Curve::LocalContinuity(const Standard_Real U1,
Index2--;
Standard_Integer MultMax;
// attention aux courbes peridiques.
if ( (myBspl->IsPeriodic()) && (Index1 == Nb) )
if ( (aBspl->IsPeriodic()) && (Index1 == Nb) )
Index1 = 1;
if ( Index2 - Index1 <= 0) {
@ -120,7 +120,7 @@ GeomAbs_Shape Geom2dAdaptor_Curve::LocalContinuity(const Standard_Real U1,
for(Standard_Integer i = Index1+1;i<=Index2;i++) {
if ( TM(i)>MultMax) MultMax=TM(i);
}
MultMax = myBspl->Degree() - MultMax;
MultMax = aBspl->Degree() - MultMax;
}
if ( MultMax <= 0) {
return GeomAbs_C0;
@ -195,6 +195,7 @@ void Geom2dAdaptor_Curve::load(const Handle(Geom2d_Curve)& C,
if ( myCurve != C) {
myCurve = C;
myCurveCache = Handle(BSplCLib_Cache)();
Handle(Standard_Type) TheType = C->DynamicType();
if ( TheType == STANDARD_TYPE(Geom2d_TrimmedCurve)) {
@ -218,12 +219,19 @@ void Geom2dAdaptor_Curve::load(const Handle(Geom2d_Curve)& C,
}
else if ( TheType == STANDARD_TYPE(Geom2d_BezierCurve)) {
myTypeCurve = GeomAbs_BezierCurve;
// Create cache for Bezier
Handle(Geom2d_BezierCurve) aBezier = Handle(Geom2d_BezierCurve)::DownCast(myCurve);
Standard_Integer aDeg = aBezier->Degree();
TColStd_Array1OfReal aFlatKnots(BSplCLib::FlatBezierKnots(aDeg), 1, 2 * (aDeg + 1));
myCurveCache = new BSplCLib_Cache(aDeg, aBezier->IsPeriodic(), aFlatKnots,
aBezier->Poles(), aBezier->Weights());
}
else if ( TheType == STANDARD_TYPE(Geom2d_BSplineCurve)) {
myTypeCurve = GeomAbs_BSplineCurve;
// Create cache for B-spline
myCurveCache = new BSplCLib_Cache(myBspl->Degree(), myBspl->IsPeriodic(),
myBspl->KnotSequence(), myBspl->Poles(), myBspl->Weights());
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
myCurveCache = new BSplCLib_Cache(aBspl->Degree(), aBspl->IsPeriodic(),
aBspl->KnotSequence(), aBspl->Poles(), aBspl->Weights());
}
else if ( TheType == STANDARD_TYPE(Geom2d_OffsetCurve))
{
@ -289,8 +297,9 @@ Standard_Integer Geom2dAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
Standard_Integer myNbIntervals = 1;
Standard_Integer NbSplit;
if (myTypeCurve == GeomAbs_BSplineCurve) {
Standard_Integer FirstIndex = myBspl->FirstUKnotIndex();
Standard_Integer LastIndex = myBspl->LastUKnotIndex();
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer FirstIndex = aBspl->FirstUKnotIndex();
Standard_Integer LastIndex = aBspl->LastUKnotIndex();
TColStd_Array1OfInteger Inter (1, LastIndex-FirstIndex+1);
if ( S > Continuity()) {
Standard_Integer Cont;
@ -310,11 +319,11 @@ Standard_Integer Geom2dAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
if ( S == GeomAbs_C1) Cont = 1;
else if ( S == GeomAbs_C2) Cont = 2;
else if ( S == GeomAbs_C3) Cont = 3;
else Cont = myBspl->Degree();
Standard_Integer Degree = myBspl->Degree();
Standard_Integer NbKnots = myBspl->NbKnots();
else Cont = aBspl->Degree();
Standard_Integer Degree = aBspl->Degree();
Standard_Integer NbKnots = aBspl->NbKnots();
TColStd_Array1OfInteger Mults (1, NbKnots);
myBspl->Multiplicities (Mults);
aBspl->Multiplicities (Mults);
NbSplit = 1;
Standard_Integer Index = FirstIndex;
Inter (NbSplit) = Index;
@ -333,19 +342,19 @@ Standard_Integer Geom2dAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
Standard_Integer NbInt = NbSplit-1;
Standard_Integer Nb = myBspl->NbKnots();
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myFirst,
myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myFirst,
aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myLast,
myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myLast,
aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
// On decale eventuellement les indices
@ -393,8 +402,9 @@ void Geom2dAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
Standard_Integer myNbIntervals = 1;
Standard_Integer NbSplit;
if (myTypeCurve == GeomAbs_BSplineCurve) {
Standard_Integer FirstIndex = myBspl->FirstUKnotIndex();
Standard_Integer LastIndex = myBspl->LastUKnotIndex();
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer FirstIndex = aBspl->FirstUKnotIndex();
Standard_Integer LastIndex = aBspl->LastUKnotIndex();
TColStd_Array1OfInteger Inter (1, LastIndex-FirstIndex+1);
if ( S > Continuity()) {
Standard_Integer Cont;
@ -414,11 +424,11 @@ void Geom2dAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
if ( S == GeomAbs_C1) Cont = 1;
else if ( S == GeomAbs_C2) Cont = 2;
else if ( S == GeomAbs_C3) Cont = 3;
else Cont = myBspl->Degree();
Standard_Integer Degree = myBspl->Degree();
Standard_Integer NbKnots = myBspl->NbKnots();
else Cont = aBspl->Degree();
Standard_Integer Degree = aBspl->Degree();
Standard_Integer NbKnots = aBspl->NbKnots();
TColStd_Array1OfInteger Mults (1, NbKnots);
myBspl->Multiplicities (Mults);
aBspl->Multiplicities (Mults);
NbSplit = 1;
Standard_Integer Index = FirstIndex;
Inter (NbSplit) = Index;
@ -436,19 +446,19 @@ void Geom2dAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
Inter (NbSplit) = Index;
Standard_Integer NbInt = NbSplit-1;
Standard_Integer Nb = myBspl->NbKnots();
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myFirst,
myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myFirst,
aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myLast,
myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myLast,
aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
@ -561,9 +571,21 @@ Standard_Real Geom2dAdaptor_Curve::Period() const
//=======================================================================
void Geom2dAdaptor_Curve::RebuildCache(const Standard_Real theParameter) const
{
myCurveCache->BuildCache(theParameter, myBspl->Degree(),
myBspl->IsPeriodic(), myBspl->KnotSequence(),
myBspl->Poles(), myBspl->Weights());
if (myTypeCurve == GeomAbs_BezierCurve)
{
Handle(Geom2d_BezierCurve) aBezier = Handle(Geom2d_BezierCurve)::DownCast(myCurve);
Standard_Integer aDeg = aBezier->Degree();
TColStd_Array1OfReal aFlatKnots(BSplCLib::FlatBezierKnots(aDeg), 1, 2 * (aDeg + 1));
myCurveCache->BuildCache(theParameter, aDeg, aBezier->IsPeriodic(), aFlatKnots,
aBezier->Poles(), aBezier->Weights());
}
else if (myTypeCurve == GeomAbs_BSplineCurve)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
myCurveCache->BuildCache(theParameter, aBspl->Degree(),
aBspl->IsPeriodic(), aBspl->KnotSequence(),
aBspl->Poles(), aBspl->Weights());
}
}
//=======================================================================
@ -577,6 +599,12 @@ gp_Pnt2d Geom2dAdaptor_Curve::Value(const Standard_Real U) const
return ValueBSpline(U);
else if (myTypeCurve == GeomAbs_OffsetCurve)
return ValueOffset(U);
else if (myTypeCurve == GeomAbs_BezierCurve)
{ // use cached data
gp_Pnt2d aRes;
myCurveCache->D0(U, aRes);
return aRes;
}
return myCurve->Value(U);
}
@ -589,20 +617,21 @@ gp_Pnt2d Geom2dAdaptor_Curve::ValueBSpline(const Standard_Real theU) const
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst)
{
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast)
{
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
return myBspl->LocalValue(theU, Ideb, Ifin);
return aBspl->LocalValue(theU, Ideb, Ifin);
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
{
@ -642,17 +671,13 @@ gp_Pnt2d Geom2dAdaptor_Curve::ValueOffset(const Standard_Real theU) const
void Geom2dAdaptor_Curve::D0(const Standard_Real U, gp_Pnt2d& P) const
{
if (myTypeCurve == GeomAbs_BSplineCurve)
{
D0BSpline(U, P);
return;
}
else if (myTypeCurve == GeomAbs_OffsetCurve)
{
D0Offset(U, P);
return;
}
myCurve->D0(U, P);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D0(U, P);
else
myCurve->D0(U, P);
}
//=======================================================================
@ -663,18 +688,19 @@ void Geom2dAdaptor_Curve::D0BSpline(const Standard_Real theU, gp_Pnt2d& theP) co
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD0(theU, Ideb, Ifin, theP);
aBspl->LocalD0(theU, Ideb, Ifin, theP);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -705,17 +731,13 @@ void Geom2dAdaptor_Curve::D1(const Standard_Real U,
gp_Pnt2d& P, gp_Vec2d& V) const
{
if (myTypeCurve == GeomAbs_BSplineCurve)
{
D1BSpline(U, P, V);
return;
}
else if (myTypeCurve == GeomAbs_OffsetCurve)
{
D1Offset(U, P, V);
return;
}
myCurve->D1(U, P, V);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D1(U, P, V);
else
myCurve->D1(U, P, V);
}
//=======================================================================
@ -726,18 +748,19 @@ void Geom2dAdaptor_Curve::D1BSpline(const Standard_Real theU, gp_Pnt2d& theP, gp
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD1(theU, Ideb, Ifin, theP, theV);
aBspl->LocalD1(theU, Ideb, Ifin, theP, theV);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -781,17 +804,13 @@ void Geom2dAdaptor_Curve::D2(const Standard_Real U,
gp_Pnt2d& P, gp_Vec2d& V1, gp_Vec2d& V2) const
{
if (myTypeCurve == GeomAbs_BSplineCurve)
{
D2BSpline(U, P, V1, V2);
return;
}
else if (myTypeCurve == GeomAbs_OffsetCurve)
{
D2Offset(U, P, V1, V2);
return;
}
myCurve->D2(U, P, V1, V2);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D2(U, P, V1, V2);
else
myCurve->D2(U, P, V1, V2);
}
//=======================================================================
@ -803,18 +822,19 @@ void Geom2dAdaptor_Curve::D2BSpline(const Standard_Real theU, gp_Pnt2d& theP,
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD2(theU, Ideb, Ifin, theP, theV1, theV2);
aBspl->LocalD2(theU, Ideb, Ifin, theP, theV1, theV2);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -862,17 +882,13 @@ void Geom2dAdaptor_Curve::D3(const Standard_Real U,
gp_Vec2d& V2, gp_Vec2d& V3) const
{
if (myTypeCurve == GeomAbs_BSplineCurve)
{
D3BSpline(U, P, V1, V2, V3);
return;
}
else if (myTypeCurve == GeomAbs_OffsetCurve)
{
D3Offset(U, P, V1, V2, V3);
return;
}
myCurve->D3(U, P, V1, V2, V3);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D3(U, P, V1, V2, V3);
else
myCurve->D3(U, P, V1, V2, V3);
}
//=======================================================================
@ -884,18 +900,19 @@ void Geom2dAdaptor_Curve::D3BSpline(const Standard_Real theU, gp_Pnt2d& theP,
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD3(theU, Ideb, Ifin, theP, theV1, theV2, theV3);
aBspl->LocalD3(theU, Ideb, Ifin, theP, theV1, theV2, theV3);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -961,18 +978,19 @@ gp_Vec2d Geom2dAdaptor_Curve::DNBSpline(const Standard_Real U,
{
if (U==myFirst || U==myLast)
{
Handle(Geom2d_BSplineCurve) aBspl = Handle(Geom2d_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (U==myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (U==myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
return myBspl->LocalDN( U, Ideb, Ifin, N);
return aBspl->LocalDN( U, Ideb, Ifin, N);
}
return myCurve->DN( U, N);

View File

@ -60,10 +60,9 @@
#include <TColStd_HArray1OfInteger.hxx>
//#include <GeomConvert_BSplineCurveKnotSplitting.hxx>
#define myBspl Handle(Geom_BSplineCurve)::DownCast (myCurve)
#define PosTol Precision::PConfusion()/2
static const Standard_Real PosTol = Precision::PConfusion() / 2;
static const int maxDerivOrder = 3;
static const Standard_Integer maxDerivOrder = 3;
static const Standard_Real MinStep = 1e-7;
static gp_Vec dummyDerivative; // used as empty value for unused derivatives in AdjustDerivative
@ -88,17 +87,18 @@ GeomAbs_Shape GeomAdaptor_Curve::LocalContinuity(const Standard_Real U1,
const
{
Standard_NoSuchObject_Raise_if(myTypeCurve!=GeomAbs_BSplineCurve," ");
Standard_Integer Nb = myBspl->NbKnots();
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,U1,myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,U1,aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,U2,myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,U2,aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
if ( Abs(newFirst-TK(Index1+1))<Precision::PConfusion()) {
if (Index1 < Nb) Index1++;
@ -107,7 +107,7 @@ GeomAbs_Shape GeomAdaptor_Curve::LocalContinuity(const Standard_Real U1,
Index2--;
Standard_Integer MultMax;
// attention aux courbes peridiques.
if ( (myBspl->IsPeriodic()) && (Index1 == Nb) )
if ( (aBspl->IsPeriodic()) && (Index1 == Nb) )
Index1 = 1;
if ( Index2 - Index1 <= 0) {
@ -118,7 +118,7 @@ GeomAbs_Shape GeomAdaptor_Curve::LocalContinuity(const Standard_Real U1,
for(Standard_Integer i = Index1+1;i<=Index2;i++) {
if ( TM(i)>MultMax) MultMax=TM(i);
}
MultMax = myBspl->Degree() - MultMax;
MultMax = aBspl->Degree() - MultMax;
}
if ( MultMax <= 0) {
return GeomAbs_C0;
@ -152,7 +152,8 @@ void GeomAdaptor_Curve::load(const Handle(Geom_Curve)& C,
if ( myCurve != C) {
myCurve = C;
myCurveCache = Handle(BSplCLib_Cache)();
const Handle(Standard_Type)& TheType = C->DynamicType();
if ( TheType == STANDARD_TYPE(Geom_TrimmedCurve)) {
Load(Handle(Geom_TrimmedCurve)::DownCast (C)->BasisCurve(),UFirst,ULast);
@ -174,12 +175,19 @@ void GeomAdaptor_Curve::load(const Handle(Geom_Curve)& C,
}
else if ( TheType == STANDARD_TYPE(Geom_BezierCurve)) {
myTypeCurve = GeomAbs_BezierCurve;
// Create cache for Bezier
Handle(Geom_BezierCurve) aBezier = Handle(Geom_BezierCurve)::DownCast(myCurve);
Standard_Integer aDeg = aBezier->Degree();
TColStd_Array1OfReal aFlatKnots(BSplCLib::FlatBezierKnots(aDeg), 1, 2 * (aDeg + 1));
myCurveCache = new BSplCLib_Cache(aDeg, aBezier->IsPeriodic(), aFlatKnots,
aBezier->Poles(), aBezier->Weights());
}
else if ( TheType == STANDARD_TYPE(Geom_BSplineCurve)) {
myTypeCurve = GeomAbs_BSplineCurve;
// Create cache for B-spline
myCurveCache = new BSplCLib_Cache(myBspl->Degree(), myBspl->IsPeriodic(),
myBspl->KnotSequence(), myBspl->Poles(), myBspl->Weights());
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
myCurveCache = new BSplCLib_Cache(aBspl->Degree(), aBspl->IsPeriodic(),
aBspl->KnotSequence(), aBspl->Poles(), aBspl->Weights());
}
else if ( TheType == STANDARD_TYPE(Geom_OffsetCurve)) {
myTypeCurve = GeomAbs_OffsetCurve;
@ -240,8 +248,9 @@ Standard_Integer GeomAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
Standard_Integer myNbIntervals = 1;
Standard_Integer NbSplit;
if (myTypeCurve == GeomAbs_BSplineCurve) {
Standard_Integer FirstIndex = myBspl->FirstUKnotIndex();
Standard_Integer LastIndex = myBspl->LastUKnotIndex();
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer FirstIndex = aBspl->FirstUKnotIndex();
Standard_Integer LastIndex = aBspl->LastUKnotIndex();
TColStd_Array1OfInteger Inter (1, LastIndex-FirstIndex+1);
if ( S > Continuity()) {
Standard_Integer Cont;
@ -261,11 +270,11 @@ Standard_Integer GeomAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
if ( S == GeomAbs_C1) Cont = 1;
else if ( S == GeomAbs_C2) Cont = 2;
else if ( S == GeomAbs_C3) Cont = 3;
else Cont = myBspl->Degree();
Standard_Integer Degree = myBspl->Degree();
Standard_Integer NbKnots = myBspl->NbKnots();
else Cont = aBspl->Degree();
Standard_Integer Degree = aBspl->Degree();
Standard_Integer NbKnots = aBspl->NbKnots();
TColStd_Array1OfInteger Mults (1, NbKnots);
myBspl->Multiplicities (Mults);
aBspl->Multiplicities (Mults);
NbSplit = 1;
Standard_Integer Index = FirstIndex;
Inter (NbSplit) = Index;
@ -284,19 +293,19 @@ Standard_Integer GeomAdaptor_Curve::NbIntervals(const GeomAbs_Shape S) const
Standard_Integer NbInt = NbSplit-1;
Standard_Integer Nb = myBspl->NbKnots();
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myFirst,
myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myFirst,
aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myLast,
myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myLast,
aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
// On decale eventuellement les indices
@ -361,8 +370,9 @@ void GeomAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
if (myTypeCurve == GeomAbs_BSplineCurve)
{
Standard_Integer FirstIndex = myBspl->FirstUKnotIndex();
Standard_Integer LastIndex = myBspl->LastUKnotIndex();
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer FirstIndex = aBspl->FirstUKnotIndex();
Standard_Integer LastIndex = aBspl->LastUKnotIndex();
TColStd_Array1OfInteger Inter (1, LastIndex-FirstIndex+1);
if ( S > Continuity()) {
@ -383,11 +393,11 @@ void GeomAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
if ( S == GeomAbs_C1) Cont = 1;
else if ( S == GeomAbs_C2) Cont = 2;
else if ( S == GeomAbs_C3) Cont = 3;
else Cont = myBspl->Degree();
Standard_Integer Degree = myBspl->Degree();
Standard_Integer NbKnots = myBspl->NbKnots();
else Cont = aBspl->Degree();
Standard_Integer Degree = aBspl->Degree();
Standard_Integer NbKnots = aBspl->NbKnots();
TColStd_Array1OfInteger Mults (1, NbKnots);
myBspl->Multiplicities (Mults);
aBspl->Multiplicities (Mults);
NbSplit = 1;
Standard_Integer Index = FirstIndex;
Inter (NbSplit) = Index;
@ -409,19 +419,19 @@ void GeomAdaptor_Curve::Intervals(TColStd_Array1OfReal& T,
// TColStd_Array1OfInteger Inter(1,NbInt+1);
// Convector.Splitting( Inter);
Standard_Integer Nb = myBspl->NbKnots();
Standard_Integer Nb = aBspl->NbKnots();
Standard_Integer Index1 = 0;
Standard_Integer Index2 = 0;
Standard_Real newFirst, newLast;
TColStd_Array1OfReal TK(1,Nb);
TColStd_Array1OfInteger TM(1,Nb);
myBspl->Knots(TK);
myBspl->Multiplicities(TM);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myFirst,
myBspl->IsPeriodic(),
aBspl->Knots(TK);
aBspl->Multiplicities(TM);
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myFirst,
aBspl->IsPeriodic(),
1,Nb,Index1,newFirst);
BSplCLib::LocateParameter(myBspl->Degree(),TK,TM,myLast,
myBspl->IsPeriodic(),
BSplCLib::LocateParameter(aBspl->Degree(),TK,TM,myLast,
aBspl->IsPeriodic(),
1,Nb,Index2,newLast);
FirstParam = newFirst;
LastParam = newLast;
@ -543,9 +553,21 @@ Standard_Real GeomAdaptor_Curve::Period() const
//=======================================================================
void GeomAdaptor_Curve::RebuildCache(const Standard_Real theParameter) const
{
myCurveCache->BuildCache(theParameter, myBspl->Degree(),
myBspl->IsPeriodic(), myBspl->KnotSequence(),
myBspl->Poles(), myBspl->Weights());
if (myTypeCurve == GeomAbs_BezierCurve)
{
Handle(Geom_BezierCurve) aBezier = Handle(Geom_BezierCurve)::DownCast(myCurve);
Standard_Integer aDeg = aBezier->Degree();
TColStd_Array1OfReal aFlatKnots(BSplCLib::FlatBezierKnots(aDeg), 1, 2 * (aDeg + 1));
myCurveCache->BuildCache(theParameter, aDeg, aBezier->IsPeriodic(), aFlatKnots,
aBezier->Poles(), aBezier->Weights());
}
else if (myTypeCurve == GeomAbs_BSplineCurve)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
myCurveCache->BuildCache(theParameter, aBspl->Degree(),
aBspl->IsPeriodic(), aBspl->KnotSequence(),
aBspl->Poles(), aBspl->Weights());
}
}
//=======================================================================
@ -559,6 +581,12 @@ gp_Pnt GeomAdaptor_Curve::Value(const Standard_Real U) const
return ValueBSpline(U);
else if (myTypeCurve == GeomAbs_OffsetCurve)
return ValueOffset(U);
else if (myTypeCurve == GeomAbs_BezierCurve)
{ // use cached data
gp_Pnt aRes;
myCurveCache->D0(U, aRes);
return aRes;
}
return myCurve->Value(U);
}
@ -570,18 +598,19 @@ gp_Pnt GeomAdaptor_Curve::ValueBSpline(const Standard_Real theU) const
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
return myBspl->LocalValue(theU, Ideb, Ifin);
return aBspl->LocalValue(theU, Ideb, Ifin);
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
{
@ -626,6 +655,8 @@ void GeomAdaptor_Curve::D0(const Standard_Real U, gp_Pnt& P) const
D0BSpline(U, P);
else if (myTypeCurve == GeomAbs_OffsetCurve)
D0Offset(U, P);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D0(U, P);
else
myCurve->D0(U, P);
}
@ -638,18 +669,19 @@ void GeomAdaptor_Curve::D0BSpline(const Standard_Real theU, gp_Pnt& theP) const
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD0(theU, Ideb, Ifin, theP);
aBspl->LocalD0(theU, Ideb, Ifin, theP);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -682,6 +714,8 @@ void GeomAdaptor_Curve::D1(const Standard_Real U, gp_Pnt& P, gp_Vec& V) const
D1BSpline(U, P, V);
else if (myTypeCurve == GeomAbs_OffsetCurve)
D1Offset(U, P, V);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D1(U, P, V);
else
myCurve->D1(U, P, V);
}
@ -694,18 +728,19 @@ void GeomAdaptor_Curve::D1BSpline(const Standard_Real theU, gp_Pnt& theP, gp_Vec
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD1(theU, Ideb, Ifin, theP, theV);
aBspl->LocalD1(theU, Ideb, Ifin, theP, theV);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -750,6 +785,8 @@ void GeomAdaptor_Curve::D2(const Standard_Real U,
D2BSpline(U, P, V1, V2);
else if (myTypeCurve == GeomAbs_OffsetCurve)
D2Offset(U, P, V1, V2);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D2(U, P, V1, V2);
else
myCurve->D2(U, P, V1, V2);
}
@ -763,18 +800,19 @@ void GeomAdaptor_Curve::D2BSpline(const Standard_Real theU, gp_Pnt& theP,
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD2(theU, Ideb, Ifin, theP, theV1, theV2);
aBspl->LocalD2(theU, Ideb, Ifin, theP, theV1, theV2);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -820,6 +858,8 @@ void GeomAdaptor_Curve::D3(const Standard_Real U,
D3BSpline(U, P, V1, V2, V3);
else if (myTypeCurve == GeomAbs_OffsetCurve)
D3Offset(U, P, V1, V2, V3);
else if (myTypeCurve == GeomAbs_BezierCurve) // use cached data
myCurveCache->D3(U, P, V1, V2, V3);
else
myCurve->D3(U, P, V1, V2, V3);
}
@ -834,18 +874,19 @@ void GeomAdaptor_Curve::D3BSpline(const Standard_Real theU,
{
if (theU == myFirst || theU == myLast)
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (theU == myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (theU == myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
myBspl->LocalD3(theU, Ideb, Ifin, theP, theV1, theV2, theV3);
aBspl->LocalD3(theU, Ideb, Ifin, theP, theV1, theV2, theV3);
return;
}
else if (!myCurveCache.IsNull()) // use cached B-spline data
@ -901,18 +942,19 @@ gp_Vec GeomAdaptor_Curve::DNBSpline(const Standard_Real U,
{
if ((U==myFirst || U==myLast))
{
Handle(Geom_BSplineCurve) aBspl = Handle(Geom_BSplineCurve)::DownCast(myCurve);
Standard_Integer Ideb = 0, Ifin = 0;
if (U==myFirst) {
myBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
aBspl->LocateU(myFirst, PosTol, Ideb, Ifin);
if (Ideb<1) Ideb=1;
if (Ideb>=Ifin) Ifin = Ideb+1;
}
if (U==myLast) {
myBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin>myBspl->NbKnots()) Ifin = myBspl->NbKnots();
aBspl->LocateU(myLast, PosTol, Ideb, Ifin);
if (Ifin > aBspl->NbKnots()) Ifin = aBspl->NbKnots();
if (Ideb>=Ifin) Ideb = Ifin-1;
}
return myBspl->LocalDN( U, Ideb, Ifin, N);
return aBspl->LocalDN( U, Ideb, Ifin, N);
}
return myCurve->DN( U, N);
}

View File

@ -20,8 +20,6 @@
#define No_Standard_RangeError
#define No_Standard_OutOfRange
#define PosTol (Precision::PConfusion()*0.5)
#include <Adaptor3d_HCurve.hxx>
#include <Adaptor3d_HSurface.hxx>
@ -70,6 +68,8 @@
#include <TColStd_Array1OfReal.hxx>
#include <TColStd_HArray1OfInteger.hxx>
static const Standard_Real PosTol = Precision::PConfusion()*0.5;
//=======================================================================
//function : LocalContinuity
//purpose :
@ -139,9 +139,7 @@ void GeomAdaptor_Surface::load(const Handle(Geom_Surface)& S,
myNestedEvaluator = Handle(GeomEvaluator_Surface)();
const Handle(Standard_Type)& TheType = S->DynamicType();
if ( TheType == STANDARD_TYPE(Geom_BezierSurface))
mySurfaceType = GeomAbs_BezierSurface;
else if (TheType == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
if (TheType == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
Load(Handle(Geom_RectangularTrimmedSurface)::DownCast (S)->BasisSurface(),
UFirst,ULast,VFirst,VLast);
}
@ -155,7 +153,7 @@ void GeomAdaptor_Surface::load(const Handle(Geom_Surface)& S,
mySurfaceType = GeomAbs_Sphere;
else if ( TheType == STANDARD_TYPE(Geom_ToroidalSurface))
mySurfaceType = GeomAbs_Torus;
else if (TheType == STANDARD_TYPE(Geom_SurfaceOfRevolution))
else if ( TheType == STANDARD_TYPE(Geom_SurfaceOfRevolution))
{
mySurfaceType = GeomAbs_SurfaceOfRevolution;
Handle(Geom_SurfaceOfRevolution) myRevSurf =
@ -168,7 +166,7 @@ void GeomAdaptor_Surface::load(const Handle(Geom_Surface)& S,
Handle(Adaptor3d_HCurve)::DownCast(aBaseAdaptor),
myRevSurf->Direction(), myRevSurf->Location());
}
else if (TheType == STANDARD_TYPE(Geom_SurfaceOfLinearExtrusion))
else if ( TheType == STANDARD_TYPE(Geom_SurfaceOfLinearExtrusion))
{
mySurfaceType = GeomAbs_SurfaceOfExtrusion;
Handle(Geom_SurfaceOfLinearExtrusion) myExtSurf =
@ -180,16 +178,30 @@ void GeomAdaptor_Surface::load(const Handle(Geom_Surface)& S,
myNestedEvaluator = new GeomEvaluator_SurfaceOfExtrusion(
Handle(Adaptor3d_HCurve)::DownCast(aBaseAdaptor), myExtSurf->Direction());
}
else if ( TheType == STANDARD_TYPE(Geom_BSplineSurface)) {
else if (TheType == STANDARD_TYPE(Geom_BezierSurface))
{
mySurfaceType = GeomAbs_BezierSurface;
// Create cache for Bezier
Handle(Geom_BezierSurface) aBezier = Handle(Geom_BezierSurface)::DownCast(mySurface);
Standard_Integer aDegU = aBezier->UDegree();
Standard_Integer aDegV = aBezier->VDegree();
TColStd_Array1OfReal aFlatKnotsU(BSplCLib::FlatBezierKnots(aDegU), 1, 2 * (aDegU + 1));
TColStd_Array1OfReal aFlatKnotsV(BSplCLib::FlatBezierKnots(aDegV), 1, 2 * (aDegV + 1));
mySurfaceCache = new BSplSLib_Cache(
aDegU, aBezier->IsUPeriodic(), aFlatKnotsU,
aDegV, aBezier->IsVPeriodic(), aFlatKnotsV,
aBezier->Poles(), aBezier->Weights());
}
else if (TheType == STANDARD_TYPE(Geom_BSplineSurface)) {
mySurfaceType = GeomAbs_BSplineSurface;
Handle(Geom_BSplineSurface) myBspl = Handle(Geom_BSplineSurface)::DownCast(mySurface);
// Create cache for B-spline
mySurfaceCache = new BSplSLib_Cache(
mySurfaceCache = new BSplSLib_Cache(
myBspl->UDegree(), myBspl->IsUPeriodic(), myBspl->UKnotSequence(),
myBspl->VDegree(), myBspl->IsVPeriodic(), myBspl->VKnotSequence(),
myBspl->Poles(), myBspl->Weights());
}
else if (TheType == STANDARD_TYPE(Geom_OffsetSurface))
else if ( TheType == STANDARD_TYPE(Geom_OffsetSurface))
{
mySurfaceType = GeomAbs_OffsetSurface;
Handle(Geom_OffsetSurface) myOffSurf = Handle(Geom_OffsetSurface)::DownCast(mySurface);
@ -668,11 +680,26 @@ Standard_Real GeomAdaptor_Surface::VPeriod() const
void GeomAdaptor_Surface::RebuildCache(const Standard_Real theU,
const Standard_Real theV) const
{
Handle(Geom_BSplineSurface) myBspl = Handle(Geom_BSplineSurface)::DownCast(mySurface);
mySurfaceCache->BuildCache(theU, theV,
myBspl->UDegree(), myBspl->IsUPeriodic(), myBspl->UKnotSequence(),
myBspl->VDegree(), myBspl->IsVPeriodic(), myBspl->VKnotSequence(),
myBspl->Poles(), myBspl->Weights());
if (mySurfaceType == GeomAbs_BezierSurface)
{
Handle(Geom_BezierSurface) aBezier = Handle(Geom_BezierSurface)::DownCast(mySurface);
Standard_Integer aDegU = aBezier->UDegree();
Standard_Integer aDegV = aBezier->VDegree();
TColStd_Array1OfReal aFlatKnotsU(BSplCLib::FlatBezierKnots(aDegU), 1, 2 * (aDegU + 1));
TColStd_Array1OfReal aFlatKnotsV(BSplCLib::FlatBezierKnots(aDegV), 1, 2 * (aDegV + 1));
mySurfaceCache->BuildCache(theU, theV,
aDegU, aBezier->IsUPeriodic(), aFlatKnotsU,
aDegV, aBezier->IsVPeriodic(), aFlatKnotsV,
aBezier->Poles(), aBezier->Weights());
}
else if (mySurfaceType == GeomAbs_BSplineSurface)
{
Handle(Geom_BSplineSurface) myBspl = Handle(Geom_BSplineSurface)::DownCast(mySurface);
mySurfaceCache->BuildCache(theU, theV,
myBspl->UDegree(), myBspl->IsUPeriodic(), myBspl->UKnotSequence(),
myBspl->VDegree(), myBspl->IsVPeriodic(), myBspl->VKnotSequence(),
myBspl->Poles(), myBspl->Weights());
}
}
//=======================================================================
@ -698,6 +725,7 @@ void GeomAdaptor_Surface::D0(const Standard_Real U,
{
switch (mySurfaceType)
{
case GeomAbs_BezierSurface:
case GeomAbs_BSplineSurface:
if (!mySurfaceCache.IsNull())
{
@ -742,9 +770,11 @@ void GeomAdaptor_Surface::D1(const Standard_Real U,
else if (Abs(V-myVLast) <= myTolV) {VSide= -1; v = myVLast;}
switch(mySurfaceType) {
case GeomAbs_BezierSurface:
case GeomAbs_BSplineSurface: {
Handle(Geom_BSplineSurface) myBspl = Handle(Geom_BSplineSurface)::DownCast(mySurface);
if ((USide != 0 || VSide != 0) &&
if (!myBspl.IsNull() &&
(USide != 0 || VSide != 0) &&
IfUVBound(u, v, Ideb, Ifin, IVdeb, IVfin, USide, VSide))
myBspl->LocalD1(u, v, Ideb, Ifin, IVdeb, IVfin, P, D1U, D1V);
else if (!mySurfaceCache.IsNull())
@ -754,7 +784,7 @@ void GeomAdaptor_Surface::D1(const Standard_Real U,
mySurfaceCache->D1(U, V, P, D1U, D1V);
}
else
myBspl->D1(u, v, P, D1U, D1V);
mySurface->D1(u, v, P, D1U, D1V);
break;
}
@ -793,9 +823,11 @@ void GeomAdaptor_Surface::D2(const Standard_Real U,
else if (Abs(V-myVLast) <= myTolV) {VSide= -1; v = myVLast;}
switch(mySurfaceType) {
case GeomAbs_BezierSurface:
case GeomAbs_BSplineSurface: {
Handle(Geom_BSplineSurface) myBspl = Handle(Geom_BSplineSurface)::DownCast(mySurface);
if ((USide != 0 || VSide != 0) &&
if (!myBspl.IsNull() &&
(USide != 0 || VSide != 0) &&
IfUVBound(u, v, Ideb, Ifin, IVdeb, IVfin, USide, VSide))
myBspl->LocalD2(u, v, Ideb, Ifin, IVdeb, IVfin, P, D1U, D1V, D2U, D2V, D2UV);
else if (!mySurfaceCache.IsNull())
@ -804,12 +836,13 @@ void GeomAdaptor_Surface::D2(const Standard_Real U,
RebuildCache(U, V);
mySurfaceCache->D2(U, V, P, D1U, D1V, D2U, D2V, D2UV);
}
else myBspl->D2(u, v, P, D1U, D1V, D2U, D2V, D2UV);
else
mySurface->D2(u, v, P, D1U, D1V, D2U, D2V, D2UV);
break;
}
case GeomAbs_SurfaceOfExtrusion:
case GeomAbs_SurfaceOfRevolution:
case GeomAbs_SurfaceOfExtrusion :
case GeomAbs_SurfaceOfRevolution :
case GeomAbs_OffsetSurface :
Standard_NoSuchObject_Raise_if(myNestedEvaluator.IsNull(),
"GeomAdaptor_Surface::D2: evaluator is not initialized");
@ -854,9 +887,9 @@ void GeomAdaptor_Surface::D3(const Standard_Real U, const Standard_Real V,
}
break;
}
case GeomAbs_SurfaceOfExtrusion:
case GeomAbs_SurfaceOfRevolution:
case GeomAbs_SurfaceOfExtrusion :
case GeomAbs_SurfaceOfRevolution :
case GeomAbs_OffsetSurface:
Standard_NoSuchObject_Raise_if(myNestedEvaluator.IsNull(),
"GeomAdaptor_Surface::D3: evaluator is not initialized");
@ -864,7 +897,7 @@ void GeomAdaptor_Surface::D3(const Standard_Real U, const Standard_Real V,
break;
default : { mySurface->D3(u,v,P,D1U,D1V,D2U,D2V,D2UV,D3U,D3V,D3UUV,D3UVV);
break;}
break;}
}
}