parallel/broadphase.hpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2022-2024 INRIA
3 //
4 
5 #ifndef __pinocchio_collision_parallel_broadphase_hpp__
6 #define __pinocchio_collision_parallel_broadphase_hpp__
7 
12 
13 #include <cstdint>
14 
15 namespace pinocchio
16 {
17 
18  template<
19  typename BroadPhaseManagerDerived,
20  typename Scalar,
21  int Options,
22  template<typename, int>
23  class JointCollectionTpl,
24  typename ConfigVectorPool,
25  typename CollisionVectorResult>
27  const size_t num_threads,
29  const Eigen::MatrixBase<ConfigVectorPool> & q,
30  const Eigen::MatrixBase<CollisionVectorResult> & res,
31  const bool stopAtFirstCollisionInConfiguration = false,
32  const bool stopAtFirstCollisionInBatch = false)
33  {
35  Pool;
36  typedef typename Pool::Model Model;
37  typedef typename Pool::Data Data;
38  typedef typename Pool::ModelVector ModelVector;
39  typedef typename Pool::DataVector DataVector;
40  typedef typename Pool::BroadPhaseManager BroadPhaseManager;
41  typedef typename Pool::BroadPhaseManagerVector BroadPhaseManagerVector;
42 
43  const ModelVector & models = pool.getModels();
44  const Model & model_check = models[0];
45  DataVector & datas = pool.getDatas();
46  BroadPhaseManagerVector & broadphase_managers = pool.getBroadPhaseManagers();
47  CollisionVectorResult & res_ = res.const_cast_derived();
48 
49  PINOCCHIO_CHECK_INPUT_ARGUMENT(num_threads <= pool.size(), "The pool is too small");
50  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.rows(), model_check.nq);
51  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.cols(), res.size());
52  res_.fill(false);
53 
55  const Eigen::DenseIndex batch_size = res.size();
56 
57  if (stopAtFirstCollisionInBatch)
58  {
59  bool is_colliding = false;
60  Eigen::DenseIndex i = 0;
61 #pragma omp parallel for schedule(static)
62  for (i = 0; i < batch_size; i++)
63  {
64  if (is_colliding)
65  continue;
66 
67  const int thread_id = omp_get_thread_num();
68  const Model & model = models[(size_t)thread_id];
69  Data & data = datas[(size_t)thread_id];
70  BroadPhaseManager & manager = broadphase_managers[(size_t)thread_id];
71  res_[i] =
72  computeCollisions(model, data, manager, q.col(i), stopAtFirstCollisionInConfiguration);
73 
74  if (res_[i])
75  {
76  is_colliding = true;
77  }
78  }
79  }
80  else
81  {
82  Eigen::DenseIndex i = 0;
83 #pragma omp parallel for schedule(static)
84  for (i = 0; i < batch_size; i++)
85  {
86  const int thread_id = omp_get_thread_num();
87  const Model & model = models[(size_t)thread_id];
88  Data & data = datas[(size_t)thread_id];
89  BroadPhaseManager & manager = broadphase_managers[(size_t)thread_id];
90  res_[i] =
91  computeCollisions(model, data, manager, q.col(i), stopAtFirstCollisionInConfiguration);
92  }
93  }
94  }
95 
100  template<
101  typename BroadPhaseManagerDerived,
102  typename Scalar,
103  int Options,
104  template<typename, int>
105  class JointCollectionTpl>
107  const size_t num_threads,
109  const std::vector<Eigen::MatrixXd> & trajectories,
110  std::vector<VectorXb> & res,
111  const bool stopAtFirstCollisionInTrajectory = false)
112  {
114  Pool;
115  typedef typename Pool::Model Model;
116  typedef typename Pool::Data Data;
117  typedef typename Pool::ModelVector ModelVector;
118  typedef typename Pool::DataVector DataVector;
119  typedef typename Pool::BroadPhaseManager BroadPhaseManager;
120  typedef typename Pool::BroadPhaseManagerVector BroadPhaseManagerVector;
121 
122  const ModelVector & models = pool.getModels();
123  const Model & model_check = models[0];
124  DataVector & datas = pool.getDatas();
125  BroadPhaseManagerVector & broadphase_managers = pool.getBroadPhaseManagers();
126 
127  PINOCCHIO_CHECK_INPUT_ARGUMENT(num_threads <= pool.size(), "The pool is too small");
128  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories.size(), res.size());
129 
130  for (size_t k = 0; k < trajectories.size(); ++k)
131  {
132  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories[k].cols(), res[k].size());
133  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories[k].rows(), model_check.nq);
134  }
135 
137  const int64_t batch_size = static_cast<int64_t>(trajectories.size());
138 
139  int64_t omp_i = 0;
140  // Visual studio support OpenMP 2 that only support signed indexes in for loops
141  // See
142  // https://stackoverflow.com/questions/2820621/why-arent-unsigned-openmp-index-variables-allowed
143  // -openmp:llvm could solve this:
144  // https://learn.microsoft.com/en-us/cpp/build/reference/openmp-enable-openmp-2-0-support?view=msvc-160
145 #pragma omp parallel for schedule(static)
146  for (omp_i = 0; omp_i < batch_size; omp_i++)
147  {
148  size_t i = static_cast<size_t>(omp_i);
149  const int thread_id = omp_get_thread_num();
150  const Model & model = models[size_t(thread_id)];
151  Data & data = datas[(size_t)thread_id];
152  const Eigen::MatrixXd & current_traj = trajectories[i];
153  VectorXb & res_current_traj = res[i];
154  res_current_traj.fill(false);
155  BroadPhaseManager & manager = broadphase_managers[size_t(thread_id)];
156 
157  for (Eigen::DenseIndex col_id = 0; col_id < current_traj.cols(); ++col_id)
158  {
159  res_current_traj[col_id] =
160  computeCollisions(model, data, manager, current_traj.col(col_id), true);
161  if (res_current_traj[col_id] && stopAtFirstCollisionInTrajectory)
162  break;
163  }
164  }
165  }
166 } // namespace pinocchio
167 
168 #endif // ifndef __pinocchio_collision_parallel_broadphase_hpp__
PINOCCHIO_CHECK_ARGUMENT_SIZE
#define PINOCCHIO_CHECK_ARGUMENT_SIZE(...)
Macro to check if the size of an element is equal to the expected size.
Definition: include/pinocchio/macros.hpp:217
pinocchio::DataTpl
Definition: context/generic.hpp:25
PINOCCHIO_CHECK_INPUT_ARGUMENT
#define PINOCCHIO_CHECK_INPUT_ARGUMENT(...)
Macro to check an assert-like condition and throw a std::invalid_argument exception (with a message) ...
Definition: include/pinocchio/macros.hpp:193
pinocchio::computeCollisions
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
pinocchio::BroadPhaseManagerPoolBase
Definition: collision/pool/broadphase-manager.hpp:22
pinocchio::Options
Options
Definition: joint-configuration.hpp:1116
omp.hpp
inverse-kinematics.i
int i
Definition: inverse-kinematics.py:20
run-algo-in-parallel.num_threads
num_threads
Definition: run-algo-in-parallel.py:10
rows
int rows
setup.data
data
Definition: cmake/cython/setup.in.py:48
pinocchio::python::Scalar
context::Scalar Scalar
Definition: admm-solver.cpp:29
pinocchio::res
ReturnType res
Definition: spatial/classic-acceleration.hpp:57
run-algo-in-parallel.batch_size
int batch_size
Definition: run-algo-in-parallel.py:11
pinocchio::python::context::Model
ModelTpl< Scalar, Options > Model
Definition: bindings/python/context/generic.hpp:61
broadphase.hpp
geometry.hpp
pinocchio::python::context::Data
DataTpl< Scalar, Options > Data
Definition: bindings/python/context/generic.hpp:62
pinocchio::Data
DataTpl< context::Scalar, context::Options > Data
Definition: multibody/fwd.hpp:34
size
FCL_REAL size() const
run-algo-in-parallel.pool
pool
Definition: run-algo-in-parallel.py:8
int64_t
int64_t
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1117
broadphase-manager.hpp
pinocchio::computeCollisionsInParallel
void computeCollisionsInParallel(const size_t num_threads, BroadPhaseManagerPoolBase< BroadPhaseManagerDerived, Scalar, Options, JointCollectionTpl > &pool, const Eigen::MatrixBase< ConfigVectorPool > &q, const Eigen::MatrixBase< CollisionVectorResult > &res, const bool stopAtFirstCollisionInConfiguration=false, const bool stopAtFirstCollisionInBatch=false)
Definition: parallel/broadphase.hpp:26
cols
int cols
pinocchio::set_default_omp_options
void set_default_omp_options(const size_t num_threads=(size_t) omp_get_max_threads())
Definition: omp.hpp:12
pinocchio::VectorXb
Eigen::Matrix< bool, Eigen::Dynamic, 1 > VectorXb
Definition: fwd.hpp:144
pinocchio::ModelTpl
Definition: context/generic.hpp:20
pinocchio::ModelTpl::nq
int nq
Dimension of the configuration vector representation.
Definition: multibody/model.hpp:99
pinocchio::Model
ModelTpl< context::Scalar, context::Options > Model
Definition: multibody/fwd.hpp:33
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1116
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:27


pinocchio
Author(s):
autogenerated on Sat Jun 22 2024 02:41:43