Go to the documentation of this file.
24 int main(
int argc,
const char ** argv)
26 using namespace Eigen;
31 const int NBT = 1000 * 100;
34 std::cout <<
"(the time score in debug mode is not relevant) " << std::endl;
47 const std::string ff_option = argv[2];
48 if (ff_option ==
"-no-ff")
60 const std::string RA =
"RARM_LINK6";
62 const std::string LA =
"LARM_LINK6";
64 const std::string RF =
"RLEG_LINK6";
66 const std::string LF =
"LLEG_LINK6";
84 contact_models_6D.push_back(ci_RF_6D);
90 contact_models_6D6D.push_back(ci_RF_6D);
91 contact_models_6D6D.push_back(ci_LF_6D);
101 std::cout <<
"nq = " <<
model.nq << std::endl;
102 std::cout <<
"nv = " <<
model.nv << std::endl;
103 std::cout <<
"--" << std::endl;
106 VectorXd qmax = Eigen::VectorXd::Ones(
model.nq);
115 for (
size_t i = 0;
i < NBT; ++
i)
118 qdots[
i] = Eigen::VectorXd::Random(
model.nv);
119 qddots[
i] = Eigen::VectorXd::Random(
model.nv);
120 taus[
i] = Eigen::VectorXd::Random(
model.nv);
128 std::cout <<
"ABA = \t\t";
129 timer.
toc(std::cout, NBT);
135 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_empty,
138 std::cout <<
"contact ABA = \t\t";
139 timer.
toc(std::cout, NBT);
146 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_empty,
149 std::cout <<
"PV = \t\t";
150 timer.
toc(std::cout, NBT);
156 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_empty,
159 std::cout <<
"constrainedABA = \t\t";
160 timer.
toc(std::cout, NBT);
162 double total_time = 0;
170 std::cout <<
"Sparse Cholesky = \t\t" << (total_time / NBT) <<
" "
178 contact_chol_empty.compute(
model,
data, contact_models_empty, contact_data_empty);
181 std::cout <<
"contactCholesky {} = \t\t" << (total_time / NBT) <<
" "
185 MatrixXd H_inverse(contact_chol_empty.size(), contact_chol_empty.size());
189 contact_chol_empty.compute(
model,
data, contact_models_empty, contact_data_empty);
191 contact_chol_empty.inverse(H_inverse);
194 std::cout <<
"contactCholeskyInverse {} = \t\t" << (total_time / NBT) <<
" "
202 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_empty,
205 std::cout <<
"constraintDynamics {} = \t\t";
206 timer.
toc(std::cout, NBT);
208 std::cout <<
"--" << std::endl;
214 contact_chol_6D.compute(
model,
data, contact_models_6D, contact_data_6D);
217 std::cout <<
"contactCholesky {6D} = \t\t" << (total_time / NBT) <<
" "
221 H_inverse.resize(contact_chol_6D.size(), contact_chol_6D.size());
225 contact_chol_6D.compute(
model,
data, contact_models_6D, contact_data_6D);
227 contact_chol_6D.inverse(H_inverse);
230 std::cout <<
"contactCholeskyInverse {6D} = \t\t" << (total_time / NBT) <<
" "
233 MatrixXd
J(contact_chol_6D.constraintDim(),
model.nv);
236 model.nv + contact_chol_6D.constraintDim(),
model.nv + contact_chol_6D.constraintDim());
239 VectorXd gamma(contact_chol_6D.constraintDim());
257 std::cout <<
"KKTContactDynamicMatrixInverse {6D} = \t\t" << (total_time / NBT) <<
" "
265 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D, contact_data_6D);
267 std::cout <<
"constraintDynamics {6D} = \t\t";
268 timer.
toc(std::cout, NBT);
273 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D, contact_data_6D,
276 std::cout <<
"contact ABA {6D} = \t\t";
277 timer.
toc(std::cout, NBT);
284 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D, contact_data_6D,
287 std::cout <<
"PV = \t\t";
288 timer.
toc(std::cout, NBT);
294 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D, contact_data_6D,
297 std::cout <<
"constrainedABA = \t\t";
298 timer.
toc(std::cout, NBT);
299 std::cout <<
"--" << std::endl;
306 contact_chol_6D6D.compute(
model,
data, contact_models_6D6D, contact_data_6D6D);
309 std::cout <<
"contactCholesky {6D,6D} = \t\t" << (total_time / NBT) <<
" "
313 H_inverse.resize(contact_chol_6D6D.size(), contact_chol_6D6D.size());
317 contact_chol_6D6D.compute(
model,
data, contact_models_6D6D, contact_data_6D6D);
319 contact_chol_6D6D.inverse(H_inverse);
322 std::cout <<
"contactCholeskyInverse {6D,6D} = \t\t" << (total_time / NBT) <<
" "
325 J.resize(contact_chol_6D6D.constraintDim(),
model.nv);
328 model.nv + contact_chol_6D6D.constraintDim(),
model.nv + contact_chol_6D6D.constraintDim());
331 gamma.resize(contact_chol_6D6D.constraintDim());
348 std::cout <<
"KKTContactDynamicMatrixInverse {6D,6D} = \t\t" << (total_time / NBT) <<
" "
356 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D6D,
359 std::cout <<
"constraintDynamics {6D,6D} = \t\t";
360 timer.
toc(std::cout, NBT);
365 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D6D,
368 std::cout <<
"contact ABA {6D,6D} = \t\t";
369 timer.
toc(std::cout, NBT);
376 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D6D,
379 std::cout <<
"PV = \t\t";
380 timer.
toc(std::cout, NBT);
386 model,
data,
qs[_smooth], qdots[_smooth], taus[_smooth], contact_models_6D6D,
389 std::cout <<
"constrainedABA = \t\t";
390 timer.
toc(std::cout, NBT);
401 std::cout <<
"forwardDynamics {6D,6D} = \t\t";
402 timer.
toc(std::cout, NBT);
404 std::cout <<
"--" << std::endl;
RigidConstraintDataTpl< context::Scalar, context::Options > RigidConstraintData
PINOCCHIO_DEPRECATED void getKKTContactDynamicMatrixInverse(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< KKTMatrixType > &KKTMatrix_inv)
Computes the inverse of the KKT matrix for dynamics with contact constraints.
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & constrainedABA(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ContactModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ContactDataAllocator > &contact_datas, ProximalSettingsTpl< Scalar > &settings)
The constrained Articulated Body Algorithm (constrainedABA). It computes constrained forward dynamics...
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & pv(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ContactModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ContactDataAllocator > &contact_datas, ProximalSettingsTpl< Scalar > &settings)
The Popov-Vereshchagin algorithm. It computes constrained forward dynamics, aka the joint acceleratio...
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
JointModelFreeFlyerTpl< Scalar > JointModelFreeFlyer
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & crba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Convention convention=Convention::LOCAL)
Computes the upper triangular part of the joint space inertia matrix M by using the Composite Rigid B...
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & aba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Convention convention=Convention::LOCAL)
The Articulated-Body algorithm. It computes the forward dynamics, aka the joint accelerations given t...
static std::string unitName(Unit u)
void initPvSolver(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const std::vector< RigidConstraintModelTpl< Scalar, Options >, Allocator > &contact_models)
Init the data according to the contact information contained in contact_models.
void randomConfiguration(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &lowerLimits, const Eigen::MatrixBase< ConfigVectorIn2 > &upperLimits, const Eigen::MatrixBase< ReturnType > &qout)
Generate a configuration vector uniformly sampled among provided limits.
@ CONTACT_6D
Point contact model.
Structure containing all the settings parameters for the proximal algorithms.
void getJointJacobian(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const JointIndex joint_id, const ReferenceFrame reference_frame, const Eigen::MatrixBase< Matrix6Like > &J)
Computes the Jacobian of a specific joint frame expressed in one of the pinocchio::ReferenceFrame opt...
const PINOCCHIO_DEPRECATED DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & forwardDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< DriftVectorType > &gamma, const Scalar inv_damping=0.)
Compute the forward dynamics with contact constraints. Internally, pinocchio::computeAllTerms is call...
void initConstraintDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const std::vector< RigidConstraintModelTpl< Scalar, Options >, Allocator > &contact_models)
Init the forward dynamics data according to the contact information contained in contact_models.
ModelTpl< Scalar, Options, JointCollectionTpl > & buildModel(const std::string &filename, const typename ModelTpl< Scalar, Options, JointCollectionTpl >::JointModel &rootJoint, const std::string &rootJointName, ModelTpl< Scalar, Options, JointCollectionTpl > &model, const bool verbose=false)
Build the model from a URDF file with a particular joint as root of the model tree inside the model g...
VectorXd constraintDynamics(const Model &model, Data &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Model::JointIndex id)
#define PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(T)
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & decompose(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Compute the Cholesky decomposition of the joint space inertia matrix M contained in data.
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & contactABA(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, DataAllocator > &contact_data)
void computeAllTerms(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes efficiently all the terms needed for dynamic simulation. It is equivalent to the call at the...
Main pinocchio namespace.
#define PINOCCHIO_MODEL_DIR
pinocchio
Author(s):
autogenerated on Sun Dec 22 2024 03:41:12