exact_cubic.h
Go to the documentation of this file.
1 
19 #ifndef _CLASS_EXACTCUBIC
20 #define _CLASS_EXACTCUBIC
21 
22 #include <functional>
23 #include <vector>
24 
25 #include "MathDefs.h"
26 #include "curve_abc.h"
27 #include "curve_constraint.h"
28 #include "piecewise_curve.h"
29 #include "polynomial.h"
30 
31 namespace ndcurves {
36 template <
37  typename Time = double, typename Numeric = Time, bool Safe = false,
38  typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>,
39  typename T_Point = std::vector<Point, Eigen::aligned_allocator<Point> >,
40  typename SplineBase = polynomial<Time, Numeric, Safe, Point, T_Point> >
41 struct exact_cubic : public piecewise_curve<Time, Numeric, Safe, Point> {
42  typedef Point point_t;
43  typedef const Eigen::Ref<const point_t> point_ref_t;
44  typedef T_Point t_point_t;
45  typedef Eigen::Matrix<Numeric, Eigen::Dynamic, Eigen::Dynamic> MatrixX;
46  typedef Eigen::Matrix<Numeric, 3, 3> Matrix3;
47  typedef Time time_t;
48  typedef Numeric num_t;
49  typedef SplineBase spline_t;
50  typedef typename std::vector<spline_t> t_spline_t;
51  typedef typename t_spline_t::iterator it_spline_t;
52  typedef typename t_spline_t::const_iterator cit_spline_t;
54 
61 
62  /* Constructors - destructors */
63  public:
68 
74  template <typename In>
75  exact_cubic(In wayPointsBegin, In wayPointsEnd) : piecewise_curve_t() {
76  t_spline_t subSplines = computeWayPoints<In>(wayPointsBegin, wayPointsEnd);
77  for (cit_spline_t it = subSplines.begin(); it != subSplines.end(); ++it) {
78  this->add_curve(*it);
79  }
80  }
81 
88  template <typename In>
89  exact_cubic(In wayPointsBegin, In wayPointsEnd,
90  const spline_constraints& constraints)
91  : piecewise_curve_t() {
92  t_spline_t subSplines =
93  computeWayPoints<In>(wayPointsBegin, wayPointsEnd, constraints);
94  for (cit_spline_t it = subSplines.begin(); it != subSplines.end(); ++it) {
95  this->add_curve(*it);
96  }
97  }
98 
101  exact_cubic(const t_spline_t& subSplines) : piecewise_curve_t() {
102  for (cit_spline_t it = subSplines.begin(); it != subSplines.end(); ++it) {
103  this->add_curve(*it);
104  }
105  }
106 
107  exact_cubic(const t_curve_ptr_t& subSplines)
108  : piecewise_curve_t(subSplines) {}
109 
111  exact_cubic(const exact_cubic& other) : piecewise_curve_t(other) {}
112 
114  virtual ~exact_cubic() {}
115 
116  std::size_t getNumberSplines() { return this->getNumberCurves(); }
117 
118  spline_t getSplineAt(std::size_t index) {
119  boost::shared_ptr<spline_t> s_ptr =
120  boost::dynamic_pointer_cast<spline_t>(this->curves_.at(index));
121  if (s_ptr)
122  return *s_ptr;
123  else
124  throw std::runtime_error(
125  "Parent piecewise curve do not contain only curves created from "
126  "exact_cubic class methods");
127  }
128 
129  private:
130  static polynomial_t create_cubic(point_ref_t a, point_ref_t b, point_ref_t c,
131  point_ref_t d, const time_t t_min,
132  const time_t t_max) {
133  typename polynomial_t::t_point_t coeffs;
134  coeffs.push_back(a);
135  coeffs.push_back(b);
136  coeffs.push_back(c);
137  coeffs.push_back(d);
138  return polynomial_t(coeffs.begin(), coeffs.end(), t_min, t_max);
139  }
140  static polynomial_t create_quintic(point_ref_t a, point_ref_t b,
143  const time_t t_min, const time_t t_max) {
144  typename polynomial_t::t_point_t coeffs;
145  coeffs.push_back(a);
146  coeffs.push_back(b);
147  coeffs.push_back(c);
148  coeffs.push_back(d);
149  coeffs.push_back(e);
150  coeffs.push_back(f);
151  return polynomial_t(coeffs.begin(), coeffs.end(), t_min, t_max);
152  }
153 
161  template <typename In>
162  t_spline_t computeWayPoints(In wayPointsBegin, In wayPointsEnd) const {
163  const std::size_t dim = wayPointsBegin->second.size();
164  const std::size_t size = std::distance(wayPointsBegin, wayPointsEnd);
165  if (Safe && size < 1) {
166  throw std::length_error(
167  "size of waypoints must be superior to 0"); // TODO
168  }
169  t_spline_t subSplines;
170  subSplines.reserve(size);
171  // refer to the paper to understand all this.
172  MatrixX h1 = MatrixX::Zero(size, size);
173  MatrixX h2 = MatrixX::Zero(size, size);
174  MatrixX h3 = MatrixX::Zero(size, size);
175  MatrixX h4 = MatrixX::Zero(size, size);
176  MatrixX h5 = MatrixX::Zero(size, size);
177  MatrixX h6 = MatrixX::Zero(size, size);
178  MatrixX a = MatrixX::Zero(size, dim);
179  MatrixX b = MatrixX::Zero(size, dim);
180  MatrixX c = MatrixX::Zero(size, dim);
181  MatrixX d = MatrixX::Zero(size, dim);
182  MatrixX x = MatrixX::Zero(size, dim);
183  In it(wayPointsBegin), next(wayPointsBegin);
184  ++next;
185  // Fill the matrices H as specified in the paper.
186  for (std::size_t i(0); next != wayPointsEnd; ++next, ++it, ++i) {
187  num_t const dTi((*next).first - (*it).first);
188  num_t const dTi_sqr(dTi * dTi);
189  num_t const dTi_cube(dTi_sqr * dTi);
190  // filling matrices values
191  h3(i, i) = -3 / dTi_sqr;
192  h3(i, i + 1) = 3 / dTi_sqr;
193  h4(i, i) = -2 / dTi;
194  h4(i, i + 1) = -1 / dTi;
195  h5(i, i) = 2 / dTi_cube;
196  h5(i, i + 1) = -2 / dTi_cube;
197  h6(i, i) = 1 / dTi_sqr;
198  h6(i, i + 1) = 1 / dTi_sqr;
199  if (i + 2 < size) {
200  In it2(next);
201  ++it2;
202  num_t const dTi_1((*it2).first - (*next).first);
203  num_t const dTi_1sqr(dTi_1 * dTi_1);
204  // this can be optimized but let's focus on clarity as long as not
205  // needed
206  h1(i + 1, i) = 2 / dTi;
207  h1(i + 1, i + 1) = 4 / dTi + 4 / dTi_1;
208  h1(i + 1, i + 2) = 2 / dTi_1;
209  h2(i + 1, i) = -6 / dTi_sqr;
210  h2(i + 1, i + 1) = (6 / dTi_1sqr) - (6 / dTi_sqr);
211  h2(i + 1, i + 2) = 6 / dTi_1sqr;
212  }
213  x.row(i) = (*it).second.transpose();
214  }
215  // adding last x
216  x.row(size - 1) = (*it).second.transpose();
217  // Compute coefficients of polynom.
218  a = x;
219  PseudoInverse(h1);
220  b = h1 * h2 * x; // h1 * b = h2 * x => b = (h1)^-1 * h2 * x
221  c = h3 * x + h4 * b;
222  d = h5 * x + h6 * b;
223  // create splines along waypoints.
224  it = wayPointsBegin, next = wayPointsBegin;
225  ++next;
226  for (int i = 0; next != wayPointsEnd; ++i, ++it, ++next) {
227  subSplines.push_back(create_cubic(a.row(i), b.row(i), c.row(i), d.row(i),
228  (*it).first, (*next).first));
229  }
230  return subSplines;
231  }
232 
233  template <typename In>
234  t_spline_t computeWayPoints(In wayPointsBegin, In wayPointsEnd,
235  const spline_constraints& constraints) const {
236  std::size_t const size(std::distance(wayPointsBegin, wayPointsEnd));
237  if (Safe && size < 1) {
238  throw std::length_error(
239  "number of waypoints should be superior to one"); // TODO
240  }
241  t_spline_t subSplines;
242  subSplines.reserve(size - 1);
243  spline_constraints cons = constraints;
244  In it(wayPointsBegin), next(wayPointsBegin), end(wayPointsEnd - 1);
245  ++next;
246  for (std::size_t i(0); next != end; ++next, ++it, ++i) {
247  compute_one_spline<In>(it, next, cons, subSplines);
248  }
249  compute_end_spline<In>(it, next, cons, subSplines);
250  return subSplines;
251  }
252 
253  template <typename In>
254  void compute_one_spline(In wayPointsBegin, In wayPointsNext,
255  spline_constraints& constraints,
256  t_spline_t& subSplines) const {
257  const point_t &a0 = wayPointsBegin->second, a1 = wayPointsNext->second;
258  const point_t &b0 = constraints.init_vel, c0 = constraints.init_acc / 2.;
259  const num_t &init_t = wayPointsBegin->first, end_t = wayPointsNext->first;
260  const num_t dt = end_t - init_t, dt_2 = dt * dt, dt_3 = dt_2 * dt;
261  const point_t d0 = (a1 - a0 - b0 * dt - c0 * dt_2) / dt_3;
262  subSplines.push_back(create_cubic(a0, b0, c0, d0, init_t, end_t));
263  constraints.init_vel = subSplines.back().derivate(end_t, 1);
264  constraints.init_acc = subSplines.back().derivate(end_t, 2);
265  }
266 
267  template <typename In>
268  void compute_end_spline(In wayPointsBegin, In wayPointsNext,
269  spline_constraints& constraints,
270  t_spline_t& subSplines) const {
271  const std::size_t dim = wayPointsBegin->second.size();
272  const point_t &a0 = wayPointsBegin->second, a1 = wayPointsNext->second;
273  const point_t &b0 = constraints.init_vel, b1 = constraints.end_vel,
274  c0 = constraints.init_acc / 2., c1 = constraints.end_acc;
275  const num_t &init_t = wayPointsBegin->first, end_t = wayPointsNext->first;
276  const num_t dt = end_t - init_t, dt_2 = dt * dt, dt_3 = dt_2 * dt,
277  dt_4 = dt_3 * dt, dt_5 = dt_4 * dt;
278  // solving a system of four linear eq with 4 unknows: d0, e0
279  const point_t alpha_0 = a1 - a0 - b0 * dt - c0 * dt_2;
280  const point_t alpha_1 = b1 - b0 - 2 * c0 * dt;
281  const point_t alpha_2 = c1 - 2 * c0;
282  const num_t x_d_0 = dt_3, x_d_1 = 3 * dt_2, x_d_2 = 6 * dt;
283  const num_t x_e_0 = dt_4, x_e_1 = 4 * dt_3, x_e_2 = 12 * dt_2;
284  const num_t x_f_0 = dt_5, x_f_1 = 5 * dt_4, x_f_2 = 20 * dt_3;
285  point_t d, e, f;
286  MatrixX rhs = MatrixX::Zero(3, dim);
287  rhs.row(0) = alpha_0;
288  rhs.row(1) = alpha_1;
289  rhs.row(2) = alpha_2;
290  Matrix3 eq = Matrix3::Zero(3, 3);
291  eq(0, 0) = x_d_0;
292  eq(0, 1) = x_e_0;
293  eq(0, 2) = x_f_0;
294  eq(1, 0) = x_d_1;
295  eq(1, 1) = x_e_1;
296  eq(1, 2) = x_f_1;
297  eq(2, 0) = x_d_2;
298  eq(2, 1) = x_e_2;
299  eq(2, 2) = x_f_2;
300  rhs = eq.inverse().eval() * rhs;
301  d = rhs.row(0);
302  e = rhs.row(1);
303  f = rhs.row(2);
304  subSplines.push_back(create_quintic(a0, b0, c0, d, e, f, init_t, end_t));
305  }
306 
307  public:
308  // Serialization of the class
310  template <class Archive>
311  void serialize(Archive& ar, const unsigned int version) {
312  if (version) {
313  // Do something depending on version ?
314  }
315  ar& BOOST_SERIALIZATION_BASE_OBJECT_NVP(piecewise_curve_t);
316  }
317 };
318 } // namespace ndcurves
319 
321  SINGLE_ARG(typename Time, typename Numeric, bool Safe, typename Point,
322  typename T_Point, typename SplineBase),
323  SINGLE_ARG(
325 #endif //_CLASS_EXACTCUBIC
#define DEFINE_CLASS_TEMPLATE_VERSION(Template, Type)
Definition: archive.hpp:27
#define SINGLE_ARG(...)
Definition: archive.hpp:23
interface for a Curve of arbitrary dimension.
struct to define constraints on start / end velocities and acceleration on a curve
Eigen::Matrix< Numeric, Eigen::Dynamic, 1 > Point
Definition: effector_spline.h:28
double Numeric
Definition: effector_spline.h:26
double Time
Definition: effector_spline.h:27
std::vector< Point, Eigen::aligned_allocator< Point > > T_Point
Definition: effector_spline.h:29
Definition: bernstein.h:20
void PseudoInverse(_Matrix_Type_ &pinvmat)
An inverse kinematics architecture enforcing an arbitrary number of strict priority levels (Reference...
Definition: MathDefs.h:26
class allowing to create a piecewise curve.
Definition of a cubic spline.
Represents a curve of dimension Dim. If value of parameter Safe is false, no verification is made on ...
Definition: curve_abc.h:37
Definition: curve_constraint.h:20
Definition: exact_cubic.h:41
void serialize(Archive &ar, const unsigned int version)
Definition: exact_cubic.h:311
curve_abc< Time, Numeric, Safe, point_t > curve_abc_t
Definition: exact_cubic.h:57
virtual ~exact_cubic()
Destructor.
Definition: exact_cubic.h:114
t_spline_t::const_iterator cit_spline_t
Definition: exact_cubic.h:52
exact_cubic()
Empty constructor. Add at least one curve to call other class functions.
Definition: exact_cubic.h:67
Numeric num_t
Definition: exact_cubic.h:48
spline_t getSplineAt(std::size_t index)
Definition: exact_cubic.h:118
exact_cubic< Time, Numeric, Safe, point_t, T_Point, SplineBase > exact_cubic_t
Definition: exact_cubic.h:56
Time time_t
Definition: exact_cubic.h:47
SplineBase spline_t
Definition: exact_cubic.h:49
Eigen::Matrix< Numeric, 3, 3 > Matrix3
Definition: exact_cubic.h:46
t_spline_t::iterator it_spline_t
Definition: exact_cubic.h:51
T_Point t_point_t
Definition: exact_cubic.h:44
piecewise_curve_t::t_curve_ptr_t t_curve_ptr_t
Definition: exact_cubic.h:60
exact_cubic(const t_curve_ptr_t &subSplines)
Definition: exact_cubic.h:107
exact_cubic(const t_spline_t &subSplines)
Constructor.
Definition: exact_cubic.h:101
piecewise_curve< Time, Numeric, Safe, point_t > piecewise_curve_t
Definition: exact_cubic.h:58
polynomial< Time, Numeric, Safe, point_t > polynomial_t
Definition: exact_cubic.h:59
exact_cubic(const exact_cubic &other)
Copy Constructor.
Definition: exact_cubic.h:111
const Eigen::Ref< const point_t > point_ref_t
Definition: exact_cubic.h:43
curve_constraints< Point > spline_constraints
Definition: exact_cubic.h:53
Eigen::Matrix< Numeric, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Definition: exact_cubic.h:45
friend class boost::serialization::access
Definition: exact_cubic.h:309
exact_cubic(In wayPointsBegin, In wayPointsEnd)
Constructor.
Definition: exact_cubic.h:75
std::size_t getNumberSplines()
Definition: exact_cubic.h:116
std::vector< spline_t > t_spline_t
Definition: exact_cubic.h:50
Point point_t
Definition: exact_cubic.h:42
exact_cubic(In wayPointsBegin, In wayPointsEnd, const spline_constraints &constraints)
Constructor.
Definition: exact_cubic.h:89
Definition: piecewise_curve.h:37
virtual std::size_t dim() const
Get dimension of curve.
Definition: piecewise_curve.h:641
std::vector< curve_ptr_t > t_curve_ptr_t
Definition: piecewise_curve.h:50
Represents a polynomial of an arbitrary order defined on the interval . It follows the equation : ...
Definition: polynomial.h:35
T_Point t_point_t
Definition: polynomial.h:37