10 #ifndef EIGEN_SPLINE_FITTING_H 11 #define EIGEN_SPLINE_FITTING_H 44 template <
typename KnotVectorType>
47 knots.resize(parameters.size()+degree+1);
50 knots(
j+degree) = parameters.segment(
j,degree).mean();
52 knots.segment(0,degree+1) = KnotVectorType::Zero(degree+1);
53 knots.segment(knots.size()-degree-1,degree+1) = KnotVectorType::Ones(degree+1);
77 template <
typename KnotVectorType,
typename ParameterVectorType,
typename IndexArray>
80 const IndexArray& derivativeIndices,
81 KnotVectorType& knots)
86 DenseIndex numDerivatives = derivativeIndices.size();
88 if (numDerivatives < 1)
97 DenseIndex numInternalDerivatives = numDerivatives;
99 if (derivativeIndices[0] == 0)
102 --numInternalDerivatives;
108 if (derivativeIndices[numDerivatives - 1] == numParameters - 1)
110 endIndex = numParameters -
degree;
111 --numInternalDerivatives;
115 endIndex = numParameters - degree - 1;
120 DenseIndex numAverageKnots = endIndex - startIndex + 3;
121 KnotVectorType averageKnots(numAverageKnots);
122 averageKnots[0] = parameters[0];
124 int newKnotIndex = 0;
126 averageKnots[++newKnotIndex] = parameters.segment(
i, degree).mean();
127 averageKnots[++newKnotIndex] = parameters[numParameters - 1];
131 ParameterVectorType temporaryParameters(numParameters + 1);
132 KnotVectorType derivativeKnots(numInternalDerivatives);
135 temporaryParameters[0] = averageKnots[
i];
136 ParameterVectorType parameterIndices(numParameters);
137 int temporaryParameterIndex = 1;
140 Scalar parameter = parameters[
j];
141 if (parameter >= averageKnots[
i] && parameter < averageKnots[
i + 1])
143 parameterIndices[temporaryParameterIndex] =
j;
144 temporaryParameters[temporaryParameterIndex++] = parameter;
147 temporaryParameters[temporaryParameterIndex] = averageKnots[
i + 1];
149 for (
int j = 0;
j <= temporaryParameterIndex - 2; ++
j)
151 for (
DenseIndex k = 0; k < derivativeIndices.size(); ++k)
153 if (parameterIndices[
j + 1] == derivativeIndices[k]
154 && parameterIndices[
j + 1] != 0
155 && parameterIndices[
j + 1] != numParameters - 1)
157 derivativeKnots[++newKnotIndex] = temporaryParameters.segment(
j, 3).mean();
164 KnotVectorType temporaryKnots(averageKnots.size() + derivativeKnots.size());
166 std::merge(averageKnots.data(), averageKnots.data() + averageKnots.size(),
167 derivativeKnots.data(), derivativeKnots.data() + derivativeKnots.size(),
168 temporaryKnots.data());
171 DenseIndex numKnots = numParameters + numDerivatives + degree + 1;
172 knots.resize(numKnots);
174 knots.head(degree).fill(temporaryKnots[0]);
175 knots.tail(degree).fill(temporaryKnots.template tail<1>()[0]);
176 knots.segment(degree, temporaryKnots.size()) = temporaryKnots;
188 template <
typename Po
intArrayType,
typename KnotVectorType>
189 void ChordLengths(
const PointArrayType& pts, KnotVectorType& chord_lengths)
196 chord_lengths.resize(pts.cols());
197 chord_lengths[0] = 0;
198 chord_lengths.rightCols(n-1) = (pts.array().leftCols(n-1) - pts.array().rightCols(n-1)).
matrix().colwise().norm();
201 std::partial_sum(chord_lengths.data(), chord_lengths.data()+
n, chord_lengths.data());
204 chord_lengths /= chord_lengths(n-1);
205 chord_lengths(n-1) =
Scalar(1);
212 template <
typename SplineType>
226 template <
typename Po
intArrayType>
238 template <
typename Po
intArrayType>
258 template <
typename Po
intArrayType,
typename IndexArray>
260 const PointArrayType& derivatives,
261 const IndexArray& derivativeIndices,
262 const unsigned int degree);
280 template <
typename Po
intArrayType,
typename IndexArray>
282 const PointArrayType& derivatives,
283 const IndexArray& derivativeIndices,
284 const unsigned int degree,
288 template <
typename SplineType>
289 template <
typename Po
intArrayType>
293 typedef typename SplineType::ControlPointVectorType ControlPointVectorType;
301 MatrixType
A = MatrixType::Zero(n,n);
304 const DenseIndex span = SplineType::Span(knot_parameters[
i], degree, knots);
307 A.row(i).segment(span-degree, degree+1) = SplineType::BasisFunctions(knot_parameters[i], degree, knots);
315 ControlPointVectorType ctrls = qr.
solve(
MatrixType(pts.transpose())).transpose();
317 return SplineType(knots, ctrls);
320 template <
typename SplineType>
321 template <
typename Po
intArrayType>
329 template <
typename SplineType>
330 template <
typename Po
intArrayType,
typename IndexArray>
333 const PointArrayType& derivatives,
334 const IndexArray& derivativeIndices,
335 const unsigned int degree,
339 typedef typename SplineType::ControlPointVectorType ControlPointVectorType;
343 const DenseIndex n = points.cols() + derivatives.cols();
350 MatrixType
A = MatrixType::Zero(n, n);
353 MatrixType
b(points.rows(),
n);
359 if (derivativeIndices[0] == 0)
361 A.template block<1, 2>(1, 0) << -1, 1;
363 Scalar
y = (knots(degree + 1) - knots(0)) / degree;
364 b.col(1) = y*derivatives.col(0);
374 if (derivativeIndices[derivatives.cols() - 1] == points.cols() - 1)
376 A.template block<1, 2>(n - 2, n - 2) << -1, 1;
378 Scalar
y = (knots(knots.size() - 1) - knots(knots.size() - (degree + 2))) / degree;
379 b.col(
b.cols() - 2) = y*derivatives.col(derivatives.cols() - 1);
386 const DenseIndex span = SplineType::Span(parameters[
i], degree, knots);
388 if (derivativeIndices[derivativeIndex] == i)
390 A.block(row, span - degree, 2, degree + 1)
391 = SplineType::BasisFunctionDerivatives(parameters[i], 1, degree, knots);
393 b.col(row++) = points.col(i);
394 b.col(row++) = derivatives.col(derivativeIndex++);
398 A.row(row++).segment(span - degree, degree + 1)
399 = SplineType::BasisFunctions(parameters[i], degree, knots);
402 b.col(0) = points.col(0);
403 b.col(
b.cols() - 1) = points.col(points.cols() - 1);
409 ControlPointVectorType controlPoints = lu.
solve(
MatrixType(
b.transpose())).transpose();
411 SplineType spline(knots, controlPoints);
416 template <
typename SplineType>
417 template <
typename Po
intArrayType,
typename IndexArray>
420 const PointArrayType& derivatives,
421 const IndexArray& derivativeIndices,
422 const unsigned int degree)
430 #endif // EIGEN_SPLINE_FITTING_H
Namespace containing all symbols from the Eigen library.
void ChordLengths(const PointArrayType &pts, KnotVectorType &chord_lengths)
Computes chord length parameters which are required for spline interpolation.
HouseholderQR< MatrixXf > qr(A)
void KnotAveragingWithDerivatives(const ParameterVectorType ¶meters, const unsigned int degree, const IndexArray &derivativeIndices, KnotVectorType &knots)
Computes knot averages when derivative constraints are present. Note that this is a technical interpr...
static SplineType InterpolateWithDerivatives(const PointArrayType &points, const PointArrayType &derivatives, const IndexArray &derivativeIndices, const unsigned int degree)
Fits an interpolating spline to the given data points and derivatives.
void KnotAveraging(const KnotVectorType ¶meters, DenseIndex degree, KnotVectorType &knots)
Computes knot averages.The knots are computed as where is the degree and the number knots of the d...
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
static ConjugateGradientParameters parameters
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
static SplineType Interpolate(const PointArrayType &pts, DenseIndex degree)
Fits an interpolating Spline to the given data points.
SplineType::ParameterVectorType ParameterVectorType
LU decomposition of a matrix with complete pivoting, and related features.
SplineType::KnotVectorType KnotVectorType
const Solve< HouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const