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> class JointCollectionTpl,
23  typename ConfigVectorPool,
24  typename CollisionVectorResult>
26  const size_t num_threads,
28  const Eigen::MatrixBase<ConfigVectorPool> & q,
29  const Eigen::MatrixBase<CollisionVectorResult> & res,
30  const bool stopAtFirstCollisionInConfiguration = false,
31  const bool stopAtFirstCollisionInBatch = false)
32  {
34  Pool;
35  typedef typename Pool::Model Model;
36  typedef typename Pool::Data Data;
37  typedef typename Pool::ModelVector ModelVector;
38  typedef typename Pool::DataVector DataVector;
39  typedef typename Pool::BroadPhaseManager BroadPhaseManager;
40  typedef typename Pool::BroadPhaseManagerVector BroadPhaseManagerVector;
41 
42  const ModelVector & models = pool.getModels();
43  const Model & model_check = models[0];
44  DataVector & datas = pool.getDatas();
45  BroadPhaseManagerVector & broadphase_managers = pool.getBroadPhaseManagers();
46  CollisionVectorResult & res_ = res.const_cast_derived();
47 
48  PINOCCHIO_CHECK_INPUT_ARGUMENT(num_threads <= pool.size(), "The pool is too small");
49  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.rows(), model_check.nq);
50  PINOCCHIO_CHECK_ARGUMENT_SIZE(q.cols(), res.size());
51  res_.fill(false);
52 
54  const Eigen::DenseIndex batch_size = res.size();
55 
56  if (stopAtFirstCollisionInBatch)
57  {
58  bool is_colliding = false;
59  Eigen::DenseIndex i = 0;
60 #pragma omp parallel for schedule(static)
61  for (i = 0; i < batch_size; i++)
62  {
63  if (is_colliding)
64  continue;
65 
66  const int thread_id = omp_get_thread_num();
67  const Model & model = models[(size_t)thread_id];
68  Data & data = datas[(size_t)thread_id];
69  BroadPhaseManager & manager = broadphase_managers[(size_t)thread_id];
70  res_[i] =
71  computeCollisions(model, data, manager, q.col(i), stopAtFirstCollisionInConfiguration);
72 
73  if (res_[i])
74  {
75  is_colliding = true;
76  }
77  }
78  }
79  else
80  {
81  Eigen::DenseIndex i = 0;
82 #pragma omp parallel for schedule(static)
83  for (i = 0; i < batch_size; i++)
84  {
85  const int thread_id = omp_get_thread_num();
86  const Model & model = models[(size_t)thread_id];
87  Data & data = datas[(size_t)thread_id];
88  BroadPhaseManager & manager = broadphase_managers[(size_t)thread_id];
89  res_[i] =
90  computeCollisions(model, data, manager, q.col(i), stopAtFirstCollisionInConfiguration);
91  }
92  }
93  }
94 
99  template<
100  typename BroadPhaseManagerDerived,
101  typename Scalar,
102  int Options,
103  template<typename, int> class JointCollectionTpl>
105  const size_t num_threads,
107  const std::vector<Eigen::MatrixXd> & trajectories,
108  std::vector<VectorXb> & res,
109  const bool stopAtFirstCollisionInTrajectory = false)
110  {
112  Pool;
113  typedef typename Pool::Model Model;
114  typedef typename Pool::Data Data;
115  typedef typename Pool::ModelVector ModelVector;
116  typedef typename Pool::DataVector DataVector;
117  typedef typename Pool::BroadPhaseManager BroadPhaseManager;
118  typedef typename Pool::BroadPhaseManagerVector BroadPhaseManagerVector;
119 
120  const ModelVector & models = pool.getModels();
121  const Model & model_check = models[0];
122  DataVector & datas = pool.getDatas();
123  BroadPhaseManagerVector & broadphase_managers = pool.getBroadPhaseManagers();
124 
125  PINOCCHIO_CHECK_INPUT_ARGUMENT(num_threads <= pool.size(), "The pool is too small");
126  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories.size(), res.size());
127 
128  for (size_t k = 0; k < trajectories.size(); ++k)
129  {
130  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories[k].cols(), res[k].size());
131  PINOCCHIO_CHECK_ARGUMENT_SIZE(trajectories[k].rows(), model_check.nq);
132  }
133 
135  const int64_t batch_size = static_cast<int64_t>(trajectories.size());
136 
137  int64_t omp_i = 0;
138  // Visual studio support OpenMP 2 that only support signed indexes in for loops
139  // See
140  // https://stackoverflow.com/questions/2820621/why-arent-unsigned-openmp-index-variables-allowed
141  // -openmp:llvm could solve this:
142  // https://learn.microsoft.com/en-us/cpp/build/reference/openmp-enable-openmp-2-0-support?view=msvc-160
143 #pragma omp parallel for schedule(static)
144  for (omp_i = 0; omp_i < batch_size; omp_i++)
145  {
146  size_t i = static_cast<size_t>(omp_i);
147  const int thread_id = omp_get_thread_num();
148  const Model & model = models[size_t(thread_id)];
149  Data & data = datas[(size_t)thread_id];
150  const Eigen::MatrixXd & current_traj = trajectories[i];
151  VectorXb & res_current_traj = res[i];
152  res_current_traj.fill(false);
153  BroadPhaseManager & manager = broadphase_managers[size_t(thread_id)];
154 
155  for (Eigen::DenseIndex col_id = 0; col_id < current_traj.cols(); ++col_id)
156  {
157  res_current_traj[col_id] =
158  computeCollisions(model, data, manager, current_traj.col(col_id), true);
159  if (res_current_traj[col_id] && stopAtFirstCollisionInTrajectory)
160  break;
161  }
162  }
163  }
164 } // namespace pinocchio
165 
166 #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:216
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:192
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:21
pinocchio::Options
Options
Definition: joint-configuration.hpp:1082
omp.hpp
inverse-kinematics.i
int i
Definition: inverse-kinematics.py:17
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:63
broadphase.hpp
geometry.hpp
pinocchio::python::context::Data
DataTpl< Scalar, Options > Data
Definition: bindings/python/context/generic.hpp:64
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:1083
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:25
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:98
pinocchio::Model
ModelTpl< context::Scalar, context::Options > Model
Definition: multibody/fwd.hpp:33
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1082
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:27


pinocchio
Author(s):
autogenerated on Fri Nov 1 2024 02:41:39