Skip to content

Commit 7fd15d0

Browse files
authored
Merge pull request #2428 from jcarpent/topic/examples
Add inverse-kinematics for translation error only
2 parents 500e6a9 + c0284f3 commit 7fd15d0

File tree

4 files changed

+121
-3
lines changed

4 files changed

+121
-3
lines changed

CHANGELOG.md

+1
Original file line numberDiff line numberDiff line change
@@ -9,6 +9,7 @@ The format is based on [Keep a Changelog](https://keepachangelog.com/en/1.0.0/).
99
### Added
1010

1111
- Default visualizer can be changed with `PINOCCHIO_VIEWER` environment variable ([#2419](https://github.com/stack-of-tasks/pinocchio/pull/2419))
12+
- Add more Python and C++ examples related to inverse kinematics with 3d tasks ([#2428](https://github.com/stack-of-tasks/pinocchio/pull/2428))
1213

1314
### Fixed
1415
- Fix linkage of Boost.Serialization on Windows ([#2400](https://github.com/stack-of-tasks/pinocchio/pull/2400))

examples/CMakeLists.txt

+5-3
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,5 @@
11
#
2-
# Copyright (c) 2015-2022 CNRS INRIA
2+
# Copyright (c) 2015-2024 CNRS INRIA
33
#
44

55
# Compute flags outside the macro to avoid recomputing it for each tests
@@ -46,6 +46,7 @@ function(ADD_PINOCCHIO_CPP_EXAMPLE EXAMPLE)
4646
endfunction()
4747

4848
add_pinocchio_cpp_example(inverse-kinematics)
49+
add_pinocchio_cpp_example(inverse-kinematics-3d)
4950
add_pinocchio_cpp_example(overview-simple)
5051
add_pinocchio_cpp_example(overview-lie)
5152
add_pinocchio_cpp_example(overview-SE3)
@@ -70,8 +71,9 @@ if(BUILD_WITH_COLLISION_SUPPORT)
7071
endif()
7172

7273
if(BUILD_PYTHON_INTERFACE)
73-
set(${PROJECT_NAME}_PYTHON_EXAMPLES inverse-kinematics overview-simple kinematics-derivatives
74-
forward-dynamics-derivatives inverse-dynamics-derivatives)
74+
set(${PROJECT_NAME}_PYTHON_EXAMPLES
75+
inverse-kinematics inverse-kinematics-3d overview-simple kinematics-derivatives
76+
forward-dynamics-derivatives inverse-dynamics-derivatives)
7577

7678
if(BUILD_WITH_URDF_SUPPORT)
7779
list(

examples/inverse-kinematics-3d.cpp

+68
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,68 @@
1+
#include <iostream>
2+
3+
#include "pinocchio/multibody/sample-models.hpp"
4+
#include "pinocchio/spatial/explog.hpp"
5+
#include "pinocchio/algorithm/kinematics.hpp"
6+
#include "pinocchio/algorithm/jacobian.hpp"
7+
#include "pinocchio/algorithm/joint-configuration.hpp"
8+
9+
int main(int /* argc */, char ** /* argv */)
10+
{
11+
pinocchio::Model model;
12+
pinocchio::buildModels::manipulator(model);
13+
pinocchio::Data data(model);
14+
15+
const int JOINT_ID = 6;
16+
const pinocchio::SE3 oMdes(Eigen::Matrix3d::Identity(), Eigen::Vector3d(1., 0., 1.));
17+
18+
Eigen::VectorXd q = pinocchio::neutral(model);
19+
const double eps = 1e-4;
20+
const int IT_MAX = 1000;
21+
const double DT = 1e-1;
22+
const double damp = 1e-12;
23+
24+
pinocchio::Data::Matrix6x joint_jacobian(6, model.nv);
25+
joint_jacobian.setZero();
26+
27+
bool success = false;
28+
Eigen::Vector3d err;
29+
Eigen::VectorXd v(model.nv);
30+
for (int i = 0;; i++)
31+
{
32+
pinocchio::forwardKinematics(model, data, q);
33+
const pinocchio::SE3 iMd = data.oMi[JOINT_ID].actInv(oMdes);
34+
err = iMd.translation(); // in joint frame
35+
if (err.norm() < eps)
36+
{
37+
success = true;
38+
break;
39+
}
40+
if (i >= IT_MAX)
41+
{
42+
success = false;
43+
break;
44+
}
45+
pinocchio::computeJointJacobian(
46+
model, data, q, JOINT_ID, joint_jacobian); // joint_jacobian expressed in the joint frame
47+
const auto J = -joint_jacobian.topRows<3>(); // Jacobian associated with the error
48+
const Eigen::Matrix3d JJt = J * J.transpose() + damp * Eigen::Matrix3d::Identity();
49+
v.noalias() = -J.transpose() * JJt.ldlt().solve(err);
50+
q = pinocchio::integrate(model, q, v * DT);
51+
if (!(i % 10))
52+
std::cout << i << ": error = " << err.transpose() << std::endl;
53+
}
54+
55+
if (success)
56+
{
57+
std::cout << "Convergence achieved!" << std::endl;
58+
}
59+
else
60+
{
61+
std::cout
62+
<< "\nWarning: the iterative algorithm has not reached convergence to the desired precision"
63+
<< std::endl;
64+
}
65+
66+
std::cout << "\nresult: " << q.transpose() << std::endl;
67+
std::cout << "\nfinal error: " << err.transpose() << std::endl;
68+
}

examples/inverse-kinematics-3d.py

+47
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,47 @@
1+
from __future__ import print_function
2+
3+
import numpy as np
4+
from numpy.linalg import norm, solve
5+
6+
import pinocchio
7+
8+
model = pinocchio.buildSampleModelManipulator()
9+
data = model.createData()
10+
11+
JOINT_ID = 6
12+
oMdes = pinocchio.SE3(np.eye(3), np.array([1.0, 0.0, 1.0]))
13+
14+
q = pinocchio.neutral(model)
15+
eps = 1e-4
16+
IT_MAX = 1000
17+
DT = 1e-1
18+
damp = 1e-12
19+
20+
it = 0
21+
while True:
22+
pinocchio.forwardKinematics(model, data, q)
23+
iMd = data.oMi[JOINT_ID].actInv(oMdes)
24+
err = iMd.translation
25+
if norm(err) < eps:
26+
success = True
27+
break
28+
if it >= IT_MAX:
29+
success = False
30+
break
31+
J = pinocchio.computeJointJacobian(model, data, q, JOINT_ID) # in joint frame
32+
J = -J[:3, :] # linear part of the Jacobian
33+
v = -J.T.dot(solve(J.dot(J.T) + damp * np.eye(3), err))
34+
q = pinocchio.integrate(model, q, v * DT)
35+
if not it % 10:
36+
print("%d: error = %s" % (it, err.T))
37+
it += 1
38+
39+
if success:
40+
print("Convergence achieved!")
41+
else:
42+
print(
43+
"\nWarning: the iterative algorithm has not reached convergence to the desired precision"
44+
)
45+
46+
print("\nresult: %s" % q.flatten().tolist())
47+
print("\nfinal error: %s" % err.T)

0 commit comments

Comments
 (0)