38 #ifndef PCL_TRACKING_IMPL_PYRAMIDAL_KLT_HPP
39 #define PCL_TRACKING_IMPL_PYRAMIDAL_KLT_HPP
42 #include <pcl/common/utils.h>
43 #include <pcl/common/io.h>
44 #include <pcl/common/utils.h>
47 template <
typename Po
intInT,
typename IntensityT>
inline void
51 track_height_ = height;
55 template <
typename Po
intInT,
typename IntensityT>
inline void
58 if (keypoints->
size () <= keypoints_nbr_)
59 keypoints_ = keypoints;
64 for (std::size_t i = 0; i < keypoints_nbr_; ++i)
70 keypoints_status_->indices.
resize (keypoints_->size (), 0);
74 template <
typename Po
intInT,
typename IntensityT>
inline void
77 assert ((input_ || ref_) &&
"[pcl::tracking::PyramidalKLTTracker] CALL setInputCloud FIRST!");
80 keypoints->
reserve (keypoints_nbr_);
81 for (std::size_t i = 0; i < keypoints_nbr_; ++i)
84 uv.
u = points->indices[i] % input_->width;
85 uv.
v = points->indices[i] / input_->width;
88 setPointsToTrack (keypoints);
92 template <
typename Po
intInT,
typename IntensityT>
bool
98 PCL_ERROR (
"[pcl::tracking::%s::initCompute] PCLBase::Init failed.\n",
99 tracker_name_.c_str ());
103 if (!input_->isOrganized ())
105 PCL_ERROR (
"[pcl::tracking::%s::initCompute] Need an organized point cloud to proceed!",
106 tracker_name_.c_str ());
110 if (!keypoints_ || keypoints_->empty ())
112 PCL_ERROR (
"[pcl::tracking::%s::initCompute] No keypoints aborting!",
113 tracker_name_.c_str ());
123 if ((track_height_ * track_width_)%2 == 0)
125 PCL_ERROR (
"[pcl::tracking::%s::initCompute] Tracking window (%dx%d) must be odd!\n",
126 tracker_name_.c_str (), track_width_, track_height_);
130 if (track_height_ < 3 || track_width_ < 3)
132 PCL_ERROR (
"[pcl::tracking::%s::initCompute] Tracking window (%dx%d) must be >= 3x3!\n",
133 tracker_name_.c_str (), track_width_, track_height_);
137 track_width_2_ = track_width_ / 2;
138 track_height_2_ = track_height_ / 2;
142 PCL_ERROR (
"[pcl::tracking::%s::initCompute] Number of pyramid levels should be at least 2!",
143 tracker_name_.c_str ());
149 PCL_ERROR (
"[pcl::tracking::%s::initCompute] Number of pyramid levels should not exceed 5!",
150 tracker_name_.c_str ());
163 template <
typename Po
intInT,
typename IntensityT>
void
182 float *row0 =
new float [src.
width + 2];
183 float *row1 =
new float [src.
width + 2];
184 float *trow0 = row0; ++trow0;
185 float *trow1 = row1; ++trow1;
186 const float* src_ptr = &(src.
points[0]);
188 for (
int y = 0; y < height; y++)
190 const float* srow0 = src_ptr + (y > 0 ? y-1 : height > 1 ? 1 : 0) * width;
191 const float* srow1 = src_ptr + y * width;
192 const float* srow2 = src_ptr + (y < height-1 ? y+1 : height > 1 ? height-2 : 0) * width;
193 float* grad_x_row = &(grad_x.
points[y * width]);
194 float* grad_y_row = &(grad_y.
points[y * width]);
197 for (
int x = 0; x < width; x++)
199 trow0[x] = (srow0[x] + srow2[x])*3 + srow1[x]*10;
200 trow1[x] = srow2[x] - srow0[x];
204 int x0 = width > 1 ? 1 : 0, x1 = width > 1 ? width-2 : 0;
205 trow0[-1] = trow0[x0]; trow0[width] = trow0[x1];
206 trow1[-1] = trow1[x0]; trow1[width] = trow1[x1];
209 for (
int x = 0; x < width; x++)
211 grad_x_row[x] = trow0[x+1] - trow0[x-1];
212 grad_y_row[x] = (trow1[x+1] + trow1[x-1])*3 + trow1[x]*10;
218 template <
typename Po
intInT,
typename IntensityT>
void
222 FloatImage smoothed (input->width, input->height);
223 convolve (input, smoothed);
225 int width = (smoothed.
width +1) / 2;
226 int height = (smoothed.
height +1) / 2;
227 std::vector<int> ii (width);
228 for (
int i = 0; i < width; ++i)
233 #pragma omp parallel for shared (output) firstprivate (ii) num_threads (threads_)
235 for (
int j = 0; j < height; ++j)
238 for (
int i = 0; i < width; ++i)
239 (*down) (i,j) = smoothed (ii[i],jj);
246 template <
typename Po
intInT,
typename IntensityT>
void
252 downsample (input, output);
255 derivatives (*output, *grad_x, *grad_y);
256 output_grad_x = grad_x;
257 output_grad_y = grad_y;
261 template <
typename Po
intInT,
typename IntensityT>
void
265 convolveRows (input, *tmp);
266 convolveCols (tmp, output);
270 template <
typename Po
intInT,
typename IntensityT>
void
273 int width = input->width;
274 int height = input->height;
275 int last = input->width - kernel_size_2_;
279 #pragma omp parallel for shared (output) num_threads (threads_)
281 for (
int j = 0; j < height; ++j)
283 for (
int i = kernel_size_2_; i < last; ++i)
286 for (
int k = kernel_last_, l = i - kernel_size_2_; k > -1; --k, ++l)
287 result+= kernel_[k] * (*input) (l,j);
289 output (i,j) = static_cast<float> (result);
292 for (
int i = last; i < width; ++i)
293 output (i,j) = output (w, j);
295 for (
int i = 0; i < kernel_size_2_; ++i)
296 output (i,j) = output (kernel_size_2_, j);
301 template <
typename Po
intInT,
typename IntensityT>
void
304 output =
FloatImage (input->width, input->height);
306 int width = input->width;
307 int height = input->height;
308 int last = input->height - kernel_size_2_;
312 #pragma omp parallel for shared (output) num_threads (threads_)
314 for (
int i = 0; i < width; ++i)
316 for (
int j = kernel_size_2_; j < last; ++j)
319 for (
int k = kernel_last_, l = j - kernel_size_2_; k > -1; --k, ++l)
320 result += kernel_[k] * (*input) (i,l);
321 output (i,j) = static_cast<float> (result);
324 for (
int j = last; j < height; ++j)
325 output (i,j) = output (i,h);
327 for (
int j = 0; j < kernel_size_2_; ++j)
328 output (i,j) = output (i, kernel_size_2_);
333 template <
typename Po
intInT,
typename IntensityT>
void
335 std::vector<FloatImageConstPtr>& pyramid,
339 pyramid.resize (step * nb_levels_);
344 #pragma omp parallel for num_threads (threads_)
346 for (
int i = 0; i < static_cast<int> (input->size ()); ++i)
347 tmp->points[i] = intensity_ (input->points[i]);
351 previous->height + 2*track_height_));
360 derivatives (*img, *g_x, *g_y);
363 previous->height + 2*track_height_));
364 pcl::copyPointCloud (*g_x, *grad_x, track_height_, track_height_, track_width_, track_width_,
369 previous->height + 2*track_height_));
370 pcl::copyPointCloud (*g_y, *grad_y, track_height_, track_height_, track_width_, track_width_,
374 for (
int level = 1; level < nb_levels_; ++level)
380 downsample (previous, current, g_x, g_y);
383 current->height + 2*track_height_));
384 pcl::copyPointCloud (*current, *image, track_height_, track_height_, track_width_, track_width_,
386 pyramid[level*step] = image;
388 pcl::copyPointCloud (*g_x, *gradx, track_height_, track_height_, track_width_, track_width_,
390 pyramid[level*step + 1] = gradx;
392 pcl::copyPointCloud (*g_y, *grady, track_height_, track_height_, track_width_, track_width_,
394 pyramid[level*step + 2] = grady;
401 template <
typename Po
intInT,
typename IntensityT>
void
405 const Eigen::Array2i& location,
406 const Eigen::Array4f& weight,
407 Eigen::ArrayXXf& win,
408 Eigen::ArrayXXf& grad_x_win,
409 Eigen::ArrayXXf& grad_y_win,
410 Eigen::Array3f &covariance)
const
412 const int step = img.
width;
413 covariance.setZero ();
415 for (
int y = 0; y < track_height_; y++)
417 const float* img_ptr = &(img.
points[0]) + (y + location[1])*step + location[0];
418 const float* grad_x_ptr = &(grad_x.
points[0]) + (y + location[1])*step + location[0];
419 const float* grad_y_ptr = &(grad_y.
points[0]) + (y + location[1])*step + location[0];
421 float* win_ptr = win.data () + y*win.cols ();
422 float* grad_x_win_ptr = grad_x_win.data () + y*grad_x_win.cols ();
423 float* grad_y_win_ptr = grad_y_win.data () + y*grad_y_win.cols ();
425 for (
int x =0; x < track_width_; ++x, ++grad_x_ptr, ++grad_y_ptr)
427 *win_ptr++ = img_ptr[x]*weight[0] + img_ptr[x+1]*weight[1] + img_ptr[x+step]*weight[2] + img_ptr[x+step+1]*weight[3];
428 float ixval = grad_x_ptr[0]*weight[0] + grad_x_ptr[1]*weight[1] + grad_x_ptr[step]*weight[2] + grad_x_ptr[step+1]*weight[3];
429 float iyval = grad_y_ptr[0]*weight[0] + grad_y_ptr[1]*weight[1] + grad_y_ptr[step]*weight[2] + grad_y_ptr[step+1]*weight[3];
431 *grad_x_win_ptr++ = ixval;
432 *grad_y_win_ptr++ = iyval;
434 covariance[0] += ixval*ixval;
435 covariance[1] += ixval*iyval;
436 covariance[2] += iyval*iyval;
442 template <
typename Po
intInT,
typename IntensityT>
void
444 const Eigen::ArrayXXf& prev_grad_x,
445 const Eigen::ArrayXXf& prev_grad_y,
447 const Eigen::Array2i& location,
448 const Eigen::Array4f& weight,
449 Eigen::Array2f &b)
const
451 const int step = next.
width;
453 for (
int y = 0; y < track_height_; y++)
455 const float* next_ptr = &(next.
points[0]) + (y + location[1])*step + location[0];
456 const float* prev_ptr = prev.data () + y*prev.cols ();
457 const float* prev_grad_x_ptr = prev_grad_x.data () + y*prev_grad_x.cols ();
458 const float* prev_grad_y_ptr = prev_grad_y.data () + y*prev_grad_y.cols ();
460 for (
int x = 0; x < track_width_; ++x, ++prev_grad_y_ptr, ++prev_grad_x_ptr)
462 float diff = next_ptr[x]*weight[0] + next_ptr[x+1]*weight[1]
463 + next_ptr[x+step]*weight[2] + next_ptr[x+step+1]*weight[3] - prev_ptr[x];
464 b[0] += *prev_grad_x_ptr * diff;
465 b[1] += *prev_grad_y_ptr * diff;
471 template <
typename Po
intInT,
typename IntensityT>
void
474 const std::vector<FloatImageConstPtr>& prev_pyramid,
475 const std::vector<FloatImageConstPtr>& pyramid,
478 std::vector<int>& status,
479 Eigen::Affine3f& motion)
const
481 std::vector<Eigen::Array2f, Eigen::aligned_allocator<Eigen::Array2f> > next_pts (prev_keypoints->
size ());
482 Eigen::Array2f half_win ((track_width_-1)*0.5f, (track_height_-1)*0.5f);
484 const int nb_points = prev_keypoints->
size ();
485 for (
int level = nb_levels_ - 1; level >= 0; --level)
487 const FloatImage& prev = *(prev_pyramid[level*3]);
489 const FloatImage& grad_x = *(prev_pyramid[level*3+1]);
490 const FloatImage& grad_y = *(prev_pyramid[level*3+2]);
492 Eigen::ArrayXXf prev_win (track_height_, track_width_);
493 Eigen::ArrayXXf grad_x_win (track_height_, track_width_);
494 Eigen::ArrayXXf grad_y_win (track_height_, track_width_);
495 float ratio (1./(1 << level));
496 for (
int ptidx = 0; ptidx < nb_points; ptidx++)
498 Eigen::Array2f prev_pt (prev_keypoints->
points[ptidx].u * ratio,
499 prev_keypoints->
points[ptidx].v * ratio);
500 Eigen::Array2f next_pt;
501 if (level == nb_levels_ -1)
504 next_pt = next_pts[ptidx]*2.f;
506 next_pts[ptidx] = next_pt;
508 Eigen::Array2i iprev_point;
510 iprev_point[0] = std::floor (prev_pt[0]);
511 iprev_point[1] = std::floor (prev_pt[1]);
513 if (iprev_point[0] < -track_width_ || (std::uint32_t) iprev_point[0] >= grad_x.
width ||
514 iprev_point[1] < -track_height_ || (std::uint32_t) iprev_point[1] >= grad_y.
height)
521 float a = prev_pt[0] - iprev_point[0];
522 float b = prev_pt[1] - iprev_point[1];
523 Eigen::Array4f weight;
524 weight[0] = (1.f - a)*(1.f - b);
525 weight[1] = a*(1.f - b);
526 weight[2] = (1.f - a)*b;
527 weight[3] = 1 - weight[0] - weight[1] - weight[2];
529 Eigen::Array3f covar = Eigen::Array3f::Zero ();
530 spatialGradient (prev, grad_x, grad_y, iprev_point, weight, prev_win, grad_x_win, grad_y_win, covar);
532 float det = covar[0]*covar[2] - covar[1]*covar[1];
533 float min_eigenvalue = (covar[2] + covar[0] - std::sqrt ((covar[0]-covar[2])*(covar[0]-covar[2]) + 4.f*covar[1]*covar[1]))/2.f;
535 if (min_eigenvalue < min_eigenvalue_threshold_ || det < std::numeric_limits<float>::epsilon ())
544 Eigen::Array2f prev_delta (0, 0);
545 for (
unsigned int j = 0; j < max_iterations_; j++)
547 Eigen::Array2i inext_pt = next_pt.floor ().cast<
int> ();
549 if (inext_pt[0] < -track_width_ || (std::uint32_t) inext_pt[0] >= next.
width ||
550 inext_pt[1] < -track_height_ || (std::uint32_t) inext_pt[1] >= next.
height)
557 a = next_pt[0] - inext_pt[0];
558 b = next_pt[1] - inext_pt[1];
559 weight[0] = (1.f - a)*(1.f - b);
560 weight[1] = a*(1.f - b);
561 weight[2] = (1.f - a)*b;
562 weight[3] = 1 - weight[0] - weight[1] - weight[2];
564 Eigen::Array2f beta = Eigen::Array2f::Zero ();
565 mismatchVector (prev_win, grad_x_win, grad_y_win, next, inext_pt, weight, beta);
567 Eigen::Vector2f delta ((covar[1]*beta[1] - covar[2]*beta[0])*det, (covar[1]*beta[0] - covar[0]*beta[1])*det);
569 next_pt[0] += delta[0]; next_pt[1] += delta[1];
570 next_pts[ptidx] = next_pt + half_win;
572 if (delta.squaredNorm () <= epsilon_)
575 if (j > 0 && std::abs (delta[0] + prev_delta[0]) < 0.01 &&
576 std::abs (delta[1] + prev_delta[1]) < 0.01 )
578 next_pts[ptidx][0] -= delta[0]*0.5f;
579 next_pts[ptidx][1] -= delta[1]*0.5f;
587 if (level == 0 && !status[ptidx])
589 Eigen::Array2f next_point = next_pts[ptidx] - half_win;
590 Eigen::Array2i inext_point;
592 inext_point[0] = std::floor (next_point[0]);
593 inext_point[1] = std::floor (next_point[1]);
595 if (inext_point[0] < -track_width_ || (std::uint32_t) inext_point[0] >= next.
width ||
596 inext_point[1] < -track_height_ || (std::uint32_t) inext_point[1] >= next.
height)
603 n.
u = next_pts[ptidx][0];
604 n.
v = next_pts[ptidx][1];
607 inext_point[0] = std::floor (next_pts[ptidx][0]);
608 inext_point[1] = std::floor (next_pts[ptidx][1]);
609 iprev_point[0] = std::floor (prev_keypoints->
points[ptidx].u);
610 iprev_point[1] = std::floor (prev_keypoints->
points[ptidx].v);
611 const PointInT& prev_pt = prev_input->points[iprev_point[1]*prev_input->width + iprev_point[0]];
612 const PointInT& next_pt = input->points[inext_point[1]*input->width + inext_point[0]];
613 transformation_computer.
add (prev_pt.getVector3fMap (), next_pt.getVector3fMap (), 1.0);
621 template <
typename Po
intInT,
typename IntensityT>
void
627 std::vector<FloatImageConstPtr> pyramid;
630 keypoints->
reserve (keypoints_->size ());
631 std::vector<int> status (keypoints_->size (), 0);
632 track (ref_, input_, ref_pyramid_, pyramid, keypoints_, keypoints, status, motion_);
635 ref_pyramid_ = pyramid;
636 keypoints_ = keypoints;
637 keypoints_status_->indices = status;