18 #ifndef __se3_check_hxx__ 19 #define __se3_check_hxx__ 21 #include <boost/fusion/algorithm.hpp> 22 #include <boost/foreach.hpp> 30 struct AlgoFusionChecker
32 typedef bool result_type;
35 AlgoFusionChecker(
const Model&model) : model(model) {}
37 inline bool operator()(
const bool& accumul,
const boost::fusion::void_ &)
const 41 inline bool operator()(
const bool& accumul,
const AlgorithmCheckerBase<T> & t)
const 42 {
return accumul && t.checkModel(model); }
47 inline bool ParentChecker::checkModel_impl(
const Model& model )
const 49 for( JointIndex j=1;(int)j<model.njoints;++j )
50 if( model.parents[j]>=j )
return false;
55 #if !defined(BOOST_FUSION_HAS_VARIADIC_LIST) 56 template<BOOST_PP_ENUM_PARAMS(PINOCCHIO_ALGO_CHECKER_LIST_MAX_LIST_SIZE,
class T)>
57 bool AlgorithmCheckerList<BOOST_PP_ENUM_PARAMS(PINOCCHIO_ALGO_CHECKER_LIST_MAX_LIST_SIZE,T)>::checkModel_impl(
const Model& model)
const 59 return boost::fusion::accumulate(checkerList,
true,internal::AlgoFusionChecker(model));
63 bool AlgorithmCheckerList<T...>::checkModel_impl(
const Model& model)
const 65 return boost::fusion::accumulate(checkerList,
true,internal::AlgoFusionChecker(model));
71 #define CHECK_DATA(a) if(!(a)) return false; 78 CHECK_DATA( (
int)data.
a.size() == model.
njoints );
79 CHECK_DATA( (
int)data.
a_gf.size() == model.
njoints );
80 CHECK_DATA( (
int)data.
v.size() == model.
njoints );
81 CHECK_DATA( (
int)data.
f.size() == model.
njoints );
82 CHECK_DATA( (
int)data.
oMi.size() == model.
njoints );
83 CHECK_DATA( (
int)data.
liMi.size() == model.
njoints );
84 CHECK_DATA( (
int)data.
Ycrb.size() == model.
njoints );
85 CHECK_DATA( (
int)data.
Yaba.size() == model.
njoints );
86 CHECK_DATA( (
int)data.
Fcrb.size() == model.
njoints );
88 CHECK_DATA( (
int)data.
iMf.size() == model.
njoints );
89 CHECK_DATA( (
int)data.
iMf.size() == model.
njoints );
90 CHECK_DATA( (
int)data.
com.size() == model.
njoints );
91 CHECK_DATA( (
int)data.
vcom.size() == model.
njoints );
92 CHECK_DATA( (
int)data.
acom.size() == model.
njoints );
93 CHECK_DATA( (
int)data.
mass.size() == model.
njoints );
95 CHECK_DATA( data.
tau.size() == model.
nv );
96 CHECK_DATA( data.
nle.size() == model.
nv );
97 CHECK_DATA( data.
ddq.size() == model.
nv );
98 CHECK_DATA( data.
u.size() == model.
nv );
99 CHECK_DATA( data.
M.rows() == model.
nv );
100 CHECK_DATA( data.
M.cols() == model.
nv );
101 CHECK_DATA( data.
Ag.cols() == model.
nv );
102 CHECK_DATA( data.
U.cols() == model.
nv );
103 CHECK_DATA( data.
U.rows() == model.
nv );
104 CHECK_DATA( data.
D.size() == model.
nv );
105 CHECK_DATA( data.
tmp.size() == model.
nv );
106 CHECK_DATA( data.
J.cols() == model.
nv );
107 CHECK_DATA( data.
Jcom.cols() == model.
nv );
109 CHECK_DATA( data.
dq_after.size() == model.
nv );
112 CHECK_DATA( (
int)data.
oMf.size() == model.
nframes );
119 for( JointIndex j=1;int(j)<model.
njoints;++j )
121 JointIndex c = (JointIndex)data.
lastChild[j];
122 CHECK_DATA((
int)c<model.
njoints);
124 for( JointIndex d=j+1;d<=c;++d )
126 CHECK_DATA( model.
parents[d]>=j );
131 for( JointIndex d=c+1;(int)d<model.
njoints;++d)
134 int row = model.
joints[j].idx_v();
139 else { CHECK_DATA(jparent.idx_v()+jparent.nv()-1 == data.
parents_fromRow[(size_t)row]); }
151 #endif // ifndef __se3_check_hxx__ Eigen::VectorXd tau
Vector of joint torques (dim model.nv).
Matrix6x Ag
Centroidal Momentum Matrix.
container::aligned_vector< Motion > a
Vector of joint accelerations expressed at the centers of the joints.
Matrix3x Jcom
Jacobien of center of mass.
container::aligned_vector< Eigen::Vector3d > vcom
Vector of subtree center of mass linear velocities expressed in the root joint of the subtree...
int nv
Dimension of the velocity vector space.
std::vector< int > nvSubtree
Dimension of the subtree motion space (for CRBA)
container::aligned_vector< SE3 > iMf
Vector of joint placements wrt to algorithm end effector.
Eigen::VectorXd u
Intermediate quantity corresponding to apparent torque [ABA].
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6x
The 6d jacobian type (temporary)
JointModelVector joints
Model of joint i, encapsulated in a JointModelAccessor.
container::aligned_vector< Inertia > Ycrb
Vector of sub-tree composite rigid body inertias, i.e. the apparent inertia of the subtree supported ...
Matrix6x J
Jacobian of joint placements.
bool checkData(const Model &model, const Data &data)
Eigen::VectorXd D
Diagonal of the joint space intertia matrix obtained by a Cholesky Decomposition. ...
container::aligned_vector< Eigen::Vector3d > acom
Vector of subtree center of mass linear accelerations expressed in the root joint of the subtree...
Eigen::VectorXd ddq
The joint accelerations computed from ABA.
Eigen::VectorXd torque_residual
Temporary corresponding to the residual torque .
std::vector< int > nvSubtree_fromRow
Subtree of the current row index (used in Cholesky Decomposition).
std::vector< JointIndex > parents
Joint parent of joint i, denoted li (li==parents[i]).
int nv(const JointModelVariant &jmodel)
Visit a JointModelVariant through JointNvVisitor to get the dimension of the joint tangent space...
container::aligned_vector< Inertia::Matrix6 > Yaba
Inertia matrix of the subtree expressed as dense matrix [ABA].
bool check() const
Run check(fusion::list) with DEFAULT_CHECKERS as argument.
container::aligned_vector< Matrix6x > Fcrb
Spatial forces set, used in CRBA and CCRBA.
Eigen::VectorXd nle
Vector of Non Linear Effects (dim model.nv). It corresponds to concatenation of the Coriolis...
container::aligned_vector< Motion > a_gf
Vector of joint accelerations due to the gravity field.
container::aligned_vector< Eigen::Vector3d > com
Vector of subtree center of mass positions expressed in the root joint of the subtree. In other words, com[j] is the CoM position of the subtree supported by joint and expressed in the joint frame . The element com[0] corresponds to the center of mass position of the whole model and expressed in the global frame.
JointDataVector joints
Vector of se3::JointData associated to the se3::JointModel stored in model, encapsulated in JointData...
Eigen::VectorXd dq_after
Generalized velocity after impact.
Eigen::MatrixXd M
The joint space inertia matrix (a square matrix of dim model.nv).
Eigen::MatrixXd U
Joint space intertia matrix square root (upper trianglular part) computed with a Cholesky Decompositi...
container::aligned_vector< SE3 > oMf
Vector of absolute operationnel frame placements (wrt the world).
container::aligned_vector< SE3 > liMi
Vector of relative joint placements (wrt the body parent).
container::aligned_vector< SE3 > oMi
Vector of absolute joint placements (wrt the world).
container::aligned_vector< Force > f
Vector of body forces expressed in the local frame of the joint. For each body, the force represents ...
std::vector< int > parents_fromRow
First previous non-zero row in M (used in Cholesky Decomposition).
Eigen::VectorXd tmp
Temporary of size NV used in Cholesky Decomposition.
container::aligned_vector< Motion > v
Vector of joint velocities expressed at the centers of the joints.
int njoints
Number of joints.
std::vector< double > mass
Vector of subtree mass. In other words, mass[j] is the mass of the subtree supported by joint ...
std::vector< int > lastChild
Index of the last child (for CRBA)
int nframes
Number of operational frames.