STK++ 0.9.13
STK::BSplineCoefficients< Data > Class Template Reference

Compute the regression B-splines coefficients. More...

#include <STK_BSplineCoefficients.h>

Inheritance diagram for STK::BSplineCoefficients< Data >:
Inheritance graph

Public Types

typedef IBasis< Data, ArrayXXBase
 
- Public Types inherited from STK::IBasis< Data, ArrayXX >
typedef Data::Type Type
 

Public Member Functions

 BSplineCoefficients (Data const *p_data=0, int nbControlPoints=1, int degree=3, Regress::KnotsPosition position=Regress::uniformKnotsPositions_, bool useDataValues=true)
 Default constructor : initialize the data members with default values.
 
 BSplineCoefficients (Data const &data, int nbControlPoints, int degree=3, Regress::KnotsPosition position=Regress::uniformKnotsPositions_, bool useDataValues=true)
 Constructor : initialize the data members.
 
 BSplineCoefficients (BSplineCoefficients const &coefs)
 copy constructor.
 
virtual ~BSplineCoefficients ()
 Destructor.
 
BSplineCoefficientsclone () const
 clone pattern implementation
 
virtual bool run ()
 run the computations.
 
int degree () const
 
int nbKnots () const
 
int nbControlPoints () const
 
VectorX constknots () const
 
void setNbControlPoints (int nbControlPoints)
 Set the number of control point (the number of BSpline)
 
void setDim (int dim)
 Alias for setting the number of control point (the number of BSpline)
 
void setDegree (int degree=3)
 Set the degree of the BSpline basis.
 
void setPosition (Regress::KnotsPosition position)
 Set the kind of position to use for the knots.
 
void setParameters (Data const &data, int nbControlPoints, int degree=3, Regress::KnotsPosition position=Regress::uniformKnotsPositions_)
 Set the whole set of parameters for computing the coefficients.
 
void setData (Data const &data, int nbControlPoints, int degree=3, Regress::KnotsPosition position=Regress::uniformKnotsPositions_)
 [DEPRECATED] Set the whole set of parameters for computing the coefficients
 
template<class OtherVector >
ArrayXX extrapolate (OtherVector const &x) const
 
- Public Member Functions inherited from STK::IBasis< Data, ArrayXX >
 IBasis ()
 default constructor
 
 IBasis (Data const *p_data, int dim, bool useDataValues=true)
 constructor
 
 IBasis (Data const &data, int dim, bool useDataValues=true)
 constructor
 
 IBasis (IBasis const &basis)
 copy constructor.
 
virtual ~IBasis ()
 destructor
 
int dim () const
 
ArrayXX constcoefficients () const
 
Type minValue () const
 
Type maxValue () const
 
void setData (Data const &data)
 Set the data set.
 
void setDim (int dim)
 
void setMinValue (Type const &minValue)
 
void setMaxValue (Type const &maxValue)
 
bool initializeStep ()
 Initialize the parameters minValue_ and maxValue_ using data set.
 
- Public Member Functions inherited from STK::IRunnerBase
String consterror () const
 get the last error message.
 

Protected Member Functions

bool computeKnots ()
 compute the position of the knots of the B-spline curves.
 
void computeCoefficients ()
 Compute the coefficients of the B-spline curves.
 
- Protected Member Functions inherited from STK::IBasis< Data, ArrayXX >
virtual void update ()
 update IBasis if a parameter or a new data set is set, update the state of this runner.
 
- Protected Member Functions inherited from STK::IRunnerBase
 IRunnerBase ()
 default constructor
 
 IRunnerBase (IRunnerBase const &runner)
 copy constructor
 
virtual ~IRunnerBase ()
 destructor
 

Protected Attributes

int degree_
 degree of the B-splines curves.
 
Regress::KnotsPosition position_
 Method used in order to position the knots.
 
int nbKnots_
 number of knots of the B-spline curves.
 
int lastControlPoint_
 Index of the last control point of the B-spline curves.
 
VectorX knots_
 Data of the knots.
 
- Protected Attributes inherited from STK::IBasis< Data, ArrayXX >
Data constp_data_
 the input data set
 
int dim_
 number of dimension to build
 
bool useDataValues_
 
Type minValue_
 Minimal value of the data.
 
Type maxValue_
 Maximal value of the data.
 
ArrayXX coefficients_
 Array2D of the coefficients.
 
- Protected Attributes inherited from STK::IRunnerBase
String msg_error_
 String with the last error message.
 
bool hasRun_
 true if run has been used, false otherwise
 

Private Member Functions

void computeUniformKnots ()
 compute the position of the uniform knots.
 
void computePeriodicKnots ()
 compute the position of the periodic knots.
 
void computeDensityKnots (bool isSorted)
 compute the position of the density knots
 
void computeCoefficientsRow (int irow, Real const &value)
 Compute a row of the coefficients matrix for a given value.
 

Detailed Description

template<class Data>
class STK::BSplineCoefficients< Data >

Compute the regression B-splines coefficients.

The BSplineCoefficients class computes the coefficients of a BSpline curve using the de Boor's algorithm. The knots can be uniform (the default), periodic or density placed.

Given the number of control points and the degree of the B-splines, the number of knots is given by the formula nbKnots = nbControlPoints + degree + 1.

Note
If the input data set is a vector of size n the output matrix of the coefficients Coefficients() is a matrix of size (n, nbControlPoints).

Definition at line 62 of file STK_BSplineCoefficients.h.

Member Typedef Documentation

◆ Base

template<class Data >
typedef IBasis<Data, ArrayXX> STK::BSplineCoefficients< Data >::Base

Definition at line 65 of file STK_BSplineCoefficients.h.

Constructor & Destructor Documentation

◆ BSplineCoefficients() [1/3]

template<class Data >
STK::BSplineCoefficients< Data >::BSplineCoefficients ( Data const p_data = 0,
int  nbControlPoints = 1,
int  degree = 3,
Regress::KnotsPosition  position = Regress::uniformKnotsPositions_,
bool  useDataValues = true 
)

Default constructor : initialize the data members with default values.

The number of knots is given by the formula nbKnots = nbControlPoints + degree +1.

Parameters
p_datathe input data values
nbControlPointsnumber of control points
degreedegree of the B-spline curves
positionmethod to use for positioning the knots
useDataValuesif true use data in order to find the minValue_ and maxValue_, use values set otherwise.

Definition at line 209 of file STK_BSplineCoefficients.h.

215 : Base(p_data, nbControlPoints, useDataValues)
216 , degree_(degree)
217 , position_(position)
220 , knots_( Range(0, nbKnots_) )
221{ }
int nbKnots_
number of knots of the B-spline curves.
int degree_
degree of the B-splines curves.
IBasis< Data, ArrayXX > Base
Regress::KnotsPosition position_
Method used in order to position the knots.
VectorX knots_
Data of the knots.
int lastControlPoint_
Index of the last control point of the B-spline curves.
int dim_
number of dimension to build
Definition STK_IBasis.h:112
TRange< UnknownSize > Range
Definition STK_Range.h:59

◆ BSplineCoefficients() [2/3]

template<class Data >
STK::BSplineCoefficients< Data >::BSplineCoefficients ( Data const data,
int  nbControlPoints,
int  degree = 3,
Regress::KnotsPosition  position = Regress::uniformKnotsPositions_,
bool  useDataValues = true 
)

Constructor : initialize the data members.

The number of knots is given by the formula nbKnots = nbControlPoints + degree +1.

Parameters
datathe input data values
nbControlPointsnumber of control points
degreedegree of the B-spline curves
positionmethod to use for positioning the knots
useDataValuesif true use data in order to find the minValue_ and maxValue_, use values set otherwise.

Definition at line 225 of file STK_BSplineCoefficients.h.

231 : Base(data, nbControlPoints, useDataValues)
232 , degree_(degree)
233 , position_(position)
236 , knots_( Range(0, nbKnots_) )
237{}

◆ BSplineCoefficients() [3/3]

template<class Data >
STK::BSplineCoefficients< Data >::BSplineCoefficients ( BSplineCoefficients< Data > const coefs)

copy constructor.

Parameters
coefsthe coefficients to copy

Definition at line 242 of file STK_BSplineCoefficients.h.

243 : Base(coefs)
244 , degree_(coefs.degree_)
245 , position_(coefs.position_)
246 , nbKnots_(coefs.nbKnots_)
247 , lastControlPoint_(coefs.lastControlPoint_)
248 , knots_(coefs.knots_)
249{}

◆ ~BSplineCoefficients()

template<class Data >
virtual STK::BSplineCoefficients< Data >::~BSplineCoefficients ( )
inlinevirtual

Destructor.

Definition at line 111 of file STK_BSplineCoefficients.h.

111{}

Member Function Documentation

◆ clone()

template<class Data >
BSplineCoefficients * STK::BSplineCoefficients< Data >::clone ( ) const
inline

clone pattern implementation

Definition at line 113 of file STK_BSplineCoefficients.h.

113{ return new BSplineCoefficients(*this);}
BSplineCoefficients(Data const *p_data=0, int nbControlPoints=1, int degree=3, Regress::KnotsPosition position=Regress::uniformKnotsPositions_, bool useDataValues=true)
Default constructor : initialize the data members with default values.

◆ computeCoefficients()

template<class Data >
void STK::BSplineCoefficients< Data >::computeCoefficients ( )
protected

Compute the coefficients of the B-spline curves.

Definition at line 442 of file STK_BSplineCoefficients.h.

443{
444#ifdef STK_REGRESS_VERBOSE
445 stk_cout << _T("BSplineCoefficients::computeCoefficients()\n");
446#endif
447
448 // compute the coefficients
449 for (int i=p_data_->begin(); i< p_data_->end(); i++)
450 { computeCoefficientsRow(i, (*p_data_)[i]);}
451
452#ifdef STK_REGRESS_VERBOSE
453 stk_cout << _T("BSplineCoefficients::computeCoefficients() done\n");
454#endif
455}
#define stk_cout
Standard stk output stream.
#define _T(x)
Let x unmodified.
void computeCoefficientsRow(int irow, Real const &value)
Compute a row of the coefficients matrix for a given value.
Data const * p_data_
the input data set
Definition STK_IBasis.h:110

References _T, and stk_cout.

◆ computeCoefficientsRow()

template<class Data >
void STK::BSplineCoefficients< Data >::computeCoefficientsRow ( int  irow,
Real const value 
)
private

Compute a row of the coefficients matrix for a given value.

Parameters
irowindex of the row
valuethe value to which we want to compute the row

Definition at line 515 of file STK_BSplineCoefficients.h.

516{
517 // value outside the range of the knots case
518 if (value <= minValue_)
519 {
520 coefficients_(irow, 0) = 1.0;
521 return;
522 }
523 if (value >= maxValue_)
524 {
525 coefficients_(irow, lastControlPoint_) = 1.0;
526 return;
527 }
528 // find interval
529 int k, k1;
530 for (k=0, k1=1; k<lastControlPoint_; k++, k1++)
531 {
532 if (value < knots_[k1]) break;
533 }
534 // begin recursion
535 coefficients_(irow, k) = 1.0;
536 for (int d=1; d<=degree_; d++)
537 {
538 // right (south-west corner) term only
539 coefficients_(irow, k-d) = ( (knots_[k1] - value)/(knots_[k1] - knots_[k1-d]) )
540 * coefficients_(irow, k1-d);
541 // compute internal terms
542 for (int i = k1-d; i<k; i++)
543 {
544 const Real knots_i = knots_[i], knots_id1 = knots_[i+d+1];
545 coefficients_(irow, i) = ( (value - knots_i)/(knots_[i+d] - knots_i) )
546 * coefficients_(irow, i)
547 + ( (knots_id1 - value)/(knots_id1 - knots_[i+1]) )
548 * coefficients_(irow, i+1);
549 }
550 // left (north-west corner) term only
551 coefficients_(irow, k) *= (value - knots_[k])/(knots_[k+d] - knots_[k]);
552 }
553}
Type minValue_
Minimal value of the data.
Definition STK_IBasis.h:116
Type maxValue_
Maximal value of the data.
Definition STK_IBasis.h:118
ArrayXX coefficients_
Array2D of the coefficients.
Definition STK_IBasis.h:121
double Real
STK fundamental type of Real values.

◆ computeDensityKnots()

template<class Data >
void STK::BSplineCoefficients< Data >::computeDensityKnots ( bool  isSorted)
private

compute the position of the density knots

Parameters
isSortedtrue if the data is sorted, false otherwise

Definition at line 486 of file STK_BSplineCoefficients.h.

487{
488 // sorted data
489 Data xtri(*p_data_, true);
490 // sort the data
491 if (!isSorted) heapSort< Data >(xtri);
492
493 // compute step
494 Real step = xtri.size()/(Real)(lastControlPoint_-degree_+1);
495
496 // set internal knots
497 int first = xtri.begin();
498 for (int k = degree_ + 1, kcell =1; k <= lastControlPoint_; k++, kcell++)
499 {
500 int cell = first + int(kcell* step);
501 knots_[k] = (xtri[cell] + xtri[cell+1])/2.;
502 }
503 // set external knots
504 for ( int k=0, j = lastControlPoint_+1; k <= degree_; j++, k++)
505 {
506 knots_[k] = xtri.front();
507 knots_[j] = xtri.back();
508 }
509}

◆ computeKnots()

template<class Data >
bool STK::BSplineCoefficients< Data >::computeKnots ( )
protected

compute the position of the knots of the B-spline curves.

Definition at line 406 of file STK_BSplineCoefficients.h.

407{
408 // resize and initialize knots
410 // set knots values
411 switch (position_)
412 {
413 // uniform position
416 break;
417 // periodic position
420 break;
421 // density position
424 break;
425 default:
426
428 return false;
429 break;
430 }
431 // shift knots
433 {
434 Real range = (maxValue_ - minValue_);
435 for (int k = 0; k < nbKnots_; k++) knots_[k] = minValue_ + range * knots_[k];
436 }
437 return true;
438}
#define STKERROR_NO_ARG(Where, Error)
Definition STK_Macros.h:49
void computeDensityKnots(bool isSorted)
compute the position of the density knots
void computePeriodicKnots()
compute the position of the periodic knots.
void computeUniformKnots()
compute the position of the uniform knots.
bool computeKnots()
compute the position of the knots of the B-spline curves.
Derived & resize(Range const &I, Range const &J)
resize the array.
String msg_error_
String with the last error message.
Definition STK_IRunner.h:96
@ uniformKnotsPositions_
uniform knots
@ densityKnotsPositions_
knots using density of the data
@ periodicKnotsPositions_
periodic knots

References STK::BSplineCoefficients< Data >::computeKnots(), STK::Regress::densityKnotsPositions_, STK::IRunnerBase::msg_error_, STK::Regress::periodicKnotsPositions_, STKERROR_NO_ARG, and STK::Regress::uniformKnotsPositions_.

Referenced by STK::BSplineCoefficients< Data >::computeKnots().

◆ computePeriodicKnots()

template<class Data >
void STK::BSplineCoefficients< Data >::computePeriodicKnots ( )
private

compute the position of the periodic knots.

Definition at line 475 of file STK_BSplineCoefficients.h.

476{
477 // compute step
478 Real step = 1.0/(dim_ - degree_);
479 // set knots
480 for (int k = 0, j = -degree_; k < nbKnots_; j++, k++)
481 knots_[k] = j * step;
482;
483}

◆ computeUniformKnots()

template<class Data >
void STK::BSplineCoefficients< Data >::computeUniformKnots ( )
private

compute the position of the uniform knots.

Definition at line 459 of file STK_BSplineCoefficients.h.

460{
461 // compute step
462 Real step = 1.0/(dim_ - degree_);
463 // set internal knots
464 const int first = degree_ + 1;
465 for (int k = first, j = 1; k <= lastControlPoint_; j++, k++) knots_[k] = j * step;
466 // set external knots
467 for ( int k=0, j = lastControlPoint_+1; k < first; j++, k++)
468 {
469 knots_[k] = 0;
470 knots_[j] = 1;
471 }
472}

◆ degree()

template<class Data >
int STK::BSplineCoefficients< Data >::degree ( ) const
inline
Returns
the degree of the B-Splines

Definition at line 118 of file STK_BSplineCoefficients.h.

118{ return degree_;}

References STK::BSplineCoefficients< Data >::degree_.

◆ extrapolate()

template<class Data >
template<class OtherVector >
ArrayXX STK::BSplineCoefficients< Data >::extrapolate ( OtherVector const x) const
Returns
the extrapolated matrix of coefficients for a given set of x-values.
Parameters
xthe values to extrapolate

Definition at line 354 of file STK_BSplineCoefficients.h.

355{
356 // check if knots exists
357 if (!this->hasRun_)
359 // resize coeficients
360 ArrayXX coefs(x.range(), Range(0, lastControlPoint_, 0), 0.0);
361 // check if the original data set was not reduced to a single point
362 if (minValue_ == maxValue_) return coefs;
363 // compute the coefficients
364 for (int irow=x.begin(); irow< x.end(); irow++)
365 {
366 const Real value = x[irow];
367 // value outside the range of the knots case
368 if (value <= minValue_)
369 {
370 coefs(irow, 0) = 1.0;
371 continue;
372 }
373 if (value >= maxValue_)
374 {
375 coefs(irow, lastControlPoint_) = 1.0;
376 continue;
377 }
378 // find interval
379 int k, k1;
380 for (k=0, k1=1; k<lastControlPoint_; k++, k1++)
381 {
382 if (value < knots_[k1]) break;
383 }
384 // begin recursion
385 coefs(irow, k) = 1.0;
386 for (int d=1; d<=degree_; d++)
387 {
388 // right (south-west corner) term only
389 coefs(irow, k-d) = ( (knots_[k1] - value)/(knots_[k1] - knots_[k1-d]) ) * coefs(irow, k1-d);
390 // compute internal terms
391 for (int i = k1-d; i<k; i++)
392 {
393 const Real knots_i = knots_[i], knots_id1 = knots_[i+d+1];
394 coefs(irow, i) = ( (value - knots_i)/(knots_[i+d] - knots_i) ) * coefs(irow, i)
395 + ( (knots_id1 - value)/(knots_id1 - knots_[i+1]) ) * coefs(irow, i+1);
396 }
397 // left (north-west corner) term only
398 coefs(irow, k) *= (value - knots_[k])/(knots_[k+d] - knots_[k]);
399 }
400 }
401 return coefs;
402}
#define STKRUNTIME_ERROR_NO_ARG(Where, Error)
Definition STK_Macros.h:138
virtual bool run()
run the computations.
ArrayXX extrapolate(OtherVector const &x) const
bool hasRun_
true if run has been used, false otherwise
Definition STK_IRunner.h:98
Array2D< Real > ArrayXX
Specialization of the Array2D class for Real values.
Definition STK_Array2D.h:53

References STK::MultidimRegression< Array, Weight >::coefs(), STK::BSplineCoefficients< Data >::extrapolate(), STK::IRunnerBase::hasRun_, STK::IRegression< Array, Array, Weight >::run(), and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::BSplineCoefficients< Data >::extrapolate().

◆ knots()

template<class Data >
VectorX const & STK::BSplineCoefficients< Data >::knots ( ) const
inline
Returns
the vector of knots of the B-spline curve

Definition at line 124 of file STK_BSplineCoefficients.h.

124{ return knots_;}

References STK::BSplineCoefficients< Data >::knots_.

Referenced by STK::BSplineRegression< YArray, XVector, Weights >::knots().

◆ nbControlPoints()

template<class Data >
int STK::BSplineCoefficients< Data >::nbControlPoints ( ) const
inline
Returns
the number of control points of the curve (as dim())

Definition at line 122 of file STK_BSplineCoefficients.h.

122{ return dim_;}

References STK::IBasis< Data, ArrayXX >::dim_.

◆ nbKnots()

template<class Data >
int STK::BSplineCoefficients< Data >::nbKnots ( ) const
inline
Returns
the number of knots of the B-spline curve

Definition at line 120 of file STK_BSplineCoefficients.h.

120{ return nbKnots_;}

References STK::BSplineCoefficients< Data >::nbKnots_.

◆ run()

template<class Data >
bool STK::BSplineCoefficients< Data >::run ( )
virtual

run the computations.

Implements STK::IRunnerBase.

Definition at line 331 of file STK_BSplineCoefficients.h.

332{
333 // check if data exists
334 if (!p_data_)
335 {
337 return false;
338 }
339 if (!this->initializeStep()) return false;
341 // compute the knots and coefficients
343 else { return false;}
344 this->hasRun_ = true;
345 return true;
346}
void computeCoefficients()
Compute the coefficients of the B-spline curves.
bool initializeStep()
Initialize the parameters minValue_ and maxValue_ using data set.
Definition STK_IBasis.h:181

References STK::IRunnerBase::hasRun_, STK::IRegression< Array, Array, Weight >::initializeStep(), STK::IRunnerBase::msg_error_, STK::BSplineCoefficients< Data >::run(), and STKERROR_NO_ARG.

Referenced by STK::BSplineCoefficients< Data >::run().

◆ setData()

template<class Data >
void STK::BSplineCoefficients< Data >::setData ( Data const data,
int  nbControlPoints,
int  degree = 3,
Regress::KnotsPosition  position = Regress::uniformKnotsPositions_ 
)

[DEPRECATED] Set the whole set of parameters for computing the coefficients

Parameters
datathe input data values
nbControlPointsnumber of control points
degreedegree of the B-spline curves (default is 3)
positionmethod to use for positioning the knots (default is uniform)

Definition at line 257 of file STK_BSplineCoefficients.h.

262{ // set data
263 p_data_ = &data;
265 degree_ = degree;
266 position_ = position;
267 nbKnots_ = dim_ + degree_ +1;
269 Base::update();
270}
virtual void update()
update IBasis if a parameter or a new data set is set, update the state of this runner.
Definition STK_IBasis.h:170

References STK::IRunnerBase::update().

◆ setDegree()

template<class Data >
void STK::BSplineCoefficients< Data >::setDegree ( int  degree = 3)

Set the degree of the BSpline basis.

Parameters
degreedegree of the B-spline curves (default is 3)

Definition at line 309 of file STK_BSplineCoefficients.h.

310{
311 degree_ = degree;
312 nbKnots_ = dim_ + degree_ +1;
315 Base::update();
316}

References STK::IRunnerBase::update().

◆ setDim()

template<class Data >
void STK::BSplineCoefficients< Data >::setDim ( int  dim)
inline

Alias for setting the number of control point (the number of BSpline)

Parameters
dimnumber of control points

Definition at line 133 of file STK_BSplineCoefficients.h.

void setNbControlPoints(int nbControlPoints)
Set the number of control point (the number of BSpline)

References STK::IBasis< Data, ArrayXX >::dim(), and STK::BSplineCoefficients< Data >::setNbControlPoints().

◆ setNbControlPoints()

template<class Data >
void STK::BSplineCoefficients< Data >::setNbControlPoints ( int  nbControlPoints)

Set the number of control point (the number of BSpline)

Parameters
nbControlPointsnumber of control points

Definition at line 296 of file STK_BSplineCoefficients.h.

297{
299 nbKnots_ = dim_ + degree_ +1;
302 Base::update();
303}

References STK::IRunnerBase::update().

Referenced by STK::BSplineCoefficients< Data >::setDim().

◆ setParameters()

template<class Data >
void STK::BSplineCoefficients< Data >::setParameters ( Data const data,
int  nbControlPoints,
int  degree = 3,
Regress::KnotsPosition  position = Regress::uniformKnotsPositions_ 
)

Set the whole set of parameters for computing the coefficients.

Parameters
datathe input data values
nbControlPointsnumber of control points
degreedegree of the B-spline curves (default is 3)
positionmethod to use for positioning the knots (default is uniform)

Definition at line 276 of file STK_BSplineCoefficients.h.

281{ // set data
282 p_data_ = &data;
284 degree_ = degree;
285 position_ = position;
286 nbKnots_ = dim_ + degree_ +1;
289 Base::update();
290}

References STK::IRunnerBase::update().

◆ setPosition()

template<class Data >
void STK::BSplineCoefficients< Data >::setPosition ( Regress::KnotsPosition  position)

Set the kind of position to use for the knots.

Parameters
positionmethod to use for positioning the knots (default is uniform)

Definition at line 322 of file STK_BSplineCoefficients.h.

323{
324 position_ = position;
325 Base::update();
326}

References STK::IRunnerBase::update().

Member Data Documentation

◆ degree_

template<class Data >
int STK::BSplineCoefficients< Data >::degree_
protected

degree of the B-splines curves.

Definition at line 174 of file STK_BSplineCoefficients.h.

Referenced by STK::BSplineCoefficients< Data >::degree().

◆ knots_

template<class Data >
VectorX STK::BSplineCoefficients< Data >::knots_
protected

Data of the knots.

Definition at line 184 of file STK_BSplineCoefficients.h.

Referenced by STK::BSplineCoefficients< Data >::knots().

◆ lastControlPoint_

template<class Data >
int STK::BSplineCoefficients< Data >::lastControlPoint_
protected

Index of the last control point of the B-spline curves.

This is dim_ - 1.

Definition at line 182 of file STK_BSplineCoefficients.h.

◆ nbKnots_

template<class Data >
int STK::BSplineCoefficients< Data >::nbKnots_
protected

number of knots of the B-spline curves.

Definition at line 178 of file STK_BSplineCoefficients.h.

Referenced by STK::BSplineCoefficients< Data >::nbKnots().

◆ position_

template<class Data >
Regress::KnotsPosition STK::BSplineCoefficients< Data >::position_
protected

Method used in order to position the knots.

Definition at line 176 of file STK_BSplineCoefficients.h.


The documentation for this class was generated from the following file: