40 #ifndef PCL_SURFACE_IMPL_MLS_H_
41 #define PCL_SURFACE_IMPL_MLS_H_
43 #include <pcl/type_traits.h>
44 #include <pcl/surface/mls.h>
45 #include <pcl/common/io.h>
46 #include <pcl/common/copy_point.h>
47 #include <pcl/common/centroid.h>
48 #include <pcl/common/eigen.h>
49 #include <pcl/common/geometry.h>
56 template <
typename Po
intInT,
typename Po
intOutT>
void
67 normals_->header = input_->header;
69 normals_->width = normals_->height = 0;
70 normals_->points.clear ();
74 output.
header = input_->header;
78 if (search_radius_ <= 0 || sqr_gauss_param_ <= 0)
80 PCL_ERROR (
"[pcl::%s::process] Invalid search radius (%f) or Gaussian parameter (%f)!\n", getClassName ().c_str (), search_radius_, sqr_gauss_param_);
85 if (upsample_method_ == DISTINCT_CLOUD && !distinct_cloud_)
87 PCL_ERROR (
"[pcl::%s::process] Upsample method was set to DISTINCT_CLOUD, but no distinct cloud was specified.\n", getClassName ().c_str ());
98 if (input_->isOrganized ())
102 setSearchMethod (tree);
106 tree_->setInputCloud (input_);
108 switch (upsample_method_)
111 case (RANDOM_UNIFORM_DENSITY):
113 std::random_device rd;
115 const double tmp = search_radius_ / 2.0;
116 rng_uniform_distribution_.reset (
new std::uniform_real_distribution<> (-tmp, tmp));
120 case (VOXEL_GRID_DILATION):
121 case (DISTINCT_CLOUD):
123 if (!cache_mls_results_)
124 PCL_WARN (
"The cache mls results is forced when using upsampling method VOXEL_GRID_DILATION or DISTINCT_CLOUD.\n");
126 cache_mls_results_ =
true;
133 if (cache_mls_results_)
135 mls_results_.resize (input_->size ());
139 mls_results_.resize (1);
143 performProcessing (output);
145 if (compute_normals_)
147 normals_->height = 1;
148 normals_->width =
static_cast<std::uint32_t
> (normals_->size ());
150 for (std::size_t i = 0; i < output.
size (); ++i)
163 output.
width =
static_cast<std::uint32_t
> (output.
size ());
169 template <
typename Po
intInT,
typename Po
intOutT>
void
171 const std::vector<int> &nn_indices,
180 mls_result.
computeMLSSurface<PointInT> (*input_, index, nn_indices, search_radius_, order_);
182 switch (upsample_method_)
187 addProjectedPointNormal (index, proj.
point, proj.
normal, mls_result.
curvature, projected_points, projected_points_normals, corresponding_input_indices);
191 case (SAMPLE_LOCAL_PLANE):
194 for (
float u_disp = -static_cast<float> (upsampling_radius_); u_disp <= upsampling_radius_; u_disp += static_cast<float> (upsampling_step_))
195 for (
float v_disp = -static_cast<float> (upsampling_radius_); v_disp <= upsampling_radius_; v_disp += static_cast<float> (upsampling_step_))
196 if (u_disp * u_disp + v_disp * v_disp < upsampling_radius_ * upsampling_radius_)
199 addProjectedPointNormal (index, proj.
point, proj.
normal, mls_result.
curvature, projected_points, projected_points_normals, corresponding_input_indices);
204 case (RANDOM_UNIFORM_DENSITY):
207 const int num_points_to_add =
static_cast<int> (std::floor (desired_num_points_in_radius_ / 2.0 / static_cast<double> (nn_indices.size ())));
210 if (num_points_to_add <= 0)
214 addProjectedPointNormal (index, proj.
point, proj.
normal, mls_result.
curvature, projected_points, projected_points_normals, corresponding_input_indices);
219 for (
int num_added = 0; num_added < num_points_to_add;)
221 const double u = (*rng_uniform_distribution_) (rng_);
222 const double v = (*rng_uniform_distribution_) (rng_);
225 if (u * u + v * v > search_radius_ * search_radius_ / 4)
234 addProjectedPointNormal (index, proj.
point, proj.
normal, mls_result.
curvature, projected_points, projected_points_normals, corresponding_input_indices);
247 template <
typename Po
intInT,
typename Po
intOutT>
void
249 const Eigen::Vector3d &point,
250 const Eigen::Vector3d &normal,
257 aux.x =
static_cast<float> (point[0]);
258 aux.y =
static_cast<float> (point[1]);
259 aux.z =
static_cast<float> (point[2]);
262 copyMissingFields (input_->points[index], aux);
265 corresponding_input_indices.
indices.push_back (index);
267 if (compute_normals_)
270 aux_normal.normal_x =
static_cast<float> (normal[0]);
271 aux_normal.normal_y =
static_cast<float> (normal[1]);
272 aux_normal.normal_z =
static_cast<float> (normal[2]);
274 projected_points_normals.
push_back (aux_normal);
279 template <
typename Po
intInT,
typename Po
intOutT>
void
283 nr_coeff_ = (order_ + 1) * (order_ + 2) / 2;
287 const unsigned int threads = threads_ == 0 ? 1 : threads_;
291 std::vector<PointIndices> corresponding_input_indices (threads);
295 #pragma omp parallel for \
297 shared(corresponding_input_indices, projected_points, projected_points_normals) \
298 schedule(dynamic,1000) \
300 for (
int cp = 0; cp < static_cast<int> (indices_->size ()); ++cp)
304 std::vector<int> nn_indices;
305 std::vector<float> nn_sqr_dists;
308 if (searchForNeighbors ((*indices_)[cp], nn_indices, nn_sqr_dists))
311 if (nn_indices.size () >= 3)
315 const int tn = omp_get_thread_num ();
317 std::size_t pp_size = projected_points[tn].size ();
324 const int index = (*indices_)[cp];
326 std::size_t mls_result_index = 0;
327 if (cache_mls_results_)
328 mls_result_index = index;
331 computeMLSPointNormal (index, nn_indices, projected_points[tn], projected_points_normals[tn], corresponding_input_indices[tn], mls_results_[mls_result_index]);
334 for (std::size_t pp = pp_size; pp < projected_points[tn].
size (); ++pp)
335 copyMissingFields (input_->points[(*indices_)[cp]], projected_points[tn][pp]);
337 computeMLSPointNormal (index, nn_indices, projected_points, projected_points_normals, *corresponding_input_indices_, mls_results_[mls_result_index]);
340 output.
insert (output.
end (), projected_points.
begin (), projected_points.
end ());
341 if (compute_normals_)
342 normals_->insert (normals_->end (), projected_points_normals.
begin (), projected_points_normals.
end ());
350 for (
unsigned int tn = 0; tn < threads; ++tn)
352 output.
insert (output.
end (), projected_points[tn].begin (), projected_points[tn].end ());
353 corresponding_input_indices_->indices.insert (corresponding_input_indices_->indices.end (),
354 corresponding_input_indices[tn].indices.begin (), corresponding_input_indices[tn].indices.end ());
355 if (compute_normals_)
356 normals_->insert (normals_->end (), projected_points_normals[tn].begin (), projected_points_normals[tn].end ());
361 performUpsampling (output);
365 template <
typename Po
intInT,
typename Po
intOutT>
void
369 if (upsample_method_ == DISTINCT_CLOUD)
372 for (std::size_t dp_i = 0; dp_i < distinct_cloud_->size (); ++dp_i)
375 if (!std::isfinite (distinct_cloud_->points[dp_i].x))
380 std::vector<int> nn_indices;
381 std::vector<float> nn_dists;
382 tree_->nearestKSearch (distinct_cloud_->points[dp_i], 1, nn_indices, nn_dists);
383 int input_index = nn_indices.front ();
387 if (mls_results_[input_index].valid ==
false)
390 Eigen::Vector3d add_point = distinct_cloud_->points[dp_i].getVector3fMap ().template cast<double> ();
392 addProjectedPointNormal (input_index, proj.point, proj.normal, mls_results_[input_index].curvature, output, *normals_, *corresponding_input_indices_);
398 if (upsample_method_ == VOXEL_GRID_DILATION)
402 MLSVoxelGrid voxel_grid (input_, indices_, voxel_size_);
403 for (
int iteration = 0; iteration < dilation_iteration_num_; ++iteration)
406 for (
typename MLSVoxelGrid::HashMap::iterator m_it = voxel_grid.
voxel_grid_.begin (); m_it != voxel_grid.
voxel_grid_.end (); ++m_it)
417 std::vector<int> nn_indices;
418 std::vector<float> nn_dists;
419 tree_->nearestKSearch (p, 1, nn_indices, nn_dists);
420 int input_index = nn_indices.front ();
424 if (mls_results_[input_index].valid ==
false)
427 Eigen::Vector3d add_point = p.getVector3fMap ().template cast<double> ();
429 addProjectedPointNormal (input_index, proj.point, proj.normal, mls_results_[input_index].curvature, output, *normals_, *corresponding_input_indices_);
436 const Eigen::Vector3d &a_mean,
437 const Eigen::Vector3d &a_plane_normal,
438 const Eigen::Vector3d &a_u,
439 const Eigen::Vector3d &a_v,
440 const Eigen::VectorXd &a_c_vec,
441 const int a_num_neighbors,
442 const float a_curvature,
444 query_point (a_query_point), mean (a_mean), plane_normal (a_plane_normal), u_axis (a_u), v_axis (a_v), c_vec (a_c_vec), num_neighbors (a_num_neighbors),
445 curvature (a_curvature), order (a_order), valid (true)
451 Eigen::Vector3d delta = pt - mean;
452 u = delta.dot (u_axis);
453 v = delta.dot (v_axis);
454 w = delta.dot (plane_normal);
460 Eigen::Vector3d delta = pt - mean;
461 u = delta.dot (u_axis);
462 v = delta.dot (v_axis);
473 for (
int ui = 0; ui <= order; ++ui)
476 for (
int vi = 0; vi <= order - ui; ++vi)
478 result += c_vec[j++] * u_pow * v_pow;
493 Eigen::VectorXd u_pow (order + 2), v_pow (order + 2);
496 d.z = d.z_u = d.z_v = d.z_uu = d.z_vv = d.z_uv = 0;
497 u_pow (0) = v_pow (0) = 1;
498 for (
int ui = 0; ui <= order; ++ui)
500 for (
int vi = 0; vi <= order - ui; ++vi)
503 d.z += u_pow (ui) * v_pow (vi) * c_vec[j];
507 d.z_u += c_vec[j] * ui * u_pow (ui - 1) * v_pow (vi);
510 d.z_v += c_vec[j] * vi * u_pow (ui) * v_pow (vi - 1);
512 if (ui >= 1 && vi >= 1)
513 d.z_uv += c_vec[j] * ui * u_pow (ui - 1) * vi * v_pow (vi - 1);
516 d.z_uu += c_vec[j] * ui * (ui - 1) * u_pow (ui - 2) * v_pow (vi);
519 d.z_vv += c_vec[j] * vi * (vi - 1) * u_pow (ui) * v_pow (vi - 2);
522 v_pow (vi + 1) = v_pow (vi) * v;
526 u_pow (ui + 1) = u_pow (ui) * u;
535 Eigen::Vector2f k (1e-5, 1e-5);
541 if (order > 1 && c_vec.size () >= (order + 1) * (order + 2) / 2 && std::isfinite (c_vec[0]))
545 const double Zlen = std::sqrt (Z);
548 const double disc2 = H * H -
K;
549 assert (disc2 >= 0.0);
550 const double disc = std::sqrt (disc2);
554 if (std::abs (k[0]) > std::abs (k[1])) std::swap (k[0], k[1]);
558 PCL_ERROR (
"No Polynomial fit data, unable to calculate the principle curvatures!\n");
572 result.
normal = plane_normal;
573 if (order > 1 && c_vec.size () >= (order + 1) * (order + 2) / 2 && std::isfinite (c_vec[0]))
578 const double dist1 = std::abs (gw - w);
582 double e1 = (gu - u) + d.
z_u * gw - d.
z_u * w;
583 double e2 = (gv - v) + d.
z_v * gw - d.
z_v * w;
591 Eigen::MatrixXd J (2, 2);
597 Eigen::Vector2d err (e1, e2);
598 Eigen::Vector2d update = J.inverse () * err;
602 d = getPolynomialPartialDerivative (gu, gv);
604 dist2 = std::sqrt ((gu - u) * (gu - u) + (gv - v) * (gv - v) + (gw - w) * (gw - w));
606 err_total = std::sqrt (e1 * e1 + e2 * e2);
608 }
while (err_total > 1e-8 && dist2 < dist1);
614 d = getPolynomialPartialDerivative (u, v);
621 result.
normal.normalize ();
624 result.
point = mean + gu * u_axis + gv * v_axis + gw * plane_normal;
635 result.
normal = plane_normal;
636 result.
point = mean + u * u_axis + v * v_axis;
649 result.
normal = plane_normal;
651 if (order > 1 && c_vec.size () >= (order + 1) * (order + 2) / 2 && std::isfinite (c_vec[0]))
656 result.
normal.normalize ();
659 result.
point = mean + u * u_axis + v * v_axis + w * plane_normal;
668 getMLSCoordinates (pt, u, v, w);
671 if (order > 1 && num_neighbors >= required_neighbors && std::isfinite (c_vec[0]) && method != NONE)
673 if (method == ORTHOGONAL)
674 proj = projectPointOrthogonalToPolynomialSurface (u, v, w);
676 proj = projectPointSimpleToPolynomialSurface (u, v);
680 proj = projectPointToMLSPlane (u, v);
690 if (order > 1 && num_neighbors >= required_neighbors && std::isfinite (c_vec[0]) && method != NONE)
692 if (method == ORTHOGONAL)
695 getMLSCoordinates (query_point, u, v, w);
696 proj = projectPointOrthogonalToPolynomialSurface (u, v, w);
701 proj.
point = mean + (c_vec[0] * plane_normal);
704 proj.
normal = plane_normal - c_vec[order + 1] * u_axis - c_vec[1] * v_axis;
710 proj.
normal = plane_normal;
717 template <
typename Po
intT>
void
720 const std::vector<int> &nn_indices,
721 double search_radius,
722 int polynomial_order,
723 std::function<
double(
const double)> weight_func)
726 EIGEN_ALIGN16 Eigen::Matrix3d covariance_matrix;
727 Eigen::Vector4d xyz_centroid;
734 EIGEN_ALIGN16 Eigen::Vector3d::Scalar eigen_value;
735 EIGEN_ALIGN16 Eigen::Vector3d eigen_vector;
736 Eigen::Vector4d model_coefficients (0, 0, 0, 0);
737 pcl::eigen33 (covariance_matrix, eigen_value, eigen_vector);
738 model_coefficients.head<3> ().matrix () = eigen_vector;
739 model_coefficients[3] = -1 * model_coefficients.dot (xyz_centroid);
741 query_point = cloud.
points[index].getVector3fMap ().template cast<double> ();
743 if (!std::isfinite(eigen_vector[0]) || !std::isfinite(eigen_vector[1]) || !std::isfinite(eigen_vector[2]))
754 const double distance = query_point.dot (model_coefficients.head<3> ()) + model_coefficients[3];
755 mean = query_point - distance * model_coefficients.head<3> ();
757 curvature = covariance_matrix.trace ();
760 curvature = std::abs (eigen_value / curvature);
763 plane_normal = model_coefficients.head<3> ();
766 v_axis = plane_normal.unitOrthogonal ();
767 u_axis = plane_normal.cross (v_axis);
771 num_neighbors =
static_cast<int> (nn_indices.size ());
772 order = polynomial_order;
775 const int nr_coeff = (order + 1) * (order + 2) / 2;
777 if (num_neighbors >= nr_coeff)
780 weight_func = [=] (
const double sq_dist) {
return this->computeMLSWeight (sq_dist, search_radius * search_radius); };
783 Eigen::VectorXd weight_vec (num_neighbors);
784 Eigen::MatrixXd P (nr_coeff, num_neighbors);
785 Eigen::VectorXd f_vec (num_neighbors);
786 Eigen::MatrixXd P_weight_Pt (nr_coeff, nr_coeff);
790 std::vector<Eigen::Vector3d, Eigen::aligned_allocator<Eigen::Vector3d> > de_meaned (num_neighbors);
791 for (std::size_t ni = 0; ni < static_cast<std::size_t>(num_neighbors); ++ni)
793 de_meaned[ni][0] = cloud.
points[nn_indices[ni]].x - mean[0];
794 de_meaned[ni][1] = cloud.
points[nn_indices[ni]].y - mean[1];
795 de_meaned[ni][2] = cloud.
points[nn_indices[ni]].z - mean[2];
796 weight_vec (ni) = weight_func (de_meaned[ni].dot (de_meaned[ni]));
801 for (std::size_t ni = 0; ni < static_cast<std::size_t>(num_neighbors); ++ni)
804 const double u_coord = de_meaned[ni].dot(u_axis);
805 const double v_coord = de_meaned[ni].dot(v_axis);
806 f_vec (ni) = de_meaned[ni].dot (plane_normal);
811 for (
int ui = 0; ui <= order; ++ui)
814 for (
int vi = 0; vi <= order - ui; ++vi)
816 P (j++, ni) = u_pow * v_pow;
824 const Eigen::MatrixXd P_weight = P * weight_vec.asDiagonal();
825 P_weight_Pt = P_weight * P.transpose ();
826 c_vec = P_weight * f_vec;
827 P_weight_Pt.llt ().solveInPlace (c_vec);
833 template <
typename Po
intInT,
typename Po
intOutT>
837 voxel_grid_ (), data_size_ (), voxel_size_ (voxel_size)
842 const double max_size = (std::max) ((std::max)(bounding_box_size.x (), bounding_box_size.y ()), bounding_box_size.z ());
845 for (std::size_t i = 0; i < indices->size (); ++i)
846 if (std::isfinite (cloud->points[(*indices)[i]].x))
849 getCellIndex (cloud->points[(*indices)[i]].getVector3fMap (), pos);
851 std::uint64_t index_1d;
859 template <
typename Po
intInT,
typename Po
intOutT>
void
862 HashMap new_voxel_grid = voxel_grid_;
863 for (
typename MLSVoxelGrid::HashMap::iterator m_it = voxel_grid_.begin (); m_it != voxel_grid_.end (); ++m_it)
865 Eigen::Vector3i index;
866 getIndexIn3D (m_it->first, index);
869 for (
int x = -1; x <= 1; ++x)
870 for (
int y = -1; y <= 1; ++y)
871 for (
int z = -1; z <= 1; ++z)
872 if (x != 0 || y != 0 || z != 0)
874 Eigen::Vector3i new_index;
875 new_index = index + Eigen::Vector3i (x, y, z);
877 std::uint64_t index_1d;
878 getIndexIn1D (new_index, index_1d);
880 new_voxel_grid[index_1d] = leaf;
883 voxel_grid_ = new_voxel_grid;
888 template <
typename Po
intInT,
typename Po
intOutT>
void
890 PointOutT &point_out)
const
892 PointOutT temp = point_out;
894 point_out.x = temp.x;
895 point_out.y = temp.y;
896 point_out.z = temp.z;
899 #define PCL_INSTANTIATE_MovingLeastSquares(T,OutT) template class PCL_EXPORTS pcl::MovingLeastSquares<T,OutT>;
900 #define PCL_INSTANTIATE_MovingLeastSquaresOMP(T,OutT) template class PCL_EXPORTS pcl::MovingLeastSquaresOMP<T,OutT>;
902 #endif // PCL_SURFACE_IMPL_MLS_H_
A point structure representing normal coordinates and the surface curvature estimate.
Data structure used to store the MLS polynomial partial derivatives.
double z_u
The partial derivative dz/du.
typename PointCloudIn::ConstPtr PointCloudInConstPtr
search::KdTree is a wrapper class which inherits the pcl::KdTree class for performing search function...
A helper functor that can set a specific value in a field if the field exists.
std::vector< PointT, Eigen::aligned_allocator< PointT > > points
The point data.
double u
The u-coordinate of the projected point in local MLS frame.
MLSProjectionResults projectPointOrthogonalToPolynomialSurface(const double u, const double v, const double w) const
Project a point orthogonal to the polynomial surface.
MLSProjectionResults projectPointToMLSPlane(const double u, const double v) const
Project a point onto the MLS plane.
shared_ptr< Indices > IndicesPtr
double z_vv
The partial derivative d^2z/dv^2.
Eigen::Vector4f bounding_max_
void push_back(const PointT &pt)
Insert a new point in the cloud, at the end of the container.
double z_uv
The partial derivative d^2z/dudv.
void process(PointCloudOut &output) override
Base method for surface reconstruction for all points given in ...
std::map< std::uint64_t, Leaf > HashMap
Eigen::Vector4f bounding_min_
double z_uu
The partial derivative d^2z/du^2.
std::uint32_t width
The point cloud width (if organized as an image-structure).
double getPolynomialValue(const double u, const double v) const
Calculate the polynomial.
Data structure used to store the MLS projection results.
Eigen::Vector3d point
The projected point.
double v
The u-coordinate of the projected point in local MLS frame.
double z_v
The partial derivative dz/dv.
Eigen::Vector3d normal
The projected point's normal.
void getIndexIn1D(const Eigen::Vector3i &index, std::uint64_t &index_1d) const
MLSProjectionResults projectPoint(const Eigen::Vector3d &pt, ProjectionMethod method, int required_neighbors=0) const
Project a point using the specified method.
PolynomialPartialDerivative getPolynomialPartialDerivative(const double u, const double v) const
Calculate the polynomial's first and second partial derivatives.
unsigned int computeCovarianceMatrix(const pcl::PointCloud< PointT > &cloud, const Eigen::Matrix< Scalar, 4, 1 > ¢roid, Eigen::Matrix< Scalar, 3, 3 > &covariance_matrix)
Compute the 3x3 covariance matrix of a given set of points.
Data structure used to store the results of the MLS fitting.
int num_neighbors
The number of neighbors used to create the mls surface.
void performUpsampling(PointCloudOut &output)
Perform upsampling for the distinct-cloud and voxel-grid methods.
std::vector< PointCloud< PointOutT >, Eigen::aligned_allocator< PointCloud< PointOutT > > > CloudVectorType
void getMinMax3D(const pcl::PointCloud< PointT > &cloud, PointT &min_pt, PointT &max_pt)
Get the minimum and maximum values on each of the 3 (x-y-z) dimensions in a given pointcloud...
std::uint32_t height
The point cloud height (if organized as an image-structure).
pcl::PCLHeader header
The point cloud header.
void eigen33(const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector)
determines the eigenvector and eigenvalue of the smallest eigenvalue of the symmetric positive semi d...
void performProcessing(PointCloudOut &output) override
Abstract surface reconstruction method.
A minimalistic implementation of a voxel grid, necessary for the point cloud upsampling.
MLSVoxelGrid(PointCloudInConstPtr &cloud, IndicesPtr &indices, float voxel_size)
void getMLSCoordinates(const Eigen::Vector3d &pt, double &u, double &v, double &w) const
Given a point calculate it's 3D location in the MLS frame.
void computeMLSSurface(const pcl::PointCloud< PointT > &cloud, int index, const std::vector< int > &nn_indices, double search_radius, int polynomial_order=2, std::function< double(const double)> weight_func={})
Smooth a given point and its neighborghood using Moving Least Squares.
OrganizedNeighbor is a class for optimized nearest neigbhor search in organized point clouds...
float curvature
The curvature at the query point.
MLSProjectionResults projectPointSimpleToPolynomialSurface(const double u, const double v) const
Project a point along the MLS plane normal to the polynomial surface.
MLSProjectionResults projectQueryPoint(ProjectionMethod method, int required_neighbors=0) const
Project the query point used to generate the mls surface about using the specified method...
double z
The z component of the polynomial evaluated at z(u, v).
unsigned int compute3DCentroid(ConstCloudIterator< PointT > &cloud_iterator, Eigen::Matrix< Scalar, 4, 1 > ¢roid)
Compute the 3D (X-Y-Z) centroid of a set of points and return it as a 3D vector.
void getCellIndex(const Eigen::Vector3f &p, Eigen::Vector3i &index) const
iterator insert(iterator position, const PointT &pt)
Insert a new point in the cloud, given an iterator.
void copyMissingFields(const PointInT &point_in, PointOutT &point_out) const
void addProjectedPointNormal(int index, const Eigen::Vector3d &point, const Eigen::Vector3d &normal, double curvature, PointCloudOut &projected_points, NormalCloud &projected_points_normals, PointIndices &corresponding_input_indices) const
This is a helper function for add projected points.
void copyPoint(const PointInT &point_in, PointOutT &point_out)
Copy the fields of a source point into a target point.
void getPosition(const std::uint64_t &index_1d, Eigen::Vector3f &point) const
typename KdTree::Ptr KdTreePtr
void computeMLSPointNormal(int index, const std::vector< int > &nn_indices, PointCloudOut &projected_points, NormalCloud &projected_points_normals, PointIndices &corresponding_input_indices, MLSResult &mls_result) const
Smooth a given point and its neighborghood using Moving Least Squares.
Eigen::Vector2f calculatePrincipleCurvatures(const double u, const double v) const
Calculate the principle curvatures using the polynomial surface.