tinyrobotics
Loading...
Searching...
No Matches
kinematics.hpp
Go to the documentation of this file.
1#ifndef TR_KINEMATICS_HPP
2#define TR_KINEMATICS_HPP
3
4#include <Eigen/Core>
5#include <Eigen/Geometry>
6
7#include "math.hpp"
8#include "model.hpp"
9
13namespace tinyrobotics {
14
25 template <typename Scalar, int nq, typename TargetLink>
26 int get_link_idx(const Model<Scalar, nq>& model, const TargetLink& target_link) {
27 if constexpr (std::is_integral<TargetLink>::value) {
28 return static_cast<int>(target_link);
29 }
30 else if constexpr (std::is_same<TargetLink, std::string>::value) {
31 return model.get_link(target_link).idx;
32 }
33 else {
34 throw std::invalid_argument("Invalid target_link type. Must be int or std::string.");
35 }
36 }
37
46 template <typename Scalar, int nq>
47 std::vector<Eigen::Transform<Scalar, 3, Eigen::Isometry>> forward_kinematics(
48 Model<Scalar, nq>& model,
49 const Eigen::Matrix<Scalar, nq, 1>& q) {
50 model.forward_kinematics.resize(model.links.size(), Eigen::Transform<Scalar, 3, Eigen::Isometry>::Identity());
51 for (auto const link : model.links) {
52 model.forward_kinematics[link.idx] = link.joint.parent_transform;
53 if (link.joint.idx != -1) {
54 model.forward_kinematics[link.idx] =
55 model.forward_kinematics[link.idx] * link.joint.get_joint_transform(q[link.joint.idx]);
56 }
57 if (link.parent != -1) {
58 model.forward_kinematics[link.idx] =
59 model.forward_kinematics[link.parent] * model.forward_kinematics[link.idx];
60 }
61 }
62 return model.forward_kinematics;
63 }
64
76 template <typename Scalar, int nq, typename TargetLink>
77 Eigen::Transform<Scalar, 3, Eigen::Isometry> forward_kinematics(const Model<Scalar, nq>& model,
78 const Eigen::Matrix<Scalar, nq, 1>& q,
79 const TargetLink& target_link) {
80 Link<Scalar> current_link = model.links[get_link_idx(model, target_link)];
81
82 // Build kinematic tree from base {b} to target {t} link
83 Eigen::Transform<Scalar, 3, Eigen::Isometry> Htb = Eigen::Transform<Scalar, 3, Eigen::Isometry>::Identity();
84 while (current_link.idx != model.base_link_idx) {
85 // Apply joint transform if joint is actuatable
86 if (current_link.joint.idx != -1) {
87 Htb = Htb * current_link.joint.get_joint_transform(q[current_link.joint.idx]).inverse();
88 }
89 Htb = Htb * current_link.joint.parent_transform.inverse();
90 current_link = model.links[current_link.parent];
91 }
92 return Htb.inverse();
93 }
94
108 template <typename Scalar, int nq, typename TargetLink, typename SourceLink>
109 Eigen::Transform<Scalar, 3, Eigen::Isometry> forward_kinematics(const Model<Scalar, nq>& model,
110 const Eigen::Matrix<Scalar, nq, 1>& q,
111 const TargetLink& target_link,
112 const SourceLink& source_link) {
113 return forward_kinematics(model, q, get_link_idx(model, source_link)).inverse()
114 * forward_kinematics(model, q, get_link_idx(model, target_link));
115 }
116
125 template <typename Scalar, int nq>
126 std::vector<Eigen::Transform<Scalar, 3, Eigen::Isometry>> forward_kinematics_com(
127 Model<Scalar, nq>& model,
128 const Eigen::Matrix<Scalar, nq, 1>& q) {
129 // Compute forward kinematics for all the links
130 forward_kinematics(model, q);
131 model.forward_kinematics_com.resize(model.links.size(),
132 Eigen::Transform<Scalar, 3, Eigen::Isometry>::Identity());
133 // Apply center of mass transform for each link
134 for (auto const link : model.links) {
135 model.forward_kinematics_com[link.idx] =
136 model.forward_kinematics[link.idx] * model.links[link.idx].center_of_mass;
137 }
138 return model.forward_kinematics_com;
139 }
140
154 template <typename Scalar, int nq, typename TargetLink, typename SourceLink = int>
155 Eigen::Transform<Scalar, 3, Eigen::Isometry> forward_kinematics_com(const Model<Scalar, nq>& model,
156 const Eigen::Matrix<Scalar, nq, 1>& q,
157 const TargetLink& target_link,
158 const SourceLink& source_link = 0) {
159 // Apply center of mass transform to forward kinematics between source and target link
160 return forward_kinematics(model, q, target_link, source_link)
161 * model.links[get_link_idx(model, target_link)].center_of_mass;
162 }
163
174 template <typename Scalar, int nq, typename TargetLink>
175 Eigen::Matrix<Scalar, 6, nq> jacobian(Model<Scalar, nq>& model,
176 const Eigen::Matrix<Scalar, nq, 1>& q,
177 const TargetLink& target_link) {
178 // Compute forward kinematics for all the links
179 forward_kinematics(model, q);
180
181 model.J.setZero();
182 Eigen::Matrix<Scalar, 3, 1> zIBb;
183 Eigen::Matrix<Scalar, 3, 1> rIBb;
184 Link<Scalar> current_link = model.links[get_link_idx(model, target_link)];
185 Eigen::Matrix<Scalar, 3, 1> rTBb = model.forward_kinematics[current_link.idx].translation();
186
187 // Compute the jacobian
188 while (current_link.idx != model.base_link_idx) {
189 if (current_link.joint.idx != -1) {
190 zIBb = model.forward_kinematics[current_link.idx].linear() * current_link.joint.axis;
191 rIBb = model.forward_kinematics[current_link.idx].translation();
192 if (current_link.joint.type == JointType::PRISMATIC) {
193 model.J.block(0, current_link.joint.idx, 3, 1) = zIBb;
194 model.J.block(3, current_link.joint.idx, 3, 1) = Eigen::Matrix<Scalar, 3, 1>::Zero();
195 }
196 else if (current_link.joint.type == JointType::REVOLUTE) {
197 model.J.block(0, current_link.joint.idx, 3, 1) = (zIBb).cross(rTBb - rIBb);
198 model.J.block(3, current_link.joint.idx, 3, 1) = zIBb;
199 }
200 }
201 // Move up the tree to parent towards the base
202 current_link = model.links[current_link.parent];
203 }
204 return model.J;
205 }
206
218 template <typename Scalar, int nq, typename SourceLink = int>
219 Eigen::Matrix<Scalar, 3, 1> center_of_mass(Model<Scalar, nq>& model,
220 const Eigen::Matrix<Scalar, nq, 1>& q,
221 const SourceLink& source_link = 0) {
222 forward_kinematics_com(model, q);
223 model.center_of_mass.setZero();
224 for (auto link : model.links) {
225 model.center_of_mass += model.forward_kinematics_com[link.idx].translation() * link.mass;
226 }
227 model.center_of_mass /= model.mass;
228 return model.center_of_mass;
229 }
230
231} // namespace tinyrobotics
232
233#endif
std::vector< Eigen::Transform< Scalar, 3, Eigen::Isometry > > forward_kinematics_com(Model< Scalar, nq > &model, const Eigen::Matrix< Scalar, nq, 1 > &q)
Computes the transform between center of mass of each link and the source link.
Definition: kinematics.hpp:126
int get_link_idx(const Model< Scalar, nq > &model, const TargetLink &target_link)
Retrieves the index of the target link in the tinyrobotics model.
Definition: kinematics.hpp:26
Eigen::Matrix< Scalar, 6, nq > jacobian(Model< Scalar, nq > &model, const Eigen::Matrix< Scalar, nq, 1 > &q, const TargetLink &target_link)
Computes the geometric jacobian of the target link from the base link, in the base link frame.
Definition: kinematics.hpp:175
std::vector< Eigen::Transform< Scalar, 3, Eigen::Isometry > > forward_kinematics(Model< Scalar, nq > &model, const Eigen::Matrix< Scalar, nq, 1 > &q)
Computes the transform to all the links in the tinyrobotics model.
Definition: kinematics.hpp:47
Eigen::Matrix< Scalar, 3, 1 > center_of_mass(Model< Scalar, nq > &model, const Eigen::Matrix< Scalar, nq, 1 > &q, const SourceLink &source_link=0)
Computes the center of mass expressed in source link frame.
Definition: kinematics.hpp:219
Contains various math related functions.
Contains struct for representing a tinyrobotics model.
A tinyrobotics model.
Definition: model.hpp:25
Eigen::Matrix< Scalar, 6, nq > J
Jacobian.
Definition: model.hpp:131
Scalar mass
Total mass of the model.
Definition: model.hpp:48
std::vector< Eigen::Transform< Scalar, 3, Eigen::Isometry > > forward_kinematics_com
Vector of forward kinematics com data.
Definition: model.hpp:65
Link< Scalar > get_link(const std::string &name) const
Get a link in the model by name.
Definition: model.hpp:144
Eigen::Matrix< Scalar, 3, 1 > center_of_mass
center of mass position
Definition: model.hpp:68
int base_link_idx
Index of the base link in the models links vector.
Definition: model.hpp:36
std::vector< Eigen::Transform< Scalar, 3, Eigen::Isometry > > forward_kinematics
Vector of forward kinematics data.
Definition: model.hpp:62
std::vector< Link< Scalar > > links
Vector of links in the model.
Definition: model.hpp:51