Rigid Body Dynamics Library
rbdl_utils.cc
Go to the documentation of this file.
1/*
2 * RBDL - Rigid Body Dynamics Library
3 * Copyright (c) 2011-2018 Martin Felis <martin@fysx.org>
4 *
5 * Licensed under the zlib license. See LICENSE for more details.
6 */
7
8#include "rbdl/rbdl_utils.h"
9#include "rbdl/rbdl_errors.h"
10
11#include "rbdl/rbdl_math.h"
12#include "rbdl/Model.h"
13#include "rbdl/Kinematics.h"
14
15#include <sstream>
16#include <iomanip>
17
18namespace RigidBodyDynamics
19{
20
21namespace Utils
22{
23
24using namespace std;
25using namespace Math;
26
27#ifndef RBDL_USE_CASADI_MATH
28string get_dof_name (const SpatialVector &joint_dof)
29{
30 if (joint_dof == SpatialVector (1., 0., 0., 0., 0., 0.)) {
31 return "RX";
32 } else if (joint_dof == SpatialVector (0., 1., 0., 0., 0., 0.)) {
33 return "RY";
34 } else if (joint_dof == SpatialVector (0., 0., 1., 0., 0., 0.)) {
35 return "RZ";
36 } else if (joint_dof == SpatialVector (0., 0., 0., 1., 0., 0.)) {
37 return "TX";
38 } else if (joint_dof == SpatialVector (0., 0., 0., 0., 1., 0.)) {
39 return "TY";
40 } else if (joint_dof == SpatialVector (0., 0., 0., 0., 0., 1.)) {
41 return "TZ";
42 }
43
44 ostringstream dof_stream(ostringstream::out);
45 dof_stream << "custom_axis (" << joint_dof.transpose() << ")";
46 return dof_stream.str();
47}
48#endif
49
50string get_body_name (const RigidBodyDynamics::Model &model,
51 unsigned int body_id)
52{
53 if (model.mBodies[body_id].mIsVirtual) {
54 // if there is not a unique child we do not know what to do...
55 if (model.mu[body_id].size() != 1) {
56 return "";
57 }
58
59 return get_body_name (model, model.mu[body_id][0]);
60 }
61
62 return model.GetBodyName(body_id);
63}
64
65#ifndef RBDL_USE_CASADI_MATH
66RBDL_DLLAPI std::string GetModelDOFOverview (const Model &model) {
67 stringstream result ("");
68
69 unsigned int q_index = 0;
70 for (unsigned int i = 1; i < model.mBodies.size(); i++) {
71 if (model.mJoints[i].mDoFCount == 1) {
72 result << setfill(' ') << setw(3) << q_index << ": " << get_body_name(model,
73 i) << "_" << get_dof_name (model.S[i]) << endl;
74 q_index++;
75 } else {
76 for (unsigned int j = 0; j < model.mJoints[i].mDoFCount; j++) {
77 result << setfill(' ') << setw(3) << q_index << ": " << get_body_name(model,
78 i) << "_" << get_dof_name (model.mJoints[i].mJointAxes[j]) << endl;
79 q_index++;
80 }
81 }
82 }
83
84 return result.str();
85}
86#endif
87
88#ifndef RBDL_USE_CASADI_MATH
89std::string print_hierarchy (const RigidBodyDynamics::Model &model,
90 unsigned int body_index = 0, int indent = 0)
91{
92 stringstream result ("");
93
94 for (int j = 0; j < indent; j++) {
95 result << " ";
96 }
97
98 result << get_body_name (model, body_index);
99
100 if (body_index > 0) {
101 result << " [ ";
102 }
103
104 while (model.mBodies[body_index].mIsVirtual) {
105 if (model.mu[body_index].size() == 0) {
106 result << " end";
107 break;
108 } else if (model.mu[body_index].size() > 1) {
109 std::ostringstream errormsg;
110 errormsg << endl <<
111 "Error: Cannot determine multi-dof joint as massless body with id " <<
112 body_index << " (name: " << model.GetBodyName(body_index) <<
113 ") has more than one child:" << endl;
114 for (unsigned int ci = 0; ci < model.mu[body_index].size(); ci++) {
115 errormsg << " id: " << model.mu[body_index][ci] << " name: " <<
116 model.GetBodyName(model.mu[body_index][ci]) << endl;
117 }
118 throw Errors::RBDLError(errormsg.str());
119 }
120
121 result << get_dof_name(model.S[body_index]) << ", ";
122
123 body_index = model.mu[body_index][0];
124 }
125
126 if (body_index > 0) {
127 result << get_dof_name(model.S[body_index]) << " ]";
128 }
129 result << endl;
130
131 unsigned int child_index = 0;
132 for (child_index = 0; child_index < model.mu[body_index].size();
133 child_index++) {
134 result << print_hierarchy (model, model.mu[body_index][child_index],
135 indent + 1);
136 }
137
138 // print fixed children
139 for (unsigned int fbody_index = 0; fbody_index < model.mFixedBodies.size();
140 fbody_index++) {
141 if (model.mFixedBodies[fbody_index].mMovableParent == body_index) {
142 for (int j = 0; j < indent + 1; j++) {
143 result << " ";
144 }
145
146 result << model.GetBodyName(model.fixed_body_discriminator + fbody_index) <<
147 " [fixed]" << endl;
148 }
149 }
150
151
152 return result.str();
153}
154#endif
155
156#ifndef RBDL_USE_CASADI_MATH
157RBDL_DLLAPI std::string GetModelHierarchy (const Model &model)
158{
159 stringstream result ("");
160
161 result << print_hierarchy (model);
162
163 return result.str();
164}
165#endif
166
167RBDL_DLLAPI std::string GetNamedBodyOriginsOverview (Model &model)
168{
169 stringstream result ("");
170
171 VectorNd Q (VectorNd::Zero(model.dof_count));
172 UpdateKinematicsCustom (model, &Q, NULL, NULL);
173
174 for (unsigned int body_id = 0; body_id < model.mBodies.size(); body_id++) {
175 std::string body_name = model.GetBodyName (body_id);
176
177 if (body_name.size() == 0) {
178 continue;
179 }
180
181 Vector3d position = CalcBodyToBaseCoordinates (model, Q, body_id, Vector3d (0.,
182 0., 0.), false);
183
184 result << body_name << ": " << position.transpose() << endl;
185 }
186
187 return result.str();
188}
189
190RBDL_DLLAPI void CalcCenterOfMass (
191 Model &model,
192 const Math::VectorNd &q,
193 const Math::VectorNd &qdot,
194 const Math::VectorNd *qddot,
195 Scalar &mass,
196 Math::Vector3d &com,
197 Math::Vector3d *com_velocity,
198 Math::Vector3d *com_acceleration,
199 Math::Vector3d *angular_momentum,
200 Math::Vector3d *change_of_angular_momentum,
201 bool update_kinematics)
202{
203 // If we want to compute com_acceleration or change of angular momentum
204 // we must have qddot provided.
205 assert( (com_acceleration == NULL && change_of_angular_momentum == NULL)
206 || (qddot != NULL) );
207
208 if (update_kinematics) {
209 UpdateKinematicsCustom (model, &q, &qdot, qddot);
210 }
211
212 for (size_t i = 1; i < model.mBodies.size(); i++) {
213 model.Ic[i] = model.I[i];
214 model.hc[i] = model.Ic[i].toMatrix() * model.v[i];
215 model.hdotc[i] = model.Ic[i] * model.a[i] + crossf(model.v[i],
216 model.Ic[i] * model.v[i]);
217 }
218
219 if (qddot && (com_acceleration || change_of_angular_momentum)) {
220 for (size_t i = 1; i < model.mBodies.size(); i++) {
221 model.hdotc[i] = model.Ic[i] * model.a[i] + crossf(model.v[i],
222 model.Ic[i] * model.v[i]);
223 }
224 }
225
226 SpatialRigidBodyInertia Itot (0., Vector3d (0., 0., 0.), Matrix3d::Zero());
227 SpatialVector htot (SpatialVector::Zero());
228 SpatialVector hdot_tot (SpatialVector::Zero());
229
230 for (size_t i = model.mBodies.size() - 1; i > 0; i--) {
231 unsigned int lambda = model.lambda[i];
232
233 if (lambda != 0) {
234 model.Ic[lambda] = model.Ic[lambda] + model.X_lambda[i].applyTranspose (
235 model.Ic[i]);
236 model.hc[lambda] = model.hc[lambda] + model.X_lambda[i].applyTranspose (
237 model.hc[i]);
238 } else {
239 Itot = Itot + model.X_lambda[i].applyTranspose (model.Ic[i]);
240 htot = htot + model.X_lambda[i].applyTranspose (model.hc[i]);
241 }
242 }
243
244 if (qddot && (com_acceleration || change_of_angular_momentum)) {
245 for (size_t i = model.mBodies.size() - 1; i > 0; i--) {
246 unsigned int lambda = model.lambda[i];
247
248 if (lambda != 0) {
249 model.hdotc[lambda] = model.hdotc[lambda] + model.X_lambda[i].applyTranspose (
250 model.hdotc[i]);
251 } else {
252 hdot_tot = hdot_tot + model.X_lambda[i].applyTranspose (model.hdotc[i]);
253 }
254 }
255 }
256
257 mass = Itot.m;
258 com = Itot.h / mass;
259 LOG << "mass = " << mass << " com = " << com.transpose() << " htot = " <<
260 htot.transpose() << std::endl;
261
262 if (com_velocity) {
263 *com_velocity = Vector3d (htot[3] / mass, htot[4] / mass, htot[5] / mass);
264 }
265
266 if (angular_momentum) {
267 htot = Xtrans (com).applyAdjoint (htot);
268 angular_momentum->set (htot[0], htot[1], htot[2]);
269 }
270
271 if (com_acceleration) {
272 *com_acceleration = Vector3d (hdot_tot[3] / mass, hdot_tot[4] / mass,
273 hdot_tot[5] / mass);
274 }
275
276 if (change_of_angular_momentum) {
277 hdot_tot = Xtrans (com).applyAdjoint (hdot_tot);
278 change_of_angular_momentum->set (hdot_tot[0], hdot_tot[1], hdot_tot[2]);
279 }
280}
281
282RBDL_DLLAPI void CalcZeroMomentPoint (
283 Model &model,
284 const Math::VectorNd &q,
285 const Math::VectorNd &qdot,
286 const Math::VectorNd &qddot,
287 Vector3d* zmp,
288 const Math::Vector3d &normal,
289 const Math::Vector3d &point,
290 bool update_kinematics
291)
292{
293 if (zmp == NULL) {
294 throw Errors::RBDLError("ZMP (output) is 'NULL'!\n");
295 }
296
297 // update kinematics if required
298 // NOTE UpdateKinematics computes model.a[i] and model.v[i] required for
299 // change of momentum
300 if (update_kinematics) {
301 UpdateKinematicsCustom (model, &q, &qdot, &qddot);
302 }
303
304 // compute change of momentum of each single body (same as in RNEA/InverseDynamics)
305 for (size_t i = 1; i < model.mBodies.size(); i++) {
306 model.Ic[i] = model.I[i];
307 model.hdotc[i] = model.Ic[i] * model.a[i] + crossf(model.v[i],
308 model.Ic[i] * model.v[i]);
309 }
310
311 SpatialRigidBodyInertia I_tot (0., Vector3d (0., 0., 0.), Matrix3d::Zero());
312 SpatialVector h_tot (SpatialVector::Zero());
313 SpatialVector hdot_tot (SpatialVector::Zero());
314
315 // compute total change of momentum and CoM wrt to root body (idx = 0)
316 // by recursively summing up local change of momentum
317 for (size_t i = model.mBodies.size() - 1; i > 0; i--) {
318 unsigned int lambda = model.lambda[i];
319
320 if (lambda != 0) {
321 model.Ic[lambda] = model.Ic[lambda] + model.X_lambda[i].applyTranspose (
322 model.Ic[i]);
323 model.hc[lambda] = model.hc[lambda] + model.X_lambda[i].applyTranspose (
324 model.hc[i]);
325 model.hdotc[lambda] = model.hdotc[lambda] + model.X_lambda[i].applyTranspose (
326 model.hdotc[i]);
327 } else {
328 I_tot = I_tot + model.X_lambda[i].applyTranspose (model.Ic[i]);
329 h_tot = h_tot + model.X_lambda[i].applyTranspose (model.hc[i]);
330 hdot_tot = hdot_tot + model.X_lambda[i].applyTranspose (model.hdotc[i]);
331 }
332 }
333
334 // compute CoM from mass and total inertia
335 const Scalar mass = I_tot.m;
336 const Vector3d com = I_tot.h / mass;
337
338 // project angular momentum onto CoM
339 SpatialTransform Xcom = Xtrans (com);
340 hdot_tot = Xcom.applyAdjoint (hdot_tot);
341
342 // compute net external force at CoM by removing effects due to gravity
343 hdot_tot = hdot_tot - mass * SpatialVector (0., 0., 0., model.gravity[0],
344 model.gravity[1], model.gravity[2]);
345
346 // express total change of momentum in world coordinates
347 hdot_tot = Xcom.inverse().applyAdjoint (hdot_tot);
348
349 // project purified change of momentum onto surface
350 // z = n x n_0
351 // -------
352 // n * f
353 Vector3d n_0 = hdot_tot.block<3,1>(0,0);
354 Vector3d f = hdot_tot.block<3,1>(3,0);
355 *zmp = normal.cross(n_0) / normal.dot(f);
356
357 // double distance = (hdot_tot - point).dot(normal);
358 // zmp = hdot_tot - distance * normal;
359
360 return;
361}
362
364 Model &model,
365 const Math::VectorNd &q,
366 bool update_kinematics)
367{
368 Scalar mass;
369 Vector3d com;
371 model,
372 q,
373 VectorNd::Zero (model.qdot_size),
374 NULL,
375 mass,
376 com,
377 NULL,
378 NULL,
379 NULL,
380 NULL,
381 update_kinematics);
382
383 Vector3d g = - Vector3d (model.gravity[0], model.gravity[1], model.gravity[2]);
384 LOG << "pot_energy: " << " mass = " << mass << " com = " << com.transpose() <<
385 std::endl;
386
387 return mass * com.dot(g);
388}
389
391 Model &model,
392 const Math::VectorNd &q,
393 const Math::VectorNd &qdot,
394 bool update_kinematics)
395{
396 if (update_kinematics) {
397 UpdateKinematicsCustom (model, &q, &qdot, NULL);
398 }
399
400 Scalar result = 0.;
401
402 for (size_t i = 1; i < model.mBodies.size(); i++) {
403 result += 0.5 * model.v[i].transpose() * (model.I[i] * model.v[i]);
404 }
405 return result;
406}
407
408}
409}
static MX_Xd_dynamic Zero(unsigned int nrows, unsigned int ncols=1)
Definition: MX_Xd_dynamic.h:52
Base class for all RBDL exceptions.
Definition: rbdl_errors.h:35
void set(const double &v0, const double &v1, const double &v2)
RBDL_DLLAPI Vector3d CalcBodyToBaseCoordinates(Model &model, const VectorNd &Q, unsigned int body_id, const Vector3d &point_body_coordinates, bool update_kinematics)
Returns the base coordinates of a point given in body coordinates.
Definition: Kinematics.cc:157
RBDL_DLLAPI void UpdateKinematicsCustom(Model &model, const VectorNd *Q, const VectorNd *QDot, const VectorNd *QDDot)
Selectively updates model internal states of body positions, velocities and/or accelerations.
Definition: Kinematics.cc:78
SpatialMatrix crossf(const SpatialVector &v)
SpatialVector_t SpatialVector
Definition: rbdl_math.h:61
SpatialTransform Xtrans(const Vector3d &r)
RBDL_DLLAPI std::string GetModelHierarchy(const Model &model)
Creates a human readable overview of the model.
Definition: rbdl_utils.cc:157
RBDL_DLLAPI Scalar CalcKineticEnergy(Model &model, const Math::VectorNd &q, const Math::VectorNd &qdot, bool update_kinematics)
Computes the kinetic energy of the full model.
Definition: rbdl_utils.cc:390
RBDL_DLLAPI std::string GetModelDOFOverview(const Model &model)
Creates a human readable overview of the Degrees of Freedom.
Definition: rbdl_utils.cc:66
RBDL_DLLAPI void CalcCenterOfMass(Model &model, const Math::VectorNd &q, const Math::VectorNd &qdot, const Math::VectorNd *qddot, Scalar &mass, Math::Vector3d &com, Math::Vector3d *com_velocity, Math::Vector3d *com_acceleration, Math::Vector3d *angular_momentum, Math::Vector3d *change_of_angular_momentum, bool update_kinematics)
Computes the Center of Mass (COM) and optionally its linear velocity.
Definition: rbdl_utils.cc:190
RBDL_DLLAPI void CalcZeroMomentPoint(Model &model, const Math::VectorNd &q, const Math::VectorNd &qdot, const Math::VectorNd &qddot, Vector3d *zmp, const Math::Vector3d &normal, const Math::Vector3d &point, bool update_kinematics)
Computes the Zero-Moment-Point (ZMP) on a given contact surface.
Definition: rbdl_utils.cc:282
RBDL_DLLAPI Scalar CalcPotentialEnergy(Model &model, const Math::VectorNd &q, bool update_kinematics)
Computes the potential energy of the full model.
Definition: rbdl_utils.cc:363
RBDL_DLLAPI std::string GetNamedBodyOriginsOverview(Model &model)
Creates a human readable overview of the locations of all bodies that have names.
Definition: rbdl_utils.cc:167
std::string print_hierarchy(const RigidBodyDynamics::Model &model, unsigned int body_index=0, int indent=0)
Definition: rbdl_utils.cc:89
string get_body_name(const RigidBodyDynamics::Model &model, unsigned int body_id)
Definition: rbdl_utils.cc:50
string get_dof_name(const SpatialVector &joint_dof)
Definition: rbdl_utils.cc:28
STL namespace.
Compact representation for Spatial Inertia.
Vector3d h
Coordinates of the center of mass.
Compact representation of spatial transformations.
SpatialVector applyAdjoint(const SpatialVector &f_sp)