92 #include <pcl/cuda/point_cloud.h>
93 #include <pcl/cuda/cutil_math.h>
106 std::numeric_limits<float>::epsilon() * std::numeric_limits<float>::epsilon();
107 return x * x <= prec_sqr * y * y;
123 float invnm = 1.0f / sqrtf (src.x*src.x + src.y*src.y);
124 perp.x = -src.y*invnm;
125 perp.y = src.x*invnm;
134 float invnm = 1.0f / sqrtf (src.z*src.z + src.y*src.y);
136 perp.y = -src.z*invnm;
137 perp.z = src.y*invnm;
143 inline __host__ __device__
void computeRoots2 (
const float& b,
const float& c, float3& roots)
146 float d = b * b - 4.0f * c;
152 roots.z = 0.5f * (b + sd);
153 roots.y = 0.5f * (b - sd);
156 inline __host__ __device__
void swap (
float& a,
float& b)
158 float c(a); a=b; b=c;
163 inline __host__ __device__
void
174 float c1 = m.
data[0].x*m.
data[1].y -
183 if (std::abs(c0) < std::numeric_limits<float>::epsilon())
187 const float s_inv3 = 1.0f/3.0f;
188 const float s_sqrt3 = sqrtf (3.0f);
191 float c2_over_3 = c2 * s_inv3;
192 float a_over_3 = (c1 - c2 * c2_over_3) * s_inv3;
196 float half_b = 0.5f * (c0 + c2_over_3 * (2.0f * c2_over_3 * c2_over_3 - c1));
198 float q = half_b * half_b + a_over_3 * a_over_3 * a_over_3;
203 float rho = sqrtf (-a_over_3);
204 float theta = std::atan2 (sqrtf (-q), half_b) * s_inv3;
205 float cos_theta = std::cos (theta);
206 float sin_theta = sin (theta);
207 roots.x = c2_over_3 + 2.f * rho * cos_theta;
208 roots.y = c2_over_3 - rho * (cos_theta + s_sqrt3 * sin_theta);
209 roots.z = c2_over_3 - rho * (cos_theta - s_sqrt3 * sin_theta);
212 if (roots.x >= roots.y)
213 swap (roots.x, roots.y);
214 if (roots.y >= roots.z)
216 swap (roots.y, roots.z);
217 if (roots.x >= roots.y)
218 swap (roots.x, roots.y);
226 inline __host__ __device__
void
229 evals = evecs.
data[0] = evecs.
data[1] = evecs.
data[2] = make_float3 (0.0f, 0.0f, 0.0f);
235 float3 scale_tmp = fmaxf (fmaxf (fabs (mat.
data[0]), fabs (mat.
data[1])), fabs (mat.
data[2]));
236 float scale = fmaxf (fmaxf (scale_tmp.x, scale_tmp.y), scale_tmp.z);
237 if (scale <= std::numeric_limits<float>::min())
241 scaledMat.
data[0] = mat.
data[0] / scale;
242 scaledMat.
data[1] = mat.
data[1] / scale;
243 scaledMat.
data[2] = mat.
data[2] / scale;
248 if ((evals.z - evals.x) <= std::numeric_limits<float>::epsilon())
251 evecs.
data[0] = make_float3 (1.0f, 0.0f, 0.0f);
252 evecs.
data[1] = make_float3 (0.0f, 1.0f, 0.0f);
253 evecs.
data[2] = make_float3 (0.0f, 0.0f, 1.0f);
255 else if ((evals.y - evals.x) <= std::numeric_limits<float>::epsilon())
263 tmp.
data[0].x -= evals.z;
264 tmp.
data[1].y -= evals.z;
265 tmp.
data[2].z -= evals.z;
267 float3 vec1 = cross (tmp.
data[0], tmp.
data[1]);
268 float3 vec2 = cross (tmp.
data[0], tmp.
data[2]);
269 float3 vec3 = cross (tmp.
data[1], tmp.
data[2]);
271 float len1 = dot (vec1, vec1);
272 float len2 = dot (vec2, vec2);
273 float len3 = dot (vec3, vec3);
275 if (len1 >= len2 && len1 >= len3)
276 evecs.
data[2] = vec1 / sqrtf (len1);
277 else if (len2 >= len1 && len2 >= len3)
278 evecs.
data[2] = vec2 / sqrtf (len2);
280 evecs.
data[2] = vec3 / sqrtf (len3);
285 else if ((evals.z - evals.y) <= std::numeric_limits<float>::epsilon())
292 tmp.
data[0].x -= evals.x;
293 tmp.
data[1].y -= evals.x;
294 tmp.
data[2].z -= evals.x;
296 float3 vec1 = cross (tmp.
data[0], tmp.
data[1]);
297 float3 vec2 = cross (tmp.
data[0], tmp.
data[2]);
298 float3 vec3 = cross (tmp.
data[1], tmp.
data[2]);
300 float len1 = dot (vec1, vec1);
301 float len2 = dot (vec2, vec2);
302 float len3 = dot (vec3, vec3);
304 if (len1 >= len2 && len1 >= len3)
305 evecs.
data[0] = vec1 / sqrtf (len1);
306 else if (len2 >= len1 && len2 >= len3)
307 evecs.
data[0] = vec2 / sqrtf (len2);
309 evecs.
data[0] = vec3 / sqrtf (len3);
320 tmp.
data[0].x -= evals.z;
321 tmp.
data[1].y -= evals.z;
322 tmp.
data[2].z -= evals.z;
324 float3 vec1 = cross (tmp.
data[0], tmp.
data[1]);
325 float3 vec2 = cross (tmp.
data[0], tmp.
data[2]);
326 float3 vec3 = cross (tmp.
data[1], tmp.
data[2]);
328 float len1 = dot (vec1, vec1);
329 float len2 = dot (vec2, vec2);
330 float len3 = dot (vec3, vec3);
333 unsigned int min_el = 2;
334 unsigned int max_el = 2;
335 if (len1 >= len2 && len1 >= len3)
338 evecs.
data[2] = vec1 / sqrtf (len1);
340 else if (len2 >= len1 && len2 >= len3)
343 evecs.
data[2] = vec2 / sqrtf (len2);
348 evecs.
data[2] = vec3 / sqrtf (len3);
354 tmp.
data[0].x -= evals.y;
355 tmp.
data[1].y -= evals.y;
356 tmp.
data[2].z -= evals.y;
358 vec1 = cross (tmp.
data[0], tmp.
data[1]);
359 vec2 = cross (tmp.
data[0], tmp.
data[2]);
360 vec3 = cross (tmp.
data[1], tmp.
data[2]);
362 len1 = dot (vec1, vec1);
363 len2 = dot (vec2, vec2);
364 len3 = dot (vec3, vec3);
365 if (len1 >= len2 && len1 >= len3)
368 evecs.
data[1] = vec1 / sqrtf (len1);
369 min_el = len1 <= mmax[min_el]? 1: min_el;
370 max_el = len1 > mmax[max_el]? 1: max_el;
372 else if (len2 >= len1 && len2 >= len3)
375 evecs.
data[1] = vec2 / sqrtf (len2);
376 min_el = len2 <= mmax[min_el]? 1: min_el;
377 max_el = len2 > mmax[max_el]? 1: max_el;
382 evecs.
data[1] = vec3 / sqrtf (len3);
383 min_el = len3 <= mmax[min_el]? 1: min_el;
384 max_el = len3 > mmax[max_el]? 1: max_el;
390 tmp.
data[0].x -= evals.x;
391 tmp.
data[1].y -= evals.x;
392 tmp.
data[2].z -= evals.x;
394 vec1 = cross (tmp.
data[0], tmp.
data[1]);
395 vec2 = cross (tmp.
data[0], tmp.
data[2]);
396 vec3 = cross (tmp.
data[1], tmp.
data[2]);
398 len1 = dot (vec1, vec1);
399 len2 = dot (vec2, vec2);
400 len3 = dot (vec3, vec3);
401 if (len1 >= len2 && len1 >= len3)
404 evecs.
data[0] = vec1 / sqrtf (len1);
405 min_el = len3 <= mmax[min_el]? 0: min_el;
406 max_el = len3 > mmax[max_el]? 0: max_el;
408 else if (len2 >= len1 && len2 >= len3)
411 evecs.
data[0] = vec2 / sqrtf (len2);
412 min_el = len3 <= mmax[min_el]? 0: min_el;
413 max_el = len3 > mmax[max_el]? 0: max_el;
418 evecs.
data[0] = vec3 / sqrtf (len3);
419 min_el = len3 <= mmax[min_el]? 0: min_el;
420 max_el = len3 > mmax[max_el]? 0: max_el;
423 unsigned mid_el = 3 - min_el - max_el;
424 evecs.
data[min_el] = normalize (cross (evecs.
data[(min_el+1)%3], evecs.
data[(min_el+2)%3] ));
425 evecs.
data[mid_el] = normalize (cross (evecs.
data[(mid_el+1)%3], evecs.
data[(mid_el+2)%3] ));
434 __inline__ __host__ __device__ float3
444 __inline__ __host__ __device__
459 __inline__ __host__ __device__ float3
467 __inline__ __host__ __device__
475 cov.
data[1].y = pt.y * pt.y;
476 cov.
data[1].z = pt.y * pt.z;
477 cov.
data[2].z = pt.z * pt.z;
480 cov.
data[0].x = pt.x;
481 cov.
data[0].y = pt.y;
482 cov.
data[0].z = pt.z;
488 template <
class IteratorT>
492 centroid.x = centroid.y = centroid.z = 0;
495 centroid /= (float) (end - begin);
499 template <
class IteratorT>
502 cov.
data[0] = make_float3 (0.0f, 0.0f, 0.0f);
503 cov.
data[1] = make_float3 (0.0f, 0.0f, 0.0f);
504 cov.
data[2] = make_float3 (0.0f, 0.0f, 0.0f);
506 cov = transform_reduce (begin, end,
517 cov.
data[0] /= (float) (end - begin);
518 cov.
data[1] /= (float) (end - begin);
519 cov.
data[2] /= (float) (end - begin);
523 template <
typename CloudPtr>
528 :
points_(thrust::raw_pointer_cast (&input->points[0]))
538 inline __host__ __device__
543 float r_quadr, z_sqr;
544 float sqrt_term_y, sqrt_term_x, norm;
545 float x_times_z, y_times_z;
551 z_sqr = point_arg.z * point_arg.z;
559 x_times_z = point_arg.x * point_arg.z;
560 y_times_z = point_arg.y * point_arg.z;
563 bounds.x = (x_times_z - sqrt_term_x) * norm;
564 bounds.y = (x_times_z + sqrt_term_x) * norm;
565 bounds.z = (y_times_z - sqrt_term_y) * norm;
566 bounds.w = (y_times_z + sqrt_term_y) * norm;
570 bounds.x +=
width_ / 2.0f;
571 bounds.y +=
width_ / 2.0f;
575 res.x = (int)std::floor (bounds.x);
576 res.y = (int)std::ceil (bounds.y);
577 res.z = (int)std::floor (bounds.z);
578 res.w = (int)std::ceil (bounds.w);
581 res.x = clamp (res.x, 0,
width_-1);
582 res.y = clamp (res.y, 0,
width_-1);
583 res.z = clamp (res.z, 0,
height_-1);
584 res.w = clamp (res.w, 0,
height_-1);
589 inline __host__ __device__
590 int radiusSearch (
const float3 &query_pt,
int k_indices[],
int max_nnn)
597 for (
int x = bounds.x; (x <= bounds.y) & (nnn < max_nnn); ++x)
599 for (
int y = bounds.z; (y <= bounds.w) & (nnn < max_nnn); ++y)
606 float3 point_dif =
points_[idx] - query_pt;
611 k_indices[nnn] = idx;
621 inline __host__ __device__
641 bounds.x = clamp (bounds.x, 0,
width_-1);
642 bounds.y = clamp (bounds.y, 0,
width_-1);
643 bounds.z = clamp (bounds.z, 0,
height_-1);
644 bounds.w = clamp (bounds.w, 0,
height_-1);
650 float skipX = max (sqrtf ((
float)bounds.y-bounds.x) / sqrt_desired_nr_neighbors, 1.0f);
651 float skipY = max (sqrtf ((
float)bounds.w-bounds.z) / sqrt_desired_nr_neighbors, 1.0f);
655 cov.
data[0] = make_float3(0,0,0);
656 cov.
data[1] = make_float3(0,0,0);
657 cov.
data[2] = make_float3(0,0,0);
658 float3 centroid = make_float3(0,0,0);
661 for (
float y = (
float) bounds.z; y <= bounds.w; y += skipY)
663 for (
float x = (
float) bounds.x; x <= bounds.y; x += skipX)
666 int idx = ((int)y) *
width_ + ((int)x);
672 float3 point_dif =
points_[idx] - query_pt;
678 float3 demean_old =
points_[idx] - centroid;
679 centroid += demean_old / (float) nnn;
680 float3 demean_new =
points_[idx] - centroid;
682 cov.
data[1].y += demean_new.y * demean_old.y;
683 cov.
data[1].z += demean_new.y * demean_old.z;
684 cov.
data[2].z += demean_new.z * demean_old.z;
686 demean_old *= demean_new.x;
687 cov.
data[0].x += demean_old.x;
688 cov.
data[0].y += demean_old.y;
689 cov.
data[0].z += demean_old.z;
697 cov.
data[0] /= (float) nnn;
698 cov.
data[1] /= (float) nnn;
699 cov.
data[2] /= (float) nnn;
704 inline __host__ __device__
724 bounds.x = clamp (bounds.x, 0,
width_-1);
725 bounds.y = clamp (bounds.y, 0,
width_-1);
726 bounds.z = clamp (bounds.z, 0,
height_-1);
727 bounds.w = clamp (bounds.w, 0,
height_-1);
732 float skipX = max (sqrtf ((
float)bounds.y-bounds.x) / sqrt_desired_nr_neighbors, 1.0f);
733 float skipY = max (sqrtf ((
float)bounds.w-bounds.z) / sqrt_desired_nr_neighbors, 1.0f);
737 float3 centroid = make_float3(0,0,0);
740 for (
float y = (
float) bounds.z; y <= bounds.w; y += skipY)
742 for (
float x = (
float) bounds.x; x <= bounds.y; x += skipX)
745 int idx = ((int)y) *
width_ + ((int)x);
751 float3 point_dif =
points_[idx] - query_pt;
762 return centroid / (float) nnn;
Kernel to compute a radius neighborhood given a organized point cloud (aka range image cloud)
__host__ __device__ int radiusSearch(const float3 &query_pt, int k_indices[], int max_nnn)
__host__ __device__ float3 computeCentroid(const float3 &query_pt, CovarianceMatrix &cov, float sqrt_desired_nr_neighbors)
OrganizedRadiusSearch(const CloudPtr &input, float focalLength, float sqr_radius)
__host__ __device__ int computeCovarianceOnline(const float3 &query_pt, CovarianceMatrix &cov, float sqrt_desired_nr_neighbors)
__host__ __device__ int4 getProjectedRadiusSearchBox(const float3 &point_arg)
const PointXYZRGB * points_
__host__ __device__ float3 unitOrthogonal(const float3 &src)
void computeCovariance(IteratorT begin, IteratorT end, CovarianceMatrix &cov, float3 centroid)
Computes a covariance matrix for a given range of points.
__host__ __device__ void eigen33(const CovarianceMatrix &mat, CovarianceMatrix &evecs, float3 &evals)
__host__ __device__ void computeRoots(const CovarianceMatrix &m, float3 &roots)
__host__ __device__ void swap(float &a, float &b)
__host__ __device__ bool isMuchSmallerThan(float x, float y)
__host__ __device__ void computeRoots2(const float &b, const float &c, float3 &roots)
void compute3DCentroid(IteratorT begin, IteratorT end, float3 ¢roid)
Computes a centroid for a given range of points.
Adds two matrices element-wise.
__inline__ __host__ __device__ CovarianceMatrix operator()(CovarianceMatrix lhs, CovarianceMatrix rhs)
Simple kernel to add two points.
__inline__ __host__ __device__ float3 operator()(float3 lhs, float3 rhs)
Kernel to compute a `‘covariance matrix’' for a single point.
__inline__ __host__ __device__ CovarianceMatrix operator()(const PointXYZRGB &point)
__inline__ __host__ __device__ ComputeCovarianceForPoint(const float3 ¢roid)
misnamed class holding a 3x3 matrix
Default point xyz-rgb structure.
Simple kernel to convert a PointXYZRGB to float3.
__inline__ __host__ __device__ float3 operator()(const PointXYZRGB &pt)