5 #ifndef __pinocchio_collision_parallel_geometry_hpp__
6 #define __pinocchio_collision_parallel_geometry_hpp__
22 const bool stopAtFirstCollision =
false)
24 bool is_colliding =
false;
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)
51 if (!is_colliding &&
res)
65 template<
typename,
int>
class JointCollectionTpl,
66 typename ConfigVectorType>
73 const Eigen::MatrixBase<ConfigVectorType> &
q,
74 const bool stopAtFirstCollision =
false)
83 template<
typename,
int>
class JointCollectionTpl,
84 typename ConfigVectorPool,
85 typename CollisionVectorResult>
89 const Eigen::MatrixBase<ConfigVectorPool> &
q,
90 const Eigen::MatrixBase<CollisionVectorResult> &
res,
91 const bool stopAtFirstCollisionInConfiguration =
false,
92 const bool stopAtFirstCollisionInBatch =
false)
99 typedef typename Pool::ModelVector ModelVector;
100 typedef typename Pool::DataVector DataVector;
101 typedef typename Pool::GeometryModelVector GeometryModelVector;
102 typedef typename Pool::GeometryDataVector GeometryDataVector;
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();
121 if (stopAtFirstCollisionInBatch)
123 bool is_colliding =
false;
124 Eigen::DenseIndex
i = 0;
125 #pragma omp parallel for schedule(static)
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)
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);
167 #endif // ifndef __pinocchio_collision_parallel_geometry_hpp__