coal  3.0.1
Coal, The Collision Detection Library. Previously known as HPP-FCL, fork of FCL -- The Flexible Collision Library
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 
39 #ifndef COAL_BVH_MODEL_H
40 #define COAL_BVH_MODEL_H
41 
42 #include "coal/fwd.hh"
43 #include "coal/collision_object.h"
44 #include "coal/BVH/BVH_internal.h"
45 #include "coal/BV/BV_node.h"
46 
47 #include <vector>
48 #include <memory>
49 #include <iostream>
50 
51 namespace coal {
52 
55 
56 template <typename _IndexType>
57 class ConvexBaseTpl;
58 
59 template <typename BV>
60 class BVFitter;
61 template <typename BV>
62 class BVSplitter;
63 
67  public:
68  typedef typename Triangle32::IndexType IndexType;
70 
72  std::shared_ptr<std::vector<Vec3s>> vertices;
73 
75  std::shared_ptr<std::vector<Triangle32>> tri_indices;
76 
78  std::shared_ptr<std::vector<Vec3s>> prev_vertices;
79 
81  unsigned int num_tris;
82 
84  unsigned int num_vertices;
85 
88 
90  // TODO: deprecated
91  shared_ptr<ConvexType> convex;
92 
95  if (num_tris && num_vertices)
96  return BVH_MODEL_TRIANGLES;
97  else if (num_vertices)
98  return BVH_MODEL_POINTCLOUD;
99  else
100  return BVH_MODEL_UNKNOWN;
101  }
102 
105 
107  BVHModelBase(const BVHModelBase& other);
108 
110  virtual ~BVHModelBase() {}
111 
113  OBJECT_TYPE getObjectType() const { return OT_BVH; }
114 
117 
119  int beginModel(unsigned int num_tris = 0, unsigned int num_vertices = 0);
120 
122  int addVertex(const Vec3s& p);
123 
125  int addVertices(const MatrixX3s& points);
126 
128  int addTriangles(const Matrixx3i& triangles);
129 
131  int addTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
132 
134  int addSubModel(const std::vector<Vec3s>& ps,
135  const std::vector<Triangle32>& ts);
136 
138  int addSubModel(const std::vector<Vec3s>& ps);
139 
142  int endModel();
143 
147 
149  int replaceVertex(const Vec3s& p);
150 
152  int replaceTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
153 
155  int replaceSubModel(const std::vector<Vec3s>& ps);
156 
159  int endReplaceModel(bool refit = true, bool bottomup = true);
160 
165 
167  int updateVertex(const Vec3s& p);
168 
170  int updateTriangle(const Vec3s& p1, const Vec3s& p2, const Vec3s& p3);
171 
173  int updateSubModel(const std::vector<Vec3s>& ps);
174 
177  int endUpdateModel(bool refit = true, bool bottomup = true);
178 
183  void buildConvexRepresentation(bool share_memory);
184 
196  bool buildConvexHull(bool keepTriangle, const char* qhullCommand = NULL);
197 
198  virtual int memUsage(const bool msg = false) const = 0;
199 
204  virtual void makeParentRelative() = 0;
205 
206  Vec3s computeCOM() const {
207  Scalar vol = 0;
208  Vec3s com(0, 0, 0);
209  if (!(vertices.get())) {
210  std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
211  "vertices."
212  << std::endl;
213  return com;
214  }
215  const std::vector<Vec3s>& vertices_ = *vertices;
216  if (!(tri_indices.get())) {
217  std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
218  "triangles."
219  << std::endl;
220  return com;
221  }
222  const std::vector<Triangle32>& tri_indices_ = *tri_indices;
223 
224  for (unsigned int i = 0; i < num_tris; ++i) {
225  const Triangle32& tri = tri_indices_[i];
226  Scalar d_six_vol =
227  (vertices_[tri[0]].cross(vertices_[tri[1]])).dot(vertices_[tri[2]]);
228  vol += d_six_vol;
229  com += (vertices_[tri[0]] + vertices_[tri[1]] + vertices_[tri[2]]) *
230  d_six_vol;
231  }
232 
233  return com / (vol * 4);
234  }
235 
237  Scalar vol = 0;
238  if (!(vertices.get())) {
239  std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
240  "vertices."
241  << std::endl;
242  return vol;
243  }
244  const std::vector<Vec3s>& vertices_ = *vertices;
245  if (!(tri_indices.get())) {
246  std::cerr << "BVH Error in `computeCOM`! The BVHModel does not contain "
247  "triangles."
248  << std::endl;
249  return vol;
250  }
251  const std::vector<Triangle32>& tri_indices_ = *tri_indices;
252  for (unsigned int i = 0; i < num_tris; ++i) {
253  const Triangle32& tri = tri_indices_[i];
254  Scalar d_six_vol =
255  (vertices_[tri[0]].cross(vertices_[tri[1]])).dot(vertices_[tri[2]]);
256  vol += d_six_vol;
257  }
258 
259  return vol / 6;
260  }
261 
263  Matrix3s C = Matrix3s::Zero();
264 
265  Matrix3s C_canonical;
266  C_canonical << Scalar(1 / 60.0), //
267  Scalar(1 / 120.0), //
268  Scalar(1 / 120.0), //
269  Scalar(1 / 120.0), //
270  Scalar(1 / 60.0), //
271  Scalar(1 / 120.0), //
272  Scalar(1 / 120.0), //
273  Scalar(1 / 120.0), //
274  Scalar(1 / 60.0);
275 
276  if (!(vertices.get())) {
277  std::cerr << "BVH Error in `computeMomentofInertia`! The BVHModel does "
278  "not contain vertices."
279  << std::endl;
280  return C;
281  }
282  const std::vector<Vec3s>& vertices_ = *vertices;
283  if (!(vertices.get())) {
284  std::cerr << "BVH Error in `computeMomentofInertia`! The BVHModel does "
285  "not contain vertices."
286  << std::endl;
287  return C;
288  }
289  const std::vector<Triangle32>& tri_indices_ = *tri_indices;
290  for (unsigned int i = 0; i < num_tris; ++i) {
291  const Triangle32& tri = tri_indices_[i];
292  const Vec3s& v1 = vertices_[tri[0]];
293  const Vec3s& v2 = vertices_[tri[1]];
294  const Vec3s& v3 = vertices_[tri[2]];
295  Matrix3s A;
296  A << v1.transpose(), v2.transpose(), v3.transpose();
297  C += A.derived().transpose() * C_canonical * A * (v1.cross(v2)).dot(v3);
298  }
299 
300  return C.trace() * Matrix3s::Identity() - C;
301  }
302 
303  protected:
304  virtual void deleteBVs() = 0;
305  virtual bool allocateBVs() = 0;
306 
308  virtual int buildTree() = 0;
309 
311  virtual int refitTree(bool bottomup) = 0;
312 
313  unsigned int num_tris_allocated;
315  unsigned int num_vertex_updated;
316 
317  protected:
319  virtual bool isEqual(const CollisionGeometry& other) const;
320 };
321 
325 template <typename BV>
327  typedef BVHModelBase Base;
328 
329  public:
331  std::vector<BVNode<BV>, Eigen::aligned_allocator<BVNode<BV>>>;
332 
334  shared_ptr<BVSplitter<BV>> bv_splitter;
335 
337  shared_ptr<BVFitter<BV>> bv_fitter;
338 
341 
346  BVHModel(const BVHModel& other);
347 
349  virtual BVHModel<BV>* clone() const { return new BVHModel(*this); }
350 
353 
356 
358  const BVNode<BV>& getBV(unsigned int i) const {
359  assert(i < num_bvs);
360  return (*bvs)[i];
361  }
362 
364  BVNode<BV>& getBV(unsigned int i) {
365  assert(i < num_bvs);
366  return (*bvs)[i];
367  }
368 
370  unsigned int getNumBVs() const { return num_bvs; }
371 
373  NODE_TYPE getNodeType() const { return BV_UNKNOWN; }
374 
376  int memUsage(const bool msg) const;
377 
383  Matrix3s I(Matrix3s::Identity());
384  makeParentRelativeRecurse(0, I, Vec3s::Zero());
385  }
386 
387  protected:
388  void deleteBVs();
389  bool allocateBVs();
390 
391  unsigned int num_bvs_allocated;
392  std::shared_ptr<std::vector<unsigned int>> primitive_indices;
393 
395  std::shared_ptr<bv_node_vector_t> bvs;
396 
398  unsigned int num_bvs;
399 
401  int buildTree();
402 
404  int refitTree(bool bottomup);
405 
409 
413 
415  int recursiveBuildTree(int bv_id, unsigned int first_primitive,
416  unsigned int num_primitives);
417 
420 
424  void makeParentRelativeRecurse(int bv_id, Matrix3s& parent_axes,
425  const Vec3s& parent_c) {
426  bv_node_vector_t& bvs_ = *bvs;
427  if (!bvs_[static_cast<size_t>(bv_id)].isLeaf()) {
428  makeParentRelativeRecurse(bvs_[static_cast<size_t>(bv_id)].first_child,
429  parent_axes,
430  bvs_[static_cast<size_t>(bv_id)].getCenter());
431 
432  makeParentRelativeRecurse(
433  bvs_[static_cast<size_t>(bv_id)].first_child + 1, parent_axes,
434  bvs_[static_cast<size_t>(bv_id)].getCenter());
435  }
436 
437  bvs_[static_cast<size_t>(bv_id)].bv =
438  translate(bvs_[static_cast<size_t>(bv_id)].bv, -parent_c);
439  }
440 
441  private:
442  virtual bool isEqual(const CollisionGeometry& _other) const {
443  const BVHModel* other_ptr = dynamic_cast<const BVHModel*>(&_other);
444  if (other_ptr == nullptr) return false;
445  const BVHModel& other = *other_ptr;
446 
447  bool res = Base::isEqual(other);
448  if (!res) return false;
449 
450  // unsigned int other_num_primitives = 0;
451  // if(other.primitive_indices)
452  // {
453 
454  // switch(other.getModelType())
455  // {
456  // case BVH_MODEL_TRIANGLES:
457  // other_num_primitives = num_tris;
458  // break;
459  // case BVH_MODEL_POINTCLOUD:
460  // other_num_primitives = num_vertices;
461  // break;
462  // default:
463  // ;
464  // }
465  // }
466 
467  // unsigned int num_primitives = 0;
468  // if(primitive_indices)
469  // {
470  //
471  // switch(other.getModelType())
472  // {
473  // case BVH_MODEL_TRIANGLES:
474  // num_primitives = num_tris;
475  // break;
476  // case BVH_MODEL_POINTCLOUD:
477  // num_primitives = num_vertices;
478  // break;
479  // default:
480  // ;
481  // }
482  // }
483  //
484  // if(num_primitives != other_num_primitives)
485  // return false;
486  //
487  // for(int k = 0; k < num_primitives; ++k)
488  // {
489  // if(primitive_indices[k] != other.primitive_indices[k])
490  // return false;
491  // }
492 
493  if (num_bvs != other.num_bvs) return false;
494 
495  if ((!(bvs.get()) && other.bvs.get()) || (bvs.get() && !(other.bvs.get())))
496  return false;
497  if (bvs.get() && other.bvs.get()) {
498  const bv_node_vector_t& bvs_ = *bvs;
499  const bv_node_vector_t& other_bvs_ = *(other.bvs);
500  for (unsigned int k = 0; k < num_bvs; ++k) {
501  if (bvs_[k] != other_bvs_[k]) return false;
502  }
503  }
504 
505  return true;
506  }
507 };
508 
510 
511 template <>
513  const Vec3s& parent_c);
514 
515 template <>
517  const Vec3s& parent_c);
518 
519 template <>
521  Matrix3s& parent_axes,
522  const Vec3s& parent_c);
523 
525 template <>
527 
528 template <>
530 
531 template <>
533 
534 template <>
536 
537 template <>
539 
540 template <>
541 NODE_TYPE BVHModel<KDOP<16>>::getNodeType() const;
542 
543 template <>
544 NODE_TYPE BVHModel<KDOP<18>>::getNodeType() const;
545 
546 template <>
547 NODE_TYPE BVHModel<KDOP<24>>::getNodeType() const;
548 
549 } // namespace coal
550 
551 #endif
A base class describing the bounding hierarchy of a mesh model or a point cloud model (which is viewe...
Definition: BVH_model.h:66
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.
ConvexBaseTpl< IndexType > ConvexType
Definition: BVH_model.h:69
OBJECT_TYPE getObjectType() const
Get the object type: it is a BVH.
Definition: BVH_model.h:113
shared_ptr< ConvexType > convex
Convex<Triangle> representation of this object.
Definition: BVH_model.h:91
unsigned int num_vertices
Number of points.
Definition: BVH_model.h:84
virtual int refitTree(bool bottomup)=0
Refit the bounding volume hierarchy.
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.
Triangle32::IndexType IndexType
Definition: BVH_model.h:68
virtual ~BVHModelBase()
deconstruction, delete mesh data related.
Definition: BVH_model.h:110
Scalar computeVolume() const
compute the volume
Definition: BVH_model.h:236
int addSubModel(const std::vector< Vec3s > &ps, const std::vector< Triangle32 > &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:78
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:313
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:262
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:87
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:81
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.
std::shared_ptr< std::vector< Triangle32 > > tri_indices
Geometry triangle index data, will be NULL for point clouds.
Definition: BVH_model.h:75
unsigned int num_vertices_allocated
Definition: BVH_model.h:314
BVHModelType getModelType() const
Model type described by the instance.
Definition: BVH_model.h:94
std::shared_ptr< std::vector< Vec3s > > vertices
Geometry point data.
Definition: BVH_model.h:72
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:206
unsigned int num_vertex_updated
Definition: BVH_model.h:315
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:326
std::shared_ptr< bv_node_vector_t > bvs
Bounding volume hierarchy.
Definition: BVH_model.h:395
unsigned int num_bvs
Number of BV nodes in bounding volume hierarchy.
Definition: BVH_model.h:398
std::shared_ptr< std::vector< unsigned int > > primitive_indices
Definition: BVH_model.h:392
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:352
std::vector< BVNode< BV >, Eigen::aligned_allocator< BVNode< BV > >> bv_node_vector_t
Definition: BVH_model.h:331
int refitTree(bool bottomup)
Refit the bounding volume hierarchy.
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:358
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:424
BVNode< BV > & getBV(unsigned int i)
Access the bv giving the its index.
Definition: BVH_model.h:364
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:373
void makeParentRelative()
This is a special acceleration: BVH_model default stores the BV's transform in world coordinate....
Definition: BVH_model.h:382
BVHModel()
Default constructor to build an empty BVH.
shared_ptr< BVSplitter< BV > > bv_splitter
Split rule to split one BV node into two children.
Definition: BVH_model.h:334
virtual BVHModel< BV > * clone() const
Clone *this into a new BVHModel.
Definition: BVH_model.h:349
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:337
unsigned int num_bvs_allocated
Definition: BVH_model.h:391
unsigned int getNumBVs() const
Get the number of bv in the BVH.
Definition: BVH_model.h:370
BVHModel(const BVHModel &other)
Copy constructor from another BVH.
The geometry for the object for collision or distance computation.
Definition: collision_object.h:96
Triangle with 3 indices for points.
Definition: data_types.h:122
_IndexType IndexType
Definition: data_types.h:127
#define COAL_DLLAPI
Definition: config.hh:88
@ BV_UNKNOWN
Definition: collision_object.h:65
@ OT_BVH
Definition: collision_object.h:54
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
BVHBuildState
States for BVH construction empty->begun->processed ->replace_begun->processed -> ....
Definition: BVH_internal.h:49
Eigen::Matrix< Scalar, Eigen::Dynamic, 3, Eigen::RowMajor > MatrixX3s
Definition: data_types.h:75
NODE_TYPE
traversal node type: bounding volume (AABB, OBB, RSS, kIOS, OBBRSS, KDOP16, KDOP18,...
Definition: collision_object.h:64
KDOP< N > translate(const KDOP< N > &bv, const Vec3s &t)
translate the KDOP BV
Eigen::Matrix< Scalar, 3, 1 > Vec3s
Definition: data_types.h:70
double Scalar
Definition: data_types.h:68
Eigen::Matrix< Eigen::DenseIndex, Eigen::Dynamic, 3, Eigen::RowMajor > Matrixx3i
Definition: data_types.h:78
OBJECT_TYPE
object type: BVH (mesh, points), basic geometry, octree
Definition: collision_object.h:52
Eigen::Matrix< Scalar, 3, 3 > Matrix3s
Definition: data_types.h:74
bool isEqual(const Eigen::MatrixBase< Derived > &lhs, const Eigen::MatrixBase< OtherDerived > &rhs, const Scalar tol=std::numeric_limits< Scalar >::epsilon() *100)
Definition: tools.h:204
A class describing a bounding volume node. It includes the tree structure providing in BVNodeBase and...
Definition: BV_node.h:106