Point Cloud Library (PCL)  1.14.0-dev
List of all members | Public Types | Public Member Functions | Protected Attributes
pcl::RangeImagePlanar Class Reference

RangeImagePlanar is derived from the original range image and differs from it because it's not a spherical projection, but using a projection plane (as normal cameras do), therefore being better applicable for range sensors that already provide a range image by themselves (stereo cameras, ToF-cameras), so that a conversion to point cloud and then to a spherical range image becomes unnecessary. More...

#include <pcl/range_image/range_image_planar.h>

+ Inheritance diagram for pcl::RangeImagePlanar:
+ Collaboration diagram for pcl::RangeImagePlanar:

Public Types

using BaseClass = RangeImage
 
using Ptr = shared_ptr< RangeImagePlanar >
 
using ConstPtr = shared_ptr< const RangeImagePlanar >
 
- Public Types inherited from pcl::RangeImage
enum  CoordinateFrame { CAMERA_FRAME = 0 , LASER_FRAME = 1 }
 
using BaseClass = pcl::PointCloud< PointWithRange >
 
using VectorOfEigenVector3f = std::vector< Eigen::Vector3f, Eigen::aligned_allocator< Eigen::Vector3f > >
 
using Ptr = shared_ptr< RangeImage >
 
using ConstPtr = shared_ptr< const RangeImage >
 
- Public Types inherited from pcl::PointCloud< PointWithRange >
using PointType = PointWithRange
 
using VectorType = std::vector< PointWithRange, Eigen::aligned_allocator< PointWithRange > >
 
using CloudVectorType = std::vector< PointCloud< PointWithRange >, Eigen::aligned_allocator< PointCloud< PointWithRange > > >
 
using Ptr = shared_ptr< PointCloud< PointWithRange > >
 
using ConstPtr = shared_ptr< const PointCloud< PointWithRange > >
 
using value_type = PointWithRange
 
using reference = PointWithRange &
 
using const_reference = const PointWithRange &
 
using difference_type = typename VectorType::difference_type
 
using size_type = typename VectorType::size_type
 
using iterator = typename VectorType::iterator
 
using const_iterator = typename VectorType::const_iterator
 
using reverse_iterator = typename VectorType::reverse_iterator
 
using const_reverse_iterator = typename VectorType::const_reverse_iterator
 

Public Member Functions

PCL_EXPORTS RangeImagePlanar ()
 Constructor. More...
 
PCL_EXPORTS ~RangeImagePlanar () override
 Destructor. More...
 
RangeImagegetNew () const override
 Return a newly created RangeImagePlanar. More...
 
PCL_EXPORTS void copyTo (RangeImage &other) const override
 Copy *this to other. More...
 
Ptr makeShared ()
 Get a boost shared pointer of a copy of this. More...
 
PCL_EXPORTS void setDisparityImage (const float *disparity_image, int di_width, int di_height, float focal_length, float base_line, float desired_angular_resolution=-1)
 Create the image from an existing disparity image. More...
 
PCL_EXPORTS void setDepthImage (const float *depth_image, int di_width, int di_height, float di_center_x, float di_center_y, float di_focal_length_x, float di_focal_length_y, float desired_angular_resolution=-1)
 Create the image from an existing depth image. More...
 
PCL_EXPORTS void setDepthImage (const unsigned short *depth_image, int di_width, int di_height, float di_center_x, float di_center_y, float di_focal_length_x, float di_focal_length_y, float desired_angular_resolution=-1)
 Create the image from an existing depth image. More...
 
template<typename PointCloudType >
void createFromPointCloudWithFixedSize (const PointCloudType &point_cloud, int di_width, int di_height, float di_center_x, float di_center_y, float di_focal_length_x, float di_focal_length_y, const Eigen::Affine3f &sensor_pose, CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f)
 Create the image from an existing point cloud. More...
 
void calculate3DPoint (float image_x, float image_y, float range, Eigen::Vector3f &point) const override
 Calculate the 3D point according to the given image point and range. More...
 
void getImagePoint (const Eigen::Vector3f &point, float &image_x, float &image_y, float &range) const override
 Calculate the image point and range from the given 3D point. More...
 
PCL_EXPORTS void getSubImage (int sub_image_image_offset_x, int sub_image_image_offset_y, int sub_image_width, int sub_image_height, int combine_pixels, RangeImage &sub_image) const override
 Get a sub part of the complete image as a new range image. More...
 
PCL_EXPORTS void getHalfImage (RangeImage &half_image) const override
 Get a range image with half the resolution. More...
 
float getFocalLengthX () const
 Getter for the focal length in X. More...
 
float getFocalLengthY () const
 Getter for the focal length in Y. More...
 
float getCenterX () const
 Getter for the principal point in X. More...
 
float getCenterY () const
 Getter for the principal point in Y. More...
 
void calculate3DPoint (float image_x, float image_y, float range, PointWithRange &point) const
 Calculate the 3D point according to the given image point and range. More...
 
void calculate3DPoint (float image_x, float image_y, PointWithRange &point) const
 Calculate the 3D point according to the given image point and the range value at the closest pixel. More...
 
virtual void calculate3DPoint (float image_x, float image_y, float range, Eigen::Vector3f &point) const
 Calculate the 3D point according to the given image point and range. More...
 
void calculate3DPoint (float image_x, float image_y, Eigen::Vector3f &point) const
 Calculate the 3D point according to the given image point and the range value at the closest pixel. More...
 
virtual void getImagePoint (const Eigen::Vector3f &point, float &image_x, float &image_y, float &range) const
 Get imagePoint from 3D point in world coordinates. More...
 
void getImagePoint (const Eigen::Vector3f &point, int &image_x, int &image_y, float &range) const
 Same as above. More...
 
void getImagePoint (const Eigen::Vector3f &point, float &image_x, float &image_y) const
 Same as above. More...
 
void getImagePoint (const Eigen::Vector3f &point, int &image_x, int &image_y) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, float &image_x, float &image_y, float &range) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, float &image_x, float &image_y) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, int &image_x, int &image_y) const
 Same as above. More...
 
- Public Member Functions inherited from pcl::RangeImage
PCL_EXPORTS RangeImage ()
 Constructor. More...
 
virtual PCL_EXPORTS ~RangeImage ()=default
 Destructor. More...
 
Ptr makeShared ()
 Get a boost shared pointer of a copy of this. More...
 
PCL_EXPORTS void reset ()
 Reset all values to an empty range image. More...
 
template<typename PointCloudType >
void createFromPointCloud (const PointCloudType &point_cloud, float angular_resolution=pcl::deg2rad(0.5f), float max_angle_width=pcl::deg2rad(360.0f), float max_angle_height=pcl::deg2rad(180.0f), const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud. More...
 
template<typename PointCloudType >
void createFromPointCloud (const PointCloudType &point_cloud, float angular_resolution_x=pcl::deg2rad(0.5f), float angular_resolution_y=pcl::deg2rad(0.5f), float max_angle_width=pcl::deg2rad(360.0f), float max_angle_height=pcl::deg2rad(180.0f), const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud. More...
 
template<typename PointCloudType >
void createFromPointCloudWithKnownSize (const PointCloudType &point_cloud, float angular_resolution, const Eigen::Vector3f &point_cloud_center, float point_cloud_radius, const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud, getting a hint about the size of the scene for faster calculation. More...
 
template<typename PointCloudType >
void createFromPointCloudWithKnownSize (const PointCloudType &point_cloud, float angular_resolution_x, float angular_resolution_y, const Eigen::Vector3f &point_cloud_center, float point_cloud_radius, const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud, getting a hint about the size of the scene for faster calculation. More...
 
template<typename PointCloudTypeWithViewpoints >
void createFromPointCloudWithViewpoints (const PointCloudTypeWithViewpoints &point_cloud, float angular_resolution, float max_angle_width, float max_angle_height, CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud, using the average viewpoint of the points (vp_x,vp_y,vp_z in the point type) in the point cloud as sensor pose (assuming a rotation of (0,0,0)). More...
 
template<typename PointCloudTypeWithViewpoints >
void createFromPointCloudWithViewpoints (const PointCloudTypeWithViewpoints &point_cloud, float angular_resolution_x, float angular_resolution_y, float max_angle_width, float max_angle_height, CoordinateFrame coordinate_frame=CAMERA_FRAME, float noise_level=0.0f, float min_range=0.0f, int border_size=0)
 Create the depth image from a point cloud, using the average viewpoint of the points (vp_x,vp_y,vp_z in the point type) in the point cloud as sensor pose (assuming a rotation of (0,0,0)). More...
 
void createEmpty (float angular_resolution, const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), RangeImage::CoordinateFrame coordinate_frame=CAMERA_FRAME, float angle_width=pcl::deg2rad(360.0f), float angle_height=pcl::deg2rad(180.0f))
 Create an empty depth image (filled with unobserved points) More...
 
void createEmpty (float angular_resolution_x, float angular_resolution_y, const Eigen::Affine3f &sensor_pose=Eigen::Affine3f::Identity(), RangeImage::CoordinateFrame coordinate_frame=CAMERA_FRAME, float angle_width=pcl::deg2rad(360.0f), float angle_height=pcl::deg2rad(180.0f))
 Create an empty depth image (filled with unobserved points) More...
 
template<typename PointCloudType >
void doZBuffer (const PointCloudType &point_cloud, float noise_level, float min_range, int &top, int &right, int &bottom, int &left)
 Integrate the given point cloud into the current range image using a z-buffer. More...
 
template<typename PointCloudType >
void integrateFarRanges (const PointCloudType &far_ranges)
 Integrates the given far range measurements into the range image. More...
 
PCL_EXPORTS void cropImage (int border_size=0, int top=-1, int right=-1, int bottom=-1, int left=-1)
 Cut the range image to the minimal size so that it still contains all actual range readings. More...
 
PCL_EXPORTS float * getRangesArray () const
 Get all the range values in one float array of size width*height. More...
 
const Eigen::Affine3f & getTransformationToRangeImageSystem () const
 Getter for the transformation from the world system into the range image system (the sensor coordinate frame) More...
 
void setTransformationToRangeImageSystem (const Eigen::Affine3f &to_range_image_system)
 Setter for the transformation from the range image system (the sensor coordinate frame) into the world system. More...
 
const Eigen::Affine3f & getTransformationToWorldSystem () const
 Getter for the transformation from the range image system (the sensor coordinate frame) into the world system. More...
 
float getAngularResolution () const
 Getter for the angular resolution of the range image in x direction in radians per pixel. More...
 
float getAngularResolutionX () const
 Getter for the angular resolution of the range image in x direction in radians per pixel. More...
 
float getAngularResolutionY () const
 Getter for the angular resolution of the range image in y direction in radians per pixel. More...
 
void getAngularResolution (float &angular_resolution_x, float &angular_resolution_y) const
 Getter for the angular resolution of the range image in x and y direction (in radians). More...
 
void setAngularResolution (float angular_resolution)
 Set the angular resolution of the range image. More...
 
void setAngularResolution (float angular_resolution_x, float angular_resolution_y)
 Set the angular resolution of the range image. More...
 
const PointWithRangegetPoint (int image_x, int image_y) const
 Return the 3D point with range at the given image position. More...
 
PointWithRangegetPoint (int image_x, int image_y)
 Non-const-version of getPoint. More...
 
const PointWithRangegetPoint (float image_x, float image_y) const
 Return the 3d point with range at the given image position. More...
 
PointWithRangegetPoint (float image_x, float image_y)
 Non-const-version of the above. More...
 
const PointWithRangegetPointNoCheck (int image_x, int image_y) const
 Return the 3D point with range at the given image position. More...
 
PointWithRangegetPointNoCheck (int image_x, int image_y)
 Non-const-version of getPointNoCheck. More...
 
void getPoint (int image_x, int image_y, Eigen::Vector3f &point) const
 Same as above. More...
 
void getPoint (int index, Eigen::Vector3f &point) const
 Same as above. More...
 
const Eigen::Map< const Eigen::Vector3f > getEigenVector3f (int x, int y) const
 Same as above. More...
 
const Eigen::Map< const Eigen::Vector3f > getEigenVector3f (int index) const
 Same as above. More...
 
const PointWithRangegetPoint (int index) const
 Return the 3d point with range at the given index (whereas index=y*width+x) More...
 
void calculate3DPoint (float image_x, float image_y, float range, PointWithRange &point) const
 Calculate the 3D point according to the given image point and range. More...
 
void calculate3DPoint (float image_x, float image_y, PointWithRange &point) const
 Calculate the 3D point according to the given image point and the range value at the closest pixel. More...
 
void calculate3DPoint (float image_x, float image_y, Eigen::Vector3f &point) const
 Calculate the 3D point according to the given image point and the range value at the closest pixel. More...
 
PCL_EXPORTS void recalculate3DPointPositions ()
 Recalculate all 3D point positions according to their pixel position and range. More...
 
void getImagePoint (const Eigen::Vector3f &point, int &image_x, int &image_y, float &range) const
 Same as above. More...
 
void getImagePoint (const Eigen::Vector3f &point, float &image_x, float &image_y) const
 Same as above. More...
 
void getImagePoint (const Eigen::Vector3f &point, int &image_x, int &image_y) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, float &image_x, float &image_y, float &range) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, float &image_x, float &image_y) const
 Same as above. More...
 
void getImagePoint (float x, float y, float z, int &image_x, int &image_y) const
 Same as above. More...
 
float checkPoint (const Eigen::Vector3f &point, PointWithRange &point_in_image) const
 point_in_image will be the point in the image at the position the given point would be. More...
 
float getRangeDifference (const Eigen::Vector3f &point) const
 Returns the difference in range between the given point and the range of the point in the image at the position the given point would be. More...
 
void getImagePointFromAngles (float angle_x, float angle_y, float &image_x, float &image_y) const
 Get the image point corresponding to the given angles. More...
 
void getAnglesFromImagePoint (float image_x, float image_y, float &angle_x, float &angle_y) const
 Get the angles corresponding to the given image point. More...
 
void real2DToInt2D (float x, float y, int &xInt, int &yInt) const
 Transforms an image point in float values to an image point in int values. More...
 
bool isInImage (int x, int y) const
 Check if a point is inside of the image. More...
 
bool isValid (int x, int y) const
 Check if a point is inside of the image and has a finite range. More...
 
bool isValid (int index) const
 Check if a point has a finite range. More...
 
bool isObserved (int x, int y) const
 Check if a point is inside of the image and has either a finite range or a max reading (range=INFINITY) More...
 
bool isMaxRange (int x, int y) const
 Check if a point is a max range (range=INFINITY) - please check isInImage or isObserved first! More...
 
bool getNormal (int x, int y, int radius, Eigen::Vector3f &normal, int step_size=1) const
 Calculate the normal of an image point using the neighbors with a maximum pixel distance of radius. More...
 
bool getNormalForClosestNeighbors (int x, int y, int radius, const PointWithRange &point, int no_of_nearest_neighbors, Eigen::Vector3f &normal, int step_size=1) const
 Same as above, but only the no_of_nearest_neighbors points closest to the given point are considered. More...
 
bool getNormalForClosestNeighbors (int x, int y, int radius, const Eigen::Vector3f &point, int no_of_nearest_neighbors, Eigen::Vector3f &normal, Eigen::Vector3f *point_on_plane=nullptr, int step_size=1) const
 Same as above. More...
 
bool getNormalForClosestNeighbors (int x, int y, Eigen::Vector3f &normal, int radius=2) const
 Same as above, using default values. More...
 
bool getSurfaceInformation (int x, int y, int radius, const Eigen::Vector3f &point, int no_of_closest_neighbors, int step_size, float &max_closest_neighbor_distance_squared, Eigen::Vector3f &normal, Eigen::Vector3f &mean, Eigen::Vector3f &eigen_values, Eigen::Vector3f *normal_all_neighbors=nullptr, Eigen::Vector3f *mean_all_neighbors=nullptr, Eigen::Vector3f *eigen_values_all_neighbors=nullptr) const
 Same as above but extracts some more data and can also return the extracted information for all neighbors in radius if normal_all_neighbors is not NULL. More...
 
float getSquaredDistanceOfNthNeighbor (int x, int y, int radius, int n, int step_size) const
 
float getImpactAngle (const PointWithRange &point1, const PointWithRange &point2) const
 Calculate the impact angle based on the sensor position and the two given points - will return -INFINITY if one of the points is unobserved. More...
 
float getImpactAngle (int x1, int y1, int x2, int y2) const
 Same as above. More...
 
float getImpactAngleBasedOnLocalNormal (int x, int y, int radius) const
 Extract a local normal (with a heuristic not to include background points) and calculate the impact angle based on this. More...
 
PCL_EXPORTS float * getImpactAngleImageBasedOnLocalNormals (int radius) const
 Uses the above function for every point in the image. More...
 
float getNormalBasedAcutenessValue (int x, int y, int radius) const
 Calculate a score [0,1] that tells how acute the impact angle is (1.0f - getImpactAngle/90deg) This uses getImpactAngleBasedOnLocalNormal Will return -INFINITY if no normal could be calculated. More...
 
float getAcutenessValue (const PointWithRange &point1, const PointWithRange &point2) const
 Calculate a score [0,1] that tells how acute the impact angle is (1.0f - getImpactAngle/90deg) will return -INFINITY if one of the points is unobserved. More...
 
float getAcutenessValue (int x1, int y1, int x2, int y2) const
 Same as above. More...
 
PCL_EXPORTS void getAcutenessValueImages (int pixel_distance, float *&acuteness_value_image_x, float *&acuteness_value_image_y) const
 Calculate getAcutenessValue for every point. More...
 
PCL_EXPORTS float getSurfaceChange (int x, int y, int radius) const
 Calculates, how much the surface changes at a point. More...
 
PCL_EXPORTS float * getSurfaceChangeImage (int radius) const
 Uses the above function for every point in the image. More...
 
void getSurfaceAngleChange (int x, int y, int radius, float &angle_change_x, float &angle_change_y) const
 Calculates, how much the surface changes at a point. More...
 
PCL_EXPORTS void getSurfaceAngleChangeImages (int radius, float *&angle_change_image_x, float *&angle_change_image_y) const
 Uses the above function for every point in the image. More...
 
float getCurvature (int x, int y, int radius, int step_size) const
 Calculates the curvature in a point using pca. More...
 
const Eigen::Vector3f getSensorPos () const
 Get the sensor position. More...
 
PCL_EXPORTS void setUnseenToMaxRange ()
 Sets all -INFINITY values to INFINITY. More...
 
int getImageOffsetX () const
 Getter for image_offset_x_. More...
 
int getImageOffsetY () const
 Getter for image_offset_y_. More...
 
void setImageOffsets (int offset_x, int offset_y)
 Setter for image offsets. More...
 
PCL_EXPORTS void getMinMaxRanges (float &min_range, float &max_range) const
 Find the minimum and maximum range in the image. More...
 
PCL_EXPORTS void change3dPointsToLocalCoordinateFrame ()
 This function sets the sensor pose to 0 and transforms all point positions to this local coordinate frame. More...
 
PCL_EXPORTS float * getInterpolatedSurfaceProjection (const Eigen::Affine3f &pose, int pixel_size, float world_size) const
 Calculate a range patch as the z values of the coordinate frame given by pose. More...
 
PCL_EXPORTS float * getInterpolatedSurfaceProjection (const Eigen::Vector3f &point, int pixel_size, float world_size) const
 Same as above, but using the local coordinate frame defined by point and the viewing direction. More...
 
Eigen::Affine3f getTransformationToViewerCoordinateFrame (const Eigen::Vector3f &point) const
 Get the local coordinate frame with 0,0,0 in point, upright and Z as the viewing direction. More...
 
void getTransformationToViewerCoordinateFrame (const Eigen::Vector3f &point, Eigen::Affine3f &transformation) const
 Same as above, using a reference for the retrurn value. More...
 
void getRotationToViewerCoordinateFrame (const Eigen::Vector3f &point, Eigen::Affine3f &transformation) const
 Same as above, but only returning the rotation. More...
 
PCL_EXPORTS bool getNormalBasedUprightTransformation (const Eigen::Vector3f &point, float max_dist, Eigen::Affine3f &transformation) const
 Get a local coordinate frame at the given point based on the normal. More...
 
PCL_EXPORTS void getIntegralImage (float *&integral_image, int *&valid_points_num_image) const
 Get the integral image of the range values (used for fast blur operations). More...
 
PCL_EXPORTS void getBlurredImageUsingIntegralImage (int blur_radius, float *integral_image, int *valid_points_num_image, RangeImage &range_image) const
 Get a blurred version of the range image using box filters on the provided integral image. More...
 
virtual PCL_EXPORTS void getBlurredImage (int blur_radius, RangeImage &range_image) const
 Get a blurred version of the range image using box filters. More...
 
float getEuclideanDistanceSquared (int x1, int y1, int x2, int y2) const
 Get the squared euclidean distance between the two image points. More...
 
float getAverageEuclideanDistance (int x, int y, int offset_x, int offset_y, int max_steps) const
 Doing the above for some steps in the given direction and averaging. More...
 
PCL_EXPORTS void getRangeImageWithSmoothedSurface (int radius, RangeImage &smoothed_range_image) const
 Project all points on the local plane approximation, thereby smoothing the surface of the scan. More...
 
void get1dPointAverage (int x, int y, int delta_x, int delta_y, int no_of_points, PointWithRange &average_point) const
 Calculates the average 3D position of the no_of_points points described by the start point x,y in the direction delta. More...
 
PCL_EXPORTS float getOverlap (const RangeImage &other_range_image, const Eigen::Affine3f &relative_transformation, int search_radius, float max_distance, int pixel_step=1) const
 Calculates the overlap of two range images given the relative transformation (from the given image to *this) More...
 
bool getViewingDirection (int x, int y, Eigen::Vector3f &viewing_direction) const
 Get the viewing direction for the given point. More...
 
void getViewingDirection (const Eigen::Vector3f &point, Eigen::Vector3f &viewing_direction) const
 Get the viewing direction for the given point. More...
 
- Public Member Functions inherited from pcl::PointCloud< PointWithRange >
 PointCloud ()=default
 Default constructor. More...
 
 PointCloud (const PointCloud< PointWithRange > &pc, const Indices &indices)
 Copy constructor from point cloud subset. More...
 
 PointCloud (std::uint32_t width_, std::uint32_t height_, const PointWithRange &value_=PointWithRange())
 Allocate constructor from point cloud subset. More...
 
PointCloudoperator+= (const PointCloud &rhs)
 Add a point cloud to the current cloud. More...
 
PointCloud operator+ (const PointCloud &rhs)
 Add a point cloud to another cloud. More...
 
const PointWithRangeat (int column, int row) const
 Obtain the point given by the (column, row) coordinates. More...
 
PointWithRangeat (int column, int row)
 Obtain the point given by the (column, row) coordinates. More...
 
const PointWithRangeat (std::size_t n) const
 
PointWithRangeat (std::size_t n)
 
const PointWithRangeoperator() (std::size_t column, std::size_t row) const
 Obtain the point given by the (column, row) coordinates. More...
 
PointWithRangeoperator() (std::size_t column, std::size_t row)
 Obtain the point given by the (column, row) coordinates. More...
 
bool isOrganized () const
 Return whether a dataset is organized (e.g., arranged in a structured grid). More...
 
Eigen::Map< Eigen::MatrixXf, Eigen::Aligned, Eigen::OuterStride<> > getMatrixXfMap (int dim, int stride, int offset)
 Return an Eigen MatrixXf (assumes float values) mapped to the specified dimensions of the PointCloud. More...
 
const Eigen::Map< const Eigen::MatrixXf, Eigen::Aligned, Eigen::OuterStride<> > getMatrixXfMap (int dim, int stride, int offset) const
 Return an Eigen MatrixXf (assumes float values) mapped to the specified dimensions of the PointCloud. More...
 
Eigen::Map< Eigen::MatrixXf, Eigen::Aligned, Eigen::OuterStride<> > getMatrixXfMap ()
 
const Eigen::Map< const Eigen::MatrixXf, Eigen::Aligned, Eigen::OuterStride<> > getMatrixXfMap () const
 
iterator begin () noexcept
 
const_iterator begin () const noexcept
 
iterator end () noexcept
 
const_iterator end () const noexcept
 
const_iterator cbegin () const noexcept
 
const_iterator cend () const noexcept
 
reverse_iterator rbegin () noexcept
 
const_reverse_iterator rbegin () const noexcept
 
reverse_iterator rend () noexcept
 
const_reverse_iterator rend () const noexcept
 
const_reverse_iterator crbegin () const noexcept
 
const_reverse_iterator crend () const noexcept
 
std::size_t size () const
 
index_t max_size () const noexcept
 
void reserve (std::size_t n)
 
bool empty () const
 
PointWithRangedata () noexcept
 
const PointWithRangedata () const noexcept
 
void resize (std::size_t count)
 Resizes the container to contain count elements. More...
 
void resize (uindex_t new_width, uindex_t new_height)
 Resizes the container to contain new_width * new_height elements. More...
 
void resize (index_t count, const PointWithRange &value)
 Resizes the container to contain count elements. More...
 
void resize (index_t new_width, index_t new_height, const PointWithRange &value)
 Resizes the container to contain count elements. More...
 
const PointWithRangeoperator[] (std::size_t n) const
 
PointWithRangeoperator[] (std::size_t n)
 
const PointWithRangefront () const
 
PointWithRangefront ()
 
const PointWithRangeback () const
 
PointWithRangeback ()
 
void assign (index_t count, const PointWithRange &value)
 Replaces the points with count copies of value More...
 
void assign (index_t new_width, index_t new_height, const PointWithRange &value)
 Replaces the points with new_width * new_height copies of value More...
 
void assign (InputIterator first, InputIterator last)
 Replaces the points with copies of those in the range [first, last) More...
 
void assign (InputIterator first, InputIterator last, index_t new_width)
 Replaces the points with copies of those in the range [first, last) More...
 
void assign (std::initializer_list< PointWithRange > ilist)
 Replaces the points with the elements from the initializer list ilist More...
 
void assign (std::initializer_list< PointWithRange > ilist, index_t new_width)
 Replaces the points with the elements from the initializer list ilist More...
 
void push_back (const PointWithRange &pt)
 Insert a new point in the cloud, at the end of the container. More...
 
void transient_push_back (const PointWithRange &pt)
 Insert a new point in the cloud, at the end of the container. More...
 
reference emplace_back (Args &&...args)
 Emplace a new point in the cloud, at the end of the container. More...
 
reference transient_emplace_back (Args &&...args)
 Emplace a new point in the cloud, at the end of the container. More...
 
iterator insert (iterator position, const PointWithRange &pt)
 Insert a new point in the cloud, given an iterator. More...
 
void insert (iterator position, std::size_t n, const PointWithRange &pt)
 Insert a new point in the cloud N times, given an iterator. More...
 
void insert (iterator position, InputIterator first, InputIterator last)
 Insert a new range of points in the cloud, at a certain position. More...
 
iterator transient_insert (iterator position, const PointWithRange &pt)
 Insert a new point in the cloud, given an iterator. More...
 
void transient_insert (iterator position, std::size_t n, const PointWithRange &pt)
 Insert a new point in the cloud N times, given an iterator. More...
 
void transient_insert (iterator position, InputIterator first, InputIterator last)
 Insert a new range of points in the cloud, at a certain position. More...
 
iterator emplace (iterator position, Args &&...args)
 Emplace a new point in the cloud, given an iterator. More...
 
iterator transient_emplace (iterator position, Args &&...args)
 Emplace a new point in the cloud, given an iterator. More...
 
iterator erase (iterator position)
 Erase a point in the cloud. More...
 
iterator erase (iterator first, iterator last)
 Erase a set of points given by a (first, last) iterator pair. More...
 
iterator transient_erase (iterator position)
 Erase a point in the cloud. More...
 
iterator transient_erase (iterator first, iterator last)
 Erase a set of points given by a (first, last) iterator pair. More...
 
void swap (PointCloud< PointWithRange > &rhs)
 Swap a point cloud with another cloud. More...
 
void clear ()
 Removes all points in a cloud and sets the width and height to 0. More...
 
Ptr makeShared () const
 Copy the cloud to the heap and return a smart pointer Note that deep copy is performed, so avoid using this function on non-empty clouds. More...
 

Protected Attributes

float focal_length_x_ {0.0f}
 
float focal_length_y_ {0.0f}
 The focal length of the image in pixels. More...
 
float focal_length_x_reciprocal_ {0.0f}
 
float focal_length_y_reciprocal_ {0.0f}
 1/focal_length -> for internal use More...
 
float center_x_ {0.0f}
 
float center_y_ {0.0f}
 The principle point of the image. More...
 
- Protected Attributes inherited from pcl::RangeImage
Eigen::Affine3f to_range_image_system_
 Inverse of to_world_system_. More...
 
Eigen::Affine3f to_world_system_
 Inverse of to_range_image_system_. More...
 
float angular_resolution_x_ {0.0f}
 Angular resolution of the range image in x direction in radians per pixel. More...
 
float angular_resolution_y_ {0.0f}
 Angular resolution of the range image in y direction in radians per pixel. More...
 
float angular_resolution_x_reciprocal_ {0.0f}
 1.0/angular_resolution_x_ - provided for better performance of multiplication compared to division More...
 
float angular_resolution_y_reciprocal_ {0.0f}
 1.0/angular_resolution_y_ - provided for better performance of multiplication compared to division More...
 
int image_offset_x_ {0}
 
int image_offset_y_ {0}
 Position of the top left corner of the range image compared to an image of full size (360x180 degrees) More...
 
PointWithRange unobserved_point
 This point is used to be able to return a reference to a non-existing point. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from pcl::RangeImage
static float getMaxAngleSize (const Eigen::Affine3f &viewer_pose, const Eigen::Vector3f &center, float radius)
 Get the size of a certain area when seen from the given pose. More...
 
static Eigen::Vector3f getEigenVector3f (const PointWithRange &point)
 Get Eigen::Vector3f from PointWithRange. More...
 
static PCL_EXPORTS void getCoordinateFrameTransformation (RangeImage::CoordinateFrame coordinate_frame, Eigen::Affine3f &transformation)
 Get the transformation that transforms the given coordinate frame into CAMERA_FRAME. More...
 
template<typename PointCloudTypeWithViewpoints >
static Eigen::Vector3f getAverageViewPoint (const PointCloudTypeWithViewpoints &point_cloud)
 Get the average viewpoint of a point cloud where each point carries viewpoint information as vp_x, vp_y, vp_z. More...
 
static PCL_EXPORTS void extractFarRanges (const pcl::PCLPointCloud2 &point_cloud_data, PointCloud< PointWithViewpoint > &far_ranges)
 Check if the provided data includes far ranges and add them to far_ranges. More...
 
- Static Public Member Functions inherited from pcl::PointCloud< PointWithRange >
static bool concatenate (pcl::PointCloud< PointWithRange > &cloud1, const pcl::PointCloud< PointWithRange > &cloud2)
 
static bool concatenate (const pcl::PointCloud< PointWithRange > &cloud1, const pcl::PointCloud< PointWithRange > &cloud2, pcl::PointCloud< PointWithRange > &cloud_out)
 
- Public Attributes inherited from pcl::PointCloud< PointWithRange >
pcl::PCLHeader header
 The point cloud header. More...
 
std::vector< PointWithRange, Eigen::aligned_allocator< PointWithRange > > points
 The point data. More...
 
std::uint32_t width
 The point cloud width (if organized as an image-structure). More...
 
std::uint32_t height
 The point cloud height (if organized as an image-structure). More...
 
bool is_dense
 True if no points are invalid (e.g., have NaN or Inf values in any of their floating point fields). More...
 
Eigen::Vector4f sensor_origin_
 Sensor acquisition pose (origin/translation). More...
 
Eigen::Quaternionf sensor_orientation_
 Sensor acquisition pose (rotation). More...
 
- Static Public Attributes inherited from pcl::RangeImage
static int max_no_of_threads
 The maximum number of openmp threads that can be used in this class. More...
 
static bool debug
 Just for... More...
 
- Static Protected Member Functions inherited from pcl::RangeImage
static void createLookupTables ()
 Create lookup tables for trigonometric functions. More...
 
static float asinLookUp (float value)
 Query the asin lookup table. More...
 
static float atan2LookUp (float y, float x)
 Query the std::atan2 lookup table. More...
 
static float cosLookUp (float value)
 Query the cos lookup table. More...
 
- Static Protected Attributes inherited from pcl::RangeImage
static const int lookup_table_size
 
static std::vector< float > asin_lookup_table
 
static std::vector< float > atan_lookup_table
 
static std::vector< float > cos_lookup_table
 

Detailed Description

RangeImagePlanar is derived from the original range image and differs from it because it's not a spherical projection, but using a projection plane (as normal cameras do), therefore being better applicable for range sensors that already provide a range image by themselves (stereo cameras, ToF-cameras), so that a conversion to point cloud and then to a spherical range image becomes unnecessary.

Author
Bastian Steder

Definition at line 51 of file range_image_planar.h.

Member Typedef Documentation

◆ BaseClass

Definition at line 55 of file range_image_planar.h.

◆ ConstPtr

Definition at line 57 of file range_image_planar.h.

◆ Ptr

Definition at line 56 of file range_image_planar.h.

Constructor & Destructor Documentation

◆ RangeImagePlanar()

PCL_EXPORTS pcl::RangeImagePlanar::RangeImagePlanar ( )

Constructor.

Referenced by getNew().

◆ ~RangeImagePlanar()

PCL_EXPORTS pcl::RangeImagePlanar::~RangeImagePlanar ( )
override

Destructor.

Member Function Documentation

◆ calculate3DPoint() [1/5]

void pcl::RangeImage::calculate3DPoint
inline

Calculate the 3D point according to the given image point and the range value at the closest pixel.

Definition at line 445 of file range_image.hpp.

◆ calculate3DPoint() [2/5]

void pcl::RangeImage::calculate3DPoint
inline

Calculate the 3D point according to the given image point and range.

Definition at line 442 of file range_image.hpp.

References pcl_lrintf.

◆ calculate3DPoint() [3/5]

void pcl::RangeImagePlanar::calculate3DPoint ( float  image_x,
float  image_y,
float  range,
Eigen::Vector3f &  point 
) const
inlineoverridevirtual

Calculate the 3D point according to the given image point and range.

Parameters
image_xthe x image position
image_ythe y image position
rangethe range
pointthe resulting 3D point
Note
Implementation according to planar range images (compared to spherical as in the original)

Reimplemented from pcl::RangeImage.

Definition at line 92 of file range_image_planar.hpp.

References center_x_, center_y_, focal_length_x_reciprocal_, focal_length_y_reciprocal_, pcl::RangeImage::image_offset_x_, pcl::RangeImage::image_offset_y_, and pcl::RangeImage::to_world_system_.

◆ calculate3DPoint() [4/5]

void pcl::RangeImage::calculate3DPoint
inline

◆ calculate3DPoint() [5/5]

void pcl::RangeImage::calculate3DPoint
inline

Calculate the 3D point according to the given image point and the range value at the closest pixel.

Definition at line 438 of file range_image.hpp.

◆ copyTo()

PCL_EXPORTS void pcl::RangeImagePlanar::copyTo ( RangeImage other) const
overridevirtual

Copy *this to other.

Derived version - also copying additional RangeImagePlanar members

Reimplemented from pcl::RangeImage.

◆ createFromPointCloudWithFixedSize()

template<typename PointCloudType >
void pcl::RangeImagePlanar::createFromPointCloudWithFixedSize ( const PointCloudType &  point_cloud,
int  di_width,
int  di_height,
float  di_center_x,
float  di_center_y,
float  di_focal_length_x,
float  di_focal_length_y,
const Eigen::Affine3f &  sensor_pose,
CoordinateFrame  coordinate_frame = CAMERA_FRAME,
float  noise_level = 0.0f,
float  min_range = 0.0f 
)

Create the image from an existing point cloud.

Parameters
point_cloudthe source point cloud
di_widththe disparity image width
di_heightthe disparity image height
di_center_xthe x-coordinate of the camera's center of projection
di_center_ythe y-coordinate of the camera's center of projection
di_focal_length_xthe camera's focal length in the horizontal direction
di_focal_length_ythe camera's focal length in the vertical direction
sensor_posethe pose of the virtual depth camera
coordinate_framethe used coordinate frame of the point cloud
noise_levelwhat is the typical noise of the sensor - is used for averaging in the z-buffer
min_rangeminimum range to consifder points

Definition at line 50 of file range_image_planar.hpp.

References center_x_, center_y_, pcl::RangeImage::doZBuffer(), focal_length_x_, focal_length_x_reciprocal_, focal_length_y_, focal_length_y_reciprocal_, pcl::RangeImage::getCoordinateFrameTransformation(), pcl::PointCloud< PointWithRange >::height, pcl::PointCloud< PointWithRange >::is_dense, pcl::PointCloud< PointWithRange >::points, pcl::RangeImage::recalculate3DPointPositions(), pcl::PointCloud< PointWithRange >::size(), pcl::RangeImage::to_range_image_system_, pcl::RangeImage::to_world_system_, pcl::RangeImage::unobserved_point, and pcl::PointCloud< PointWithRange >::width.

◆ getCenterX()

float pcl::RangeImagePlanar::getCenterX ( ) const
inline

Getter for the principal point in X.

Definition at line 201 of file range_image_planar.h.

References center_x_.

◆ getCenterY()

float pcl::RangeImagePlanar::getCenterY ( ) const
inline

Getter for the principal point in Y.

Definition at line 205 of file range_image_planar.h.

References center_y_.

◆ getFocalLengthX()

float pcl::RangeImagePlanar::getFocalLengthX ( ) const
inline

Getter for the focal length in X.

Definition at line 193 of file range_image_planar.h.

References focal_length_x_.

◆ getFocalLengthY()

float pcl::RangeImagePlanar::getFocalLengthY ( ) const
inline

Getter for the focal length in Y.

Definition at line 197 of file range_image_planar.h.

References focal_length_y_.

◆ getHalfImage()

PCL_EXPORTS void pcl::RangeImagePlanar::getHalfImage ( RangeImage half_image) const
overridevirtual

Get a range image with half the resolution.

Reimplemented from pcl::RangeImage.

◆ getImagePoint() [1/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 461 of file range_image.hpp.

◆ getImagePoint() [2/8]

void pcl::RangeImage::getImagePoint
inline

Get imagePoint from 3D point in world coordinates.

Definition at line 453 of file range_image.hpp.

◆ getImagePoint() [3/8]

void pcl::RangeImagePlanar::getImagePoint ( const Eigen::Vector3f &  point,
float &  image_x,
float &  image_y,
float &  range 
) const
inlineoverridevirtual

Calculate the image point and range from the given 3D point.

Parameters
pointthe resulting 3D point
image_xthe resulting x image position
image_ythe resulting y image position
rangethe resulting range
Note
Implementation according to planar range images (compared to spherical as in the original)

Reimplemented from pcl::RangeImage.

Definition at line 105 of file range_image_planar.hpp.

References center_x_, center_y_, focal_length_x_, focal_length_y_, pcl::RangeImage::image_offset_x_, pcl::RangeImage::image_offset_y_, and pcl::RangeImage::to_range_image_system_.

◆ getImagePoint() [4/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 465 of file range_image.hpp.

References pcl::RangeImage::getPoint().

◆ getImagePoint() [5/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 457 of file range_image.hpp.

References pcl::RangeImage::getPoint(), and pcl::RangeImage::isInImage().

◆ getImagePoint() [6/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 473 of file range_image.hpp.

◆ getImagePoint() [7/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 469 of file range_image.hpp.

◆ getImagePoint() [8/8]

void pcl::RangeImage::getImagePoint
inline

Same as above.

Definition at line 477 of file range_image.hpp.

◆ getNew()

RangeImage* pcl::RangeImagePlanar::getNew ( ) const
inlineoverridevirtual

Return a newly created RangeImagePlanar.

Reimplementation to return an image of the same type.

Reimplemented from pcl::RangeImage.

Definition at line 68 of file range_image_planar.h.

References RangeImagePlanar().

◆ getSubImage()

PCL_EXPORTS void pcl::RangeImagePlanar::getSubImage ( int  sub_image_image_offset_x,
int  sub_image_image_offset_y,
int  sub_image_width,
int  sub_image_height,
int  combine_pixels,
RangeImage sub_image 
) const
overridevirtual

Get a sub part of the complete image as a new range image.

Parameters
sub_image_image_offset_x- The x coordinate of the top left pixel of the sub image. This is always according to absolute 0,0 meaning -180°,-90° and it is already in the system of the new image, so the actual pixel used in the original image is combine_pixels* (image_offset_x-image_offset_x_)
sub_image_image_offset_y- Same as image_offset_x for the y coordinate
sub_image_width- width of the new image
sub_image_height- height of the new image
combine_pixels- shrinking factor, meaning the new angular resolution is combine_pixels times the old one
sub_image- the output image

Reimplemented from pcl::RangeImage.

◆ makeShared()

Ptr pcl::RangeImagePlanar::makeShared ( )
inline

Get a boost shared pointer of a copy of this.

Definition at line 77 of file range_image_planar.h.

◆ setDepthImage() [1/2]

PCL_EXPORTS void pcl::RangeImagePlanar::setDepthImage ( const float *  depth_image,
int  di_width,
int  di_height,
float  di_center_x,
float  di_center_y,
float  di_focal_length_x,
float  di_focal_length_y,
float  desired_angular_resolution = -1 
)

Create the image from an existing depth image.

Parameters
depth_imagethe input depth image data as float values
di_widththe disparity image width
di_heightthe disparity image height
di_center_xthe x-coordinate of the camera's center of projection
di_center_ythe y-coordinate of the camera's center of projection
di_focal_length_xthe camera's focal length in the horizontal direction
di_focal_length_ythe camera's focal length in the vertical direction
desired_angular_resolutionIf this is set, the system will skip as many pixels as necessary to get as close to this angular resolution as possible while not going over this value (the density will not be lower than this value). The value is in radians per pixel.

◆ setDepthImage() [2/2]

PCL_EXPORTS void pcl::RangeImagePlanar::setDepthImage ( const unsigned short *  depth_image,
int  di_width,
int  di_height,
float  di_center_x,
float  di_center_y,
float  di_focal_length_x,
float  di_focal_length_y,
float  desired_angular_resolution = -1 
)

Create the image from an existing depth image.

Parameters
depth_imagethe input disparity image data as short values describing millimeters
di_widththe disparity image width
di_heightthe disparity image height
di_center_xthe x-coordinate of the camera's center of projection
di_center_ythe y-coordinate of the camera's center of projection
di_focal_length_xthe camera's focal length in the horizontal direction
di_focal_length_ythe camera's focal length in the vertical direction
desired_angular_resolutionIf this is set, the system will skip as many pixels as necessary to get as close to this angular resolution as possible while not going over this value (the density will not be lower than this value). The value is in radians per pixel.

◆ setDisparityImage()

PCL_EXPORTS void pcl::RangeImagePlanar::setDisparityImage ( const float *  disparity_image,
int  di_width,
int  di_height,
float  focal_length,
float  base_line,
float  desired_angular_resolution = -1 
)

Create the image from an existing disparity image.

Parameters
disparity_imagethe input disparity image data
di_widththe disparity image width
di_heightthe disparity image height
focal_lengththe focal length of the primary camera that generated the disparity image
base_linethe baseline of the stereo pair that generated the disparity image
desired_angular_resolutionIf this is set, the system will skip as many pixels as necessary to get as close to this angular resolution as possible while not going over this value (the density will not be lower than this value). The value is in radians per pixel.

Member Data Documentation

◆ center_x_

float pcl::RangeImagePlanar::center_x_ {0.0f}
protected

◆ center_y_

float pcl::RangeImagePlanar::center_y_ {0.0f}
protected

The principle point of the image.

Definition at line 211 of file range_image_planar.h.

Referenced by calculate3DPoint(), createFromPointCloudWithFixedSize(), getCenterY(), and getImagePoint().

◆ focal_length_x_

float pcl::RangeImagePlanar::focal_length_x_ {0.0f}
protected

◆ focal_length_x_reciprocal_

float pcl::RangeImagePlanar::focal_length_x_reciprocal_ {0.0f}
protected

Definition at line 210 of file range_image_planar.h.

Referenced by calculate3DPoint(), and createFromPointCloudWithFixedSize().

◆ focal_length_y_

float pcl::RangeImagePlanar::focal_length_y_ {0.0f}
protected

The focal length of the image in pixels.

Definition at line 209 of file range_image_planar.h.

Referenced by createFromPointCloudWithFixedSize(), getFocalLengthY(), and getImagePoint().

◆ focal_length_y_reciprocal_

float pcl::RangeImagePlanar::focal_length_y_reciprocal_ {0.0f}
protected

1/focal_length -> for internal use

Definition at line 210 of file range_image_planar.h.

Referenced by calculate3DPoint(), and createFromPointCloudWithFixedSize().


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