File size: 16,505 Bytes
a42735d |
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 |
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 20010-2011 Hauke Heibel <hauke.heibel@gmail.com>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
#ifndef EIGEN_SPLINE_FITTING_H
#define EIGEN_SPLINE_FITTING_H
#include <algorithm>
#include <functional>
#include <numeric>
#include <vector>
#include "SplineFwd.h"
#include "../../../../Eigen/LU"
#include "../../../../Eigen/QR"
namespace Eigen
{
/**
* \brief Computes knot averages.
* \ingroup Splines_Module
*
* The knots are computed as
* \f{align*}
* u_0 & = \hdots = u_p = 0 \\
* u_{m-p} & = \hdots = u_{m} = 1 \\
* u_{j+p} & = \frac{1}{p}\sum_{i=j}^{j+p-1}\bar{u}_i \quad\quad j=1,\hdots,n-p
* \f}
* where \f$p\f$ is the degree and \f$m+1\f$ the number knots
* of the desired interpolating spline.
*
* \param[in] parameters The input parameters. During interpolation one for each data point.
* \param[in] degree The spline degree which is used during the interpolation.
* \param[out] knots The output knot vector.
*
* \sa Les Piegl and Wayne Tiller, The NURBS book (2nd ed.), 1997, 9.2.1 Global Curve Interpolation to Point Data
**/
template <typename KnotVectorType>
void KnotAveraging(const KnotVectorType& parameters, DenseIndex degree, KnotVectorType& knots)
{
knots.resize(parameters.size()+degree+1);
for (DenseIndex j=1; j<parameters.size()-degree; ++j)
knots(j+degree) = parameters.segment(j,degree).mean();
knots.segment(0,degree+1) = KnotVectorType::Zero(degree+1);
knots.segment(knots.size()-degree-1,degree+1) = KnotVectorType::Ones(degree+1);
}
/**
* \brief Computes knot averages when derivative constraints are present.
* Note that this is a technical interpretation of the referenced article
* since the algorithm contained therein is incorrect as written.
* \ingroup Splines_Module
*
* \param[in] parameters The parameters at which the interpolation B-Spline
* will intersect the given interpolation points. The parameters
* are assumed to be a non-decreasing sequence.
* \param[in] degree The degree of the interpolating B-Spline. This must be
* greater than zero.
* \param[in] derivativeIndices The indices corresponding to parameters at
* which there are derivative constraints. The indices are assumed
* to be a non-decreasing sequence.
* \param[out] knots The calculated knot vector. These will be returned as a
* non-decreasing sequence
*
* \sa Les A. Piegl, Khairan Rajab, Volha Smarodzinana. 2008.
* Curve interpolation with directional constraints for engineering design.
* Engineering with Computers
**/
template <typename KnotVectorType, typename ParameterVectorType, typename IndexArray>
void KnotAveragingWithDerivatives(const ParameterVectorType& parameters,
const unsigned int degree,
const IndexArray& derivativeIndices,
KnotVectorType& knots)
{
typedef typename ParameterVectorType::Scalar Scalar;
DenseIndex numParameters = parameters.size();
DenseIndex numDerivatives = derivativeIndices.size();
if (numDerivatives < 1)
{
KnotAveraging(parameters, degree, knots);
return;
}
DenseIndex startIndex;
DenseIndex endIndex;
DenseIndex numInternalDerivatives = numDerivatives;
if (derivativeIndices[0] == 0)
{
startIndex = 0;
--numInternalDerivatives;
}
else
{
startIndex = 1;
}
if (derivativeIndices[numDerivatives - 1] == numParameters - 1)
{
endIndex = numParameters - degree;
--numInternalDerivatives;
}
else
{
endIndex = numParameters - degree - 1;
}
// There are (endIndex - startIndex + 1) knots obtained from the averaging
// and 2 for the first and last parameters.
DenseIndex numAverageKnots = endIndex - startIndex + 3;
KnotVectorType averageKnots(numAverageKnots);
averageKnots[0] = parameters[0];
int newKnotIndex = 0;
for (DenseIndex i = startIndex; i <= endIndex; ++i)
averageKnots[++newKnotIndex] = parameters.segment(i, degree).mean();
averageKnots[++newKnotIndex] = parameters[numParameters - 1];
newKnotIndex = -1;
ParameterVectorType temporaryParameters(numParameters + 1);
KnotVectorType derivativeKnots(numInternalDerivatives);
for (DenseIndex i = 0; i < numAverageKnots - 1; ++i)
{
temporaryParameters[0] = averageKnots[i];
ParameterVectorType parameterIndices(numParameters);
int temporaryParameterIndex = 1;
for (DenseIndex j = 0; j < numParameters; ++j)
{
Scalar parameter = parameters[j];
if (parameter >= averageKnots[i] && parameter < averageKnots[i + 1])
{
parameterIndices[temporaryParameterIndex] = j;
temporaryParameters[temporaryParameterIndex++] = parameter;
}
}
temporaryParameters[temporaryParameterIndex] = averageKnots[i + 1];
for (int j = 0; j <= temporaryParameterIndex - 2; ++j)
{
for (DenseIndex k = 0; k < derivativeIndices.size(); ++k)
{
if (parameterIndices[j + 1] == derivativeIndices[k]
&& parameterIndices[j + 1] != 0
&& parameterIndices[j + 1] != numParameters - 1)
{
derivativeKnots[++newKnotIndex] = temporaryParameters.segment(j, 3).mean();
break;
}
}
}
}
KnotVectorType temporaryKnots(averageKnots.size() + derivativeKnots.size());
std::merge(averageKnots.data(), averageKnots.data() + averageKnots.size(),
derivativeKnots.data(), derivativeKnots.data() + derivativeKnots.size(),
temporaryKnots.data());
// Number of knots (one for each point and derivative) plus spline order.
DenseIndex numKnots = numParameters + numDerivatives + degree + 1;
knots.resize(numKnots);
knots.head(degree).fill(temporaryKnots[0]);
knots.tail(degree).fill(temporaryKnots.template tail<1>()[0]);
knots.segment(degree, temporaryKnots.size()) = temporaryKnots;
}
/**
* \brief Computes chord length parameters which are required for spline interpolation.
* \ingroup Splines_Module
*
* \param[in] pts The data points to which a spline should be fit.
* \param[out] chord_lengths The resulting chord length vector.
*
* \sa Les Piegl and Wayne Tiller, The NURBS book (2nd ed.), 1997, 9.2.1 Global Curve Interpolation to Point Data
**/
template <typename PointArrayType, typename KnotVectorType>
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);
}
/**
* \brief Spline fitting methods.
* \ingroup Splines_Module
**/
template <typename SplineType>
struct SplineFitting
{
typedef typename SplineType::KnotVectorType KnotVectorType;
typedef typename SplineType::ParameterVectorType ParameterVectorType;
/**
* \brief Fits an interpolating Spline to the given data points.
*
* \param pts The points for which an interpolating spline will be computed.
* \param degree The degree of the interpolating spline.
*
* \returns A spline interpolating the initially provided points.
**/
template <typename PointArrayType>
static SplineType Interpolate(const PointArrayType& pts, DenseIndex degree);
/**
* \brief Fits an interpolating Spline to the given data points.
*
* \param pts The points for which an interpolating spline will be computed.
* \param degree The degree of the interpolating spline.
* \param knot_parameters The knot parameters for the interpolation.
*
* \returns A spline interpolating the initially provided points.
**/
template <typename PointArrayType>
static SplineType Interpolate(const PointArrayType& pts, DenseIndex degree, const KnotVectorType& knot_parameters);
/**
* \brief Fits an interpolating spline to the given data points and
* derivatives.
*
* \param points The points for which an interpolating spline will be computed.
* \param derivatives The desired derivatives of the interpolating spline at interpolation
* points.
* \param derivativeIndices An array indicating which point each derivative belongs to. This
* must be the same size as @a derivatives.
* \param degree The degree of the interpolating spline.
*
* \returns A spline interpolating @a points with @a derivatives at those points.
*
* \sa Les A. Piegl, Khairan Rajab, Volha Smarodzinana. 2008.
* Curve interpolation with directional constraints for engineering design.
* Engineering with Computers
**/
template <typename PointArrayType, typename IndexArray>
static SplineType InterpolateWithDerivatives(const PointArrayType& points,
const PointArrayType& derivatives,
const IndexArray& derivativeIndices,
const unsigned int degree);
/**
* \brief Fits an interpolating spline to the given data points and derivatives.
*
* \param points The points for which an interpolating spline will be computed.
* \param derivatives The desired derivatives of the interpolating spline at interpolation points.
* \param derivativeIndices An array indicating which point each derivative belongs to. This
* must be the same size as @a derivatives.
* \param degree The degree of the interpolating spline.
* \param parameters The parameters corresponding to the interpolation points.
*
* \returns A spline interpolating @a points with @a derivatives at those points.
*
* \sa Les A. Piegl, Khairan Rajab, Volha Smarodzinana. 2008.
* Curve interpolation with directional constraints for engineering design.
* Engineering with Computers
*/
template <typename PointArrayType, typename IndexArray>
static SplineType InterpolateWithDerivatives(const PointArrayType& points,
const PointArrayType& derivatives,
const IndexArray& derivativeIndices,
const unsigned int degree,
const ParameterVectorType& parameters);
};
template <typename SplineType>
template <typename PointArrayType>
SplineType SplineFitting<SplineType>::Interpolate(const PointArrayType& pts, DenseIndex degree, const KnotVectorType& knot_parameters)
{
typedef typename SplineType::KnotVectorType::Scalar Scalar;
typedef typename SplineType::ControlPointVectorType ControlPointVectorType;
typedef Matrix<Scalar,Dynamic,Dynamic> MatrixType;
KnotVectorType knots;
KnotAveraging(knot_parameters, degree, knots);
DenseIndex n = pts.cols();
MatrixType A = MatrixType::Zero(n,n);
for (DenseIndex i=1; i<n-1; ++i)
{
const DenseIndex span = SplineType::Span(knot_parameters[i], degree, knots);
// The segment call should somehow be told the spline order at compile time.
A.row(i).segment(span-degree, degree+1) = SplineType::BasisFunctions(knot_parameters[i], degree, knots);
}
A(0,0) = 1.0;
A(n-1,n-1) = 1.0;
HouseholderQR<MatrixType> 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);
}
template <typename SplineType>
template <typename PointArrayType>
SplineType SplineFitting<SplineType>::Interpolate(const PointArrayType& pts, DenseIndex degree)
{
KnotVectorType chord_lengths; // knot parameters
ChordLengths(pts, chord_lengths);
return Interpolate(pts, degree, chord_lengths);
}
template <typename SplineType>
template <typename PointArrayType, typename IndexArray>
SplineType
SplineFitting<SplineType>::InterpolateWithDerivatives(const PointArrayType& points,
const PointArrayType& derivatives,
const IndexArray& derivativeIndices,
const unsigned int degree,
const ParameterVectorType& parameters)
{
typedef typename SplineType::KnotVectorType::Scalar Scalar;
typedef typename SplineType::ControlPointVectorType ControlPointVectorType;
typedef Matrix<Scalar, Dynamic, Dynamic> MatrixType;
const DenseIndex n = points.cols() + derivatives.cols();
KnotVectorType knots;
KnotAveragingWithDerivatives(parameters, degree, derivativeIndices, knots);
// fill matrix
MatrixType A = MatrixType::Zero(n, n);
// Use these dimensions for quicker populating, then transpose for solving.
MatrixType b(points.rows(), n);
DenseIndex startRow;
DenseIndex derivativeStart;
// End derivatives.
if (derivativeIndices[0] == 0)
{
A.template block<1, 2>(1, 0) << -1, 1;
Scalar y = (knots(degree + 1) - knots(0)) / degree;
b.col(1) = y*derivatives.col(0);
startRow = 2;
derivativeStart = 1;
}
else
{
startRow = 1;
derivativeStart = 0;
}
if (derivativeIndices[derivatives.cols() - 1] == points.cols() - 1)
{
A.template block<1, 2>(n - 2, n - 2) << -1, 1;
Scalar y = (knots(knots.size() - 1) - knots(knots.size() - (degree + 2))) / degree;
b.col(b.cols() - 2) = y*derivatives.col(derivatives.cols() - 1);
}
DenseIndex row = startRow;
DenseIndex derivativeIndex = derivativeStart;
for (DenseIndex i = 1; i < parameters.size() - 1; ++i)
{
const DenseIndex span = SplineType::Span(parameters[i], degree, knots);
if (derivativeIndex < derivativeIndices.size() && derivativeIndices[derivativeIndex] == i)
{
A.block(row, span - degree, 2, degree + 1)
= SplineType::BasisFunctionDerivatives(parameters[i], 1, degree, knots);
b.col(row++) = points.col(i);
b.col(row++) = derivatives.col(derivativeIndex++);
}
else
{
A.row(row).segment(span - degree, degree + 1)
= SplineType::BasisFunctions(parameters[i], degree, knots);
b.col(row++) = points.col(i);
}
}
b.col(0) = points.col(0);
b.col(b.cols() - 1) = points.col(points.cols() - 1);
A(0,0) = 1;
A(n - 1, n - 1) = 1;
// Solve
FullPivLU<MatrixType> lu(A);
ControlPointVectorType controlPoints = lu.solve(MatrixType(b.transpose())).transpose();
SplineType spline(knots, controlPoints);
return spline;
}
template <typename SplineType>
template <typename PointArrayType, typename IndexArray>
SplineType
SplineFitting<SplineType>::InterpolateWithDerivatives(const PointArrayType& points,
const PointArrayType& derivatives,
const IndexArray& derivativeIndices,
const unsigned int degree)
{
ParameterVectorType parameters;
ChordLengths(points, parameters);
return InterpolateWithDerivatives(points, derivatives, derivativeIndices, degree, parameters);
}
}
#endif // EIGEN_SPLINE_FITTING_H
|