1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-08-29 14:00:49 +03:00

0029769: Uninitialized data with BSplCLib_Cache, BSplSLib_Cache

Implementation of classes BSplCLib_Cache and BSplSLib_Cache is revised:
- Common functionality dealing with spans along one parametric direction is separated to new struct BSplCLib_CacheParams
- Empty constructors are removed; copying is prohibited
- Code reconsidering degree and other parameters on each call to BuildCache() is eliminated; curve parameters must be the same in constructor and all calls to BuildCache()
- Extra call to BuildCache() from constructor is eliminated
This commit is contained in:
abv
2018-06-10 22:40:12 +03:00
committed by bugmaster
parent 3388cf17dc
commit 0a96e0bbc4
9 changed files with 279 additions and 382 deletions

View File

@@ -31,159 +31,71 @@ static Standard_Real* ConvertArray(const Handle(TColStd_HArray2OfReal)& theHArra
return (Standard_Real*) &(anArray(anArray.LowerRow(), anArray.LowerCol()));
}
BSplCLib_Cache::BSplCLib_Cache()
BSplCLib_Cache::BSplCLib_Cache(const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt2d& /* only used to distinguish from 3d variant */,
const TColStd_Array1OfReal* theWeights)
: myIsRational(theWeights != NULL),
myParams (theDegree, thePeriodic, theFlatKnots)
{
myPolesWeights.Nullify();
myIsRational = Standard_False;
mySpanStart = 0.0;
mySpanLength = 0.0;
mySpanIndex = 0;
myDegree = 0;
myFlatKnots.Nullify();
Standard_Integer aPWColNumber = (myIsRational ? 3 : 2);
myPolesWeights = new TColStd_HArray2OfReal (1, theDegree + 1, 1, aPWColNumber);
}
BSplCLib_Cache::BSplCLib_Cache(const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt2d& thePoles2d,
const TColgp_Array1OfPnt& /* only used to distinguish from 2d variant */,
const TColStd_Array1OfReal* theWeights)
: myIsRational(theWeights != NULL),
myParams (theDegree, thePeriodic, theFlatKnots)
{
Standard_Real aCacheParam = theFlatKnots.Value(theFlatKnots.Lower() + theDegree);
BuildCache(aCacheParam, theDegree, thePeriodic,
theFlatKnots, thePoles2d, theWeights);
Standard_Integer aPWColNumber = (myIsRational ? 4 : 3);
myPolesWeights = new TColStd_HArray2OfReal (1, theDegree + 1, 1, aPWColNumber);
}
BSplCLib_Cache::BSplCLib_Cache(const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt& thePoles,
const TColStd_Array1OfReal* theWeights)
{
Standard_Real aCacheParam = theFlatKnots.Value(theFlatKnots.Lower() + theDegree);
BuildCache(aCacheParam, theDegree, thePeriodic,
theFlatKnots, thePoles, theWeights);
}
Standard_Boolean BSplCLib_Cache::IsCacheValid(Standard_Real theParameter) const
{
Standard_Real aNewParam = theParameter;
if (!myFlatKnots.IsNull())
PeriodicNormalization(myFlatKnots->Array1(), aNewParam);
Standard_Real aDelta = aNewParam - mySpanStart;
return ((aDelta >= 0.0 || mySpanIndex == mySpanIndexMin) &&
(aDelta < mySpanLength || mySpanIndex == mySpanIndexMax));
return myParams.IsCacheValid (theParameter);
}
void BSplCLib_Cache::PeriodicNormalization(const TColStd_Array1OfReal& theFlatKnots,
Standard_Real& theParameter) const
{
Standard_Real aPeriod = theFlatKnots.Value(theFlatKnots.Upper() - myDegree) -
theFlatKnots.Value(myDegree + 1) ;
if (theParameter < theFlatKnots.Value(myDegree + 1))
{
Standard_Real aScale = IntegerPart(
(theFlatKnots.Value(myDegree + 1) - theParameter) / aPeriod);
theParameter += aPeriod * (aScale + 1.0);
}
if (theParameter > theFlatKnots.Value(theFlatKnots.Upper() - myDegree))
{
Standard_Real aScale = IntegerPart(
(theParameter - theFlatKnots.Value(theFlatKnots.Upper() - myDegree)) / aPeriod);
theParameter -= aPeriod * (aScale + 1.0);
}
}
void BSplCLib_Cache::BuildCache(const Standard_Real& theParameter,
const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt2d& thePoles2d,
const TColStd_Array1OfReal* theWeights)
{
// Normalize theParameter for periodical B-splines
Standard_Real aNewParam = theParameter;
if (thePeriodic)
{
PeriodicNormalization(theFlatKnots, aNewParam);
myFlatKnots = new TColStd_HArray1OfReal(1, theFlatKnots.Length());
myFlatKnots->ChangeArray1() = theFlatKnots;
}
else if (!myFlatKnots.IsNull()) // Periodical curve became non-periodical
myFlatKnots.Nullify();
// Change the size of cached data if needed
myIsRational = (theWeights != NULL);
Standard_Integer aPWColNumber = myIsRational ? 3 : 2;
if (theDegree > myDegree)
myPolesWeights = new TColStd_HArray2OfReal(1, theDegree + 1, 1, aPWColNumber);
myDegree = theDegree;
mySpanIndex = 0;
BSplCLib::LocateParameter(theDegree, theFlatKnots, BSplCLib::NoMults(),
aNewParam, thePeriodic, mySpanIndex, aNewParam);
mySpanStart = theFlatKnots.Value(mySpanIndex);
mySpanLength = theFlatKnots.Value(mySpanIndex + 1) - mySpanStart;
mySpanIndexMin = thePeriodic ? 0 : myDegree + 1;
mySpanIndexMax = theFlatKnots.Length() - 1 - theDegree;
Standard_Real aNewParam = myParams.PeriodicNormalization (theParameter);
myParams.LocateParameter (aNewParam, theFlatKnots);
// Calculate new cache data
BSplCLib::BuildCache(mySpanStart, mySpanLength, thePeriodic, theDegree,
mySpanIndex, theFlatKnots, thePoles2d, theWeights,
myPolesWeights->ChangeArray2());
BSplCLib::BuildCache (myParams.SpanStart, myParams.SpanLength, myParams.IsPeriodic,
myParams.Degree, myParams.SpanIndex, theFlatKnots, thePoles2d,
theWeights, myPolesWeights->ChangeArray2());
}
void BSplCLib_Cache::BuildCache(const Standard_Real& theParameter,
const Standard_Integer& theDegree,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt& thePoles,
const TColStd_Array1OfReal* theWeights)
{
// Create list of knots with repetitions and normalize theParameter for periodical B-splines
Standard_Real aNewParam = theParameter;
if (thePeriodic)
{
PeriodicNormalization(theFlatKnots, aNewParam);
myFlatKnots = new TColStd_HArray1OfReal(1, theFlatKnots.Length());
myFlatKnots->ChangeArray1() = theFlatKnots;
}
else if (!myFlatKnots.IsNull()) // Periodical curve became non-periodical
myFlatKnots.Nullify();
// Change the size of cached data if needed
myIsRational = (theWeights != NULL);
Standard_Integer aPWColNumber = myIsRational ? 4 : 3;
if (theDegree > myDegree)
myPolesWeights = new TColStd_HArray2OfReal(1, theDegree + 1, 1, aPWColNumber);
myDegree = theDegree;
mySpanIndex = 0;
BSplCLib::LocateParameter(theDegree, theFlatKnots, BSplCLib::NoMults(),
aNewParam, thePeriodic, mySpanIndex, aNewParam);
mySpanStart = theFlatKnots.Value(mySpanIndex);
mySpanLength = theFlatKnots.Value(mySpanIndex + 1) - mySpanStart;
mySpanIndexMin = thePeriodic ? 0 : myDegree + 1;
mySpanIndexMax = theFlatKnots.Length() - 1 - theDegree;
Standard_Real aNewParam = myParams.PeriodicNormalization (theParameter);
myParams.LocateParameter (aNewParam, theFlatKnots);
// Calculate new cache data
BSplCLib::BuildCache(mySpanStart, mySpanLength, thePeriodic, theDegree,
mySpanIndex, theFlatKnots, thePoles, theWeights,
myPolesWeights->ChangeArray2());
BSplCLib::BuildCache (myParams.SpanStart, myParams.SpanLength, myParams.IsPeriodic,
myParams.Degree, myParams.SpanIndex, theFlatKnots, thePoles,
theWeights, myPolesWeights->ChangeArray2());
}
void BSplCLib_Cache::CalculateDerivative(const Standard_Real& theParameter,
const Standard_Integer& theDerivative,
Standard_Real& theDerivArray) const
{
Standard_Real aNewParameter = theParameter;
if (!myFlatKnots.IsNull()) // B-spline is periodical
PeriodicNormalization(myFlatKnots->Array1(), aNewParameter);
aNewParameter = (aNewParameter - mySpanStart) / mySpanLength;
Standard_Real aNewParameter = myParams.PeriodicNormalization (theParameter);
aNewParameter = (aNewParameter - myParams.SpanStart) / myParams.SpanLength;
Standard_Real* aPolesArray = ConvertArray(myPolesWeights);
Standard_Integer aDimension = myPolesWeights->RowLength(); // number of columns
@@ -199,23 +111,23 @@ void BSplCLib_Cache::CalculateDerivative(const Standard_Real& theParameter,
// When the degree of curve is lesser than the requested derivative,
// nullify array cells corresponding to greater derivatives
Standard_Integer aDerivative = theDerivative;
if (myDegree < theDerivative)
if (myParams.Degree < theDerivative)
{
aDerivative = myDegree;
for (Standard_Integer ind = myDegree * aDimension; ind < (theDerivative + 1) * aDimension; ind++)
aDerivative = myParams.Degree;
for (Standard_Integer ind = myParams.Degree * aDimension; ind < (theDerivative + 1) * aDimension; ind++)
{
aPntDeriv[ind] = 0.0;
(&theDerivArray)[ind] = 0.0; // should be cleared separately, because aPntDeriv may look to another memory area
}
}
PLib::EvalPolynomial(aNewParameter, aDerivative, myDegree, aDimension,
PLib::EvalPolynomial(aNewParameter, aDerivative, myParams.Degree, aDimension,
aPolesArray[0], aPntDeriv[0]);
// Unnormalize derivatives since those are computed normalized
Standard_Real aFactor = 1.0;
for (Standard_Integer deriv = 1; deriv <= aDerivative; deriv++)
{
aFactor /= mySpanLength;
aFactor /= myParams.SpanLength;
for (Standard_Integer ind = 0; ind < aDimension; ind++)
aPntDeriv[aDimension * deriv + ind] *= aFactor;
}
@@ -227,17 +139,15 @@ void BSplCLib_Cache::CalculateDerivative(const Standard_Real& theParameter,
void BSplCLib_Cache::D0(const Standard_Real& theParameter, gp_Pnt2d& thePoint) const
{
Standard_Real aNewParameter = theParameter;
if (!myFlatKnots.IsNull()) // B-spline is periodical
PeriodicNormalization(myFlatKnots->Array1(), aNewParameter);
aNewParameter = (aNewParameter - mySpanStart) / mySpanLength;
Standard_Real aNewParameter = myParams.PeriodicNormalization (theParameter);
aNewParameter = (aNewParameter - myParams.SpanStart) / myParams.SpanLength;
Standard_Real* aPolesArray = ConvertArray(myPolesWeights);
Standard_Real aPoint[4];
Standard_Integer aDimension = myPolesWeights->RowLength(); // number of columns
PLib::NoDerivativeEvalPolynomial(aNewParameter, myDegree,
aDimension, myDegree * aDimension,
PLib::NoDerivativeEvalPolynomial(aNewParameter, myParams.Degree,
aDimension, myParams.Degree * aDimension,
aPolesArray[0], aPoint[0]);
thePoint.SetCoord(aPoint[0], aPoint[1]);
@@ -247,17 +157,15 @@ void BSplCLib_Cache::D0(const Standard_Real& theParameter, gp_Pnt2d& thePoint) c
void BSplCLib_Cache::D0(const Standard_Real& theParameter, gp_Pnt& thePoint) const
{
Standard_Real aNewParameter = theParameter;
if (!myFlatKnots.IsNull()) // B-spline is periodical
PeriodicNormalization(myFlatKnots->Array1(), aNewParameter);
aNewParameter = (aNewParameter - mySpanStart) / mySpanLength;
Standard_Real aNewParameter = myParams.PeriodicNormalization (theParameter);
aNewParameter = (aNewParameter - myParams.SpanStart) / myParams.SpanLength;
Standard_Real* aPolesArray = ConvertArray(myPolesWeights);
Standard_Real aPoint[4];
Standard_Integer aDimension = myPolesWeights->RowLength(); // number of columns
PLib::NoDerivativeEvalPolynomial(aNewParameter, myDegree,
aDimension, myDegree * aDimension,
PLib::NoDerivativeEvalPolynomial(aNewParameter, myParams.Degree,
aDimension, myParams.Degree * aDimension,
aPolesArray[0], aPoint[0]);
thePoint.SetCoord(aPoint[0], aPoint[1], aPoint[2]);

View File

@@ -19,7 +19,6 @@
#include <Standard_Type.hxx>
#include <Standard_Transient.hxx>
#include <gp_Pnt2d.hxx>
#include <gp_Pnt.hxx>
#include <gp_Vec2d.hxx>
@@ -31,6 +30,8 @@
#include <TColgp_Array1OfPnt.hxx>
#include <TColgp_Array1OfPnt2d.hxx>
#include <BSplCLib_CacheParams.hxx>
//! \brief A cache class for Bezier and B-spline curves.
//!
//! Defines all data, that can be cached on a span of a curve.
@@ -38,11 +39,10 @@
class BSplCLib_Cache : public Standard_Transient
{
public:
//! Default constructor
Standard_EXPORT BSplCLib_Cache();
//! Constructor for caching of 2D curves
//! Constructor, prepares data structures for caching values on a 2d curve.
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param thePeriodic identify whether 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
@@ -51,9 +51,10 @@ public:
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt2d& thePoles2d,
const TColStd_Array1OfReal* theWeights = NULL);
//! Constructor for caching of 3D curves
//! Constructor, prepares data structures for caching values on a 3d curve.
//! \param theDegree degree of the curve
//! \param thePeriodic identify the curve is periodic
//! \param thePeriodic identify whether 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
@@ -69,27 +70,20 @@ 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 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,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt2d& thePoles2d,
const TColStd_Array1OfReal* theWeights = NULL);
const TColStd_Array1OfReal* theWeights);
//! 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 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,
const Standard_Boolean& thePeriodic,
const TColStd_Array1OfReal& theFlatKnots,
const TColgp_Array1OfPnt& thePoles,
const TColStd_Array1OfReal* theWeights = NULL);
@@ -142,10 +136,6 @@ public:
DEFINE_STANDARD_RTTIEXT(BSplCLib_Cache,Standard_Transient)
protected:
//! 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 curve
//! \param[in] theParameter parameter of the calculation
@@ -156,21 +146,18 @@ protected:
const Standard_Integer& theDerivative,
Standard_Real& theDerivArray) const;
// copying is prohibited
BSplCLib_Cache (const BSplCLib_Cache&);
void operator = (const BSplCLib_Cache&);
private:
Handle(TColStd_HArray2OfReal) myPolesWeights; ///< array of poles and weights of calculated cache
Standard_Boolean myIsRational; //!< identifies the rationality of Bezier/B-spline curve
BSplCLib_CacheParams myParams; //!< cache parameters
Handle(TColStd_HArray2OfReal) myPolesWeights; //!< array of poles and weights of calculated cache
// the array has following structure:
// x1 y1 [z1] [w1]
// 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 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 Bezier/B-spline curve
Standard_Integer mySpanIndexMin; ///< minimal index of 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

@@ -0,0 +1,106 @@
// Copyright (c) 2018 OPEN CASCADE SAS
//
// This file is part of Open CASCADE Technology software library.
//
// This library is free software; you can redistribute it and/or modify it under
// the terms of the GNU Lesser General Public License version 2.1 as published
// by the Free Software Foundation, with special exception defined in the file
// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
// distribution for complete text of the license and disclaimer of any warranty.
//
// Alternatively, this file may be used under the terms of Open CASCADE
// commercial license or contractual agreement.
#ifndef _BSplCLib_CacheParams_Headerfile
#define _BSplCLib_CacheParams_Headerfile
#include <Standard_Real.hxx>
#include <TColStd_Array1OfReal.hxx>
#include <BSplCLib.hxx>
//! Simple structure containing parameters describing parameterization
//! of a B-spline curve or a surface in one direction (U or V),
//! and data of the current span for its caching
struct BSplCLib_CacheParams
{
const Standard_Integer Degree; ///< degree of Bezier/B-spline
const Standard_Boolean IsPeriodic; ///< true of the B-spline is periodic
const Standard_Real FirstParameter; ///< first valid parameter
const Standard_Real LastParameter; ///< last valid parameter
const Standard_Integer SpanIndexMin; ///< minimal index of span
const Standard_Integer SpanIndexMax; ///< maximal index of span
Standard_Real SpanStart; ///< parameter for the frst point of the span
Standard_Real SpanLength; ///< length of the span
Standard_Integer SpanIndex; ///< index of the span
//! Constructor, prepares data structures for caching.
//! \param theDegree degree of the B-spline (or Bezier)
//! \param thePeriodic identify whether the B-spline is periodic
//! \param theFlatKnots knots of Bezier / B-spline parameterization
BSplCLib_CacheParams (Standard_Integer theDegree, Standard_Boolean thePeriodic,
const TColStd_Array1OfReal& theFlatKnots)
: Degree(theDegree),
IsPeriodic(thePeriodic),
FirstParameter(theFlatKnots.Value(theFlatKnots.Lower() + theDegree)),
LastParameter(theFlatKnots.Value(theFlatKnots.Upper() - theDegree)),
SpanIndexMin(theFlatKnots.Lower() + theDegree),
SpanIndexMax(theFlatKnots.Upper() - theDegree - 1),
SpanStart(0.),
SpanLength(0.),
SpanIndex(0)
{}
//! Normalizes the parameter for periodic B-splines
//! \param theParameter the value to be normalized into the knots array
Standard_Real PeriodicNormalization (Standard_Real theParameter) const
{
if (IsPeriodic)
{
if (theParameter < FirstParameter)
{
Standard_Real aPeriod = LastParameter - FirstParameter;
Standard_Real aScale = IntegerPart ((FirstParameter - theParameter) / aPeriod);
return theParameter + aPeriod * (aScale + 1.0);
}
if (theParameter > LastParameter)
{
Standard_Real aPeriod = LastParameter - FirstParameter;
Standard_Real aScale = IntegerPart ((theParameter - LastParameter) / aPeriod);
return theParameter - aPeriod * (aScale + 1.0);
}
}
return theParameter;
}
//! Verifies validity of the cache using flat parameter of the point
//! \param theParameter parameter of the point placed in the span
Standard_Boolean IsCacheValid (Standard_Real theParameter) const
{
Standard_Real aNewParam = PeriodicNormalization (theParameter);
Standard_Real aDelta = aNewParam - SpanStart;
return ((aDelta >= 0.0 || SpanIndex == SpanIndexMin) &&
(aDelta < SpanLength || SpanIndex == SpanIndexMax));
}
//! Computes span for the specified parameter
//! \param theParameter parameter of the point placed in the span
//! \param theFlatKnots knots of Bezier / B-spline parameterization
void LocateParameter (Standard_Real& theParameter, const TColStd_Array1OfReal& theFlatKnots)
{
SpanIndex = 0;
BSplCLib::LocateParameter (Degree, theFlatKnots, BSplCLib::NoMults(),
theParameter, IsPeriodic, SpanIndex, theParameter);
SpanStart = theFlatKnots.Value(SpanIndex);
SpanLength = theFlatKnots.Value(SpanIndex + 1) - SpanStart;
}
private:
// copying is prohibited
BSplCLib_CacheParams (const BSplCLib_CacheParams&);
void operator = (const BSplCLib_CacheParams&);
};
#endif

View File

@@ -7,6 +7,7 @@ BSplCLib_3.cxx
BSplCLib_BzSyntaxes.cxx
BSplCLib_Cache.cxx
BSplCLib_Cache.hxx
BSplCLib_CacheParams.hxx
BSplCLib_CurveComputation.gxx
BSplCLib_EvaluatorFunction.hxx
BSplCLib_KnotDistribution.hxx