Go to the documentation of this file.
5 #ifndef __pinocchio_extra_reachable_workspace_hpp__
6 #define __pinocchio_extra_reachable_workspace_hpp__
12 #include "pinocchio/extra/config.hpp"
14 #ifdef PINOCCHIO_WITH_HPP_FCL
16 #endif // PINOCCHIO_WITH_HPP_FCL
61 template<
typename,
int>
class JointCollectionTpl,
62 typename ConfigVectorType>
65 const Eigen::MatrixBase<ConfigVectorType> &
q0,
66 const double time_horizon,
68 Eigen::MatrixXd & vertex,
87 template<
typename,
int>
class JointCollectionTpl,
88 typename ConfigVectorType>
91 const Eigen::MatrixBase<ConfigVectorType> &
q0,
92 const double time_horizon,
97 #ifdef PINOCCHIO_WITH_HPP_FCL
115 template<
typename,
int>
class JointCollectionTpl,
116 typename ConfigVectorType>
117 void reachableWorkspaceWithCollisions(
120 const Eigen::MatrixBase<ConfigVectorType> &
q0,
121 const double time_horizon,
123 Eigen::MatrixXd & vertex,
145 template<
typename,
int>
class JointCollectionTpl,
146 typename ConfigVectorType>
147 void reachableWorkspaceWithCollisionsHull(
150 const Eigen::MatrixBase<ConfigVectorType> &
q0,
151 const double time_horizon,
155 #endif // PINOCCHIO_WITH_HPP_FCL
178 template<
typename,
int>
class JointCollectionTpl,
179 typename ConfigVectorType,
180 class FilterFunction>
183 const Eigen::MatrixBase<ConfigVectorType> &
q0,
184 const double time_horizon,
186 FilterFunction config_filter,
187 Eigen::MatrixXd & vertex,
200 const Eigen::VectorXd & res1,
201 const Eigen::VectorXd & res2,
202 const Eigen::VectorXi & comb,
203 Eigen::VectorXd & qv);
221 const Eigen::VectorXd & element,
223 Eigen::VectorXi & indices,
224 Eigen::VectorXd & combination);
229 #include "pinocchio/extra/reachable-workspace.hxx"
231 #endif // ifndef __pinocchio_extra_reachable_workspace_hpp__
void computeVertex(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q0, const double time_horizon, const int frame_id, FilterFunction config_filter, Eigen::MatrixXd &vertex, const ReachableSetParams ¶ms=ReachableSetParams())
Samples points to create the reachable workspace that will respect mechanical limits of the model as ...
void productCombination(const Eigen::VectorXd &element, const int repeat, Eigen::VectorXi &indices, Eigen::VectorXd &combination)
Cartesian product of input element with itself. Number of repetition is specified with repeat argumen...
void computeJointVel(const Eigen::VectorXd &res1, const Eigen::VectorXd &res2, const Eigen::VectorXi &comb, Eigen::VectorXd &qv)
Computes the joint configuration associated with the permutation results stored in res1 and res2.
Structure containing the return value for the reachable algorithm.
PINOCCHIO_EXTRA_DLLAPI void buildConvexHull(ReachableSetResults &res)
Computes the convex hull using qhull associated with the vertex stored in res.
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 ¶ms=ReachableSetParams())
Computes the convex Hull reachable workspace on a fixed time horizon. For more information,...
Parameters for the reachable space algorithm.
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > & q0
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 ¶ms=ReachableSetParams())
Computes the reachable workspace on a fixed time horizon. For more information, please see https://gi...
void generateCombination(const int n, const int k, Eigen::VectorXi &indices)
Return a subsequence of length k of elements from range 0 to n. Inspired by https://docs....
JointCollectionTpl & model
Main pinocchio namespace.
pinocchio
Author(s):
autogenerated on Thu Dec 19 2024 03:41:32