Eigen::Spline< Scalar_, Dim_, Degree_ > Class Template Reference

A class representing multi-dimensional spline curves. More...

Public Types

enum  { Dimension }
 
enum  { Degree }
 
typedef SplineTraits< Spline >::BasisDerivativeType BasisDerivativeType
 The data type used to store the values of the basis function derivatives. More...
 
typedef SplineTraits< Spline >::BasisVectorType BasisVectorType
 The data type used to store non-zero basis functions. More...
 
typedef SplineTraits< Spline >::ControlPointVectorType ControlPointVectorType
 The data type representing the spline's control points. More...
 
typedef SplineTraits< Spline >::KnotVectorType KnotVectorType
 The data type used to store knot vectors. More...
 
typedef SplineTraits< Spline >::ParameterVectorType ParameterVectorType
 The data type used to store parameter vectors. More...
 
typedef SplineTraits< Spline >::PointType PointType
 The point type the spline is representing. More...
 
typedef Scalar_ Scalar
 

Public Member Functions

SplineTraits< Spline >::BasisDerivativeType basisFunctionDerivatives (Scalar u, DenseIndex order) const
 Computes the non-zero spline basis function derivatives up to given order. More...
 
template<int DerivativeOrder>
SplineTraits< Spline, DerivativeOrder >::BasisDerivativeType basisFunctionDerivatives (Scalar u, DenseIndex order=DerivativeOrder) const
 Computes the non-zero spline basis function derivatives up to given order. More...
 
SplineTraits< Spline >::BasisVectorType basisFunctions (Scalar u) const
 Computes the non-zero basis functions at the given site. More...
 
const ControlPointVectorTypectrls () const
 Returns the ctrls of the underlying spline. More...
 
DenseIndex degree () const
 Returns the spline degree. More...
 
SplineTraits< Spline >::DerivativeType derivatives (Scalar u, DenseIndex order) const
 Evaluation of spline derivatives of up-to given order. More...
 
template<int DerivativeOrder>
SplineTraits< Spline, DerivativeOrder >::DerivativeType derivatives (Scalar u, DenseIndex order=DerivativeOrder) const
 Evaluation of spline derivatives of up-to given order. More...
 
const KnotVectorTypeknots () const
 Returns the knots of the underlying spline. More...
 
PointType operator() (Scalar u) const
 Returns the spline value at a given site \(u\). More...
 
DenseIndex span (Scalar u) const
 Returns the span within the knot vector in which u is falling. More...
 
 Spline ()
 Creates a (constant) zero spline. For Splines with dynamic degree, the resulting degree will be 0. More...
 
template<typename OtherVectorType , typename OtherArrayType >
 Spline (const OtherVectorType &knots, const OtherArrayType &ctrls)
 Creates a spline from a knot vector and control points. More...
 
template<int OtherDegree>
 Spline (const Spline< Scalar, Dimension, OtherDegree > &spline)
 Copy constructor for splines. More...
 

Static Public Member Functions

static BasisDerivativeType BasisFunctionDerivatives (const Scalar u, const DenseIndex order, const DenseIndex degree, const KnotVectorType &knots)
 Computes the non-zero spline basis function derivatives up to given order. More...
 
static BasisVectorType BasisFunctions (Scalar u, DenseIndex degree, const KnotVectorType &knots)
 Returns the spline's non-zero basis functions. More...
 
static DenseIndex Span (typename SplineTraits< Spline >::Scalar u, DenseIndex degree, const typename SplineTraits< Spline >::KnotVectorType &knots)
 Computes the span within the provided knot vector in which u is falling. More...
 

Static Private Member Functions

template<typename DerivativeType >
static void BasisFunctionDerivativesImpl (const typename Spline< Scalar_, Dim_, Degree_ >::Scalar u, const DenseIndex order, const DenseIndex p, const typename Spline< Scalar_, Dim_, Degree_ >::KnotVectorType &U, DerivativeType &N_)
 

Private Attributes

ControlPointVectorType m_ctrls
 
KnotVectorType m_knots
 

Detailed Description

template<typename Scalar_, int Dim_, int Degree_>
class Eigen::Spline< Scalar_, Dim_, Degree_ >

A class representing multi-dimensional spline curves.

The class represents B-splines with non-uniform knot vectors. Each control point of the B-spline is associated with a basis function

\begin{align*} C(u) & = \sum_{i=0}^{n}N_{i,p}(u)P_i \end{align*}

Template Parameters
Scalar_The underlying data type (typically float or double)
Dim_The curve dimension (e.g. 2 or 3)
Degree_Per default set to Dynamic; could be set to the actual desired degree for optimization purposes (would result in stack allocation of several temporary variables).

Definition at line 37 of file Spline.h.

Member Typedef Documentation

◆ BasisDerivativeType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::BasisDerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::BasisDerivativeType

The data type used to store the values of the basis function derivatives.

Definition at line 57 of file Spline.h.

◆ BasisVectorType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::BasisVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::BasisVectorType

The data type used to store non-zero basis functions.

Definition at line 54 of file Spline.h.

◆ ControlPointVectorType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::ControlPointVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::ControlPointVectorType

The data type representing the spline's control points.

Definition at line 60 of file Spline.h.

◆ KnotVectorType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::KnotVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::KnotVectorType

The data type used to store knot vectors.

Definition at line 48 of file Spline.h.

◆ ParameterVectorType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::ParameterVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::ParameterVectorType

The data type used to store parameter vectors.

Definition at line 51 of file Spline.h.

◆ PointType

template<typename Scalar_ , int Dim_, int Degree_>
typedef SplineTraits<Spline>::PointType Eigen::Spline< Scalar_, Dim_, Degree_ >::PointType

The point type the spline is representing.

Definition at line 45 of file Spline.h.

◆ Scalar

template<typename Scalar_ , int Dim_, int Degree_>
typedef Scalar_ Eigen::Spline< Scalar_, Dim_, Degree_ >::Scalar

The spline curve's scalar type.

Definition at line 40 of file Spline.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ , int Dim_, int Degree_>
anonymous enum
Enumerator
Dimension 

The spline curve's dimension.

Definition at line 41 of file Spline.h.

41 { Dimension = Dim_ };

◆ anonymous enum

template<typename Scalar_ , int Dim_, int Degree_>
anonymous enum
Enumerator
Degree 

The spline curve's degree.

Definition at line 42 of file Spline.h.

42 { Degree = Degree_ };

Constructor & Destructor Documentation

◆ Spline() [1/3]

template<typename Scalar_ , int Dim_, int Degree_>
Eigen::Spline< Scalar_, Dim_, Degree_ >::Spline ( )
inline

Creates a (constant) zero spline. For Splines with dynamic degree, the resulting degree will be 0.

Definition at line 66 of file Spline.h.

67  : m_knots(1, (Degree==Dynamic ? 2 : 2*Degree+2))
68  , m_ctrls(ControlPointVectorType::Zero(Dimension,(Degree==Dynamic ? 1 : Degree+1)))
69  {
70  // in theory this code can go to the initializer list but it will get pretty
71  // much unreadable ...
72  enum { MinDegree = (Degree==Dynamic ? 0 : Degree) };
73  m_knots.template segment<MinDegree+1>(0) = Array<Scalar,1,MinDegree+1>::Zero();
74  m_knots.template segment<MinDegree+1>(MinDegree+1) = Array<Scalar,1,MinDegree+1>::Ones();
75  }
static const ConstantReturnType Ones()
static const ConstantReturnType Zero()
KnotVectorType m_knots
Definition: Spline.h:223
ControlPointVectorType m_ctrls
Definition: Spline.h:224
const int Dynamic

◆ Spline() [2/3]

template<typename Scalar_ , int Dim_, int Degree_>
template<typename OtherVectorType , typename OtherArrayType >
Eigen::Spline< Scalar_, Dim_, Degree_ >::Spline ( const OtherVectorType &  knots,
const OtherArrayType &  ctrls 
)
inline

Creates a spline from a knot vector and control points.

Parameters
knotsThe spline's knot vector.
ctrlsThe spline's control point vector.

Definition at line 83 of file Spline.h.

83 : m_knots(knots), m_ctrls(ctrls) {}
const ControlPointVectorType & ctrls() const
Returns the ctrls of the underlying spline.
Definition: Spline.h:101
const KnotVectorType & knots() const
Returns the knots of the underlying spline.
Definition: Spline.h:96

◆ Spline() [3/3]

template<typename Scalar_ , int Dim_, int Degree_>
template<int OtherDegree>
Eigen::Spline< Scalar_, Dim_, Degree_ >::Spline ( const Spline< Scalar, Dimension, OtherDegree > &  spline)
inline

Copy constructor for splines.

Parameters
splineThe input spline.

Definition at line 90 of file Spline.h.

90  :
91  m_knots(spline.knots()), m_ctrls(spline.ctrls()) {}

Member Function Documentation

◆ BasisFunctionDerivatives()

template<typename Scalar_ , int Dim_, int Degree_>
SplineTraits< Spline< Scalar_, Dim_, Degree_ > >::BasisDerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::BasisFunctionDerivatives ( const Scalar  u,
const DenseIndex  order,
const DenseIndex  degree,
const KnotVectorType knots 
)
static

Computes the non-zero spline basis function derivatives up to given order.

The function computes

\begin{align*} \frac{d^i}{du^i} N_{i,p}(u), \hdots, \frac{d^i}{du^i} N_{i+p+1,p}(u) \end{align*}

with i ranging from 0 up to the specified order.

Parameters
uParameter \(u \in [0;1]\) at which the non-zero basis function derivatives are computed.
orderThe order up to which the basis function derivatives are computes.
degreeThe degree of the underlying spline
knotsThe underlying spline's knot vector.

Definition at line 497 of file Spline.h.

502  {
503  typename SplineTraits<Spline>::BasisDerivativeType der;
504  BasisFunctionDerivativesImpl(u, order, degree, knots, der);
505  return der;
506  }
DenseIndex degree() const
Returns the spline degree.
Definition: Spline.h:282
static void BasisFunctionDerivativesImpl(const typename Spline< Scalar_, Dim_, Degree_ >::Scalar u, const DenseIndex order, const DenseIndex p, const typename Spline< Scalar_, Dim_, Degree_ >::KnotVectorType &U, DerivativeType &N_)
Definition: Spline.h:373

◆ basisFunctionDerivatives() [1/2]

template<typename Scalar_ , int Dim_, int Degree_>
SplineTraits< Spline< Scalar_, Dim_, Degree_ >, DerivativeOrder >::BasisDerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::basisFunctionDerivatives ( Scalar  u,
DenseIndex  order 
) const

Computes the non-zero spline basis function derivatives up to given order.

The function computes

\begin{align*} \frac{d^i}{du^i} N_{i,p}(u), \hdots, \frac{d^i}{du^i} N_{i+p+1,p}(u) \end{align*}

with i ranging from 0 up to the specified order.

Parameters
uParameter \(u \in [0;1]\) at which the non-zero basis function derivatives are computed.
orderThe order up to which the basis function derivatives are computes.

Definition at line 478 of file Spline.h.

479  {
480  typename SplineTraits<Spline<Scalar_, Dim_, Degree_> >::BasisDerivativeType der;
481  BasisFunctionDerivativesImpl(u, order, degree(), knots(), der);
482  return der;
483  }
SplineTraits< Spline >::BasisDerivativeType BasisDerivativeType
The data type used to store the values of the basis function derivatives.
Definition: Spline.h:57

◆ basisFunctionDerivatives() [2/2]

template<typename Scalar_ , int Dim_, int Degree_>
template<int DerivativeOrder>
SplineTraits<Spline,DerivativeOrder>::BasisDerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::basisFunctionDerivatives ( Scalar  u,
DenseIndex  order = DerivativeOrder 
) const

Computes the non-zero spline basis function derivatives up to given order.

The function computes

\begin{align*} \frac{d^i}{du^i} N_{i,p}(u), \hdots, \frac{d^i}{du^i} N_{i+p+1,p}(u) \end{align*}

with i ranging from 0 up to the specified order.

Parameters
uParameter \(u \in [0;1]\) at which the non-zero basis function derivatives are computed.
orderThe order up to which the basis function derivatives are computes. Using the template version of this function is more efficieent since temporary objects are allocated on the stack whenever this is possible.

◆ BasisFunctionDerivativesImpl()

template<typename Scalar_ , int Dim_, int Degree_>
template<typename DerivativeType >
void Eigen::Spline< Scalar_, Dim_, Degree_ >::BasisFunctionDerivativesImpl ( const typename Spline< Scalar_, Dim_, Degree_ >::Scalar  u,
const DenseIndex  order,
const DenseIndex  p,
const typename Spline< Scalar_, Dim_, Degree_ >::KnotVectorType U,
DerivativeType &  N_ 
)
staticprivate

Definition at line 373 of file Spline.h.

379  {
380  typedef Spline<Scalar_, Dim_, Degree_> SplineType;
381  enum { Order = SplineTraits<SplineType>::OrderAtCompileTime };
382 
383  const DenseIndex span = SplineType::Span(u, p, U);
384 
385  const DenseIndex n = (std::min)(p, order);
386 
387  N_.resize(n+1, p+1);
388 
389  BasisVectorType left = BasisVectorType::Zero(p+1);
390  BasisVectorType right = BasisVectorType::Zero(p+1);
391 
392  Matrix<Scalar,Order,Order> ndu(p+1,p+1);
393 
394  Scalar saved, temp; // FIXME These were double instead of Scalar. Was there a reason for that?
395 
396  ndu(0,0) = 1.0;
397 
398  DenseIndex j;
399  for (j=1; j<=p; ++j)
400  {
401  left[j] = u-U[span+1-j];
402  right[j] = U[span+j]-u;
403  saved = 0.0;
404 
405  for (DenseIndex r=0; r<j; ++r)
406  {
407  /* Lower triangle */
408  ndu(j,r) = right[r+1]+left[j-r];
409  temp = ndu(r,j-1)/ndu(j,r);
410  /* Upper triangle */
411  ndu(r,j) = static_cast<Scalar>(saved+right[r+1] * temp);
412  saved = left[j-r] * temp;
413  }
414 
415  ndu(j,j) = static_cast<Scalar>(saved);
416  }
417 
418  for (j = p; j>=0; --j)
419  N_(0,j) = ndu(j,p);
420 
421  // Compute the derivatives
422  DerivativeType a(n+1,p+1);
423  DenseIndex r=0;
424  for (; r<=p; ++r)
425  {
426  DenseIndex s1,s2;
427  s1 = 0; s2 = 1; // alternate rows in array a
428  a(0,0) = 1.0;
429 
430  // Compute the k-th derivative
431  for (DenseIndex k=1; k<=static_cast<DenseIndex>(n); ++k)
432  {
433  Scalar d = 0.0;
434  DenseIndex rk,pk,j1,j2;
435  rk = r-k; pk = p-k;
436 
437  if (r>=k)
438  {
439  a(s2,0) = a(s1,0)/ndu(pk+1,rk);
440  d = a(s2,0)*ndu(rk,pk);
441  }
442 
443  if (rk>=-1) j1 = 1;
444  else j1 = -rk;
445 
446  if (r-1 <= pk) j2 = k-1;
447  else j2 = p-r;
448 
449  for (j=j1; j<=j2; ++j)
450  {
451  a(s2,j) = (a(s1,j)-a(s1,j-1))/ndu(pk+1,rk+j);
452  d += a(s2,j)*ndu(rk+j,pk);
453  }
454 
455  if (r<=pk)
456  {
457  a(s2,k) = -a(s1,k-1)/ndu(pk+1,r);
458  d += a(s2,k)*ndu(r,pk);
459  }
460 
461  N_(k,r) = static_cast<Scalar>(d);
462  j = s1; s1 = s2; s2 = j; // Switch rows
463  }
464  }
465 
466  /* Multiply through by the correct factors */
467  /* (Eq. [2.9]) */
468  r = p;
469  for (DenseIndex k=1; k<=static_cast<DenseIndex>(n); ++k)
470  {
471  for (j=p; j>=0; --j) N_(k,j) *= r;
472  r *= p-k;
473  }
474  }
ArrayXXi a
int n
float * p
DenseIndex span(Scalar u) const
Returns the span within the knot vector in which u is falling.
Definition: Spline.h:291
SplineTraits< Spline >::BasisVectorType BasisVectorType
The data type used to store non-zero basis functions.
Definition: Spline.h:54
Scalar_ Scalar
Definition: Spline.h:40
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
std::ptrdiff_t j

◆ basisFunctions()

template<typename Scalar_ , int Dim_, int Degree_>
SplineTraits< Spline< Scalar_, Dim_, Degree_ > >::BasisVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::basisFunctions ( Scalar  u) const

Computes the non-zero basis functions at the given site.

Splines have local support and a point from their image is defined by exactly \(p+1\) control points \(P_i\) where \(p\) is the spline degree.

This function computes the \(p+1\) non-zero basis function values for a given parameter value \(u\). It returns

\begin{align*} N_{i,p}(u), \hdots, N_{i+p+1,p}(u) \end{align*}

Parameters
uParameter \(u \in [0;1]\) at which the non-zero basis functions are computed.

Definition at line 363 of file Spline.h.

364  {
365  return Spline::BasisFunctions(u, degree(), knots());
366  }
static BasisVectorType BasisFunctions(Scalar u, DenseIndex degree, const KnotVectorType &knots)
Returns the spline's non-zero basis functions.
Definition: Spline.h:249

◆ BasisFunctions()

template<typename Scalar_ , int Dim_, int Degree_>
Spline< Scalar_, Dim_, Degree_ >::BasisVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::BasisFunctions ( Scalar  u,
DenseIndex  degree,
const KnotVectorType knots 
)
static

Returns the spline's non-zero basis functions.

The function computes and returns

\begin{align*} N_{i,p}(u), \hdots, N_{i+p+1,p}(u) \end{align*}

Parameters
uThe site at which the basis functions are computed.
degreeThe degree of the underlying spline.
knotsThe underlying spline's knot vector.

Definition at line 249 of file Spline.h.

253  {
254  const DenseIndex p = degree;
255  const DenseIndex i = Spline::Span(u, degree, knots);
256 
257  const KnotVectorType& U = knots;
258 
259  BasisVectorType left(p+1); left(0) = Scalar(0);
260  BasisVectorType right(p+1); right(0) = Scalar(0);
261 
262  VectorBlock<BasisVectorType,Degree>(left,1,p) = u - VectorBlock<const KnotVectorType,Degree>(U,i+1-p,p).reverse();
263  VectorBlock<BasisVectorType,Degree>(right,1,p) = VectorBlock<const KnotVectorType,Degree>(U,i+1,p) - u;
264 
265  BasisVectorType N(1,p+1);
266  N(0) = Scalar(1);
267  for (DenseIndex j=1; j<=p; ++j)
268  {
269  Scalar saved = Scalar(0);
270  for (DenseIndex r=0; r<j; r++)
271  {
272  const Scalar tmp = N(r)/(right(r+1)+left(j-r));
273  N[r] = saved + right(r+1)*tmp;
274  saved = left(j-r)*tmp;
275  }
276  N(j) = saved;
277  }
278  return N;
279  }
int i
static DenseIndex Span(typename SplineTraits< Spline >::Scalar u, DenseIndex degree, const typename SplineTraits< Spline >::KnotVectorType &knots)
Computes the span within the provided knot vector in which u is falling.
Definition: Spline.h:236
SplineTraits< Spline >::KnotVectorType KnotVectorType
The data type used to store knot vectors.
Definition: Spline.h:48

◆ ctrls()

template<typename Scalar_ , int Dim_, int Degree_>
const ControlPointVectorType& Eigen::Spline< Scalar_, Dim_, Degree_ >::ctrls ( ) const
inline

Returns the ctrls of the underlying spline.

Definition at line 101 of file Spline.h.

101 { return m_ctrls; }

◆ degree()

template<typename Scalar_ , int Dim_, int Degree_>
DenseIndex Eigen::Spline< Scalar_, Dim_, Degree_ >::degree

Returns the spline degree.

Definition at line 282 of file Spline.h.

283  {
284  if (Degree_ == Dynamic)
285  return m_knots.size() - m_ctrls.cols() - 1;
286  else
287  return Degree_;
288  }

◆ derivatives() [1/2]

template<typename Scalar_ , int Dim_, int Degree_>
SplineTraits< Spline< Scalar_, Dim_, Degree_ >, DerivativeOrder >::DerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::derivatives ( Scalar  u,
DenseIndex  order 
) const

Evaluation of spline derivatives of up-to given order.

The function returns

\begin{align*} \frac{d^i}{du^i}C(u) & = \sum_{i=0}^{n} \frac{d^i}{du^i} N_{i,p}(u)P_i \end{align*}

for i ranging between 0 and order.

Parameters
uParameter \(u \in [0;1]\) at which the spline derivative is evaluated.
orderThe order up to which the derivatives are computed.

Definition at line 344 of file Spline.h.

345  {
346  typename SplineTraits< Spline >::DerivativeType res;
347  derivativesImpl(*this, u, order, res);
348  return res;
349  }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
void derivativesImpl(const SplineType &spline, typename SplineType::Scalar u, DenseIndex order, DerivativeType &der)
Definition: Spline.h:313

◆ derivatives() [2/2]

template<typename Scalar_ , int Dim_, int Degree_>
template<int DerivativeOrder>
SplineTraits<Spline,DerivativeOrder>::DerivativeType Eigen::Spline< Scalar_, Dim_, Degree_ >::derivatives ( Scalar  u,
DenseIndex  order = DerivativeOrder 
) const

Evaluation of spline derivatives of up-to given order.

The function returns

\begin{align*} \frac{d^i}{du^i}C(u) & = \sum_{i=0}^{n} \frac{d^i}{du^i} N_{i,p}(u)P_i \end{align*}

for i ranging between 0 and order.

Parameters
uParameter \(u \in [0;1]\) at which the spline derivative is evaluated.
orderThe order up to which the derivatives are computed. Using the template version of this function is more efficieent since temporary objects are allocated on the stack whenever this is possible.

◆ knots()

template<typename Scalar_ , int Dim_, int Degree_>
const KnotVectorType& Eigen::Spline< Scalar_, Dim_, Degree_ >::knots ( ) const
inline

Returns the knots of the underlying spline.

Definition at line 96 of file Spline.h.

96 { return m_knots; }

◆ operator()()

template<typename Scalar_ , int Dim_, int Degree_>
Spline< Scalar_, Dim_, Degree_ >::PointType Eigen::Spline< Scalar_, Dim_, Degree_ >::operator() ( Scalar  u) const

Returns the spline value at a given site \(u\).

The function returns

\begin{align*} C(u) & = \sum_{i=0}^{n}N_{i,p}P_i \end{align*}

Parameters
uParameter \(u \in [0;1]\) at which the spline is evaluated.
Returns
The spline value at the given location \(u\).

Definition at line 297 of file Spline.h.

298  {
299  enum { Order = SplineTraits<Spline>::OrderAtCompileTime };
300 
301  const DenseIndex span = this->span(u);
302  const DenseIndex p = degree();
303  const BasisVectorType basis_funcs = basisFunctions(u);
304 
305  const Replicate<BasisVectorType,Dimension,1> ctrl_weights(basis_funcs);
306  const Block<const ControlPointVectorType,Dimension,Order> ctrl_pts(ctrls(),0,span-p,Dimension,p+1);
307  return (ctrl_weights * ctrl_pts).rowwise().sum();
308  }
SplineTraits< Spline >::BasisVectorType basisFunctions(Scalar u) const
Computes the non-zero basis functions at the given site.
Definition: Spline.h:363

◆ span()

template<typename Scalar_ , int Dim_, int Degree_>
DenseIndex Eigen::Spline< Scalar_, Dim_, Degree_ >::span ( Scalar  u) const

Returns the span within the knot vector in which u is falling.

Parameters
uThe site for which the span is determined.

Definition at line 291 of file Spline.h.

292  {
293  return Spline::Span(u, degree(), knots());
294  }

◆ Span()

template<typename Scalar_ , int Dim_, int Degree_>
DenseIndex Eigen::Spline< Scalar_, Dim_, Degree_ >::Span ( typename SplineTraits< Spline< Scalar_, Dim_, Degree_ > >::Scalar  u,
DenseIndex  degree,
const typename SplineTraits< Spline< Scalar_, Dim_, Degree_ > >::KnotVectorType knots 
)
static

Computes the span within the provided knot vector in which u is falling.

Definition at line 236 of file Spline.h.

240  {
241  // Piegl & Tiller, "The NURBS Book", A2.1 (p. 68)
242  if (u <= knots(0)) return degree;
243  const Scalar* pos = std::upper_bound(knots.data()+degree-1, knots.data()+knots.size()-degree-1, u);
244  return static_cast<DenseIndex>( std::distance(knots.data(), pos) - 1 );
245  }

Member Data Documentation

◆ m_ctrls

template<typename Scalar_ , int Dim_, int Degree_>
ControlPointVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::m_ctrls
private

Control points.

Definition at line 224 of file Spline.h.

◆ m_knots

template<typename Scalar_ , int Dim_, int Degree_>
KnotVectorType Eigen::Spline< Scalar_, Dim_, Degree_ >::m_knots
private

Knot vector.

Definition at line 223 of file Spline.h.


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