$darkmode
pinocchio 4.0.0
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
reachable-workspace.hpp
1 //
2 // Copyright (c) 2016-2023 CNRS INRIA
3 //
4 
5 #pragma once
6 
7 // IWYU pragma: begin_keep
8 #include <Eigen/Core>
9 #include <Eigen/SparseCore>
10 
11 // #include <vector>
12 #include <cmath>
13 
14 #include <boost/math/special_functions/factorials.hpp>
15 
16 #include "pinocchio/macros.hpp"
17 
18 #include "pinocchio/multibody.hpp"
19 #include "pinocchio/multibody/liegroup.hpp"
20 
21 #include "pinocchio/algorithm/frames.hpp"
22 
23 #include "pinocchio/extra/config.hpp"
24 
25 #ifdef PINOCCHIO_WITH_COLLISION
26  #include "pinocchio/collision/collision.hpp"
27 #endif // PINOCCHIO_WITH_COLLISION
28 // IWYU pragma: end_keep
29 
30 namespace pinocchio
31 {
36  {
37  Eigen::MatrixXd vertex;
38  Eigen::MatrixXi faces;
39  };
40 
50  {
51  int n_samples = 5;
52  int facet_dims = 3;
53  };
54 
66  template<
67  typename Scalar,
68  int Options,
69  template<typename, int> class JointCollectionTpl,
70  typename ConfigVectorType>
72  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
73  const Eigen::MatrixBase<ConfigVectorType> & q0,
74  const double time_horizon,
75  const int frame_id,
76  Eigen::MatrixXd & vertex,
77  const ReachableSetParams & params = ReachableSetParams());
78 
84 
90 
92  template<
93  typename Scalar,
94  int Options,
95  template<typename, int> class JointCollectionTpl,
96  typename ConfigVectorType>
98  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
99  const Eigen::MatrixBase<ConfigVectorType> & q0,
100  const double time_horizon,
101  const int frame_id,
102  ReachableSetResults & res,
103  const ReachableSetParams & params = ReachableSetParams());
104 
105 #ifdef PINOCCHIO_WITH_COLLISION
119  template<
120  typename Scalar,
121  int Options,
122  template<typename, int> class JointCollectionTpl,
123  typename ConfigVectorType>
125  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
126  const GeometryModel & geom_model,
127  const Eigen::MatrixBase<ConfigVectorType> & q0,
128  const double time_horizon,
129  const int frame_id,
130  Eigen::MatrixXd & vertex,
131  const ReachableSetParams & params = ReachableSetParams());
132 
147  template<
148  typename Scalar,
149  int Options,
150  template<typename, int> class JointCollectionTpl,
151  typename ConfigVectorType>
153  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
154  const GeometryModel & geom_model,
155  const Eigen::MatrixBase<ConfigVectorType> & q0,
156  const double time_horizon,
157  const int frame_id,
158  ReachableSetResults & res,
159  const ReachableSetParams & params = ReachableSetParams());
160 #endif // PINOCCHIO_WITH_COLLISION
161 
162  namespace internal
163  {
178 
180  template<
181  typename Scalar,
182  int Options,
183  template<typename, int> class JointCollectionTpl,
184  typename ConfigVectorType,
185  class FilterFunction>
186  void computeVertex(
187  const ModelTpl<Scalar, Options, JointCollectionTpl> & model,
188  const Eigen::MatrixBase<ConfigVectorType> & q0,
189  const double time_horizon,
190  const int frame_id,
191  FilterFunction config_filter,
192  Eigen::MatrixXd & vertex,
193  const ReachableSetParams & params = ReachableSetParams());
194 
197  PINOCCHIO_EXTRA_DLLAPI void buildConvexHull(ReachableSetResults & res);
198 
202 
204  inline void computeJointVel(
205  const Eigen::VectorXd & res1,
206  const Eigen::VectorXd & res2,
207  const Eigen::VectorXi & comb,
208  Eigen::VectorXd & qv);
209 
213 
215  inline void generateCombination(const int n, const int k, Eigen::VectorXi & indices);
216 
223 
225  inline void productCombination(
226  const Eigen::VectorXd & element,
227  const int repeat,
228  Eigen::VectorXi & indices,
229  Eigen::VectorXd & combination);
230  } // namespace internal
231 } // namespace pinocchio
232 
233 // IWYU pragma: begin_exports
234 #include "pinocchio/src/extra/reachable-workspace.hxx"
235 // IWYU pragma: end_exports
Main pinocchio namespace.
Definition: treeview.dox:11
void reachableWorkspaceWithCollisions(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const GeometryModel &geom_model, const Eigen::MatrixBase< ConfigVectorType > &q0, const double time_horizon, const int frame_id, Eigen::MatrixXd &vertex, const ReachableSetParams &params=ReachableSetParams())
Computes the reachable workspace with respect to a geometry model on a fixed time horizon....
void reachableWorkspaceWithCollisionsHull(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const GeometryModel &geom_model, const Eigen::MatrixBase< ConfigVectorType > &q0, const double time_horizon, const int frame_id, ReachableSetResults &res, const ReachableSetParams &params=ReachableSetParams())
Computes the convex Hull of the reachable workspace with respect to a geometry model on a fixed time ...
void reachableWorkspace(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q0, const double time_horizon, const int frame_id, Eigen::MatrixXd &vertex, const ReachableSetParams &params=ReachableSetParams())
Computes the reachable workspace on a fixed time horizon. For more information, please see https://gi...
void reachableWorkspaceHull(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q0, const double time_horizon, const int frame_id, ReachableSetResults &res, const ReachableSetParams &params=ReachableSetParams())
Computes the convex Hull reachable workspace on a fixed time horizon. For more information,...
Parameters for the reachable space algorithm.
Structure containing the return value for the reachable algorithm.