From b00a33bc70f2f94ef9922e516ca095f8dff245e1 Mon Sep 17 00:00:00 2001 From: Hauke Heibel Date: Fri, 25 Nov 2011 14:53:40 +0100 Subject: Integrated spline class and simple spline fitting --- unsupported/Eigen/CMakeLists.txt | 2 +- unsupported/Eigen/Splines | 52 ++++ unsupported/Eigen/src/CMakeLists.txt | 3 +- unsupported/Eigen/src/Splines/Spline.h | 407 ++++++++++++++++++++++++++ unsupported/Eigen/src/Splines/SplineFitting.h | 115 ++++++++ unsupported/Eigen/src/Splines/SplineFwd.h | 78 +++++ unsupported/test/CMakeLists.txt | 1 + unsupported/test/splines.cpp | 239 +++++++++++++++ 8 files changed, 895 insertions(+), 2 deletions(-) create mode 100644 unsupported/Eigen/Splines create mode 100644 unsupported/Eigen/src/Splines/Spline.h create mode 100644 unsupported/Eigen/src/Splines/SplineFitting.h create mode 100644 unsupported/Eigen/src/Splines/SplineFwd.h create mode 100644 unsupported/test/splines.cpp diff --git a/unsupported/Eigen/CMakeLists.txt b/unsupported/Eigen/CMakeLists.txt index 679c96736..e961e72c5 100644 --- a/unsupported/Eigen/CMakeLists.txt +++ b/unsupported/Eigen/CMakeLists.txt @@ -1,6 +1,6 @@ set(Eigen_HEADERS AdolcForward BVH IterativeSolvers MatrixFunctions MoreVectorization AutoDiff AlignedVector3 Polynomials FFT NonLinearOptimization SparseExtra IterativeSolvers - NumericalDiff Skyline MPRealSupport OpenGLSupport KroneckerProduct + NumericalDiff Skyline MPRealSupport OpenGLSupport KroneckerProduct Splines ) install(FILES diff --git a/unsupported/Eigen/Splines b/unsupported/Eigen/Splines new file mode 100644 index 000000000..362274157 --- /dev/null +++ b/unsupported/Eigen/Splines @@ -0,0 +1,52 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 20010-2011 Hauke Heibel +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see . + +#ifndef EIGEN_SPLINES_MODULE_H +#define EIGEN_SPLINES_MODULE_H + +namespace Eigen +{ +/** \ingroup Unsupported_modules + * \defgroup Splines_Module Spline and spline fitting module + * + * This module provides a simple multi-dimensional spline class while + * offering most basic functionality to fit a spline to point sets. + * + * \code + * #include + * \endcode + */ +//@{ +} + +#include "src/Splines/SplineFwd.h" +#include "src/Splines/Spline.h" +#include "src/Splines/SplineFitting.h" + +namespace Eigen +{ +//@} +} + +#endif // EIGEN_SPLINES_MODULE_H diff --git a/unsupported/Eigen/src/CMakeLists.txt b/unsupported/Eigen/src/CMakeLists.txt index 6fd98b310..f3180b52b 100644 --- a/unsupported/Eigen/src/CMakeLists.txt +++ b/unsupported/Eigen/src/CMakeLists.txt @@ -9,4 +9,5 @@ ADD_SUBDIRECTORY(NumericalDiff) ADD_SUBDIRECTORY(Polynomials) ADD_SUBDIRECTORY(Skyline) ADD_SUBDIRECTORY(SparseExtra) -ADD_SUBDIRECTORY(KroneckerProduct) \ No newline at end of file +ADD_SUBDIRECTORY(KroneckerProduct) +ADD_SUBDIRECTORY(Splines) diff --git a/unsupported/Eigen/src/Splines/Spline.h b/unsupported/Eigen/src/Splines/Spline.h new file mode 100644 index 000000000..bb4b73e5a --- /dev/null +++ b/unsupported/Eigen/src/Splines/Spline.h @@ -0,0 +1,407 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 20010-2011 Hauke Heibel +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see . + +#ifndef EIGEN_SPLINE_H +#define EIGEN_SPLINE_H + +#include "SplineFwd.h" + +namespace Eigen +{ + /** + * \class Spline class + * \brief A class representing N-D spline curves. + * \tparam _Scalar The underlying data type (typically float or double) + * \tparam _Dim The curve dimension (e.g. 2 or 3) + * \tparam _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). + **/ + template + class Spline + { + public: + typedef _Scalar Scalar; + enum { Dimension = _Dim }; + enum { Degree = _Degree }; + + typedef typename SplineTraits::PointType PointType; + typedef typename SplineTraits::KnotVectorType KnotVectorType; + typedef typename SplineTraits::BasisVectorType BasisVectorType; + typedef typename SplineTraits::ControlPointVectorType ControlPointVectorType; + + /** + * \brief Creates a spline from a knot vector and control points. + **/ + template + Spline(const OtherVectorType& knots, const OtherArrayType& ctrls) : m_knots(knots), m_ctrls(ctrls) {} + + template + Spline(const Spline& spline) : + m_knots(spline.knots()), m_ctrls(spline.ctrls()) {} + + /* Const access methods for knots and control points. */ + const KnotVectorType& knots() const { return m_knots; } + const ControlPointVectorType& ctrls() const { return m_ctrls; } + + /* Spline evaluation. */ + PointType operator()(Scalar u) const; + + /* Evaluation of spline derivatives of up-to given order. + * The returned matrix has dimensions Dim-by-(Order+1) containing + * the 0-th order up-to Order-th order derivatives. + */ + typename SplineTraits::DerivativeType + derivatives(Scalar u, DenseIndex order) const; + + /** + * Evaluation of spline derivatives of up-to given order. + * The function performs identically to derivatives(Scalar, int) but + * does not require any heap allocations. + * \sa derivatives(Scalar, int) + **/ + template + typename SplineTraits::DerivativeType + derivatives(Scalar u, DenseIndex order = DerivativeOrder) const; + + /* Non-zero spline basis functions. */ + typename SplineTraits::BasisVectorType + basisFunctions(Scalar u) const; + + /* Non-zero spline basis function derivatives up to given order. + * The order is different from the spline order - it is the order + * up to which derivatives will be computed. + * \sa basisFunctions(Scalar) + */ + typename SplineTraits::BasisDerivativeType + basisFunctionDerivatives(Scalar u, DenseIndex order) const; + + /** + * Computes non-zero basis function derivatives up to the given derivative order. + * As opposed to basisFunctionDerivatives(Scalar, int) this function does not perform + * any heap allocations. + * \sa basisFunctionDerivatives(Scalar, int) + **/ + template + typename SplineTraits::BasisDerivativeType + basisFunctionDerivatives(Scalar u, DenseIndex order = DerivativeOrder) const; + + /** + * \brief The current spline degree. It's a function of knot size and number + * of controls and thus does not require a dedicated member. + */ + DenseIndex degree() const; + + /** Computes the span within the knot vector in which u falls. */ + DenseIndex span(Scalar u) const; + + static DenseIndex Span(typename SplineTraits::Scalar u, DenseIndex degree, const typename SplineTraits::KnotVectorType& knots); + static BasisVectorType BasisFunctions(Scalar u, DenseIndex degree, const KnotVectorType& knots); + + + private: + KnotVectorType m_knots; /* Knot vector. */ + ControlPointVectorType m_ctrls; /* Control points. */ + }; + + template + DenseIndex Spline<_Scalar, _Dim, _Degree>::Span( + typename SplineTraits< Spline<_Scalar, _Dim, _Degree> >::Scalar u, + DenseIndex degree, + const typename SplineTraits< Spline<_Scalar, _Dim, _Degree> >::KnotVectorType& knots) + { + // Piegl & Tiller, "The NURBS Book", A2.1 (p. 68) + if (u <= knots(0)) return degree; + const Scalar* pos = std::upper_bound(knots.data()+degree-1, knots.data()+knots.size()-degree-1, u); + return static_cast( std::distance(knots.data(), pos) - 1 ); + } + + template + typename Spline<_Scalar, _Dim, _Degree>::BasisVectorType + Spline<_Scalar, _Dim, _Degree>::BasisFunctions( + typename Spline<_Scalar, _Dim, _Degree>::Scalar u, + DenseIndex degree, + const typename Spline<_Scalar, _Dim, _Degree>::KnotVectorType& knots) + { + typedef typename Spline<_Scalar, _Dim, _Degree>::BasisVectorType BasisVectorType; + + const DenseIndex p = degree; + const DenseIndex i = Spline::Span(u, degree, knots); + + const KnotVectorType& U = knots; + + BasisVectorType left(p+1); left(0) = Scalar(0); + BasisVectorType right(p+1); right(0) = Scalar(0); + + VectorBlock(left,1,p) = u - VectorBlock(U,i+1-p,p).reverse(); + VectorBlock(right,1,p) = VectorBlock(U,i+1,p) - u; + + BasisVectorType N(1,p+1); + N(0) = Scalar(1); + for (DenseIndex j=1; j<=p; ++j) + { + Scalar saved = Scalar(0); + for (DenseIndex r=0; r + DenseIndex Spline<_Scalar, _Dim, _Degree>::degree() const + { + if (_Degree == Dynamic) + return m_knots.size() - m_ctrls.cols() - 1; + else + return _Degree; + } + + template + DenseIndex Spline<_Scalar, _Dim, _Degree>::span(Scalar u) const + { + return Spline::Span(u, degree(), knots()); + } + + /** + * \brief A functor for the computation of a spline point. + * \sa Piegl & Tiller, "The NURBS Book", A4.1 (p. 124) + **/ + template + typename Spline<_Scalar, _Dim, _Degree>::PointType Spline<_Scalar, _Dim, _Degree>::operator()(Scalar u) const + { + enum { Order = SplineTraits::OrderAtCompileTime }; + + const DenseIndex span = this->span(u); + const DenseIndex p = degree(); + const BasisVectorType basis_funcs = basisFunctions(u); + + const Replicate ctrl_weights(basis_funcs); + const Block ctrl_pts(ctrls(),0,span-p,Dimension,p+1); + return (ctrl_weights * ctrl_pts).rowwise().sum(); + } + + /* --------------------------------------------------------------------------------------------- */ + + template + void derivativesImpl(const SplineType& spline, typename SplineType::Scalar u, DenseIndex order, DerivativeType& der) + { + enum { Dimension = SplineTraits::Dimension }; + enum { Order = SplineTraits::OrderAtCompileTime }; + enum { DerivativeOrder = DerivativeType::ColsAtCompileTime }; + + typedef typename SplineTraits::Scalar Scalar; + + typedef typename SplineTraits::BasisVectorType BasisVectorType; + typedef typename SplineTraits::ControlPointVectorType ControlPointVectorType; + + typedef typename SplineTraits::BasisDerivativeType BasisDerivativeType; + typedef typename BasisDerivativeType::ConstRowXpr BasisDerivativeRowXpr; + + const DenseIndex p = spline.degree(); + const DenseIndex span = spline.span(u); + + const DenseIndex n = (std::min)(p, order); + + der.resize(Dimension,n+1); + + // Retrieve the basis function derivatives up to the desired order... + const BasisDerivativeType basis_func_ders = spline.template basisFunctionDerivatives(u, n+1); + + // ... and perform the linear combinations of the control points. + for (DenseIndex der_order=0; der_order ctrl_weights( basis_func_ders.row(der_order) ); + const Block ctrl_pts(spline.ctrls(),0,span-p,Dimension,p+1); + der.col(der_order) = (ctrl_weights * ctrl_pts).rowwise().sum(); + } + } + + /** + * \brief A functor for the computation of a spline point. + * \sa Piegl & Tiller, "The NURBS Book", A4.1 (p. 124) + **/ + template + typename SplineTraits< Spline<_Scalar, _Dim, _Degree> >::DerivativeType + Spline<_Scalar, _Dim, _Degree>::derivatives(Scalar u, DenseIndex order) const + { + typename SplineTraits< Spline >::DerivativeType res; + derivativesImpl(*this, u, order, res); + return res; + } + + template + template + typename SplineTraits< Spline<_Scalar, _Dim, _Degree>, DerivativeOrder >::DerivativeType + Spline<_Scalar, _Dim, _Degree>::derivatives(Scalar u, DenseIndex order) const + { + typename SplineTraits< Spline, DerivativeOrder >::DerivativeType res; + derivativesImpl(*this, u, order, res); + return res; + } + + /** + * \brief A functor for the computation of a spline's non-zero basis functions. + * \sa Piegl & Tiller, "The NURBS Book", A2.2 (p. 70) + **/ + template + typename SplineTraits< Spline<_Scalar, _Dim, _Degree> >::BasisVectorType + Spline<_Scalar, _Dim, _Degree>::basisFunctions(Scalar u) const + { + return Spline::BasisFunctions(u, degree(), knots()); + } + + /* --------------------------------------------------------------------------------------------- */ + + template + void basisFunctionDerivativesImpl(const SplineType& spline, typename SplineType::Scalar u, DenseIndex order, DerivativeType& N_) + { + enum { Order = SplineTraits::OrderAtCompileTime }; + + typedef typename SplineTraits::Scalar Scalar; + typedef typename SplineTraits::BasisVectorType BasisVectorType; + typedef typename SplineTraits::KnotVectorType KnotVectorType; + typedef typename SplineTraits::ControlPointVectorType ControlPointVectorType; + + const KnotVectorType& U = spline.knots(); + + const DenseIndex p = spline.degree(); + const DenseIndex span = spline.span(u); + + const DenseIndex n = (std::min)(p, order); + + N_.resize(n+1, p+1); + + BasisVectorType left = BasisVectorType::Zero(p+1); + BasisVectorType right = BasisVectorType::Zero(p+1); + + Matrix ndu(p+1,p+1); + + double saved, temp; + + ndu(0,0) = 1.0; + + DenseIndex j; + for (j=1; j<=p; ++j) + { + left[j] = u-U[span+1-j]; + right[j] = U[span+j]-u; + saved = 0.0; + + for (DenseIndex r=0; r(saved+right[r+1] * temp); + saved = left[j-r] * temp; + } + + ndu(j,j) = static_cast(saved); + } + + for (j = p; j>=0; --j) + N_(0,j) = ndu(j,p); + + // Compute the derivatives + DerivativeType a(n+1,p+1); + DenseIndex r=0; + for (; r<=p; ++r) + { + DenseIndex s1,s2; + s1 = 0; s2 = 1; // alternate rows in array a + a(0,0) = 1.0; + + // Compute the k-th derivative + for (DenseIndex k=1; k<=static_cast(n); ++k) + { + double d = 0.0; + DenseIndex rk,pk,j1,j2; + rk = r-k; pk = p-k; + + if (r>=k) + { + a(s2,0) = a(s1,0)/ndu(pk+1,rk); + d = a(s2,0)*ndu(rk,pk); + } + + if (rk>=-1) j1 = 1; + else j1 = -rk; + + if (r-1 <= pk) j2 = k-1; + else j2 = p-r; + + for (j=j1; j<=j2; ++j) + { + a(s2,j) = (a(s1,j)-a(s1,j-1))/ndu(pk+1,rk+j); + d += a(s2,j)*ndu(rk+j,pk); + } + + if (r<=pk) + { + a(s2,k) = -a(s1,k-1)/ndu(pk+1,r); + d += a(s2,k)*ndu(r,pk); + } + + N_(k,r) = static_cast(d); + j = s1; s1 = s2; s2 = j; // Switch rows + } + } + + /* Multiply through by the correct factors */ + /* (Eq. [2.9]) */ + r = p; + for (DenseIndex k=1; k<=static_cast(n); ++k) + { + for (DenseIndex j=p; j>=0; --j) N_(k,j) *= r; + r *= p-k; + } + } + + template + typename SplineTraits< Spline<_Scalar, _Dim, _Degree> >::BasisDerivativeType + Spline<_Scalar, _Dim, _Degree>::basisFunctionDerivatives(Scalar u, DenseIndex order) const + { + typename SplineTraits< Spline >::BasisDerivativeType der; + basisFunctionDerivativesImpl(*this, u, order, der); + return der; + } + + template + template + typename SplineTraits< Spline<_Scalar, _Dim, _Degree>, DerivativeOrder >::BasisDerivativeType + Spline<_Scalar, _Dim, _Degree>::basisFunctionDerivatives(Scalar u, DenseIndex order) const + { + typename SplineTraits< Spline, DerivativeOrder >::BasisDerivativeType der; + basisFunctionDerivativesImpl(*this, u, order, der); + return der; + } +} + +#endif // EIGEN_SPLINE_H diff --git a/unsupported/Eigen/src/Splines/SplineFitting.h b/unsupported/Eigen/src/Splines/SplineFitting.h new file mode 100644 index 000000000..6de820af4 --- /dev/null +++ b/unsupported/Eigen/src/Splines/SplineFitting.h @@ -0,0 +1,115 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 20010-2011 Hauke Heibel +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see . + +#ifndef EIGEN_SPLINE_FITTING_H +#define EIGEN_SPLINE_FITTING_H + +#include + +#include "SplineFwd.h" + +#include + +namespace Eigen +{ + template + void KnotAveraging(const KnotVectorType& parameters, DenseIndex degree, KnotVectorType& knots) + { + typedef typename KnotVectorType::Scalar Scalar; + + knots.resize(parameters.size()+degree+1); + + for (DenseIndex j=1; j + void ChordLengths(const PointArrayType& pts, KnotVectorType& chord_lengths) + { + typedef typename KnotVectorType::Scalar Scalar; + + const DenseIndex n = pts.cols(); + + // 1. compute the column-wise norms + chord_lengths.resize(pts.cols()); + chord_lengths[0] = 0; + chord_lengths.rightCols(n-1) = (pts.array().leftCols(n-1) - pts.array().rightCols(n-1)).matrix().colwise().norm(); + + // 2. compute the partial sums + std::partial_sum(chord_lengths.data(), chord_lengths.data()+n, chord_lengths.data()); + + // 3. normalize the data + chord_lengths /= chord_lengths(n-1); + chord_lengths(n-1) = Scalar(1); + } + + template + struct SplineFitting + { + template + static SplineType Interpolate(const PointArrayType& pts, DenseIndex degree); + }; + + template + template + SplineType SplineFitting::Interpolate(const PointArrayType& pts, DenseIndex degree) + { + typedef typename SplineType::KnotVectorType::Scalar Scalar; + typedef typename SplineType::KnotVectorType KnotVectorType; + typedef typename SplineType::BasisVectorType BasisVectorType; + typedef typename SplineType::ControlPointVectorType ControlPointVectorType; + + typedef Matrix MatrixType; + + KnotVectorType chord_lengths; // knot parameters + ChordLengths(pts, chord_lengths); + + KnotVectorType knots; + KnotAveraging(chord_lengths, degree, knots); + + DenseIndex n = pts.cols(); + MatrixType A = MatrixType::Zero(n,n); + for (DenseIndex i=1; i qr(A); + + // Here, we are creating a temporary due to an Eigen issue. + ControlPointVectorType ctrls = qr.solve(MatrixType(pts.transpose())).transpose(); + + return SplineType(knots, ctrls); + } +} + +#endif // EIGEN_SPLINE_FITTING_H diff --git a/unsupported/Eigen/src/Splines/SplineFwd.h b/unsupported/Eigen/src/Splines/SplineFwd.h new file mode 100644 index 000000000..da0e75205 --- /dev/null +++ b/unsupported/Eigen/src/Splines/SplineFwd.h @@ -0,0 +1,78 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 20010-2011 Hauke Heibel +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see . + +#ifndef EIGEN_SPLINES_FWD_H +#define EIGEN_SPLINES_FWD_H + +#include + +namespace Eigen +{ + template class Spline; + + template < typename SplineType, int _DerivativeOrder = Dynamic > struct SplineTraits {}; + +// hide specializations from doxygen +#ifndef DOXYGEN_SHOULD_SKIP_THIS + + template + struct SplineTraits< Spline<_Scalar, _Dim, _Degree>, Dynamic > + { + typedef _Scalar Scalar; /* The underlying scalar value. */ + enum { Dimension = _Dim }; /* The spline curve's dimension. */ + enum { Degree = _Degree }; /* The spline curve's degree. */ + + enum { OrderAtCompileTime = _Degree==Dynamic ? Dynamic : _Degree+1 }; + enum { NumOfDerivativesAtCompileTime = OrderAtCompileTime }; + + typedef Array BasisVectorType; + + typedef Array BasisDerivativeType; + typedef Array DerivativeType; + + typedef Array PointType; + typedef Array KnotVectorType; + typedef Array ControlPointVectorType; + }; + + template < typename _Scalar, int _Dim, int _Degree, int _DerivativeOrder > + struct SplineTraits< Spline<_Scalar, _Dim, _Degree>, _DerivativeOrder > : public SplineTraits< Spline<_Scalar, _Dim, _Degree> > + { + enum { OrderAtCompileTime = _Degree==Dynamic ? Dynamic : _Degree+1 }; + enum { NumOfDerivativesAtCompileTime = _DerivativeOrder==Dynamic ? Dynamic : _DerivativeOrder+1 }; + + typedef Array<_Scalar,Dynamic,Dynamic,RowMajor,NumOfDerivativesAtCompileTime,OrderAtCompileTime> BasisDerivativeType; + typedef Array<_Scalar,_Dim,Dynamic,ColMajor,_Dim,NumOfDerivativesAtCompileTime> DerivativeType; + }; + +#endif + + typedef Spline Spline2f; + typedef Spline Spline3f; + + typedef Spline Spline2d; + typedef Spline Spline3d; +} + +#endif // EIGEN_SPLINES_FWD_H diff --git a/unsupported/test/CMakeLists.txt b/unsupported/test/CMakeLists.txt index 83e947665..8580bb061 100644 --- a/unsupported/test/CMakeLists.txt +++ b/unsupported/test/CMakeLists.txt @@ -94,3 +94,4 @@ endif(GSL_FOUND) ei_add_test(polynomialsolver " " "${GSL_LIBRARIES}" ) ei_add_test(polynomialutils) ei_add_test(kronecker_product) +ei_add_test(splines) diff --git a/unsupported/test/splines.cpp b/unsupported/test/splines.cpp new file mode 100644 index 000000000..02e170c44 --- /dev/null +++ b/unsupported/test/splines.cpp @@ -0,0 +1,239 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 2010-2011 Hauke Heibel +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see . + +#include "main.h" + +#include + +// lets do some explicit instantiations and thus +// force the compilation of all spline functions... +template class Spline; +template class Spline; + +template class Spline; +template class Spline; +template class Spline; +template class Spline; + +template class Spline; +template class Spline; + +template class Spline; +template class Spline; +template class Spline; +template class Spline; + +Spline closed_spline2d() +{ + RowVectorXd knots(12); + knots << 0, + 0, + 0, + 0, + 0.867193179093898, + 1.660330955342408, + 2.605084834823134, + 3.484154586374428, + 4.252699478956276, + 4.252699478956276, + 4.252699478956276, + 4.252699478956276; + + MatrixXd ctrls(8,2); + ctrls << -0.370967741935484, 0.236842105263158, + -0.231401860693277, 0.442245185027632, + 0.344361228532831, 0.773369994120753, + 0.828990216203802, 0.106550882647595, + 0.407270163678382, -1.043452922172848, + -0.488467813584053, -0.390098582530090, + -0.494657189446427, 0.054804824897884, + -0.370967741935484, 0.236842105263158; + ctrls.transposeInPlace(); + + return Spline(knots, ctrls); +} + +/* create a reference spline */ +Spline spline3d() +{ + RowVectorXd knots(11); + knots << 0, + 0, + 0, + 0.118997681558377, + 0.162611735194631, + 0.498364051982143, + 0.655098003973841, + 0.679702676853675, + 1.000000000000000, + 1.000000000000000, + 1.000000000000000; + + MatrixXd ctrls(8,3); + ctrls << 0.959743958516081, 0.340385726666133, 0.585267750979777, + 0.223811939491137, 0.751267059305653, 0.255095115459269, + 0.505957051665142, 0.699076722656686, 0.890903252535799, + 0.959291425205444, 0.547215529963803, 0.138624442828679, + 0.149294005559057, 0.257508254123736, 0.840717255983663, + 0.254282178971531, 0.814284826068816, 0.243524968724989, + 0.929263623187228, 0.349983765984809, 0.196595250431208, + 0.251083857976031, 0.616044676146639, 0.473288848902729; + ctrls.transposeInPlace(); + + return Spline(knots, ctrls); +} + +/* compares evaluations against known results */ +void eval_spline3d() +{ + Spline3d spline = spline3d(); + + RowVectorXd u(10); + u << 0.351659507062997, + 0.830828627896291, + 0.585264091152724, + 0.549723608291140, + 0.917193663829810, + 0.285839018820374, + 0.757200229110721, + 0.753729094278495, + 0.380445846975357, + 0.567821640725221; + + MatrixXd pts(10,3); + pts << 0.707620811535916, 0.510258911240815, 0.417485437023409, + 0.603422256426978, 0.529498282727551, 0.270351549348981, + 0.228364197569334, 0.423745615677815, 0.637687289287490, + 0.275556796335168, 0.350856706427970, 0.684295784598905, + 0.514519311047655, 0.525077224890754, 0.351628308305896, + 0.724152914315666, 0.574461155457304, 0.469860285484058, + 0.529365063753288, 0.613328702656816, 0.237837040141739, + 0.522469395136878, 0.619099658652895, 0.237139665242069, + 0.677357023849552, 0.480655768435853, 0.422227610314397, + 0.247046593173758, 0.380604672404750, 0.670065791405019; + pts.transposeInPlace(); + + for (int i=0; i::Interpolate(points,3); + + KnotVectorType chord_lengths; // knot parameters + Eigen::ChordLengths(points, chord_lengths); + + for (Eigen::DenseIndex i=0; i