5 #ifndef __pinocchio_collision_parallel_geometry_hpp__
6 #define __pinocchio_collision_parallel_geometry_hpp__
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"
18 inline bool computeCollisionsInParallel(
19 const size_t num_threads,
20 const GeometryModel & geom_model,
21 GeometryData & geom_data,
22 const bool stopAtFirstCollision =
false)
24 bool is_colliding =
false;
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;
34 #pragma omp parallel for schedule(dynamic)
35 for (omp_cp_index = 0; omp_cp_index < batch_size; ++omp_cp_index)
37 size_t cp_index =
static_cast<size_t>(omp_cp_index);
39 if (stopAtFirstCollision && is_colliding)
42 const CollisionPair & collision_pair = geom_model.collisionPairs[cp_index];
45 geom_data.activeCollisionPairs[cp_index]
47 geom_model.geometryObjects[collision_pair.first].disableCollision
48 || geom_model.geometryObjects[collision_pair.second].disableCollision))
51 if (!is_colliding && res)
54 geom_data.collisionPairIndex = cp_index;
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)
77 return computeCollisionsInParallel(num_threads, geom_model, geom_data, stopAtFirstCollision);
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)
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;
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");
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();
114 PINOCCHIO_CHECK_ARGUMENT_SIZE(q.rows(), model_check.nq);
115 PINOCCHIO_CHECK_ARGUMENT_SIZE(q.cols(), res.size());
118 set_default_omp_options(num_threads);
119 const Eigen::DenseIndex batch_size = res.size();
121 if (stopAtFirstCollisionInBatch)
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++)
131 const int thread_id = omp_get_thread_num();
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];
138 model, data, geometry_model, geometry_data, q.col(i),
139 stopAtFirstCollisionInConfiguration);
149 Eigen::DenseIndex i = 0;
150 #pragma omp parallel for schedule(static)
151 for (i = 0; i < batch_size; i++)
153 const int thread_id = omp_get_thread_num();
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];
160 model, data, geometry_model, geometry_data, q.col(i),
161 stopAtFirstCollisionInConfiguration);
Main pinocchio namespace.
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...