pinocchio  3.5.0
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
geometry.hpp
1 //
2 // Copyright (c) 2021-2024 INRIA
3 //
4 
5 #ifndef __pinocchio_collision_parallel_geometry_hpp__
6 #define __pinocchio_collision_parallel_geometry_hpp__
7 
8 #include "pinocchio/multibody/pool/geometry.hpp"
9 #include "pinocchio/algorithm/geometry.hpp"
10 #include "pinocchio/collision/collision.hpp"
11 #include "pinocchio/algorithm/parallel/omp.hpp"
12 
13 #include <cstdint>
14 
15 namespace pinocchio
16 {
17 
18  inline bool computeCollisionsInParallel(
19  const size_t num_threads,
20  const GeometryModel & geom_model,
21  GeometryData & geom_data,
22  const bool stopAtFirstCollision = false)
23  {
24  bool is_colliding = false;
25 
26  set_default_omp_options(num_threads);
27  const int64_t batch_size = static_cast<int64_t>(geom_model.collisionPairs.size());
28  int64_t omp_cp_index = 0;
29  // Visual studio support OpenMP 2 that only support signed indexes in for loops
30  // See
31  // https://stackoverflow.com/questions/2820621/why-arent-unsigned-openmp-index-variables-allowed
32  // -openmp:llvm could solve this:
33  // https://learn.microsoft.com/en-us/cpp/build/reference/openmp-enable-openmp-2-0-support?view=msvc-160
34 #pragma omp parallel for schedule(dynamic)
35  for (omp_cp_index = 0; omp_cp_index < batch_size; ++omp_cp_index)
36  {
37  size_t cp_index = static_cast<size_t>(omp_cp_index);
38 
39  if (stopAtFirstCollision && is_colliding)
40  continue;
41 
42  const CollisionPair & collision_pair = geom_model.collisionPairs[cp_index];
43 
44  if (
45  geom_data.activeCollisionPairs[cp_index]
46  && !(
47  geom_model.geometryObjects[collision_pair.first].disableCollision
48  || geom_model.geometryObjects[collision_pair.second].disableCollision))
49  {
50  bool res = computeCollision(geom_model, geom_data, cp_index);
51  if (!is_colliding && res)
52  {
53  is_colliding = true;
54  geom_data.collisionPairIndex = cp_index; // first pair to be in collision
55  }
56  }
57  }
58 
59  return is_colliding;
60  }
61 
62  template<
63  typename Scalar,
64  int Options,
65  template<typename, int> class JointCollectionTpl,
66  typename ConfigVectorType>
67  bool computeCollisionsInParallel(
68  const size_t num_threads,
69  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
70  DataTpl<Scalar, Options, JointCollectionTpl> & data,
71  const GeometryModel & geom_model,
72  GeometryData & geom_data,
73  const Eigen::MatrixBase<ConfigVectorType> & q,
74  const bool stopAtFirstCollision = false)
75  {
76  updateGeometryPlacements(model, data, geom_model, geom_data, q);
77  return computeCollisionsInParallel(num_threads, geom_model, geom_data, stopAtFirstCollision);
78  }
79 
80  template<
81  typename Scalar,
82  int Options,
83  template<typename, int> class JointCollectionTpl,
84  typename ConfigVectorPool,
85  typename CollisionVectorResult>
86  void computeCollisionsInParallel(
87  const size_t num_threads,
88  GeometryPoolTpl<Scalar, Options, JointCollectionTpl> & pool,
89  const Eigen::MatrixBase<ConfigVectorPool> & q,
90  const Eigen::MatrixBase<CollisionVectorResult> & res,
91  const bool stopAtFirstCollisionInConfiguration = false,
92  const bool stopAtFirstCollisionInBatch = false)
93  {
94  typedef GeometryPoolTpl<Scalar, Options, JointCollectionTpl> Pool;
95  typedef typename Pool::Model Model;
96  typedef typename Pool::Data Data;
97  typedef typename Pool::GeometryModel GeometryModel;
98  typedef typename Pool::GeometryData GeometryData;
99  typedef typename Pool::ModelVector ModelVector;
100  typedef typename Pool::DataVector DataVector;
101  typedef typename Pool::GeometryModelVector GeometryModelVector;
102  typedef typename Pool::GeometryDataVector GeometryDataVector;
103 
104  PINOCCHIO_CHECK_INPUT_ARGUMENT(pool.size() > 0, "The pool should have at least one element");
105  PINOCCHIO_CHECK_INPUT_ARGUMENT(num_threads <= pool.size(), "The pool is too small");
106 
107  const ModelVector & models = pool.getModels();
108  const Model & model_check = models[0];
109  const GeometryModelVector & geometry_models = pool.getGeometryModels();
110  DataVector & datas = pool.getDatas();
111  GeometryDataVector & geometry_datas = pool.getGeometryDatas();
112  CollisionVectorResult & res_ = res.const_cast_derived();
113 
114  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.rows(), model_check.nq);
115  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.cols(), res.size());
116  res_.fill(false);
117 
118  set_default_omp_options(num_threads);
119  const Eigen::DenseIndex batch_size = res.size();
120 
121  if (stopAtFirstCollisionInBatch)
122  {
123  bool is_colliding = false;
124  Eigen::DenseIndex i = 0;
125 #pragma omp parallel for schedule(static)
126  for (i = 0; i < batch_size; i++)
127  {
128  if (is_colliding)
129  continue;
130 
131  const int thread_id = omp_get_thread_num();
132 
133  const Model & model = models[(size_t)thread_id];
134  Data & data = datas[(size_t)thread_id];
135  const GeometryModel & geometry_model = geometry_models[(size_t)thread_id];
136  GeometryData & geometry_data = geometry_datas[(size_t)thread_id];
137  res_[i] = computeCollisions(
138  model, data, geometry_model, geometry_data, q.col(i),
139  stopAtFirstCollisionInConfiguration);
140 
141  if (res_[i])
142  {
143  is_colliding = true;
144  }
145  }
146  }
147  else
148  {
149  Eigen::DenseIndex i = 0;
150 #pragma omp parallel for schedule(static)
151  for (i = 0; i < batch_size; i++)
152  {
153  const int thread_id = omp_get_thread_num();
154 
155  const Model & model = models[(size_t)thread_id];
156  Data & data = datas[(size_t)thread_id];
157  const GeometryModel & geometry_model = geometry_models[(size_t)thread_id];
158  GeometryData & geometry_data = geometry_datas[(size_t)thread_id];
159  res_[i] = computeCollisions(
160  model, data, geometry_model, geometry_data, q.col(i),
161  stopAtFirstCollisionInConfiguration);
162  }
163  }
164  }
165 } // namespace pinocchio
166 
167 #endif // ifndef __pinocchio_collision_parallel_geometry_hpp__
Main pinocchio namespace.
Definition: treeview.dox:11
bool computeCollision(const GeometryModel &geom_model, GeometryData &geom_data, const PairIndex pair_id, fcl::CollisionRequest &collision_request)
Compute the collision status between a SINGLE collision pair. The result is store in the collisionRes...
void updateGeometryPlacements(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const GeometryModel &geom_model, GeometryData &geom_data, const Eigen::MatrixBase< ConfigVectorType > &q)
Apply a forward kinematics and update the placement of the geometry objects.
bool computeCollisions(BroadPhaseManagerBase< BroadPhaseManagerDerived > &broadphase_manager, CollisionCallBackBase *callback)
Calls computeCollision for every active pairs of GeometryData. This function assumes that updateGeome...
Definition: broadphase.hpp:34