Skip to content

Commit

Permalink
computeTriangleMeshArea() : calculate the area of a triangle mesh
Browse files Browse the repository at this point in the history
  • Loading branch information
mynickmynick committed Sep 9, 2023
1 parent f2492bd commit 12d02bb
Show file tree
Hide file tree
Showing 2 changed files with 98 additions and 0 deletions.
45 changes: 45 additions & 0 deletions surface/include/pcl/surface/gp3.h
Original file line number Diff line number Diff line change
Expand Up @@ -513,6 +513,51 @@ namespace pcl
}
};

/** \brief Calculate the area of a triangle mesh
* \param[in] cloud the point cloud to which the mesh is applied
* \param[in] triangleMesh a triangle mesh
* \return the mesh area
* \note example of use:
* pcl::GreedyProjectionTriangulation<PointN> gp;
* ...
* gp.setInputCloud(cloud_with_normals)
* ...
* pcl::PolygonMesh pm;
* gp3.reconstruct(pm);
* calculatePolygonArea(cloud_with_normals, pm.vertices);
* gp.reconstruct(pm);
* float functArea=computeTriangleMeshArea(cloud_with_normals, pm.polygons);
* \ingroup
* common
*/
template <typename PointT>
inline float
computeTriangleMeshArea(const shared_ptr<pcl::PointCloud<PointT>>& cloud,
std::vector<pcl::Vertices>& triangleMesh);

/** \brief Calculate the area of a triangle mesh
* \param[in] cloud the point cloud to which the mesh is applied
* \param[in] indices of the point cloud
* \param[in] triangleMesh a triangle mesh
* \return the mesh area
* \note example of use:
* pcl::GreedyProjectionTriangulation<PointN> gp;
* ...
* gp.setInputCloud(cloud_with_normals)
* ...
* pcl::PolygonMesh pm;
* gp.reconstruct(pm);
* float functArea=computeTriangleMeshArea(cloud_with_normals, pm.polygons);
* \ingroup
* common
*/

template <typename PointT>
inline float
computeTriangleMeshArea(const shared_ptr<pcl::PointCloud<PointT>>& cloud,
const shared_ptr<Indices>& indices,
std::vector<pcl::Vertices>& triangleMesh);

} // namespace pcl

#ifdef PCL_NO_PRECOMPILE
Expand Down
53 changes: 53 additions & 0 deletions surface/include/pcl/surface/impl/gp3.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1667,6 +1667,59 @@ pcl::GreedyProjectionTriangulation<PointInT>::getTriangleList (const pcl::Polygo
#define PCL_INSTANTIATE_GreedyProjectionTriangulation(T) \
template class PCL_EXPORTS pcl::GreedyProjectionTriangulation<T>;

template <typename PointT>
inline float
pcl::computeTriangleMeshArea(const shared_ptr<pcl::PointCloud<PointT>>& cloud,
std::vector<pcl::Vertices>& triangleMesh)
{
const pcl::PointCloud<PointT>::ConstPtr input_ = cloud;
double area = 0;
for (auto& triangle_ : triangleMesh) {
if (triangle_.vertices.size() == 3) {
const Eigen::Matrix<double, 3, 1> P(
(*input_)[triangle_.vertices[0]].x - (*input_)[triangle_.vertices[2]].x,
(*input_)[triangle_.vertices[0]].y - (*input_)[triangle_.vertices[2]].y,
(*input_)[triangle_.vertices[0]].z - (*input_)[triangle_.vertices[2]].z);
const Eigen::Matrix<double, 3, 1> Q(
(*input_)[triangle_.vertices[1]].x - (*input_)[triangle_.vertices[2]].x,
(*input_)[triangle_.vertices[1]].y - (*input_)[triangle_.vertices[2]].y,
(*input_)[triangle_.vertices[1]].z - (*input_)[triangle_.vertices[2]].z);
area += 0.5 * P.cross(Q).norm();
}
}
return area;
}

template <typename PointT>
inline float
pcl::computeTriangleMeshArea(const shared_ptr<pcl::PointCloud<PointT>>& cloud,
const shared_ptr<Indices>& indices,
std::vector<pcl::Vertices>& triangleMesh)
{
const pcl::PointCloud<PointT>::ConstPtr input_ = cloud;
double area = 0;
for (auto& triangle_ : triangleMesh) {
if (triangle_.vertices.size() == 3) {
const Eigen::Matrix<double, 3, 1> P(
(*input_)[(*indices_)[triangle_.vertices[0]]].x -
(*input_)[(*indices_)[triangle_.vertices[2]]].x,
(*input_)[(*indices_)[triangle_.vertices[0]]].y -
(*input_)[(*indices_)[triangle_.vertices[2]]].y,
(*input_)[(*indices_)[triangle_.vertices[0]]].z -
(*input_)[(*indices_)[triangle_.vertices[2]]].z);
const Eigen::Matrix<double, 3, 1> Q(
(*input_)[(*indices_)[triangle_.vertices[1]]].x -
(*input_)[(*indices_)[triangle_.vertices[2]]].x,
(*input_)[(*indices_)[triangle_.vertices[1]]].y -
(*input_)[(*indices_)[triangle_.vertices[2]]].y,
(*input_)[(*indices_)[triangle_.vertices[1]]].z -
(*input_)[(*indices_)[triangle_.vertices[2]]].z);
area += 0.5 * P.cross(Q).norm();
}
}
return area;
}

#endif // PCL_SURFACE_IMPL_GP3_H_


0 comments on commit 12d02bb

Please sign in to comment.