pinocchio  3.7.0
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
 
Loading...
Searching...
No Matches
broadphase.hpp
1//
2// Copyright (c) 2022-2024 INRIA
3//
4
5#ifndef __pinocchio_collision_parallel_broadphase_hpp__
6#define __pinocchio_collision_parallel_broadphase_hpp__
7
8#include "pinocchio/collision/pool/broadphase-manager.hpp"
9#include "pinocchio/collision/broadphase.hpp"
10#include "pinocchio/algorithm/geometry.hpp"
11#include "pinocchio/algorithm/parallel/omp.hpp"
12
13#include <cstdint>
14
15namespace 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>
25 void computeCollisionsInParallel(
26 const size_t num_threads,
27 BroadPhaseManagerPoolBase<BroadPhaseManagerDerived, Scalar, Options, JointCollectionTpl> & pool,
28 const Eigen::MatrixBase<ConfigVectorPool> & q,
29 const Eigen::MatrixBase<CollisionVectorResult> & res,
30 const bool stopAtFirstCollisionInConfiguration = false,
31 const bool stopAtFirstCollisionInBatch = false)
32 {
33 typedef BroadPhaseManagerPoolBase<BroadPhaseManagerDerived, Scalar, Options, JointCollectionTpl>
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
53 set_default_omp_options(num_threads);
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>
104 void computeCollisionsInParallel(
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
134 set_default_omp_options(num_threads);
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 {
158 computeCollisions(model, data, manager, current_traj.col(col_id), true);
160 break;
161 }
162 }
163 }
164} // namespace pinocchio
165
166#endif // ifndef __pinocchio_collision_parallel_broadphase_hpp__
Main pinocchio namespace.
Definition treeview.dox:11
bool computeCollisions(BroadPhaseManagerBase< BroadPhaseManagerDerived > &broadphase_manager, CollisionCallBackBase *callback)
Calls computeCollision for every active pairs of GeometryData. This function assumes that updateGeome...