coal 3.0.2
Coal, The Collision Detection Library. Previously known as HPP-FCL, fork of FCL -- The Flexible Collision Library
Loading...
Searching...
No Matches
BVH_model.h
Go to the documentation of this file.
1/*
2 * Software License Agreement (BSD License)
3 *
4 * Copyright (c) 2011-2014, Willow Garage, Inc.
5 * Copyright (c) 2014-2015, Open Source Robotics Foundation
6 * Copyright (c) 2020-2022, INRIA
7 * All rights reserved.
8 *
9 * Redistribution and use in source and binary forms, with or without
10 * modification, are permitted provided that the following conditions
11 * are met:
12 *
13 * * Redistributions of source code must retain the above copyright
14 * notice, this list of conditions and the following disclaimer.
15 * * Redistributions in binary form must reproduce the above
16 * copyright notice, this list of conditions and the following
17 * disclaimer in the documentation and/or other materials provided
18 * with the distribution.
19 * * Neither the name of Open Source Robotics Foundation nor the names of its
20 * contributors may be used to endorse or promote products derived
21 * from this software without specific prior written permission.
22 *
23 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
24 * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
25 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
26 * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
27 * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
28 * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
29 * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
30 * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
31 * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
32 * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
33 * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
34 * POSSIBILITY OF SUCH DAMAGE.
35 */
36
38
39#ifndef COAL_BVH_MODEL_H
40#define COAL_BVH_MODEL_H
41
42#include "coal/fwd.hh"
45#include "coal/BV/BV_node.h"
46
47#include <vector>
48#include <memory>
49#include <iostream>
50
51namespace coal {
52
55
56class ConvexBase;
57
58template <typename BV>
59class BVFitter;
60template <typename BV>
61class BVSplitter;
62
65class COAL_DLLAPI BVHModelBase : public CollisionGeometry {
66 public:
68 std::shared_ptr<std::vector<Vec3s>> vertices;
69
71 std::shared_ptr<std::vector<Triangle>> tri_indices;
72
74 std::shared_ptr<std::vector<Vec3s>> prev_vertices;
75
77 unsigned int num_tris;
78
80 unsigned int num_vertices;
81
84
86 shared_ptr<ConvexBase> convex;
87
92 else if (num_vertices)
94 else
96 }
97
100
103
105 virtual ~BVHModelBase() {}
106
108 OBJECT_TYPE getObjectType() const { return OT_BVH; }
109
112
114 int beginModel(unsigned int num_tris = 0, unsigned int num_vertices = 0);
115
117 int addVertex(const Vec3s& p);
118
120 int addVertices(const MatrixX3s& points);
121
123 int addTriangles(const Matrixx3i& triangles);
124
126 int addTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
127
129 int addSubModel(const std::vector<Vec3s>& ps,
130 const std::vector<Triangle>& ts);
131
133 int addSubModel(const std::vector<Vec3s>& ps);
134
137 int endModel();
138
142
144 int replaceVertex(const Vec3s& p);
145
147 int replaceTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
148
150 int replaceSubModel(const std::vector<Vec3s>& ps);
151
154 int endReplaceModel(bool refit = true, bool bottomup = true);
155
160
162 int updateVertex(const Vec3s& p);
163
165 int updateTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
166
168 int updateSubModel(const std::vector<Vec3s>& ps);
169
172 int endUpdateModel(bool refit = true, bool bottomup = true);
173
178 void buildConvexRepresentation(bool share_memory);
179
191 bool buildConvexHull(bool keepTriangle, const char* qhullCommand = NULL);
192
193 virtual int memUsage(const bool msg = false) const = 0;
194
199 virtual void makeParentRelative() = 0;
200
201 Vec3s computeCOM() const {
202 CoalScalar vol = 0;
203 Vec3s com(0, 0, 0);
204 if (!(vertices.get())) {
205 std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
206 "vertices."
207 << std::endl;
208 return com;
209 }
210 const std::vector<Vec3s>& vertices_ = *vertices;
211 if (!(tri_indices.get())) {
212 std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
213 "triangles."
214 << std::endl;
215 return com;
216 }
217 const std::vector<Triangle>& tri_indices_ = *tri_indices;
218
219 for (unsigned int i = 0; i < num_tris; ++i) {
220 const Triangle& tri = tri_indices_[i];
221 CoalScalar d_six_vol =
222 (vertices_[tri[0]].cross(vertices_[tri[1]])).dot(vertices_[tri[2]]);
223 vol += d_six_vol;
224 com += (vertices_[tri[0]] + vertices_[tri[1]] + vertices_[tri[2]]) *
225 d_six_vol;
226 }
227
228 return com / (vol * 4);
229 }
230
232 CoalScalar vol = 0;
233 if (!(vertices.get())) {
234 std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
235 "vertices."
236 << std::endl;
237 return vol;
238 }
239 const std::vector<Vec3s>& vertices_ = *vertices;
240 if (!(tri_indices.get())) {
241 std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
242 "triangles."
243 << std::endl;
244 return vol;
245 }
246 const std::vector<Triangle>& tri_indices_ = *tri_indices;
247 for (unsigned int i = 0; i < num_tris; ++i) {
248 const Triangle& tri = tri_indices_[i];
249 CoalScalar d_six_vol =
250 (vertices_[tri[0]].cross(vertices_[tri[1]])).dot(vertices_[tri[2]]);
251 vol += d_six_vol;
252 }
253
254 return vol / 6;
255 }
256
258 Matrix3s C = Matrix3s::Zero();
259
260 Matrix3s C_canonical;
261 C_canonical << 1 / 60.0, 1 / 120.0, 1 / 120.0, 1 / 120.0, 1 / 60.0,
262 1 / 120.0, 1 / 120.0, 1 / 120.0, 1 / 60.0;
263
264 if (!(vertices.get())) {
265 std::cerr << "BVH Error in `computeMomentofInertia`! The BVHModel does "
266 "not contain vertices."
267 << std::endl;
268 return C;
269 }
270 const std::vector<Vec3s>& vertices_ = *vertices;
271 if (!(vertices.get())) {
272 std::cerr << "BVH Error in `computeMomentofInertia`! The BVHModel does "
273 "not contain vertices."
274 << std::endl;
275 return C;
276 }
277 const std::vector<Triangle>& tri_indices_ = *tri_indices;
278 for (unsigned int i = 0; i < num_tris; ++i) {
279 const Triangle& tri = tri_indices_[i];
280 const Vec3s& v1 = vertices_[tri[0]];
281 const Vec3s& v2 = vertices_[tri[1]];
282 const Vec3s& v3 = vertices_[tri[2]];
283 Matrix3s A;
284 A << v1.transpose(), v2.transpose(), v3.transpose();
285 C += A.derived().transpose() * C_canonical * A * (v1.cross(v2)).dot(v3);
286 }
287
288 return C.trace() * Matrix3s::Identity() - C;
289 }
290
291 protected:
292 virtual void deleteBVs() = 0;
293 virtual bool allocateBVs() = 0;
294
296 virtual int buildTree() = 0;
297
299 virtual int refitTree(bool bottomup) = 0;
300
301 unsigned int num_tris_allocated;
303 unsigned int num_vertex_updated;
304
305 protected:
307 virtual bool isEqual(const CollisionGeometry& other) const;
308};
309
313template <typename BV>
314class COAL_DLLAPI BVHModel : public BVHModelBase {
315 typedef BVHModelBase Base;
316
317 public:
319 std::vector<BVNode<BV>, Eigen::aligned_allocator<BVNode<BV>>>;
320
322 shared_ptr<BVSplitter<BV>> bv_splitter;
323
325 shared_ptr<BVFitter<BV>> bv_fitter;
326
329
334 BVHModel(const BVHModel& other);
335
337 virtual BVHModel<BV>* clone() const { return new BVHModel(*this); }
338
341
344
346 const BVNode<BV>& getBV(unsigned int i) const {
347 assert(i < num_bvs);
348 return (*bvs)[i];
349 }
350
352 BVNode<BV>& getBV(unsigned int i) {
353 assert(i < num_bvs);
354 return (*bvs)[i];
355 }
356
358 unsigned int getNumBVs() const { return num_bvs; }
359
361 NODE_TYPE getNodeType() const { return BV_UNKNOWN; }
362
364 int memUsage(const bool msg) const;
365
371 Matrix3s I(Matrix3s::Identity());
372 makeParentRelativeRecurse(0, I, Vec3s::Zero());
373 }
374
375 protected:
376 void deleteBVs();
378
379 unsigned int num_bvs_allocated;
380 std::shared_ptr<std::vector<unsigned int>> primitive_indices;
381
383 std::shared_ptr<bv_node_vector_t> bvs;
384
386 unsigned int num_bvs;
387
389 int buildTree();
390
392 int refitTree(bool bottomup);
393
397
401
403 int recursiveBuildTree(int bv_id, unsigned int first_primitive,
404 unsigned int num_primitives);
405
408
412 void makeParentRelativeRecurse(int bv_id, Matrix3s& parent_axes,
413 const Vec3s& parent_c) {
414 bv_node_vector_t& bvs_ = *bvs;
415 if (!bvs_[static_cast<size_t>(bv_id)].isLeaf()) {
416 makeParentRelativeRecurse(bvs_[static_cast<size_t>(bv_id)].first_child,
417 parent_axes,
418 bvs_[static_cast<size_t>(bv_id)].getCenter());
419
421 bvs_[static_cast<size_t>(bv_id)].first_child + 1, parent_axes,
422 bvs_[static_cast<size_t>(bv_id)].getCenter());
423 }
424
425 bvs_[static_cast<size_t>(bv_id)].bv =
426 translate(bvs_[static_cast<size_t>(bv_id)].bv, -parent_c);
427 }
428
429 private:
430 virtual bool isEqual(const CollisionGeometry& _other) const {
431 const BVHModel* other_ptr = dynamic_cast<const BVHModel*>(&_other);
432 if (other_ptr == nullptr) return false;
433 const BVHModel& other = *other_ptr;
434
435 bool res = Base::isEqual(other);
436 if (!res) return false;
437
438 // unsigned int other_num_primitives = 0;
439 // if(other.primitive_indices)
440 // {
441
442 // switch(other.getModelType())
443 // {
444 // case BVH_MODEL_TRIANGLES:
445 // other_num_primitives = num_tris;
446 // break;
447 // case BVH_MODEL_POINTCLOUD:
448 // other_num_primitives = num_vertices;
449 // break;
450 // default:
451 // ;
452 // }
453 // }
454
455 // unsigned int num_primitives = 0;
456 // if(primitive_indices)
457 // {
458 //
459 // switch(other.getModelType())
460 // {
461 // case BVH_MODEL_TRIANGLES:
462 // num_primitives = num_tris;
463 // break;
464 // case BVH_MODEL_POINTCLOUD:
465 // num_primitives = num_vertices;
466 // break;
467 // default:
468 // ;
469 // }
470 // }
471 //
472 // if(num_primitives != other_num_primitives)
473 // return false;
474 //
475 // for(int k = 0; k < num_primitives; ++k)
476 // {
477 // if(primitive_indices[k] != other.primitive_indices[k])
478 // return false;
479 // }
480
481 if (num_bvs != other.num_bvs) return false;
482
483 if ((!(bvs.get()) && other.bvs.get()) || (bvs.get() && !(other.bvs.get())))
484 return false;
485 if (bvs.get() && other.bvs.get()) {
486 const bv_node_vector_t& bvs_ = *bvs;
487 const bv_node_vector_t& other_bvs_ = *(other.bvs);
488 for (unsigned int k = 0; k < num_bvs; ++k) {
489 if (bvs_[k] != other_bvs_[k]) return false;
490 }
491 }
492
493 return true;
494 }
495};
496
498
499template <>
501 const Vec3s& parent_c);
502
503template <>
505 const Vec3s& parent_c);
506
507template <>
509 Matrix3s& parent_axes,
510 const Vec3s& parent_c);
511
513template <>
515
516template <>
518
519template <>
521
522template <>
524
525template <>
527
528template <>
530
531template <>
533
534template <>
536
537} // namespace coal
538
539#endif
The class for the default algorithm fitting a bounding volume to a set of points.
Definition BV_fitter.h:121
virtual void deleteBVs()=0
int endUpdateModel(bool refit=true, bool bottomup=true)
End BVH model update, will also refit or rebuild the bounding volume hierarchy.
OBJECT_TYPE getObjectType() const
Get the object type: it is a BVH.
Definition BVH_model.h:108
unsigned int num_vertices
Number of points.
Definition BVH_model.h:80
virtual int refitTree(bool bottomup)=0
Refit the bounding volume hierarchy.
CoalScalar computeVolume() const
compute the volume
Definition BVH_model.h:231
void computeLocalAABB()
Compute the AABB for the BVH, used for broad-phase collision.
int replaceTriangle(const Vec3s &p1, const Vec3s &p2, const Vec3s &p3)
Replace one triangle in the old BVH model.
std::shared_ptr< std::vector< Triangle > > tri_indices
Geometry triangle index data, will be NULL for point clouds.
Definition BVH_model.h:71
virtual ~BVHModelBase()
deconstruction, delete mesh data related.
Definition BVH_model.h:105
int addSubModel(const std::vector< Vec3s > &ps, const std::vector< Triangle > &ts)
Add a set of triangles in the new BVH model.
int updateSubModel(const std::vector< Vec3s > &ps)
Update a set of points in the old BVH model.
BVHModelBase(const BVHModelBase &other)
copy from another BVH
virtual int memUsage(const bool msg=false) const =0
std::shared_ptr< std::vector< Vec3s > > prev_vertices
Geometry point data in previous frame.
Definition BVH_model.h:74
int addSubModel(const std::vector< Vec3s > &ps)
Add a set of points in the new BVH model.
int beginModel(unsigned int num_tris=0, unsigned int num_vertices=0)
Begin a new BVH model.
unsigned int num_tris_allocated
Definition BVH_model.h:301
int endModel()
End BVH model construction, will build the bounding volume hierarchy.
Matrix3s computeMomentofInertia() const
compute the inertia matrix, related to the origin
Definition BVH_model.h:257
int replaceVertex(const Vec3s &p)
Replace one point in the old BVH model.
int updateTriangle(const Vec3s &p1, const Vec3s &p2, const Vec3s &p3)
Update one triangle in the old BVH model.
virtual bool allocateBVs()=0
int beginReplaceModel()
Replace the geometry information of current frame (i.e. should have the same mesh topology with the p...
int addTriangle(const Vec3s &p1, const Vec3s &p2, const Vec3s &p3)
Add one triangle in the new BVH model.
BVHBuildState build_state
The state of BVH building process.
Definition BVH_model.h:83
virtual void makeParentRelative()=0
This is a special acceleration: BVH_model default stores the BV's transform in world coordinate....
int addTriangles(const Matrixx3i &triangles)
Add triangles in the new BVH model.
void buildConvexRepresentation(bool share_memory)
Build this Convex<Triangle> representation of this model. The result is stored in attribute convex.
virtual bool isEqual(const CollisionGeometry &other) const
for ccd vertex update
bool buildConvexHull(bool keepTriangle, const char *qhullCommand=NULL)
Build a convex hull and store it in attribute convex.
virtual int buildTree()=0
Build the bounding volume hierarchy.
int addVertex(const Vec3s &p)
Add one point 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...
unsigned int num_tris
Number of triangles.
Definition BVH_model.h:77
BVHModelBase()
Constructing an empty BVH.
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:302
BVHModelType getModelType() const
Model type described by the instance.
Definition BVH_model.h:89
std::shared_ptr< std::vector< Vec3s > > vertices
Geometry point data.
Definition BVH_model.h:68
int replaceSubModel(const std::vector< Vec3s > &ps)
Replace a set of points in the old BVH model.
Vec3s computeCOM() const
compute center of mass
Definition BVH_model.h:201
shared_ptr< ConvexBase > convex
Convex<Triangle> representation of this object.
Definition BVH_model.h:86
unsigned int num_vertex_updated
Definition BVH_model.h:303
int addVertices(const MatrixX3s &points)
Add points in the new BVH model.
int updateVertex(const Vec3s &p)
Update one point in the old BVH model.
A class describing the bounding hierarchy of a mesh model or a point cloud model (which is viewed as ...
Definition BVH_model.h:314
std::shared_ptr< bv_node_vector_t > bvs
Bounding volume hierarchy.
Definition BVH_model.h:383
unsigned int num_bvs
Number of BV nodes in bounding volume hierarchy.
Definition BVH_model.h:386
BVNode< BV > & getBV(unsigned int i)
Access the bv giving the its index.
Definition BVH_model.h:352
std::shared_ptr< std::vector< unsigned int > > primitive_indices
Definition BVH_model.h:380
int refitTree_bottomup()
Refit the bounding volume hierarchy in a bottom-up way (fast but less compact)
~BVHModel()
deconstruction, delete mesh data related.
Definition BVH_model.h:340
int refitTree(bool bottomup)
Refit the bounding volume hierarchy.
int refitTree_topdown()
Refit the bounding volume hierarchy in a top-down way (slow but more compact)
void makeParentRelativeRecurse(int bv_id, Matrix3s &parent_axes, const Vec3s &parent_c)
Definition BVH_model.h:412
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:361
void makeParentRelative()
This is a special acceleration: BVH_model default stores the BV's transform in world coordinate....
Definition BVH_model.h:370
BVHModel()
Default constructor to build an empty BVH.
virtual BVHModel< BV > * clone() const
Clone *this into a new BVHModel.
Definition BVH_model.h:337
shared_ptr< BVSplitter< BV > > bv_splitter
Split rule to split one BV node into two children.
Definition BVH_model.h:322
int recursiveBuildTree(int bv_id, unsigned int first_primitive, unsigned int num_primitives)
Recursive kernel for hierarchy construction.
bool allocateBVs()
int buildTree()
Build the bounding volume hierarchy.
int recursiveRefitTree_bottomup(int bv_id)
Recursive kernel for bottomup refitting.
shared_ptr< BVFitter< BV > > bv_fitter
Fitting rule to fit a BV node to a set of geometry primitives.
Definition BVH_model.h:325
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:346
unsigned int num_bvs_allocated
Definition BVH_model.h:379
unsigned int getNumBVs() const
Get the number of bv in the BVH.
Definition BVH_model.h:358
std::vector< BVNode< BV >, Eigen::aligned_allocator< BVNode< BV > > > bv_node_vector_t
Definition BVH_model.h:318
BVHModel(const BVHModel &other)
Copy constructor from another BVH.
A class describing the split rule that splits each BV node.
Definition BV_splitter.h:58
The geometry for the object for collision or distance computation.
Definition collision_object.h:94
Base for convex polytope.
Definition geometric_shapes.h:645
Triangle with 3 indices for points.
Definition data_types.h:111
@ BV_UNKNOWN
Definition collision_object.h:65
@ OT_BVH
Definition collision_object.h:54
CollisionGeometry()
Definition collision_object.h:96
Main namespace.
Definition broadphase_bruteforce.h:44
BVHModelType
BVH model type.
Definition BVH_internal.h:79
@ BVH_MODEL_POINTCLOUD
triangle model
Definition BVH_internal.h:82
@ BVH_MODEL_TRIANGLES
unknown model type
Definition BVH_internal.h:81
@ BVH_MODEL_UNKNOWN
Definition BVH_internal.h:80
Eigen::Matrix< CoalScalar, 3, 3 > Matrix3s
Definition data_types.h:81
bool isEqual(const Eigen::MatrixBase< Derived > &lhs, const Eigen::MatrixBase< OtherDerived > &rhs, const CoalScalar tol=std::numeric_limits< CoalScalar >::epsilon() *100)
Definition tools.h:204
BVHBuildState
States for BVH construction empty->begun->processed ->replace_begun->processed -> ....
Definition BVH_internal.h:49
NODE_TYPE
traversal node type: bounding volume (AABB, OBB, RSS, kIOS, OBBRSS, KDOP16, KDOP18,...
Definition collision_object.h:64
Eigen::Matrix< Eigen::DenseIndex, Eigen::Dynamic, 3, Eigen::RowMajor > Matrixx3i
Definition data_types.h:85
OBJECT_TYPE
object type: BVH (mesh, points), basic geometry, octree
Definition collision_object.h:52
Eigen::Matrix< CoalScalar, Eigen::Dynamic, 3, Eigen::RowMajor > MatrixX3s
Definition data_types.h:82
Eigen::Matrix< CoalScalar, 3, 1 > Vec3s
Definition data_types.h:77
double CoalScalar
Definition data_types.h:76
A class describing a bounding volume node. It includes the tree structure providing in BVNodeBase and...
Definition BV_node.h:106