Point Cloud Library (PCL)  1.14.0-dev
List of all members | Classes | Public Types | Public Member Functions | Public Attributes
pcl::MLSResult Struct Reference

Data structure used to store the results of the MLS fitting. More...

#include <pcl/surface/mls.h>

Classes

struct  MLSProjectionResults
 Data structure used to store the MLS projection results. More...
 
struct  PolynomialPartialDerivative
 Data structure used to store the MLS polynomial partial derivatives. More...
 

Public Types

enum  ProjectionMethod { NONE , SIMPLE , ORTHOGONAL }
 

Public Member Functions

 MLSResult ()
 
 MLSResult (const Eigen::Vector3d &a_query_point, const Eigen::Vector3d &a_mean, const Eigen::Vector3d &a_plane_normal, const Eigen::Vector3d &a_u, const Eigen::Vector3d &a_v, const Eigen::VectorXd &a_c_vec, const int a_num_neighbors, const float a_curvature, const int a_order)
 
void getMLSCoordinates (const Eigen::Vector3d &pt, double &u, double &v, double &w) const
 Given a point calculate its 3D location in the MLS frame. More...
 
void getMLSCoordinates (const Eigen::Vector3d &pt, double &u, double &v) const
 Given a point calculate its 2D location in the MLS frame. More...
 
double getPolynomialValue (const double u, const double v) const
 Calculate the polynomial. More...
 
PolynomialPartialDerivative getPolynomialPartialDerivative (const double u, const double v) const
 Calculate the polynomial's first and second partial derivatives. More...
 
Eigen::Vector2f calculatePrincipalCurvatures (const double u, const double v) const
 Calculate the principal curvatures using the polynomial surface. More...
 
Eigen::Vector2f calculatePrincipleCurvatures (const double u, const double v) const
 Calculate the principal curvatures using the polynomial surface. More...
 
MLSProjectionResults projectPointOrthogonalToPolynomialSurface (const double u, const double v, const double w) const
 Project a point orthogonal to the polynomial surface. More...
 
MLSProjectionResults projectPointToMLSPlane (const double u, const double v) const
 Project a point onto the MLS plane. More...
 
MLSProjectionResults projectPointSimpleToPolynomialSurface (const double u, const double v) const
 Project a point along the MLS plane normal to the polynomial surface. More...
 
MLSProjectionResults projectPoint (const Eigen::Vector3d &pt, ProjectionMethod method, int required_neighbors=0) const
 Project a point using the specified method. More...
 
MLSProjectionResults projectQueryPoint (ProjectionMethod method, int required_neighbors=0) const
 Project the query point used to generate the mls surface about using the specified method. More...
 
template<typename PointT >
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. More...
 

Public Attributes

Eigen::Vector3d query_point
 The query point about which the mls surface was generated. More...
 
Eigen::Vector3d mean
 The mean point of all the neighbors. More...
 
Eigen::Vector3d plane_normal
 The normal of the local plane of the query point. More...
 
Eigen::Vector3d u_axis
 The axis corresponding to the u-coordinates of the local plane of the query point. More...
 
Eigen::Vector3d v_axis
 The axis corresponding to the v-coordinates of the local plane of the query point. More...
 
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]*u*v + c_vec[5]*u^2. More...
 
int num_neighbors
 The number of neighbors used to create the mls surface. More...
 
float curvature
 The curvature at the query point. More...
 
int order
 The order of the polynomial. More...
 
bool valid
 If True, the mls results data is valid, otherwise False. More...
 

Detailed Description

Data structure used to store the results of the MLS fitting.

Definition at line 59 of file mls.h.

Member Enumeration Documentation

◆ ProjectionMethod

Enumerator
NONE 

Project to the mls plane.

SIMPLE 

Project along the mls plane normal to the polynomial surface.

ORTHOGONAL 

Project to the closest point on the polynonomial surface.

Definition at line 61 of file mls.h.

Constructor & Destructor Documentation

◆ MLSResult() [1/2]

pcl::MLSResult::MLSResult ( )
inline

Definition at line 92 of file mls.h.

◆ MLSResult() [2/2]

pcl::MLSResult::MLSResult ( const Eigen::Vector3d &  a_query_point,
const Eigen::Vector3d &  a_mean,
const Eigen::Vector3d &  a_plane_normal,
const Eigen::Vector3d &  a_u,
const Eigen::Vector3d &  a_v,
const Eigen::VectorXd &  a_c_vec,
const int  a_num_neighbors,
const float  a_curvature,
const int  a_order 
)
inline

Definition at line 441 of file mls.hpp.

Member Function Documentation

◆ calculatePrincipalCurvatures()

Eigen::Vector2f pcl::MLSResult::calculatePrincipalCurvatures ( const double  u,
const double  v 
) const

Calculate the principal curvatures using the polynomial surface.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The principal curvature [k1, k2] at the provided uv coordinates.
Note
If an error occurs then 1e-5 is returned.

Referenced by calculatePrincipleCurvatures().

◆ calculatePrincipleCurvatures()

Eigen::Vector2f pcl::MLSResult::calculatePrincipleCurvatures ( const double  u,
const double  v 
) const
inline

Calculate the principal curvatures using the polynomial surface.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The principal curvature [k1, k2] at the provided uv coordinates.
Note
If an error occurs then 1e-5 is returned.
Deprecated:
Scheduled for removal in version 1 . 15 : "use calculatePrincipalCurvatures() instead"

Definition at line 155 of file mls.h.

References calculatePrincipalCurvatures().

◆ computeMLSSurface()

template<typename PointT >
void pcl::MLSResult::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.

Parameters
[in]cloudthe input cloud, used together with index and nn_indices
[in]indexthe index of the query point in the input cloud
[in]nn_indicesthe set of nearest neighbors indices for pt
[in]search_radiusthe search radius used to find nearest neighbors for pt
[in]polynomial_orderthe order of the polynomial to fit to the nearest neighbors
[in]weight_funcdefines the weight function for the polynomial fit

Definition at line 692 of file mls.hpp.

References pcl::compute3DCentroid(), pcl::computeCovarianceMatrix(), pcl::geometry::distance(), and pcl::eigen33().

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

◆ getMLSCoordinates() [1/2]

void pcl::MLSResult::getMLSCoordinates ( const Eigen::Vector3d &  pt,
double &  u,
double &  v 
) const
inline

Given a point calculate its 2D location in the MLS frame.

Parameters
[in]ptThe point
[out]uThe u-coordinate of the point in local MLS frame.
[out]vThe v-coordinate of the point in local MLS frame.

Definition at line 464 of file mls.hpp.

◆ getMLSCoordinates() [2/2]

void pcl::MLSResult::getMLSCoordinates ( const Eigen::Vector3d &  pt,
double &  u,
double &  v,
double &  w 
) const
inline

Given a point calculate its 3D location in the MLS frame.

Parameters
[in]ptThe point
[out]uThe u-coordinate of the point in local MLS frame.
[out]vThe v-coordinate of the point in local MLS frame.
[out]wThe w-coordinate of the point in local MLS frame.

Definition at line 455 of file mls.hpp.

◆ getPolynomialPartialDerivative()

pcl::MLSResult::PolynomialPartialDerivative pcl::MLSResult::getPolynomialPartialDerivative ( const double  u,
const double  v 
) const
inline

Calculate the polynomial's first and second partial derivatives.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The polynomial partial derivatives at the provided uv coordinates.

Definition at line 494 of file mls.hpp.

◆ getPolynomialValue()

double pcl::MLSResult::getPolynomialValue ( const double  u,
const double  v 
) const
inline

Calculate the polynomial.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The polynomial value at the provided uv coordinates.

Definition at line 472 of file mls.hpp.

◆ projectPoint()

pcl::MLSResult::MLSProjectionResults pcl::MLSResult::projectPoint ( const Eigen::Vector3d &  pt,
ProjectionMethod  method,
int  required_neighbors = 0 
) const
inline

Project a point using the specified method.

Parameters
[in]ptThe point to be project.
[in]methodThe projection method to be used.
[in]required_neighborsThe minimum number of neighbors required.
Note
If required_neighbors is 0 then any number of neighbors is allowed.
If required_neighbors is not satisfied it projects to the mls plane.
Returns
The MLSProjectionResults for the input data.

Definition at line 639 of file mls.hpp.

◆ projectPointOrthogonalToPolynomialSurface()

pcl::MLSResult::MLSProjectionResults pcl::MLSResult::projectPointOrthogonalToPolynomialSurface ( const double  u,
const double  v,
const double  w 
) const
inline

Project a point orthogonal to the polynomial surface.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
[in]wThe w-coordinate of the point in local MLS frame.
Returns
The MLSProjectionResults for the input data.
Note
If the MLSResults does not contain polynomial data it projects the point onto the mls plane.
If the optimization diverges it performs a simple projection on to the polynomial surface.
This was implemented based on this https://math.stackexchange.com/questions/1497093/shortest-distance-between-point-and-surface

Definition at line 539 of file mls.hpp.

References pcl::MLSResult::MLSProjectionResults::normal, pcl::MLSResult::MLSProjectionResults::point, pcl::MLSResult::MLSProjectionResults::u, pcl::MLSResult::MLSProjectionResults::v, pcl::MLSResult::PolynomialPartialDerivative::z, pcl::MLSResult::PolynomialPartialDerivative::z_u, pcl::MLSResult::PolynomialPartialDerivative::z_uu, pcl::MLSResult::PolynomialPartialDerivative::z_uv, pcl::MLSResult::PolynomialPartialDerivative::z_v, and pcl::MLSResult::PolynomialPartialDerivative::z_vv.

◆ projectPointSimpleToPolynomialSurface()

pcl::MLSResult::MLSProjectionResults pcl::MLSResult::projectPointSimpleToPolynomialSurface ( const double  u,
const double  v 
) const
inline

Project a point along the MLS plane normal to the polynomial surface.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The MLSProjectionResults for the input data.
Note
If the MLSResults does not contain polynomial data it projects the point onto the mls plane.

Definition at line 616 of file mls.hpp.

References pcl::MLSResult::MLSProjectionResults::normal, pcl::MLSResult::MLSProjectionResults::point, pcl::MLSResult::MLSProjectionResults::u, pcl::MLSResult::MLSProjectionResults::v, pcl::MLSResult::PolynomialPartialDerivative::z, pcl::MLSResult::PolynomialPartialDerivative::z_u, and pcl::MLSResult::PolynomialPartialDerivative::z_v.

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

◆ projectPointToMLSPlane()

pcl::MLSResult::MLSProjectionResults pcl::MLSResult::projectPointToMLSPlane ( const double  u,
const double  v 
) const
inline

Project a point onto the MLS plane.

Parameters
[in]uThe u-coordinate of the point in local MLS frame.
[in]vThe v-coordinate of the point in local MLS frame.
Returns
The MLSProjectionResults for the input data.

Definition at line 604 of file mls.hpp.

References pcl::MLSResult::MLSProjectionResults::normal, pcl::MLSResult::MLSProjectionResults::point, pcl::MLSResult::MLSProjectionResults::u, and pcl::MLSResult::MLSProjectionResults::v.

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

◆ projectQueryPoint()

pcl::MLSResult::MLSProjectionResults pcl::MLSResult::projectQueryPoint ( ProjectionMethod  method,
int  required_neighbors = 0 
) const
inline

Project the query point used to generate the mls surface about using the specified method.

Parameters
[in]methodThe projection method to be used.
[in]required_neighborsThe minimum number of neighbors required.
Note
If required_neighbors is 0 then any number of neighbors is allowed.
If required_neighbors is not satisfied it projects to the mls plane.
Returns
The MLSProjectionResults for the input data.

Definition at line 661 of file mls.hpp.

References pcl::MLSResult::MLSProjectionResults::normal, and pcl::MLSResult::MLSProjectionResults::point.

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

Member Data Documentation

◆ c_vec

Eigen::VectorXd pcl::MLSResult::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]*u*v + c_vec[5]*u^2.

Definition at line 230 of file mls.h.

◆ curvature

float pcl::MLSResult::curvature

The curvature at the query point.

Definition at line 232 of file mls.h.

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

◆ mean

Eigen::Vector3d pcl::MLSResult::mean

The mean point of all the neighbors.

Definition at line 226 of file mls.h.

◆ num_neighbors

int pcl::MLSResult::num_neighbors

The number of neighbors used to create the mls surface.

Definition at line 231 of file mls.h.

Referenced by pcl::MovingLeastSquares< PointInT, PointOutT >::computeMLSPointNormal().

◆ order

int pcl::MLSResult::order

The order of the polynomial.

If order > 1 then use polynomial fit

Definition at line 233 of file mls.h.

◆ plane_normal

Eigen::Vector3d pcl::MLSResult::plane_normal

The normal of the local plane of the query point.

Definition at line 227 of file mls.h.

◆ query_point

Eigen::Vector3d pcl::MLSResult::query_point

The query point about which the mls surface was generated.

Definition at line 225 of file mls.h.

◆ u_axis

Eigen::Vector3d pcl::MLSResult::u_axis

The axis corresponding to the u-coordinates of the local plane of the query point.

Definition at line 228 of file mls.h.

◆ v_axis

Eigen::Vector3d pcl::MLSResult::v_axis

The axis corresponding to the v-coordinates of the local plane of the query point.

Definition at line 229 of file mls.h.

◆ valid

bool pcl::MLSResult::valid

If True, the mls results data is valid, otherwise False.

Definition at line 234 of file mls.h.


The documentation for this struct was generated from the following files: