39#ifndef HPP_FCL_BVH_MODEL_H
40#define HPP_FCL_BVH_MODEL_H
88 if (num_tris && num_vertices)
90 else if (num_vertices)
105 delete[] tri_indices;
106 delete[] prev_vertices;
116 int beginModel(
unsigned int num_tris = 0,
unsigned int num_vertices = 0);
132 const std::vector<Triangle>& ts);
195 virtual int memUsage(
const bool msg =
false)
const = 0;
206 for (
unsigned int i = 0; i < num_tris; ++i) {
207 const Triangle& tri = tri_indices[i];
209 (vertices[tri[0]].cross(vertices[tri[1]])).dot(vertices[tri[2]]);
212 (vertices[tri[0]] + vertices[tri[1]] + vertices[tri[2]]) * d_six_vol;
215 return com / (vol * 4);
220 for (
unsigned int i = 0; i < num_tris; ++i) {
221 const Triangle& tri = tri_indices[i];
223 (vertices[tri[0]].cross(vertices[tri[1]])).dot(vertices[tri[2]]);
234 C_canonical << 1 / 60.0, 1 / 120.0, 1 / 120.0, 1 / 120.0, 1 / 60.0,
235 1 / 120.0, 1 / 120.0, 1 / 120.0, 1 / 60.0;
237 for (
unsigned int i = 0; i < num_tris; ++i) {
238 const Triangle& tri = tri_indices[i];
239 const Vec3f& v1 = vertices[tri[0]];
240 const Vec3f& v2 = vertices[tri[1]];
241 const Vec3f& v3 = vertices[tri[2]];
243 A << v1.transpose(), v2.transpose(), v3.transpose();
244 C += A.derived().transpose() * C_canonical * A * (v1.cross(v2)).dot(v3);
247 return C.trace() * Matrix3f::Identity() - C;
272template <
typename BV>
298 delete[] primitive_indices;
331 makeParentRelativeRecurse(0, I, Vec3f::Zero());
363 unsigned int num_primitives);
372 const Vec3f& parent_c) {
373 if (!bvs[bv_id].isLeaf()) {
374 makeParentRelativeRecurse(bvs[bv_id].first_child, parent_axes,
375 bvs[bv_id].getCenter());
377 makeParentRelativeRecurse(bvs[bv_id].first_child + 1, parent_axes,
378 bvs[bv_id].getCenter());
381 bvs[bv_id].
bv =
translate(bvs[bv_id].bv, -parent_c);
387 if (other_ptr ==
nullptr)
return false;
390 bool res = Base::isEqual(other);
391 if (!res)
return false;
436 if (num_bvs != other.
num_bvs)
return false;
438 for (
unsigned int k = 0; k < num_bvs; ++k) {
439 if (bvs[k] != other.
bvs[k])
return false;
450 const Vec3f& parent_c);
454 const Vec3f& parent_c);
459 const Vec3f& parent_c);
A base class describing the bounding hierarchy of a mesh model or a point cloud model (which is viewe...
Definition BVH_model.h:63
Vec3f * prev_vertices
Geometry point data in previous frame.
Definition BVH_model.h:72
unsigned int num_tris
Number of triangles.
Definition BVH_model.h:75
int addTriangles(const Matrixx3i &triangles)
Add triangles in the new BVH model.
int replaceTriangle(const Vec3f &p1, const Vec3f &p2, const Vec3f &p3)
Replace one triangle in the old BVH model.
int replaceSubModel(const std::vector< Vec3f > &ps)
Replace a set of points in the old BVH model.
int updateTriangle(const Vec3f &p1, const Vec3f &p2, const Vec3f &p3)
Update one triangle in the old BVH model.
Triangle * tri_indices
Geometry triangle index data, will be NULL for point clouds.
Definition BVH_model.h:69
int endReplaceModel(bool refit=true, bool bottomup=true)
End BVH model replacement, will also refit or rebuild the bounding volume hierarchy.
unsigned int num_vertices_allocated
Definition BVH_model.h:261
int updateVertex(const Vec3f &p)
Update one point in the old BVH model.
int addSubModel(const std::vector< Vec3f > &ps)
Add a set of points in the new BVH model.
int beginUpdateModel()
Replace the geometry information of current frame (i.e. should have the same mesh topology with the p...
virtual bool isEqual(const CollisionGeometry &other) const
for ccd vertex update
BVHBuildState build_state
The state of BVH building process.
Definition BVH_model.h:81
virtual void makeParentRelative()=0
This is a special acceleration: BVH_model default stores the BV's transform in world coordinate....
virtual int refitTree(bool bottomup)=0
Refit the bounding volume hierarchy.
int addSubModel(const std::vector< Vec3f > &ps, const std::vector< Triangle > &ts)
Add a set of triangles in the new BVH model.
unsigned int num_vertex_updated
Definition BVH_model.h:262
BVHModelType getModelType() const
Model type described by the instance.
Definition BVH_model.h:87
int beginReplaceModel()
Replace the geometry information of current frame (i.e. should have the same mesh topology with the p...
int addTriangle(const Vec3f &p1, const Vec3f &p2, const Vec3f &p3)
Add one triangle in the new BVH model.
virtual ~BVHModelBase()
deconstruction, delete mesh data related.
Definition BVH_model.h:103
virtual bool allocateBVs()=0
int replaceVertex(const Vec3f &p)
Replace one point in the old BVH model.
BVHModelBase(const BVHModelBase &other)
copy from another BVH
virtual int buildTree()=0
Build the bounding volume hierarchy.
Vec3f * vertices
Geometry point data.
Definition BVH_model.h:66
unsigned int num_vertices
Number of points.
Definition BVH_model.h:78
Matrix3f computeMomentofInertia() const
compute the inertia matrix, related to the origin
Definition BVH_model.h:230
unsigned int num_tris_allocated
Definition BVH_model.h:260
int addVertices(const Matrixx3f &points)
Add points in the new BVH model.
int endModel()
End BVH model construction, will build the bounding volume hierarchy.
BVHModelBase()
Constructing an empty BVH.
void computeLocalAABB()
Compute the AABB for the BVH, used for broad-phase collision.
shared_ptr< ConvexBase > convex
Convex<Triangle> representation of this object.
Definition BVH_model.h:84
int addVertex(const Vec3f &p)
Add one point in the new BVH model.
int endUpdateModel(bool refit=true, bool bottomup=true)
End BVH model update, will also refit or rebuild the bounding volume hierarchy.
virtual void deleteBVs()=0
bool buildConvexHull(bool keepTriangle, const char *qhullCommand=NULL)
Build a convex hull and store it in attribute convex.
FCL_REAL computeVolume() const
compute the volume
Definition BVH_model.h:218
OBJECT_TYPE getObjectType() const
Get the object type: it is a BVH.
Definition BVH_model.h:110
int updateSubModel(const std::vector< Vec3f > &ps)
Update a set of points in the old BVH model.
virtual int memUsage(const bool msg=false) const =0
Vec3f computeCOM() const
compute center of mass
Definition BVH_model.h:203
int beginModel(unsigned int num_tris=0, unsigned int num_vertices=0)
Begin a new BVH model.
void buildConvexRepresentation(bool share_memory)
Build this Convex representation of this model. The result is stored in attribute convex.
A class describing the bounding hierarchy of a mesh model or a point cloud model (which is viewed as ...
Definition BVH_model.h:273
~BVHModel()
deconstruction, delete mesh data related.
Definition BVH_model.h:296
int refitTree_bottomup()
Refit the bounding volume hierarchy in a bottom-up way (fast but less compact)
const BVNode< BV > & getBV(unsigned int i) const
We provide getBV() and getNumBVs() because BVH may be compressed (in future), so we must provide some...
Definition BVH_model.h:305
BVHModel(const BVHModel &other)
Copy constructor from another BVH.
virtual BVHModel< BV > * clone() const
Clone *this into a new BVHModel.
Definition BVH_model.h:293
unsigned int getNumBVs() const
Get the number of bv in the BVH.
Definition BVH_model.h:317
void makeParentRelative()
This is a special acceleration: BVH_model default stores the BV's transform in world coordinate....
Definition BVH_model.h:329
unsigned int * primitive_indices
Definition BVH_model.h:339
BVNode< BV > * bvs
Bounding volume hierarchy.
Definition BVH_model.h:342
int buildTree()
Build the bounding volume hierarchy.
shared_ptr< BVFitter< BV > > bv_fitter
Fitting rule to fit a BV node to a set of geometry primitives.
Definition BVH_model.h:281
int recursiveRefitTree_bottomup(int bv_id)
Recursive kernel for bottomup refitting.
shared_ptr< BVSplitter< BV > > bv_splitter
Split rule to split one BV node into two children.
Definition BVH_model.h:278
BVNode< BV > & getBV(unsigned int i)
Access the bv giving the its index.
Definition BVH_model.h:311
void makeParentRelativeRecurse(int bv_id, Matrix3f &parent_axes, const Vec3f &parent_c)
Definition BVH_model.h:371
unsigned int num_bvs_allocated
Definition BVH_model.h:338
BVHModel()
Default constructor to build an empty BVH.
unsigned int num_bvs
Number of BV nodes in bounding volume hierarchy.
Definition BVH_model.h:345
int memUsage(const bool msg) const
Check the number of memory used.
NODE_TYPE getNodeType() const
Get the BV type: default is unknown.
Definition BVH_model.h:320
int refitTree_topdown()
Refit the bounding volume hierarchy in a top-down way (slow but more compact)
int recursiveBuildTree(int bv_id, unsigned int first_primitive, unsigned int num_primitives)
Recursive kernel for hierarchy construction.
int refitTree(bool bottomup)
Refit the bounding volume hierarchy.
The geometry for the object for collision or distance computation.
Definition collision_object.h:95
Triangle with 3 indices for points.
Definition data_types.h:96
#define HPP_FCL_DLLAPI
Definition config.hh:88
KDOP< N > translate(const KDOP< N > &bv, const Vec3f &t)
translate the KDOP BV
BV bv
bounding volume storing the geometry
Definition BV_node.h:113
@ OT_BVH
Definition collision_object.h:55
@ BV_UNKNOWN
Definition collision_object.h:66
Eigen::Matrix< Eigen::DenseIndex, Eigen::Dynamic, 3 > Matrixx3i
Definition data_types.h:70
Eigen::Matrix< FCL_REAL, 3, 3 > Matrix3f
Definition data_types.h:68
BVHBuildState
States for BVH construction empty->begun->processed ->replace_begun->processed -> ....
Definition BVH_internal.h:50
Eigen::Matrix< FCL_REAL, 3, 1 > Vec3f
Definition data_types.h:66
Eigen::Matrix< FCL_REAL, Eigen::Dynamic, 3 > Matrixx3f
Definition data_types.h:69
OBJECT_TYPE
object type: BVH (mesh, points), basic geometry, octree
Definition collision_object.h:53
NODE_TYPE
traversal node type: bounding volume (AABB, OBB, RSS, kIOS, OBBRSS, KDOP16, KDOP18,...
Definition collision_object.h:65
BVHModelType
BVH model type.
Definition BVH_internal.h:80
@ BVH_MODEL_POINTCLOUD
triangle model
Definition BVH_internal.h:83
@ BVH_MODEL_UNKNOWN
Definition BVH_internal.h:81
@ BVH_MODEL_TRIANGLES
unknown model type
Definition BVH_internal.h:82
bool isEqual(const Eigen::MatrixBase< Derived > &lhs, const Eigen::MatrixBase< OtherDerived > &rhs, const FCL_REAL tol=std::numeric_limits< FCL_REAL >::epsilon() *100)
Definition tools.h:209
double FCL_REAL
Definition data_types.h:65
Main namespace.
Definition broadphase_bruteforce.h:44
A class describing a bounding volume node. It includes the tree structure providing in BVNodeBase and...
Definition BV_node.h:109