40#ifndef PCL_SURFACE_IMPL_MLS_H_
41#define PCL_SURFACE_IMPL_MLS_H_
45#include <pcl/common/copy_point.h>
46#include <pcl/common/eigen.h>
47#include <pcl/search/kdtree.h>
48#include <pcl/search/organized.h>
49#include <pcl/surface/mls.h>
50#include <pcl/type_traits.h>
52#include <Eigen/Geometry>
62template <
typename Po
intInT,
typename Po
intOutT>
void
86 PCL_ERROR (
"[pcl::%s::process] Invalid search radius (%f) or Gaussian parameter (%f)!\n", getClassName ().c_str (),
search_radius_,
sqr_gauss_param_);
93 PCL_ERROR (
"[pcl::%s::process] Upsample method was set to DISTINCT_CLOUD, but no distinct cloud was specified.\n", getClassName ().c_str ());
104 if (
input_->isOrganized ())
119 std::random_device rd;
122 rng_uniform_distribution_ = std::make_unique<std::uniform_real_distribution<>> (-tmp, tmp);
130 PCL_WARN (
"The cache mls results is forced when using upsampling method VOXEL_GRID_DILATION or DISTINCT_CLOUD.\n");
156 for (std::size_t i = 0; i < output.
size (); ++i)
158 using FieldList =
typename pcl::traits::fieldList<PointOutT>::type;
175template <
typename Po
intInT,
typename Po
intOutT>
void
213 const int num_points_to_add =
static_cast<int> (std::floor (
desired_num_points_in_radius_ / 2.0 /
static_cast<double> (nn_indices.size ())));
216 if (num_points_to_add <= 0)
225 for (
int num_added = 0; num_added < num_points_to_add;)
227 const double u = (*rng_uniform_distribution_) (rng_);
228 const double v = (*rng_uniform_distribution_) (rng_);
253template <
typename Po
intInT,
typename Po
intOutT>
void
255 const Eigen::Vector3d &point,
256 const Eigen::Vector3d &normal,
263 aux.x =
static_cast<float> (point[0]);
264 aux.y =
static_cast<float> (point[1]);
265 aux.z =
static_cast<float> (point[2]);
271 corresponding_input_indices.
indices.push_back (index);
276 aux_normal.normal_x =
static_cast<float> (normal[0]);
277 aux_normal.normal_y =
static_cast<float> (normal[1]);
278 aux_normal.normal_z =
static_cast<float> (normal[2]);
280 projected_points_normals.
push_back (aux_normal);
285template <
typename Po
intInT,
typename Po
intOutT>
void
297 std::vector<PointIndices> corresponding_input_indices (threads);
301#pragma omp parallel for \
303 shared(corresponding_input_indices, projected_points, projected_points_normals) \
304 schedule(dynamic,1000) \
306 for (
int cp = 0; cp < static_cast<int> (
indices_->size ()); ++cp)
311 std::vector<float> nn_sqr_dists;
317 if (nn_indices.size () >= 3)
321 const int tn = omp_get_thread_num ();
323 std::size_t pp_size = projected_points[tn].size ();
330 const int index = (*indices_)[cp];
332 std::size_t mls_result_index = 0;
334 mls_result_index = index;
340 for (std::size_t pp = pp_size; pp < projected_points[tn].
size (); ++pp)
346 output.
insert (output.
end (), projected_points.
begin (), projected_points.
end ());
356 for (
unsigned int tn = 0; tn < threads; ++tn)
358 output.
insert (output.
end (), projected_points[tn].begin (), projected_points[tn].end ());
360 corresponding_input_indices[tn].indices.begin (), corresponding_input_indices[tn].indices.end ());
362 normals_->insert (
normals_->end (), projected_points_normals[tn].begin (), projected_points_normals[tn].end ());
371template <
typename Po
intInT,
typename Po
intOutT>
void
387 std::vector<float> nn_dists;
389 const auto input_index = nn_indices.front ();
396 Eigen::Vector3d add_point = (*distinct_cloud_)[dp_i].getVector3fMap ().template cast<double> ();
424 std::vector<float> nn_dists (1);
425 tree_->nearestKSearch (p, 1, nn_indices, nn_dists);
426 const auto input_index = nn_indices.front ();
433 Eigen::Vector3d add_point = p.getVector3fMap ().template cast<double> ();
442 const Eigen::Vector3d &a_mean,
443 const Eigen::Vector3d &a_plane_normal,
444 const Eigen::Vector3d &a_u,
445 const Eigen::Vector3d &a_v,
446 const Eigen::VectorXd &a_c_vec,
447 const int a_num_neighbors,
448 const float a_curvature,
457 Eigen::Vector3d delta = pt -
mean;
466 Eigen::Vector3d delta = pt -
mean;
479 for (
int ui = 0; ui <=
order; ++ui)
482 for (
int vi = 0; vi <=
order - ui; ++vi)
484 result +=
c_vec[j++] * u_pow * v_pow;
499 Eigen::VectorXd u_pow (
order + 2), v_pow (
order + 2);
503 u_pow (0) = v_pow (0) = 1;
504 for (
int ui = 0; ui <=
order; ++ui)
506 for (
int vi = 0; vi <=
order - ui; ++vi)
509 d.
z += u_pow (ui) * v_pow (vi) *
c_vec[j];
513 d.
z_u +=
c_vec[j] * ui * u_pow (ui - 1) * v_pow (vi);
516 d.
z_v +=
c_vec[j] * vi * u_pow (ui) * v_pow (vi - 1);
518 if (ui >= 1 && vi >= 1)
519 d.
z_uv +=
c_vec[j] * ui * u_pow (ui - 1) * vi * v_pow (vi - 1);
522 d.
z_uu +=
c_vec[j] * ui * (ui - 1) * u_pow (ui - 2) * v_pow (vi);
525 d.
z_vv +=
c_vec[j] * vi * (vi - 1) * u_pow (ui) * v_pow (vi - 2);
528 v_pow (vi + 1) = v_pow (vi) * v;
532 u_pow (ui + 1) = u_pow (ui) * u;
552 const double dist1 = std::abs (gw - w);
556 double e1 = (gu - u) + d.
z_u * gw - d.
z_u * w;
557 double e2 = (gv - v) + d.
z_v * gw - d.
z_v * w;
565 Eigen::MatrixXd J (2, 2);
571 Eigen::Vector2d err (e1, e2);
572 Eigen::Vector2d update = J.inverse () * err;
578 dist2 = std::sqrt ((gu - u) * (gu - u) + (gv - v) * (gv - v) + (gw - w) * (gw - w));
580 err_total = std::sqrt (e1 * e1 + e2 * e2);
582 }
while (err_total > 1e-8 && dist2 < dist1);
595 result.
normal.normalize ();
630 result.
normal.normalize ();
691template <
typename Po
intT>
void
695 double search_radius,
696 int polynomial_order,
697 std::function<
double(
const double)> weight_func)
700 EIGEN_ALIGN16 Eigen::Matrix3d covariance_matrix;
701 Eigen::Vector4d xyz_centroid;
708 EIGEN_ALIGN16 Eigen::Vector3d::Scalar eigen_value;
709 EIGEN_ALIGN16 Eigen::Vector3d eigen_vector;
710 Eigen::Vector4d model_coefficients (0, 0, 0, 0);
711 pcl::eigen33 (covariance_matrix, eigen_value, eigen_vector);
712 model_coefficients.head<3> ().matrix () = eigen_vector;
713 model_coefficients[3] = -1 * model_coefficients.dot (xyz_centroid);
715 query_point = cloud[index].getVector3fMap ().template cast<double> ();
717 if (!std::isfinite(eigen_vector[0]) || !std::isfinite(eigen_vector[1]) || !std::isfinite(eigen_vector[2]))
728 const double distance =
query_point.dot (model_coefficients.head<3> ()) + model_coefficients[3];
746 order = polynomial_order;
749 const int nr_coeff = (
order + 1) * (
order + 2) / 2;
754 weight_func = [
this, search_radius] (
const double sq_dist) {
return this->computeMLSWeight (sq_dist, search_radius * search_radius); };
760 Eigen::MatrixXd P_weight_Pt (nr_coeff, nr_coeff);
764 std::vector<Eigen::Vector3d, Eigen::aligned_allocator<Eigen::Vector3d> > de_meaned (
num_neighbors);
765 for (std::size_t ni = 0; ni < static_cast<std::size_t>(
num_neighbors); ++ni)
767 de_meaned[ni][0] = cloud[nn_indices[ni]].x -
mean[0];
768 de_meaned[ni][1] = cloud[nn_indices[ni]].y -
mean[1];
769 de_meaned[ni][2] = cloud[nn_indices[ni]].z -
mean[2];
770 weight_vec (ni) = weight_func (de_meaned[ni].dot (de_meaned[ni]));
775 for (std::size_t ni = 0; ni < static_cast<std::size_t>(
num_neighbors); ++ni)
778 const double u_coord = de_meaned[ni].dot(
u_axis);
779 const double v_coord = de_meaned[ni].dot(
v_axis);
785 for (
int ui = 0; ui <=
order; ++ui)
788 for (
int vi = 0; vi <=
order - ui; ++vi)
790 P (j++, ni) = u_pow * v_pow;
798 const Eigen::MatrixXd P_weight = P * weight_vec.asDiagonal();
799 P_weight_Pt.noalias() = P_weight * P.transpose ();
800 c_vec.noalias() = P_weight * f_vec;
801 P_weight_Pt.llt ().solveInPlace (
c_vec);
807template <
typename Po
intInT,
typename Po
intOutT>
811 int dilation_iteration_num) :
819 const double max_size = (std::max) ((std::max)(bounding_box_size.x (), bounding_box_size.y ()), bounding_box_size.z ());
822 for (std::size_t i = 0; i < indices->size (); ++i)
823 if (std::isfinite ((*cloud)[(*indices)[i]].x))
826 getCellIndex ((*cloud)[(*indices)[i]].getVector3fMap (), pos);
828 std::uint64_t index_1d;
836template <
typename Po
intInT,
typename Po
intOutT>
void
842 Eigen::Vector3i index;
846 for (
int x = -1; x <= 1; ++x)
847 for (
int y = -1; y <= 1; ++y)
848 for (
int z = -1; z <= 1; ++z)
849 if (x != 0 || y != 0 || z != 0)
851 Eigen::Vector3i new_index;
852 new_index = index + Eigen::Vector3i (x, y, z);
854 std::uint64_t index_1d;
857 new_voxel_grid[index_1d] = leaf;
865template <
typename Po
intInT,
typename Po
intOutT>
void
867 PointOutT &point_out)
const
869 PointOutT temp = point_out;
871 point_out.x = temp.x;
872 point_out.y = temp.y;
873 point_out.z = temp.z;
876#define PCL_INSTANTIATE_MovingLeastSquares(T,OutT) template class PCL_EXPORTS pcl::MovingLeastSquares<T,OutT>;
877#define PCL_INSTANTIATE_MovingLeastSquaresOMP(T,OutT) template class PCL_EXPORTS pcl::MovingLeastSquaresOMP<T,OutT>;
Define methods for centroid estimation and covariance matrix calculus.
A minimalistic implementation of a voxel grid, necessary for the point cloud upsampling.
void getPosition(const std::uint64_t &index_1d, Eigen::Vector3f &point) const
Eigen::Vector4f bounding_min_
void getIndexIn1D(const Eigen::Vector3i &index, std::uint64_t &index_1d) const
Eigen::Vector4f bounding_max_
void getCellIndex(const Eigen::Vector3f &p, Eigen::Vector3i &index) const
void getIndexIn3D(std::uint64_t index_1d, Eigen::Vector3i &index_3d) const
MLSVoxelGrid(PointCloudInConstPtr &cloud, IndicesPtr &indices, float voxel_size, int dilation_iteration_num)
std::map< std::uint64_t, Leaf > HashMap
unsigned int threads_
The maximum number of threads the scheduler should use.
void performUpsampling(PointCloudOut &output)
Perform upsampling for the distinct-cloud and voxel-grid methods.
int order_
The order of the polynomial to be fit.
typename KdTree::Ptr KdTreePtr
MLSResult::ProjectionMethod projection_method_
Parameter that specifies the projection method to be used.
double search_radius_
The nearest neighbors search radius for each point.
pcl::PointCloud< PointOutT > PointCloudOut
double sqr_gauss_param_
Parameter for distance based weighting of neighbors (search_radius_ * search_radius_ works fine).
typename PointCloudIn::ConstPtr PointCloudInConstPtr
KdTreePtr tree_
A pointer to the spatial search object.
void copyMissingFields(const PointInT &point_in, PointOutT &point_out) const
double upsampling_step_
Step size for the local plane sampling.
NormalCloudPtr normals_
The point cloud that will hold the estimated normals, if set.
UpsamplingMethod upsample_method_
Parameter that specifies the upsampling method to be used.
int searchForNeighbors(pcl::index_t index, pcl::Indices &indices, std::vector< float > &sqr_distances) const
Search for the nearest neighbors of a given point using a radius search.
double upsampling_radius_
Radius of the circle in the local point plane that will be sampled.
@ RANDOM_UNIFORM_DENSITY
The local plane of each input point will be sampled using an uniform random distribution such that th...
@ SAMPLE_LOCAL_PLANE
The local plane of each input point will be sampled in a circular fashion using the upsampling_radius...
@ VOXEL_GRID_DILATION
The input cloud will be inserted into a voxel grid with voxels of size voxel_size_; this voxel grid w...
@ NONE
No upsampling will be done, only the input points will be projected to their own MLS surfaces.
@ DISTINCT_CLOUD
Project the points of the distinct cloud to the MLS surface.
void setSearchMethod(const KdTreePtr &tree)
Provide a pointer to the search object.
int desired_num_points_in_radius_
Parameter that specifies the desired number of points within the search radius.
pcl::PointCloud< pcl::Normal > NormalCloud
PointIndicesPtr corresponding_input_indices_
Collects for each point in output the corresponding point in the input.
int nr_coeff_
Number of coefficients, to be computed from the requested order.
bool compute_normals_
Parameter that specifies whether the normals should be computed for the input cloud or not.
void performProcessing(PointCloudOut &output) override
Abstract surface reconstruction method.
void computeMLSPointNormal(pcl::index_t index, const pcl::Indices &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.
void process(PointCloudOut &output) override
Base method for surface reconstruction for all points given in <setInputCloud (), setIndices ()>.
void addProjectedPointNormal(pcl::index_t 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 adding projected points.
int dilation_iteration_num_
Number of dilation steps for the VOXEL_GRID_DILATION upsampling method.
bool cache_mls_results_
True if the mls results for the input cloud should be stored.
std::vector< MLSResult > mls_results_
Stores the MLS result for each point in the input cloud.
float voxel_size_
Voxel size for the VOXEL_GRID_DILATION upsampling method.
PointCloudInConstPtr distinct_cloud_
The distinct point cloud that will be projected to the MLS surface.
PointCloudConstPtr input_
PointCloud represents the base class in PCL for storing collections of 3D points.
void push_back(const PointT &pt)
Insert a new point in the cloud, at the end of the container.
std::uint32_t width
The point cloud width (if organized as an image-structure).
iterator insert(iterator position, const PointT &pt)
Insert a new point in the cloud, given an iterator.
pcl::PCLHeader header
The point cloud header.
std::vector< PointCloud< PointOutT >, Eigen::aligned_allocator< PointCloud< PointOutT > > > CloudVectorType
std::uint32_t height
The point cloud height (if organized as an image-structure).
void clear()
Removes all points in a cloud and sets the width and height to 0.
iterator begin() noexcept
search::KdTree is a wrapper class which inherits the pcl::KdTree class for performing search function...
OrganizedNeighbor is a class for optimized nearest neighbor search in organized projectable point clo...
Define standard C methods and C++ classes that are common to all methods.
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.
void copyPoint(const PointInT &point_in, PointOutT &point_out)
Copy the fields of a source point into a target point.
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.
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...
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.
detail::int_type_t< detail::index_type_size, detail::index_type_signed > index_t
Type used for an index in PCL.
shared_ptr< Indices > IndicesPtr
IndicesAllocator<> Indices
Type used for indices in PCL.
Data structure used to store the MLS projection results.
Eigen::Vector3d point
The projected point.
double v
The v-coordinate of the projected point in local MLS frame.
Eigen::Vector3d normal
The projected point's normal.
double u
The u-coordinate of the projected point in local MLS frame.
Data structure used to store the MLS polynomial partial derivatives.
double z_uv
The partial derivative d^2z/dudv.
double z_u
The partial derivative dz/du.
double z_uu
The partial derivative d^2z/du^2.
double z
The z component of the polynomial evaluated at z(u, v).
double z_vv
The partial derivative d^2z/dv^2.
double z_v
The partial derivative dz/dv.
Data structure used to store the results of the MLS fitting.
MLSProjectionResults projectPoint(const Eigen::Vector3d &pt, ProjectionMethod method, int required_neighbors=0) const
Project a point using the specified method.
Eigen::Vector3d mean
The mean point of all the neighbors.
MLSProjectionResults projectPointOrthogonalToPolynomialSurface(const double u, const double v, const double w) const
Project a point orthogonal to the polynomial surface.
Eigen::Vector3d u_axis
The axis corresponding to the u-coordinates of the local plane of the query point.
Eigen::Vector3d plane_normal
The normal of the local plane of the query point.
@ ORTHOGONAL
Project to the closest point on the polynonomial surface.
@ NONE
Project to the mls plane.
Eigen::Vector3d v_axis
The axis corresponding to the v-coordinates of the local plane of the query point.
int num_neighbors
The number of neighbors used to create the mls surface.
Eigen::VectorXd c_vec
The polynomial coefficients Example: z = c_vec[0] + c_vec[1]*v + c_vec[2]*v^2 + c_vec[3]*u + c_vec[4]...
void computeMLSSurface(const pcl::PointCloud< PointT > &cloud, pcl::index_t index, const pcl::Indices &nn_indices, double search_radius, int polynomial_order=2, std::function< double(const double)> weight_func={})
Smooth a given point and its neighborhood using Moving Least Squares.
void getMLSCoordinates(const Eigen::Vector3d &pt, double &u, double &v, double &w) const
Given a point calculate its 3D location in the MLS frame.
float curvature
The curvature at the query point.
PolynomialPartialDerivative getPolynomialPartialDerivative(const double u, const double v) const
Calculate the polynomial's first and second partial derivatives.
MLSProjectionResults projectPointSimpleToPolynomialSurface(const double u, const double v) const
Project a point along the MLS plane normal to the polynomial surface.
MLSProjectionResults projectPointToMLSPlane(const double u, const double v) const
Project a point onto the MLS plane.
double getPolynomialValue(const double u, const double v) const
Calculate the polynomial.
Eigen::Vector3d query_point
The query point about which the mls surface was generated.
MLSProjectionResults projectQueryPoint(ProjectionMethod method, int required_neighbors=0) const
Project the query point used to generate the mls surface about using the specified method.
int order
The order of the polynomial.
bool valid
If True, the mls results data is valid, otherwise False.
A point structure representing normal coordinates and the surface curvature estimate.
A helper functor that can set a specific value in a field if the field exists.