UsTK : Ultrasound ToolKit  version 2.0.1 under development (2024-11-21)
usPolynomialCurve3D.cpp
1 /****************************************************************************
2  *
3  * This file is part of the ustk software.
4  * Copyright (C) 2016 - 2017 by Inria. All rights reserved.
5  *
6  * This software is free software; you can redistribute it and/or
7  * modify it under the terms of the GNU General Public License
8  * ("GPL") version 2 as published by the Free Software Foundation.
9  * See the file LICENSE.txt at the root directory of this source
10  * distribution for additional information about the GNU GPL.
11  *
12  * For using ustk with software that can not be combined with the GNU
13  * GPL, please contact Inria about acquiring a ViSP Professional
14  * Edition License.
15  *
16  * This software was developed at:
17  * Inria Rennes - Bretagne Atlantique
18  * Campus Universitaire de Beaulieu
19  * 35042 Rennes Cedex
20  * France
21  *
22  * If you have questions regarding the use of this file, please contact
23  * Inria at ustk@inria.fr
24  *
25  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
26  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
27  *
28  * Authors:
29  * Jason Chevrie
30  *
31  *****************************************************************************/
32 
33 #include <visp3/ustk_core/usPolynomialCurve3D.h>
34 
35 #include <algorithm>
36 
37 #include <visp3/core/vpException.h>
38 
40  : m_order(0), m_startParameter(0), m_endParameter(1), m_polynomialCoefficients(3, 1, 0)
41 {
42 }
43 
45  : m_order(curve.m_order), m_startParameter(curve.m_startParameter), m_endParameter(curve.m_endParameter),
46  m_polynomialCoefficients(curve.m_polynomialCoefficients)
47 {
48 }
49 
51 {
52  m_order = curve.m_order;
56 
57  return *this;
58 }
59 
61 
63  : m_order(order), m_startParameter(0), m_endParameter(1), m_polynomialCoefficients(3, m_order + 1, 0)
64 {
65 }
66 
67 void usPolynomialCurve3D::setOrder(unsigned int order)
68 {
69  if (order == m_order)
70  return;
71  else if (order > m_order) {
72  // Keep same polynomial
73  m_order = order;
74  m_polynomialCoefficients.resize(3, order + 1, false);
75  } else {
76  // Find polynomial to have the same properties at the extremities (position, direction, curvature, ...)
77  int nb_coef = order + 1;
78  int nb_constraints_begin = nb_coef / 2 + nb_coef % 2;
79  int nb_constraints_end = nb_coef / 2;
80 
81  vpColVector startParameter_power(m_order + 1, 0);
82  vpColVector endParameter_power(m_order + 1, 0);
83  startParameter_power[0] = 1;
84  endParameter_power[0] = 1;
85  for (unsigned int i = 1; i <= m_order; i++) {
86  startParameter_power[i] = m_startParameter * startParameter_power[i - 1];
87  endParameter_power[i] = m_endParameter * endParameter_power[i - 1];
88  }
89 
90  vpMatrix A(nb_coef, nb_coef, 0);
91  for (int i = 0; i < nb_constraints_begin; i++) {
92  for (int j = i; j < nb_coef; j++) {
93  A[i][j] = startParameter_power[j - i];
94  for (int n = 0; n < i; n++)
95  A[i][j] *= j - n;
96  }
97  }
98  for (int i = 0; i < nb_constraints_end; i++) {
99  for (int j = 0; j < nb_coef; j++) {
100  A[nb_constraints_begin + i][j] = endParameter_power[j - i];
101  for (int n = 0; n < i; n++)
102  A[i][j] *= j - n;
103  }
104  }
105 
106  vpMatrix Ainvt = A.inverseByLU().t();
107 
108  vpMatrix B(3, nb_coef);
109  for (int i = 0; i < nb_constraints_begin; i++) {
110  vpColVector col(3, 0);
111  for (unsigned int j = nb_constraints_begin; j <= m_order; j++) {
112  double coef = startParameter_power[j - i];
113  for (int n = 0; n < i; n++)
114  coef *= j - n;
115  col += coef * m_polynomialCoefficients.getCol(j);
116  }
117  B.insert(col, 0, i);
118  }
119  for (int i = 0; i < nb_constraints_end; i++) {
120  vpColVector col(3, 0);
121  for (unsigned int j = nb_constraints_begin; j <= m_order; j++) {
122  double coef = endParameter_power[j - i];
123  for (int n = 0; n < i; n++)
124  coef *= j - n;
125  col += coef * m_polynomialCoefficients.getCol(j);
126  }
127  B.insert(col, 0, nb_constraints_begin + i);
128  }
129 
130  m_polynomialCoefficients.resize(3, nb_coef);
131  m_polynomialCoefficients = B * Ainvt;
132  m_order = order;
133  }
134 }
135 
136 unsigned int usPolynomialCurve3D::getOrder() const { return m_order; }
137 
138 void usPolynomialCurve3D::setStartParameter(double startParameter) { m_startParameter = startParameter; }
139 
141 
142 void usPolynomialCurve3D::setEndParameter(double endParameter) { m_endParameter = endParameter; }
143 
145 
146 void usPolynomialCurve3D::setBoundaries(double startParameter, double endParameter)
147 {
148  if (startParameter <= endParameter) {
149  m_startParameter = startParameter;
150  m_endParameter = endParameter;
151  } else {
152  m_startParameter = endParameter;
153  m_endParameter = startParameter;
154  this->reverse();
155  }
156 }
157 
159 
161 
162 void usPolynomialCurve3D::setLength(double length, double precision)
163 {
164  if (precision <= 0)
165  throw vpException(vpException::badValue, "usPolynomialCurve3D::setLength: precision should be strictly positive");
166  if (length <= 0)
167  throw vpException(vpException::badValue, "usPolynomialCurve3D::setLength: length should be strictly positive");
168 
169  double lMetric = this->getLength();
170  double lParam = this->getParametricLength();
171  double lMetricPrev = 0;
172  double lParamPrev = 0;
173  double diffMetric = length - lMetric;
174 
175  while (fabs(diffMetric) > length * precision) {
176  double slope = (lParam - lParamPrev) / (lMetric - lMetricPrev);
177  m_endParameter += diffMetric * slope;
178  lMetricPrev = lMetric;
179  lParamPrev = lParam;
180  lMetric = this->getLength();
181  lParam = this->getParametricLength();
182  diffMetric = length - lMetric;
183  }
184 }
185 
186 double usPolynomialCurve3D::getLength(int nbCountSeg) const
187 {
188  vpColVector params(nbCountSeg + 1);
189  double step = this->getParametricLength() / nbCountSeg;
190  params[0] = m_startParameter;
191  for (int i = 1; i < nbCountSeg + 1; i++)
192  params[i] = params[i - 1] + step;
193 
194  vpMatrix points = this->getPoints(params);
195 
196  double length = 0.0;
197  for (int i = 0; i < nbCountSeg; i++)
198  length += (points.getCol(i) - points.getCol(i + 1)).frobeniusNorm();
199  return length;
200 }
201 
202 void usPolynomialCurve3D::setPolynomialCoefficients(const vpMatrix &polynomialCoefficients)
203 {
204  if (polynomialCoefficients.getCols() < 2)
205  throw vpException(vpException::dimensionError,
206  "usPolynomialCurve3D::setPolynomialCoefficients: empty coefficient matrix");
207  if (polynomialCoefficients.getRows() != 3)
208  throw vpException(vpException::dimensionError,
209  "usPolynomialCurve3D::setPolynomialCoefficients: space dimension should be 3");
210 
211  m_order = polynomialCoefficients.getCols() - 1;
212  m_polynomialCoefficients = polynomialCoefficients;
213 }
214 
216 
217 vpColVector usPolynomialCurve3D::getPoint(double parameter) const
218 {
219  vpColVector T(m_order + 1);
220  T[0] = 1.0;
221  for (unsigned int i = 1; i < m_order + 1; i++)
222  T[i] = T[i - 1] * parameter;
223  return m_polynomialCoefficients * T;
224 }
225 
226 vpMatrix usPolynomialCurve3D::getPoints(const vpColVector &parameters) const
227 {
228  vpMatrix T(m_order + 1, parameters.size());
229  for (unsigned int j = 0; j < parameters.size(); j++) {
230  T[0][j] = 1.0;
231  for (unsigned int i = 1; i < m_order + 1; i++)
232  T[i][j] = T[i - 1][j] * parameters[j];
233  }
234  return m_polynomialCoefficients * T;
235 }
236 
237 vpColVector usPolynomialCurve3D::getStartPoint() const { return this->getPoint(m_startParameter); }
238 
239 vpColVector usPolynomialCurve3D::getEndPoint() const { return this->getPoint(m_endParameter); }
240 
241 vpColVector usPolynomialCurve3D::getTangent(double parameter) const
242 {
243  vpColVector T(m_order + 1);
244  double tt = 1.0;
245  T[0] = 0.0;
246  for (unsigned int i = 1; i < m_order + 1; i++) {
247  T[i] = i * tt;
248  tt *= parameter;
249  }
250  return (m_polynomialCoefficients * T).normalize();
251 }
252 
254 
255 vpColVector usPolynomialCurve3D::getEndTangent() const { return this->getTangent(m_endParameter); }
256 
257 vpColVector usPolynomialCurve3D::getDerivative(double parameter, unsigned int order) const
258 {
259  vpColVector P(3, 0);
260 
261  if (order > 0 && order <= m_order) {
262  vpColVector factor_coef(m_order + 1, 0);
263  factor_coef[order] = 1;
264  for (unsigned int i = order + 1; i <= m_order; i++)
265  factor_coef[i] = parameter * factor_coef[i - 1];
266 
267  for (unsigned int i = order; i <= m_order; i++) {
268  for (unsigned int n = 0; n < order; n++)
269  factor_coef[i] *= i - n;
270  }
271 
272  P = m_polynomialCoefficients * factor_coef;
273  }
274 
275  return P;
276 }
277 
278 void usPolynomialCurve3D::defineFromPoints(const std::vector<vpColVector> &points, const std::vector<double> &param,
279  unsigned int order)
280 {
281  unsigned int nbPoints = points.size();
282  if (nbPoints < 2)
283  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPoints(const "
284  "std::vector<vpColVector> &points, const std::vector<double> "
285  "&param, unsigned int order): need at least two points to fit");
286  if (param.size() != nbPoints)
287  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPoints(const "
288  "std::vector<vpColVector> &points, const std::vector<double> "
289  "&param, unsigned int order): mismatching number of points and "
290  "parametric values");
291 
292  if (order < 1)
293  order = m_order;
294 
295  unsigned int nbCoef = order + 1;
296 
297  vpMatrix M(nbCoef, nbCoef, 0);
298  vpMatrix B(nbCoef, 3, 0);
299 
300  unsigned int nb_pow = 2 * nbCoef - 1;
301  vpMatrix t_pow(nbPoints, nb_pow);
302  for (unsigned int i = 0; i < nbPoints; i++) {
303  double ti = param.at(i);
304  t_pow[i][0] = 1;
305  for (unsigned int j = 1; j < nb_pow; j++) {
306  t_pow[i][j] = ti * t_pow[i][j - 1];
307  }
308  }
309 
310  for (unsigned int line = 0; line < nbCoef; line++) {
311  for (unsigned int i = 0; i < nbPoints; i++) {
312  for (unsigned int j = 0; j < nbCoef; j++) {
313  M[line][j] += t_pow[i][line + j];
314  }
315  for (unsigned int dim = 0; dim < 3; dim++) {
316  B[line][dim] += t_pow[i][line] * points.at(i)[dim];
317  }
318  }
319  }
320 
321  vpMatrix A;
322  try {
323  if (nbCoef > nbPoints)
324  A = M.pseudoInverse(std::numeric_limits<double>::epsilon()) * B;
325  else
326  A = M.inverseByCholeskyLapack() * B;
327  } catch (std::exception &e) {
328  throw vpException(vpException::fatalError, "usPolynomialCurve3D::defineFromPoints(const std::vector<vpColVector> "
329  "&points, const std::vector<double> &param, unsigned int order): %s",
330  e.what());
331  }
332 
333  this->setPolynomialCoefficients(A.t());
334  this->setBoundaries(param.front(), param.back());
335 }
336 
337 void usPolynomialCurve3D::defineFromPoints(const vpMatrix points, const vpColVector &param, unsigned int order)
338 {
339  unsigned int nbPoints = points.getCols();
340  if (nbPoints < 2)
341  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPoints(const "
342  "std::vector<vpColVector> &points, const std::vector<double> "
343  "&param, unsigned int order): need at least two points to fit");
344  if (param.size() != nbPoints)
345  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPoints(const "
346  "std::vector<vpColVector> &points, const std::vector<double> "
347  "&param, unsigned int order): mismatching number of points and "
348  "parametric values");
349 
350  if (order < 1)
351  order = m_order;
352 
353  unsigned int nbCoef = order + 1;
354 
355  vpMatrix M(nbCoef, nbCoef, 0);
356  vpMatrix B(nbCoef, 3, 0);
357 
358  unsigned int nb_pow = 2 * nbCoef - 1;
359  vpMatrix t_pow(nbPoints, nb_pow);
360  for (unsigned int i = 0; i < nbPoints; i++) {
361  double ti = param[i];
362  t_pow[i][0] = 1;
363  for (unsigned int j = 1; j < nb_pow; j++) {
364  t_pow[i][j] = ti * t_pow[i][j - 1];
365  }
366  }
367 
368  for (unsigned int line = 0; line < nbCoef; line++) {
369  for (unsigned int i = 0; i < nbPoints; i++) {
370  for (unsigned int j = 0; j < nbCoef; j++) {
371  M[line][j] += t_pow[i][line + j];
372  }
373  for (unsigned int dim = 0; dim < 3; dim++) {
374  B[line][dim] += t_pow[i][line] * points[dim][i];
375  }
376  }
377  }
378 
379  vpMatrix A;
380  try {
381  if (nbCoef > nbPoints)
382  A = M.pseudoInverse(std::numeric_limits<double>::epsilon()) * B;
383  else
384  A = M.inverseByCholeskyLapack() * B;
385  } catch (std::exception &e) {
386  throw vpException(vpException::fatalError, "usPolynomialCurve3D::defineFromPoints(const std::vector<vpColVector> "
387  "&points, const std::vector<double> &param, unsigned int order): %s",
388  e.what());
389  }
390 
391  this->setPolynomialCoefficients(A.t());
392  this->setBoundaries(param[0], param[param.size() - 1]);
393 }
394 
395 void usPolynomialCurve3D::defineFromPointsAuto(const std::vector<vpColVector> &points, unsigned int order)
396 {
397  unsigned int nbPoints = points.size();
398  if (nbPoints < 2)
399  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPointsAuto(const "
400  "std::vector<vpColVector> &points, unsigned int order): need at "
401  "least two points to fit");
402 
403  if (order < 1)
404  order = m_order;
405 
406  int max_i = 0;
407  int max_j = 0;
408  double max = 0;
409  for (unsigned int i = 0; i < nbPoints; i++) {
410  for (unsigned int j = 0; j < i; j++) {
411  double v = (points.at(i) - points.at(j)).sumSquare();
412  if (v > max) {
413  max = v;
414  max_i = i;
415  max_j = j;
416  }
417  }
418  }
419  vpColVector dir((points.at(max_i) - points.at(max_j)).normalize());
420 
421  return this->defineFromPointsAuto(points, dir, order);
422 }
423 
424 void usPolynomialCurve3D::defineFromPointsAuto(const vpMatrix &points, unsigned int order)
425 {
426  unsigned int nbPoints = points.getCols();
427  if (nbPoints < 2)
428  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPointsAuto(const vpMatrix &points, "
429  "unsigned int order): need at least two points to fit");
430 
431  if (order < 1)
432  order = m_order;
433 
434  int max_i = 0;
435  int max_j = 0;
436  double max = 0;
437  for (unsigned int i = 0; i < nbPoints; i++) {
438  for (unsigned int j = 0; j < i; j++) {
439  double v = (points.getCol(i) - points.getCol(j)).sumSquare();
440  if (v > max) {
441  max = v;
442  max_i = i;
443  max_j = j;
444  }
445  }
446  }
447  vpColVector dir((points.getCol(max_i) - points.getCol(max_j)).normalize());
448 
449  return this->defineFromPointsAuto(points, dir, order);
450 }
451 
452 void usPolynomialCurve3D::defineFromPointsAuto(const std::vector<vpColVector> &points, const vpColVector &direction,
453  unsigned int order)
454 {
455  unsigned int nbPoints = points.size();
456  if (nbPoints < 2)
457  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPointsAuto(const "
458  "std::vector<vpColVector> &points, const vpColVector &direction, "
459  "unsigned int order): need at least two points to fit");
460 
461  if (order < 1)
462  order = m_order;
463 
464  std::vector<double> t(nbPoints, 0);
465 
466  double min = std::numeric_limits<double>::max();
467  for (unsigned int i = 0; i < nbPoints; i++) {
468  t.at(i) = vpColVector::dotProd(points.at(i), direction);
469  if (t.at(i) < min)
470  min = t.at(i);
471  }
472  for (unsigned int i = 0; i < nbPoints; i++)
473  t.at(i) -= min;
474 
475  std::vector<int> index(nbPoints);
476  for (unsigned int i = 0; i < nbPoints; i++)
477  index.at(i) = i;
478 
479  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
480 
481  std::vector<vpColVector> newPoints(nbPoints);
482  std::vector<double> interLength(nbPoints - 1);
483  newPoints.front() = points.at(index[0]);
484  double oldLength = this->getParametricLength();
485  double newLength = 0;
486  for (unsigned int i = 1; i < nbPoints; i++) {
487  newPoints.at(i) = points.at(index[i]);
488  interLength.at(i - 1) = sqrt(vpMath::sqr(newPoints.at(i)[0] - newPoints.at(i - 1)[0]) +
489  vpMath::sqr(newPoints.at(i)[1] - newPoints.at(i - 1)[1]) +
490  +vpMath::sqr(newPoints.at(i)[2] - newPoints.at(i - 1)[2]));
491  newLength += interLength.at(i - 1);
492  }
493  double scaling = oldLength / newLength;
494  std::vector<double> param(nbPoints);
495  param.front() = m_startParameter;
496  for (unsigned int i = 1; i < nbPoints; i++)
497  param.at(i) = param.at(i - 1) + scaling * interLength.at(i - 1);
498 
499  this->defineFromPoints(newPoints, param, order);
500 }
501 
502 void usPolynomialCurve3D::defineFromPointsAuto(const vpMatrix &points, const vpColVector &direction, unsigned int order)
503 {
504  unsigned int nbPoints = points.getCols();
505  if (nbPoints < 2)
506  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromPointsAuto(const vpMatrix &points, "
507  "const vpColVector &direction, unsigned int order): need at least "
508  "two points to fit");
509 
510  if (order < 1)
511  order = m_order;
512 
513  std::vector<double> t(nbPoints, 0);
514 
515  double min = std::numeric_limits<double>::max();
516  for (unsigned int i = 0; i < nbPoints; i++) {
517  t.at(i) = vpColVector::dotProd(points.getCol(i), direction);
518  if (t.at(i) < min)
519  min = t.at(i);
520  }
521  for (unsigned int i = 0; i < nbPoints; i++)
522  t.at(i) -= min;
523 
524  std::vector<int> index(nbPoints);
525  for (unsigned int i = 0; i < nbPoints; i++)
526  index.at(i) = i;
527 
528  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
529 
530  vpMatrix newPoints(3, nbPoints);
531  vpColVector interLength(nbPoints - 1);
532  for (unsigned int j = 0; j < 3; j++)
533  newPoints[j][0] = points[j][index[0]];
534  double oldLength = this->getParametricLength();
535  double newLength = 0;
536  for (unsigned int i = 1; i < nbPoints; i++) {
537  for (unsigned int j = 0; j < 3; j++)
538  newPoints[j][i] = points[j][index[i]];
539  interLength[i - 1] =
540  sqrt(vpMath::sqr(newPoints[0][i] - newPoints[0][i - 1]) + vpMath::sqr(newPoints[1][i] - newPoints[1][i - 1]) +
541  vpMath::sqr(newPoints[2][i] - newPoints[2][i - 1]));
542  newLength += interLength[i - 1];
543  }
544  double scaling = oldLength / newLength;
545  vpColVector param(nbPoints);
546  param[0] = m_startParameter;
547  for (unsigned int i = 1; i < nbPoints; i++)
548  param[i] = param[i - 1] + scaling * interLength[i - 1];
549 
550  this->defineFromPoints(newPoints, param, order);
551 }
552 
553 void usPolynomialCurve3D::defineFromWeightedPoints(const std::vector<vpColVector> &points,
554  const std::vector<double> &param, const std::vector<double> &weights,
555  unsigned int order)
556 {
557  unsigned int nbPoints = points.size();
558  if (nbPoints < 2)
559  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const "
560  "std::vector<vpColVector> &points, const std::vector<double> "
561  "&param, const std::vector<double> &weights, unsigned int order): "
562  "need at least two points to fit");
563  if (param.size() != nbPoints)
564  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const "
565  "std::vector<vpColVector> &points, const std::vector<double> "
566  "&param, const std::vector<double> &weights, unsigned int order): "
567  "mismatching number of points and parametric values or weights");
568  if (weights.size() != nbPoints)
569  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const "
570  "std::vector<vpColVector> &points, const std::vector<double> "
571  "&param, const std::vector<double> &weights, unsigned int order): "
572  "mismatching number of points and weights");
573 
574  if (order < 1)
575  order = m_order;
576 
577  unsigned int nbCoef = order + 1;
578 
579  vpMatrix M(nbCoef, nbCoef, 0);
580  vpMatrix B(nbCoef, 3, 0);
581 
582  unsigned int nb_pow = 2 * nbCoef - 1;
583  vpMatrix t_pow(nbPoints, nb_pow);
584  for (unsigned int i = 0; i < nbPoints; i++) {
585  double ti = param.at(i);
586  t_pow[i][0] = weights.at(i);
587  for (unsigned int j = 1; j < nb_pow; j++) {
588  t_pow[i][j] = ti * t_pow[i][j - 1];
589  }
590  }
591 
592  for (unsigned int line = 0; line < nbCoef; line++) {
593  for (unsigned int i = 0; i < nbPoints; i++) {
594  for (unsigned int j = 0; j < nbCoef; j++) {
595  M[line][j] += t_pow[i][line + j];
596  }
597  for (unsigned int dim = 0; dim < 3; dim++) {
598  B[line][dim] += t_pow[i][line] * points.at(i)[dim];
599  }
600  }
601  }
602 
603  vpMatrix A;
604  try {
605  A = M.pseudoInverse(std::numeric_limits<double>::epsilon()) * B;
606  } catch (std::exception &e) {
607  throw vpException(vpException::fatalError, "usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix &points, "
608  "const vpColVector &param, const vpColVector &weights, unsigned int "
609  "order): %s",
610  e.what());
611  }
612 
613  this->setPolynomialCoefficients(A.t());
614  this->setBoundaries(param.front(), param.back());
615 }
616 
617 void usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix &points, const vpColVector &param,
618  const vpColVector &weights, unsigned int order)
619 {
620  unsigned int nbPoints = points.getCols();
621  if (nbPoints < 2)
622  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix "
623  "&points, const vpColVector &param, const vpColVector &weights, "
624  "unsigned int order): need at least two points to fit");
625  if (param.size() != nbPoints)
626  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix "
627  "&points, const vpColVector &param, const vpColVector &weights, "
628  "unsigned int order): mismatching number of points and parametric "
629  "values or weights");
630  if (weights.size() != nbPoints)
631  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix "
632  "&points, const vpColVector &param, const vpColVector &weights, "
633  "unsigned int order): mismatching number of points and weights");
634 
635  if (order < 1)
636  order = m_order;
637 
638  unsigned int nbCoef = order + 1;
639 
640  vpMatrix M(nbCoef, nbCoef, 0);
641  vpMatrix B(nbCoef, 3, 0);
642 
643  unsigned int nb_pow = 2 * nbCoef - 1;
644  vpMatrix t_pow(nbPoints, nb_pow);
645  for (unsigned int i = 0; i < nbPoints; i++) {
646  double ti = param[i];
647  t_pow[i][0] = weights[i];
648  for (unsigned int j = 1; j < nb_pow; j++) {
649  t_pow[i][j] = ti * t_pow[i][j - 1];
650  }
651  }
652 
653  for (unsigned int line = 0; line < nbCoef; line++) {
654  for (unsigned int i = 0; i < nbPoints; i++) {
655  for (unsigned int j = 0; j < nbCoef; j++) {
656  M[line][j] += t_pow[i][line + j];
657  }
658  for (unsigned int dim = 0; dim < 3; dim++) {
659  B[line][dim] += t_pow[i][line] * points[dim][i];
660  }
661  }
662  }
663 
664  vpMatrix A;
665  try {
666  A = M.pseudoInverse(std::numeric_limits<double>::epsilon()) * B;
667  } catch (std::exception &e) {
668  throw vpException(vpException::fatalError, "usPolynomialCurve3D::defineFromWeightedPoints(const vpMatrix &points, "
669  "const vpColVector &param, const vpColVector &weights, unsigned int "
670  "order): %s",
671  e.what());
672  }
673 
674  this->setPolynomialCoefficients(A.t());
675  this->setBoundaries(param[0], param[nbPoints - 1]);
676 }
677 
678 void usPolynomialCurve3D::defineFromWeightedPointsAuto(const std::vector<vpColVector> &points,
679  const std::vector<double> &weights, unsigned int order)
680 {
681  unsigned int nbPoints = points.size();
682  if (nbPoints < 2)
683  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const "
684  "std::vector<vpColVector> &points, const std::vector<double> "
685  "&weights, unsigned int order): need at least two points to fit");
686  if (weights.size() != nbPoints)
687  throw vpException(vpException::dimensionError,
688  "usPolynomialCurve3D::defineFromWeightedPointsAuto(const std::vector<vpColVector> &points, const "
689  "std::vector<double> &weights, unsigned int order): mismatching number of points and weights");
690 
691  if (order < 1)
692  order = m_order;
693 
694  std::vector<double> t(nbPoints, 0);
695 
696  vpColVector Pmean(3, 0);
697  for (unsigned int i = 0; i < nbPoints; i++)
698  Pmean += points.at(i);
699  Pmean /= nbPoints;
700 
701  vpMatrix M(nbPoints, 3);
702  for (unsigned int i = 0; i < nbPoints; i++)
703  for (int j = 0; j < 3; j++)
704  M[i][j] = points.at(i)[j] - Pmean[j];
705 
706  // Reduction to one principal component using singular value decomposition
707  vpColVector w;
708  vpMatrix V;
709  M.svd(w, V);
710 
711  vpColVector dir(V.getCol(0));
712  double min = std::numeric_limits<double>::max();
713  for (unsigned int i = 0; i < nbPoints; i++) {
714  t.at(i) = vpColVector::dotProd(points.at(i), dir);
715  if (t.at(i) < min)
716  min = t.at(i);
717  }
718  for (unsigned int i = 0; i < nbPoints; i++)
719  t.at(i) -= min;
720 
721  std::vector<int> index(nbPoints);
722  for (unsigned int i = 0; i < nbPoints; i++)
723  index.at(i) = i;
724 
725  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
726 
727  std::vector<vpColVector> newPoints(nbPoints);
728  std::vector<double> newWeights(nbPoints);
729  std::vector<double> interLength(nbPoints - 1);
730  newPoints.front() = points.at(index[0]);
731  newWeights.front() = weights.at(index[0]);
732  double oldLength = this->getParametricLength();
733  double newLength = 0;
734  for (unsigned int i = 1; i < nbPoints; i++) {
735  newPoints.at(i) = points.at(index[i]);
736  newWeights.at(i) = weights.at(index[i]);
737  interLength.at(i - 1) = sqrt(vpMath::sqr(newPoints.at(i)[0] - newPoints.at(i - 1)[0]) +
738  vpMath::sqr(newPoints.at(i)[1] - newPoints.at(i - 1)[1]) +
739  vpMath::sqr(newPoints.at(i)[2] - newPoints.at(i - 1)[2]));
740  newLength += interLength.at(i - 1);
741  }
742  double scaling = oldLength / newLength;
743  std::vector<double> param(nbPoints);
744  param.front() = m_startParameter;
745  for (unsigned int i = 1; i < nbPoints; i++)
746  param.at(i) = param.at(i - 1) + scaling * interLength.at(i - 1);
747 
748  this->defineFromWeightedPoints(newPoints, param, newWeights, order);
749 }
750 
751 void usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix &points, const vpColVector &weights,
752  unsigned int order)
753 {
754  unsigned int nbPoints = points.getCols();
755  if (nbPoints < 2)
756  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix "
757  "&points, const vpColVector &weights, unsigned int order): need at "
758  "least two points to fit");
759  if (weights.size() != nbPoints)
760  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix "
761  "&points, const vpColVector &weights, unsigned int order): "
762  "mismatching number of points and weights");
763 
764  if (order < 1)
765  order = m_order;
766 
767  std::vector<double> t(nbPoints, 0);
768 
769  vpColVector Pmean(3, 0);
770  for (unsigned int i = 0; i < nbPoints; i++)
771  Pmean += points.getCol(i);
772  Pmean /= nbPoints;
773 
774  vpMatrix M(nbPoints, 3);
775  for (unsigned int i = 0; i < nbPoints; i++)
776  for (int j = 0; j < 3; j++)
777  M[i][j] = points[j][i] - Pmean[j];
778 
779  // Reduction to one principal component using singular value decomposition
780  vpColVector w;
781  vpMatrix V;
782  M.svd(w, V);
783 
784  vpColVector dir(V.getCol(0));
785  double min = std::numeric_limits<double>::max();
786  for (unsigned int i = 0; i < nbPoints; i++) {
787  t.at(i) = vpColVector::dotProd(points.getCol(i), dir);
788  if (t.at(i) < min)
789  min = t.at(i);
790  }
791  for (unsigned int i = 0; i < nbPoints; i++)
792  t.at(i) -= min;
793 
794  std::vector<int> index(nbPoints);
795  for (unsigned int i = 0; i < nbPoints; i++)
796  index.at(i) = i;
797 
798  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
799 
800  vpMatrix newPoints(3, nbPoints);
801  vpColVector newWeights(nbPoints);
802  vpColVector interLength(nbPoints - 1);
803  for (unsigned int j = 0; j < 3; j++)
804  newPoints[j][0] = points[j][index[0]];
805  newWeights[0] = weights[index[0]];
806  double oldLength = this->getParametricLength();
807  double newLength = 0;
808  for (unsigned int i = 1; i < nbPoints; i++) {
809  for (unsigned int j = 0; j < 3; j++)
810  newPoints[j][i] = points[j][index[i]];
811  newWeights[i] = weights[index[i]];
812  interLength[i - 1] =
813  sqrt(vpMath::sqr(newPoints[0][i] - newPoints[0][i - 1]) + vpMath::sqr(newPoints[1][i] - newPoints[1][i - 1]) +
814  vpMath::sqr(newPoints[2][i] - newPoints[2][i - 1]));
815  newLength += interLength[i - 1];
816  }
817  double scaling = oldLength / newLength;
818  vpColVector param(nbPoints);
819  param[0] = m_startParameter;
820  for (unsigned int i = 1; i < nbPoints; i++)
821  param[i] = param[i - 1] + scaling * interLength[i - 1];
822 
823  this->defineFromWeightedPoints(newPoints, param, newWeights, order);
824 }
825 
826 void usPolynomialCurve3D::defineFromWeightedPointsAuto(const std::vector<vpColVector> &points,
827  const std::vector<double> &weights, const vpColVector &direction,
828  unsigned int order)
829 {
830  unsigned int nbPoints = points.size();
831  if (nbPoints < 2)
832  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const "
833  "std::vector<vpColVector> &points, const std::vector<double> "
834  "&weights, const vpColVector &direction, unsigned int order): need "
835  "at least two points to fit");
836  if (weights.size() != nbPoints)
837  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const "
838  "std::vector<vpColVector> &points, const std::vector<double> "
839  "&weights, const vpColVector &direction, unsigned int order): "
840  "mismatching number of points and weights");
841 
842  if (order < 1)
843  order = m_order;
844 
845  std::vector<double> t(nbPoints, 0);
846 
847  double min = std::numeric_limits<double>::max();
848  for (unsigned int i = 0; i < nbPoints; i++) {
849  t.at(i) = vpColVector::dotProd(points.at(i), direction);
850  if (t.at(i) < min)
851  min = t.at(i);
852  }
853  for (unsigned int i = 0; i < nbPoints; i++)
854  t.at(i) -= min;
855 
856  std::vector<int> index(nbPoints);
857  for (unsigned int i = 0; i < nbPoints; i++)
858  index.at(i) = i;
859 
860  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
861 
862  std::vector<vpColVector> newPoints(nbPoints);
863  std::vector<double> newWeights(nbPoints);
864  std::vector<double> interLength(nbPoints - 1);
865  newPoints.front() = points.at(index[0]);
866  newWeights.front() = weights.at(index[0]);
867  double oldLength = this->getParametricLength();
868  double newLength = 0;
869  for (unsigned int i = 1; i < nbPoints; i++) {
870  newPoints.at(i) = points.at(index[i]);
871  newWeights.at(i) = weights.at(index[i]);
872  interLength.at(i - 1) = sqrt(vpMath::sqr(newPoints.at(i)[0] - newPoints.at(i - 1)[0]) +
873  vpMath::sqr(newPoints.at(i)[1] - newPoints.at(i - 1)[1]) +
874  vpMath::sqr(newPoints.at(i)[2] - newPoints.at(i - 1)[2]));
875  newLength += interLength.at(i - 1);
876  }
877  double scaling = oldLength / newLength;
878  std::vector<double> param(nbPoints);
879  param.front() = m_startParameter;
880  for (unsigned int i = 1; i < nbPoints; i++)
881  param.at(i) = param.at(i - 1) + scaling * interLength.at(i - 1);
882 
883  this->defineFromWeightedPoints(newPoints, param, newWeights, order);
884 }
885 
886 void usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix &points, const vpColVector &weights,
887  const vpColVector &direction, unsigned int order)
888 {
889  unsigned int nbPoints = points.getCols();
890  if (nbPoints < 2)
891  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix "
892  "&points, const vpColVector &weights, const vpColVector &direction, "
893  "unsigned int order): need at least two points to fit");
894  if (weights.size() != nbPoints)
895  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::defineFromWeightedPointsAuto(const vpMatrix "
896  "&points, const vpColVector &weights, const vpColVector &direction, "
897  "unsigned int order): mismatching number of points and weights");
898 
899  if (order < 1)
900  order = m_order;
901 
902  std::vector<double> t(nbPoints, 0);
903 
904  double min = std::numeric_limits<double>::max();
905  for (unsigned int i = 0; i < nbPoints; i++) {
906  t.at(i) = vpColVector::dotProd(points.getCol(i), direction);
907  if (t.at(i) < min)
908  min = t.at(i);
909  }
910  for (unsigned int i = 0; i < nbPoints; i++)
911  t.at(i) -= min;
912 
913  std::vector<int> index(nbPoints);
914  for (unsigned int i = 0; i < nbPoints; i++)
915  index.at(i) = i;
916 
917  std::sort(index.begin(), index.end(), [&t](int i, int j) { return t[i] < t[j]; });
918 
919  vpMatrix newPoints(3, nbPoints);
920  vpColVector newWeights(nbPoints);
921  vpColVector interLength(nbPoints - 1);
922  for (unsigned int j = 0; j < 3; j++)
923  newPoints[j][0] = points[j][index[0]];
924  newWeights[0] = weights[index[0]];
925  double oldLength = this->getParametricLength();
926  double newLength = 0;
927  for (unsigned int i = 1; i < nbPoints; i++) {
928  for (unsigned int j = 0; j < 3; j++)
929  newPoints[j][i] = points[j][index[i]];
930  newWeights[i] = weights[index[i]];
931  interLength[i - 1] =
932  sqrt(vpMath::sqr(newPoints[0][i] - newPoints[0][i - 1]) + vpMath::sqr(newPoints[1][i] - newPoints[1][i - 1]) +
933  vpMath::sqr(newPoints[2][i] - newPoints[2][i - 1]));
934  newLength += interLength[i - 1];
935  }
936  double scaling = oldLength / newLength;
937  vpColVector param(nbPoints);
938  param[0] = m_startParameter;
939  for (unsigned int i = 1; i < nbPoints; i++)
940  param[i] = param[i - 1] + scaling * interLength[i - 1];
941 
942  this->defineFromWeightedPoints(newPoints, param, newWeights, order);
943 }
944 
945 double usPolynomialCurve3D::getCurvature(double param) const
946 {
947  vpColVector dX_dl = this->getDerivative(param, 1);
948  vpColVector dX2_dl2 = this->getDerivative(param, 2);
949 
950  double norm = dX_dl.frobeniusNorm();
951  double curvature = vpColVector::crossProd(dX_dl, dX2_dl2).frobeniusNorm() / pow(norm, 3);
952 
953  return curvature;
954 }
955 
956 double usPolynomialCurve3D::getCurvatureFromShape(double start, double end, vpColVector &center3D,
957  vpColVector &direction3D) const
958 {
959  if (start < 0)
960  start = 0;
961  double length = this->getParametricLength();
962  if (start > length)
963  start = length;
964  if (end < 0)
965  end = 0;
966  if (end > length)
967  end = length;
968  if (start > end) {
969  double tmp = start;
970  start = end;
971  end = tmp;
972  }
973 
974  int nbPoints = 10;
975 
976  // Create data matrix with centered vectors
977  double step = (end - start) / (nbPoints - 1);
978  vpColVector params(nbPoints);
979  params[0] = start;
980  for (int i = 1; i < nbPoints; i++)
981  params[i] = params[i - 1] + step;
982 
983  vpMatrix M(this->getPoints(params));
984  vpRowVector mean(3, 0);
985 
986  for (int i = 0; i < nbPoints; i++)
987  mean += M.getRow(i);
988  mean /= nbPoints;
989 
990  for (int i = 0; i < nbPoints; i++) {
991  M[i][0] -= mean[0];
992  M[i][1] -= mean[1];
993  M[i][2] -= mean[2];
994  }
995 
996  // Reduction to two principal components using singular value decomposition
997  vpMatrix U(M);
998  vpColVector w;
999  vpMatrix V;
1000  U.svd(w, V);
1001 
1002  vpMatrix S;
1003  S.diag(w);
1004 
1005  U.resize(nbPoints, 2, false);
1006 
1007  S.resize(2, 2, false);
1008 
1009  vpMatrix P = U * S;
1010 
1011  // 2D nonlinear least square fitting (Coope93)
1012  vpColVector d(nbPoints);
1013  for (int i = 0; i < nbPoints; i++) {
1014  d[i] = pow(P.t().getCol(i).frobeniusNorm(), 2);
1015  }
1016 
1017  vpColVector x(nbPoints, 1);
1018  vpMatrix B(nbPoints, 3);
1019  B.insert(P, 0, 0);
1020  B.insert(x, 0, 2);
1021 
1022  vpColVector y = B.pseudoInverse(0) * d;
1023 
1024  vpColVector center(2);
1025  center[0] = y[0] / 2;
1026  center[1] = y[1] / 2;
1027 
1028  double r = sqrt(y[2] + pow(center.frobeniusNorm(), 2));
1029 
1030  // Check validity
1031 
1032  /*vpColVector cp = P.getRow(0).t() - center;
1033  for(int i=1 ; i<nbPoints ; i++)
1034  {
1035  double dot = vpColVector::dotProd(cp, P.getRow(i).t() - center);
1036  if(dot<0) return 0;
1037  }*/
1038 
1039  if (direction3D.size() == 3) {
1040  direction3D = V.getCol(2);
1041  } else if (direction3D.size() == 4) {
1042  direction3D.insert(0, V.getCol(2));
1043  direction3D[3] = 0;
1044  }
1045 
1046  if (center3D.size() == 3) {
1047  V.resize(3, 2, false);
1048  center3D = V * center;
1049  } else if (center3D.size() == 4) {
1050  V.resize(3, 2, false);
1051  center3D.insert(0, V * center + mean.t());
1052  center3D[3] = 1;
1053  }
1054 
1055  return 1.0 / r;
1056  /*
1057  usPolynomialCurve3D seg(m_spline.back());
1058 
1059  vpColVector dX_dl = seg.getDerivative(seg.getParametricLength(),2);
1060  if(dX_dl.frobeniusNorm()<std::numeric_limits<double>::epsilon()) return 0;
1061 
1062  vpColVector d2X_dl2 = seg.getDerivative(seg.getParametricLength(),2);
1063 
1064  double k = vpColVector::crossProd(dX_dl, d2X_dl2).frobeniusNorm() / pow(dX_dl.frobeniusNorm(),3);
1065 
1066  center3D = seg.getEndPoint() + 1/k * ( d2X_dl2 - 1/pow(dX_dl.frobeniusNorm(),2) * vpColVector::dotProd(dX_dl,
1067  d2X_dl2)*dX_dl);
1068  direction3D = vpColVector::crossProd(dX_dl, d2X_dl2).normalize();
1069  return k;*/
1070 }
1071 
1072 double usPolynomialCurve3D::getMeanAxisDeviation(int nbCountSeg) const
1073 {
1074  if (nbCountSeg < 2)
1075  throw vpException(vpException::badValue, "usPolynomialCurve3D::getMeanAxisDeviation: should use at least 2 segment "
1076  "to compute approximate deviation from axis");
1077 
1078  vpColVector params(nbCountSeg + 1);
1079  double step = this->getParametricLength() / nbCountSeg;
1080  params[0] = m_startParameter;
1081  for (int i = 1; i < nbCountSeg + 1; i++)
1082  params[i] = params[i - 1] + step;
1083 
1084  vpMatrix points(this->getPoints(params));
1085 
1086  vpColVector axis = (points.getCol(nbCountSeg) - points.getCol(0)).normalize();
1087  vpColVector origin = points.getCol(0);
1088  double meanDeviation = 0;
1089  for (int i = 0; i < nbCountSeg + 1; i++)
1090  meanDeviation += vpColVector::dotProd(points.getCol(i) - origin, axis);
1091 
1092  meanDeviation /= nbCountSeg + 1;
1093  return meanDeviation;
1094 }
1095 
1096 void usPolynomialCurve3D::setControlPoints(const vpMatrix &controlPoints)
1097 {
1098  if (controlPoints.getRows() != 3)
1099  throw vpException(vpException::dimensionError,
1100  "usPolynomialCurve3D::setControlPoints(const vpMatrix&): invalid points dimension, should be 3");
1101  if (controlPoints.getCols() < 2)
1102  throw vpException(
1103  vpException::dimensionError,
1104  "usPolynomialCurve3D::setControlPoints(const vpMatrix&): invalid number of points, should greater than 1");
1105 
1106  vpColVector direction = controlPoints.getCol(controlPoints.getCols() - 1) - controlPoints.getCol(0);
1107  this->defineFromPointsAuto(controlPoints, direction, controlPoints.getCols() - 1);
1108 }
1109 
1110 void usPolynomialCurve3D::setControlPoints(double **controlPoints)
1111 {
1112  vpMatrix M(3, m_order);
1113  memcpy(M.data, *controlPoints, M.size() * sizeof(double));
1114  this->setControlPoints(M);
1115 }
1116 
1118 {
1119  if (m_order == 0)
1120  return this->getStartPoint();
1121 
1122  vpColVector params(m_order + 1);
1123  double step = this->getParametricLength() / m_order;
1124  params[0] = m_startParameter;
1125  for (unsigned int i = 1; i < m_order + 1; i++)
1126  params[i] = params[i - 1] + step;
1127 
1128  return this->getPoints(params);
1129 }
1130 
1132 {
1133  int nbRenderingPoints = (m_order < 2) ? 2 : 10;
1134  vpColVector params(nbRenderingPoints);
1135  double step = this->getParametricLength() / (nbRenderingPoints - 1);
1136  params[0] = m_startParameter;
1137  for (int i = 1; i < nbRenderingPoints; i++)
1138  params[i] = params[i - 1] + step;
1139 
1140  return this->getPoints(params);
1141 }
1142 
1144 {
1145  unsigned int order1 = n1.getOrder();
1146  unsigned int order2 = n2.getOrder();
1147  vpMatrix coords1(order1, 50);
1148  vpMatrix coords2(order2, 50);
1149  for (unsigned int j = 0; j < 50; ++j) {
1150  coords1[0][j] = 1.0;
1151  coords2[0][j] = 1.0;
1152  double t = static_cast<double>(j) / 49.0;
1153  for (unsigned int i = 1; i < order1; ++i)
1154  coords1[i][j] = coords1[i - 1][j] * t;
1155  for (unsigned int i = 1; i < order2; ++i)
1156  coords2[i][j] = coords2[i - 1][j] * t;
1157  }
1158  vpMatrix p1 = n1.getPolynomialCoefficients() * coords1;
1159  vpMatrix p2 = n2.getPolynomialCoefficients() * coords2;
1160  double distance = 0.0;
1161  for (unsigned int i = 0; i < 50; ++i)
1162  distance += (p1.getCol(i) - p2.getCol(i)).frobeniusNorm();
1163  distance /= 50;
1164  return distance;
1165 }
1166 
1167 usPolynomialCurve3D usPolynomialCurve3D::getSubPolynomialCurve(double startParameter, double endParameter) const
1168 {
1169  usPolynomialCurve3D seg(*this);
1170 
1171  seg.setBoundaries(startParameter, endParameter);
1172 
1173  return seg;
1174 }
1175 
1177 {
1178  usPolynomialCurve3D seg(*this);
1179 
1180  seg.setOrder(order);
1181 
1182  return seg;
1183 }
1184 
1185 void usPolynomialCurve3D::changeCoefficientsToFitBoundaries(double startParameter, double endParameter)
1186 {
1187  if (startParameter == m_startParameter && endParameter == m_endParameter)
1188  return;
1189  if (startParameter == endParameter)
1190  throw vpException(vpException::dimensionError, "usPolynomialCurve3D::changeCoefficientsToFitBoundaries(double "
1191  "startParameter, double endParameter): new parametric boundaries "
1192  "should be different");
1193 
1194  double beta = (m_startParameter * endParameter - m_endParameter * startParameter) / (endParameter - startParameter);
1195  vpColVector beta_pow(m_order + 1);
1196  beta_pow[0] = 1;
1197  for (unsigned int i = 1; i <= m_order; i++)
1198  beta_pow[i] = beta * beta_pow[i - 1];
1199 
1200  double alpha = (m_endParameter - m_startParameter) / (endParameter - startParameter);
1201  vpColVector alpha_pow(m_order + 1);
1202  alpha_pow[0] = 1;
1203  for (unsigned int i = 1; i <= m_order; i++)
1204  alpha_pow[i] = alpha * alpha_pow[i - 1];
1205 
1206  vpMatrix M(m_order + 1, m_order + 1, 0);
1207  M[0][0] = 1;
1208  for (unsigned int i = 1; i <= m_order; i++) {
1209  M[i][0] = 1;
1210  for (unsigned int j = 1; j <= i; j++) {
1211  M[i][j] = M[i - 1][j - 1] + M[i - 1][j];
1212  }
1213  }
1214 
1215  for (unsigned int i = 1; i <= m_order; i++) {
1216  for (unsigned int j = 0; j <= i; j++) {
1217  M[i][j] *= alpha_pow[j] * beta_pow[i - j];
1218  }
1219  }
1220 
1222 
1223  m_startParameter = startParameter;
1224  m_endParameter = endParameter;
1225 }
1226 
1228 
1230 {
1231  double endParameter = m_endParameter;
1233  this->setBoundaries(m_startParameter, endParameter);
1234 }
1235 
1236 void usPolynomialCurve3D::move(const vpHomogeneousMatrix &H)
1237 {
1238  vpMatrix R(H.getRotationMatrix());
1239  vpMatrix M = R * m_polynomialCoefficients;
1240  for (int i = 0; i < 3; i++)
1241  M[i][0] += H[i][3];
1242 
1243  this->setPolynomialCoefficients(M);
1244 }
1245 
1246 void usPolynomialCurve3D::move(double x, double y, double z, double tx, double ty, double tz)
1247 {
1248  this->move(vpHomogeneousMatrix(x, y, z, tx, ty, tz));
1249 }
1250 
1252 
1253 std::ostream &operator<<(std::ostream &s, const usPolynomialCurve3D &seg)
1254 {
1255  s << "usPolynomialCurve3D\n";
1256  s << seg.m_order << '\n';
1257  for (int i = 0; i < 3; i++) {
1258  for (unsigned int j = 0; j < seg.m_order + 1; j++)
1259  s << seg.m_polynomialCoefficients[i][j] << " ";
1260  s << '\n';
1261  }
1262  s << seg.m_startParameter << '\n';
1263  s << seg.m_endParameter << '\n';
1264  s.flush();
1265  return s;
1266 }
1267 
1268 std::istream &operator>>(std::istream &s, usPolynomialCurve3D &seg)
1269 {
1270  std::string c;
1271  s >> c;
1272  if (c != "usPolynomialCurve3D") {
1273  vpException e(vpException::ioError,
1274  "operator>>(std::istream&, Polynomial3D&): Stream does not contain usPolynomialCurve3D data");
1275  throw e;
1276  }
1277  s >> seg.m_order;
1278  seg.m_polynomialCoefficients.resize(3, seg.m_order + 1);
1279  for (unsigned int i = 0; i < 3; i++)
1280  for (unsigned int j = 0; j < seg.m_order + 1; j++)
1281  s >> seg.m_polynomialCoefficients[i][j];
1282  s >> seg.m_startParameter;
1283  s >> seg.m_endParameter;
1284  s.get();
1285  return s;
1286 }
1287 
1288 std::ostream &operator<<=(std::ostream &s, const usPolynomialCurve3D &seg)
1289 {
1290  s.write("usPolynomialCurve3D", 20);
1291  s.write((char *)&(seg.m_order), sizeof(int));
1292  for (unsigned int i = 0; i < 3; i++)
1293  for (unsigned int j = 0; j < seg.m_order + 1; j++)
1294  s.write((char *)&(seg.m_polynomialCoefficients[i][j]), sizeof(double));
1295  s.write((char *)&(seg.m_startParameter), sizeof(double));
1296  s.write((char *)&(seg.m_endParameter), sizeof(double));
1297  s.flush();
1298  return s;
1299 }
1300 
1301 std::istream &operator>>=(std::istream &s, usPolynomialCurve3D &seg)
1302 {
1303  char c[20];
1304  s.read(c, 20);
1305  if (strcmp(c, "usPolynomialCurve3D")) {
1306  vpException e(vpException::ioError,
1307  "operator>>=(std::istream&, Polynomial3D&): Stream does not contain usPolynomialCurve3D data");
1308  throw e;
1309  }
1310  s.read((char *)&(seg.m_order), sizeof(int));
1311  seg.m_polynomialCoefficients.resize(3, seg.m_order + 1);
1312  for (unsigned int i = 0; i < 3; i++)
1313  for (unsigned int j = 0; j < seg.m_order + 1; j++)
1314  s.read((char *)&(seg.m_polynomialCoefficients[i][j]), sizeof(double));
1315  s.read((char *)&(seg.m_startParameter), sizeof(double));
1316  s.read((char *)&(seg.m_endParameter), sizeof(double));
1317  return s;
1318 }
vpColVector getDerivative(double parameter, unsigned int order) const
void defineFromPoints(const std::vector< vpColVector > &points, const std::vector< double > &param, unsigned int order=0)
vpMatrix getControlPoints() const
vpColVector getStartTangent() const
vpMatrix getRenderingPoints() const
usPolynomialCurve3D getSubPolynomialCurve(double startParameter, double endParameter) const
static double curveDistance(const usPolynomialCurve3D &n1, const usPolynomialCurve3D &n2)
usPolynomialCurve3D getNewOrderPolynomialCurve(unsigned int order) const
void setStartParameter(double startParameter)
double getLength(int nbCountSeg=50) const
void defineFromWeightedPoints(const std::vector< vpColVector > &points, const std::vector< double > &param, const std::vector< double > &weights, unsigned int order=0)
vpColVector getTangent(double parameter) const
unsigned int getOrder() const
vpColVector getPoint(double parameter) const
double getMeanAxisDeviation(int nbCountSeg=50) const
double getStartParameter() const
vpMatrix getPoints(const vpColVector &parameters) const
void setParametricLength(double length)
void setPolynomialCoefficients(const vpMatrix &polynomialCoefficients)
vpColVector getEndPoint() const
vpMatrix getPolynomialCoefficients() const
void setLength(double length, double precision=1e-4)
void setOrder(unsigned int order)
void move(const vpHomogeneousMatrix &H)
vpColVector getEndTangent() const
void defineFromWeightedPointsAuto(const std::vector< vpColVector > &points, const std::vector< double > &weights, unsigned int order=0)
void setControlPoints(const vpMatrix &controlPoints)
void defineFromPointsAuto(const std::vector< vpColVector > &points, unsigned int order=0)
double getCurvatureFromShape(double start, double end, vpColVector &center3D, vpColVector &direction3D) const
double getParametricLength() const
void changeCoefficientsToFitBoundaries(double startParameter, double endParameter)
void setEndParameter(double endParameter)
const usPolynomialCurve3D & operator=(const usPolynomialCurve3D &curve)
double getEndParameter() const
vpColVector getStartPoint() const
void setBoundaries(double startParameter, double endParamter)
double getCurvature(double param) const