Biorbd
Markers.cpp
1 #define BIORBD_API_EXPORTS
2 #include "RigidBody/Markers.h"
3 
4 #include <rbdl/Model.h>
5 #include <rbdl/Kinematics.h>
6 #include "Utils/String.h"
7 #include "Utils/Matrix.h"
8 #include "RigidBody/GeneralizedCoordinates.h"
9 #include "RigidBody/GeneralizedVelocity.h"
10 #include "RigidBody/Joints.h"
11 #include "RigidBody/NodeSegment.h"
12 #include "RigidBody/Segment.h"
13 
14 
16  m_marks(std::make_shared<std::vector<biorbd::rigidbody::NodeSegment>>())
17 {
18  //ctor
19 }
20 
22  m_marks(other.m_marks)
23 {
24 
25 }
26 
28 {
29  //dtor
30 }
31 
33 {
35  copy.DeepCopy(*this);
36  return copy;
37 }
38 
40 {
41  m_marks->resize(other.m_marks->size());
42  for (unsigned int i=0; i<other.m_marks->size(); ++i)
43  (*m_marks)[i] = (*other.m_marks)[i].DeepCopy();
44 }
45 
46 // Add a new marker to the markers pool
49  const biorbd::utils::String &name,
50  const biorbd::utils::String &parentName,
51  bool technical,
52  bool anatomical,
53  const biorbd::utils::String& axesToRemove,
54  int id)
55 {
56  biorbd::rigidbody::NodeSegment tp(pos, name, parentName, technical, anatomical, axesToRemove, id);
57  m_marks->push_back(tp);
58 }
59 
61  unsigned int idx) const
62 {
63  return (*m_marks)[idx];
64 }
65 
66 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::marker(
67  const biorbd::utils::String& name) const
68 {
69  std::vector<biorbd::rigidbody::NodeSegment> pos;
70  for (unsigned int i=0; i<nbMarkers(); ++i) // Go through all the markers and select the right ones
71  if (!marker(i).parent().compare(name))
72  pos.push_back(marker(i));
73 
74  return pos;
75 }
76 
77 // Return a marker
81  bool removeAxis,
82  bool updateKin)
83 {
84  // Assuming that this is also a joint type (via BiorbdModel)
85  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
86  if (updateKin)
87  model.UpdateKinematicsCustom(&Q);
88 
89  unsigned int id = model.GetBodyId(n.parent().c_str());
90  if (removeAxis)
91  return biorbd::rigidbody::NodeSegment(RigidBodyDynamics::CalcBodyToBaseCoordinates(model, Q, id, n.removeAxes(), false));
92  else
93  return biorbd::rigidbody::NodeSegment(RigidBodyDynamics::CalcBodyToBaseCoordinates(model, Q, id, n, false));
94 }
95 
96 // Get a marker
99  unsigned int idx,
100  bool removeAxis,
101  bool updateKin)
102 {
103  // Assuming that this is also a joint type (via BiorbdModel)
104  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
105  if (updateKin)
106  model.UpdateKinematicsCustom(&Q);
107 
108  const biorbd::rigidbody::NodeSegment& node(marker(idx));
109  unsigned int id = model.GetBodyId(node.parent().c_str());
110 
111  // Retrieve the position of the marker in the local reference
112  const biorbd::rigidbody::NodeSegment& pos = marker(idx, removeAxis);
113 
114  return biorbd::rigidbody::NodeSegment(RigidBodyDynamics::CalcBodyToBaseCoordinates(model, Q, id, pos, false));
115 }
116 // Get a marker
118  unsigned int idx,
119  bool removeAxis)
120 {
121  const biorbd::rigidbody::NodeSegment& node(marker(idx));
122  if (removeAxis)
123  return node.removeAxes();
124  else
125  return node;
126 }
127 
128 // Get all the markers
129 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::markers(
131  bool removeAxis,
132  bool updateKin)
133 {
134  // Assuming that this is also a joint type (via BiorbdModel)
135  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
136  if (updateKin)
137  model.UpdateKinematicsCustom(&Q);
138 
139  std::vector<biorbd::rigidbody::NodeSegment> pos;
140  for (unsigned int i=0; i<nbMarkers(); ++i)
141  pos.push_back(marker(Q, i, removeAxis, false));
142 
143  return pos;
144 }
145 // Get all the markers in the local reference
146 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::markers(
147  bool removeAxis)
148 {
149  std::vector<biorbd::rigidbody::NodeSegment> pos;
150  for (unsigned int i=0; i<nbMarkers(); ++i)
151  pos.push_back(marker(i, removeAxis));// Forward kinematics
152 
153  return pos;
154 }
155 
156 // Get a marker's velocity
160  unsigned int idx,
161  bool removeAxis,
162  bool updateKin)
163 {
164  // Assuming that this is also a joint type (via BiorbdModel)
165  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
166  if (updateKin) {
167  model.UpdateKinematicsCustom(&Q, &Qdot);
168  }
169 
170 
171  const biorbd::rigidbody::NodeSegment& node(marker(idx));
172  unsigned int id(model.GetBodyId(node.parent().c_str()));
173 
174  // Retrieve the position of the marker in the local reference
175  const biorbd::rigidbody::NodeSegment& pos(marker(idx, removeAxis));
176 
177  // Calculate the velocity of the point
178  return biorbd::rigidbody::NodeSegment(RigidBodyDynamics::CalcPointVelocity(model, Q, Qdot, id, pos, false));
179 }
180 
181 // Get the makers'velocities
182 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::markersVelocity(
185  bool removeAxis,
186  bool updateKin)
187 {
188  // Assuming that this is also a joint type (via BiorbdModel)
189  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
190  if (updateKin) {
191  model.UpdateKinematicsCustom(&Q);
192  }
193 
194  std::vector<biorbd::rigidbody::NodeSegment> pos;
195  for (unsigned int i=0; i<nbMarkers(); ++i)
196  pos.push_back(markerVelocity(Q, Qdot, i, removeAxis, false));
197 
198  return pos;
199 }
200 
201 // Get the technical markers
202 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::technicalMarkers(
204  bool removeAxis,
205  bool updateKin)
206 {
207  // Assuming that this is also a joint type (via BiorbdModel)
208  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
209  if (updateKin) {
210  model.UpdateKinematicsCustom(&Q);
211  }
212 
213  std::vector<biorbd::rigidbody::NodeSegment> pos;
214  for (unsigned int i=0; i<nbMarkers(); ++i) {
215  if ( marker(i).isTechnical() ) {
216  pos.push_back(marker(Q, i, removeAxis, false));
217  }
218  }
219  return pos;
220 }
221 // Get the technical markers in a local reference
222 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::technicalMarkers(
223  bool removeAxis)
224 {
225  std::vector<biorbd::rigidbody::NodeSegment> pos;
226  for (unsigned int i=0; i<nbMarkers(); ++i)
227  if ( marker(i).isTechnical() )
228  pos.push_back(marker(i, removeAxis));// Forward kinematics
229  return pos;
230 }
231 // Get the anatomical markers
232 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::anatomicalMarkers(
234  bool removeAxis,
235  bool updateKin)
236 {
237  // Assuming that this is also a joint type (via BiorbdModel)
238  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
239  if (updateKin) {
240  model.UpdateKinematicsCustom(&Q);
241  }
242 
243  std::vector<biorbd::rigidbody::NodeSegment> pos;
244  for (unsigned int i=0; i<nbMarkers(); ++i) {
245  if ( marker(i).isAnatomical() ) {
246  pos.push_back(marker(Q, i, removeAxis, false));
247  }
248  }
249  return pos;
250 }
251 // Get the anatomical markers in a local reference
252 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::anatomicalMarkers(bool removeAxis)
253 {
254  std::vector<biorbd::rigidbody::NodeSegment> pos;
255  for (unsigned int i=0; i<nbMarkers(); ++i)
256  if ( marker(i).isAnatomical() )
257  pos.push_back(marker(i, removeAxis));// Forward kinematics
258  return pos;
259 }
260 
261 
262 
263 std::vector<biorbd::rigidbody::NodeSegment> biorbd::rigidbody::Markers::segmentMarkers(
265  unsigned int idx,
266  bool removeAxis,
267  bool updateKin)
268 {
269  // Assuming that this is also a joint type (via BiorbdModel)
270  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
271  if (updateKin)
272  model.UpdateKinematicsCustom(&Q, nullptr, nullptr);
273 
274  // Name of the segment to find
275  const biorbd::utils::String& name(model.segment(idx).name());
276 
277  std::vector<biorbd::rigidbody::NodeSegment> pos;
278  for (unsigned int i=0; i<nbMarkers(); ++i) // Go through all the markers and select the right ones
279  if ((*m_marks)[i].parent().compare(name))
280  pos.push_back(marker(Q,i,removeAxis,false));
281 
282  return pos;
283 }
284 
285 
287 {
288  return static_cast<unsigned int>(m_marks->size());
289 }
290 
291 unsigned int biorbd::rigidbody::Markers::nbMarkers(unsigned int idxSegment) const
292 {
293  // Assuming that this is also a joint type (via BiorbdModel)
294  const biorbd::rigidbody::Joints &model = dynamic_cast<const biorbd::rigidbody::Joints &>(*this);
295 
296  // Name of the segment to find
297  const biorbd::utils::String& name(model.segment(idxSegment).name());
298 
299  unsigned int n = 0;
300  for (unsigned int i=0; i<nbMarkers(); ++i) // Go through all the markers and select the right ones
301  if ((*m_marks)[i].parent().compare(name))
302  ++n;
303 
304  return n;
305 }
306 
307 // Get the Jacobian of the markers
308 std::vector<biorbd::utils::Matrix> biorbd::rigidbody::Markers::markersJacobian(
310  bool removeAxis,
311  bool updateKin)
312 {
313  return markersJacobian(Q, removeAxis, updateKin, false);
314 }
315 
316 std::vector<biorbd::utils::Matrix> biorbd::rigidbody::Markers::technicalMarkersJacobian(
318  bool removeAxis,
319  bool updateKin)
320 {
321  return markersJacobian(Q, removeAxis, updateKin, true);
322 }
323 
324 // Get the Jacobian of the technical markers
327  const biorbd::utils::String& parentName,
329  bool updateKin)
330 {
331  // Assuming that this is also a joint type (via BiorbdModel)
332  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
333  if (updateKin) {
334  model.UpdateKinematicsCustom(&Q);
335  }
336 
337  biorbd::utils::Matrix G(biorbd::utils::Matrix::Zero(3, model.nbQ()));;
338 
339  // Calculate the Jacobien of this Tag
340  unsigned int id = model.GetBodyId(parentName.c_str());
341  RigidBodyDynamics::CalcPointJacobian(model, Q, id, p, G, false);
342 
343  return G;
344 }
345 
346 #ifndef BIORBD_USE_CASADI_MATH
348  const std::vector<biorbd::rigidbody::NodeSegment> &markers,
351  bool removeAxes)
352 {
353  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
354  model.UpdateKinematicsCustom(&Q); // also assert for dimensions
355 
356  // Find the technical markers only (body_point)
357  std::vector<biorbd::rigidbody::NodeSegment> body_point(
358  technicalMarkers(removeAxes));
359  std::vector<RigidBodyDynamics::Math::Vector3d> body_pointEigen;
360  for (unsigned int i=0; i<body_point.size(); ++i)
361  body_pointEigen.push_back(body_point[i]);
362 
363  std::vector<RigidBodyDynamics::Math::Vector3d> markersInRbdl;
364  for (unsigned int i = 0; i<markers.size(); ++i)
365  markersInRbdl.push_back(markers[i]);
366 
367  // Associate the body number to each technical marker (body_id)
368  std::vector<unsigned int> body_id;
369  for (unsigned int i=0; i<body_point.size(); ++i)
370  body_id.push_back( static_cast<unsigned int>((*(body_point.begin()+i)).parentId()) );
371 
372  // Call the base function
373  return RigidBodyDynamics::InverseKinematics(
374  dynamic_cast<biorbd::rigidbody::Joints &>(*this),
375  Qinit, body_id, body_pointEigen, markersInRbdl, Q);
376 }
377 #endif
378 
379 // Get the Jacobian of the technical markers
380 std::vector<biorbd::utils::Matrix> biorbd::rigidbody::Markers::markersJacobian(
382  bool removeAxis,
383  bool updateKin,
384  bool lookForTechnical)
385 {
386  // Assuming that this is also a joint type (via BiorbdModel)
387  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
388  if (updateKin) {
389  model.UpdateKinematicsCustom(&Q);
390  }
391 
392  std::vector<biorbd::utils::Matrix> G;
393 
394  for (unsigned int idx=0; idx<nbMarkers(); ++idx){
395  // Actual marker
396  const biorbd::rigidbody::NodeSegment& node(marker(idx));
397  if (lookForTechnical && !node.isTechnical())
398  continue;
399 
400  unsigned int id = model.GetBodyId(node.parent().c_str());
401  const biorbd::utils::Vector3d& pos(marker(idx, removeAxis));
402  biorbd::utils::Matrix G_tp(biorbd::utils::Matrix::Zero(3,model.nbQ()));
403 
404  // Calculate the Jacobian of this Tag
405  RigidBodyDynamics::CalcPointJacobian(model, Q, id, pos, G_tp, false);
406 
407  G.push_back(G_tp);
408  }
409 
410  return G;
411 }
412 
414 {
415  unsigned int nTechMarkers = 0;
416  if (nTechMarkers == 0) // If the function has never been called before
417  for (auto mark : *m_marks)
418  if (mark.isTechnical())
419  ++nTechMarkers;
420 
421  return nTechMarkers;
422 }
423 
424 unsigned int biorbd::rigidbody::Markers::nbTechnicalMarkers(unsigned int idxSegment)
425 {
426  // Assuming that this is also a joint type (via BiorbdModel)
427  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
428 
429  unsigned int nTechMarkers = 0;
430 
431  // Name of the segment to find
432  const biorbd::utils::String& name(model.segment(idxSegment).name());
433 
434  if (nTechMarkers == 0) // If the function has never been called before
435  for (auto mark : *m_marks)
436  if (mark.isTechnical() && !mark.parent().compare(name))
437  ++nTechMarkers;
438 
439  return nTechMarkers;
440 }
441 
442 
444 {
445  unsigned int nAnatMarkers = 0;
446  if (nAnatMarkers == 0) // If the function has never been called before
447  for (auto mark : *m_marks)
448  if (mark.isAnatomical())
449  ++nAnatMarkers;
450 
451  return nAnatMarkers;
452 }
453 
454 std::vector<biorbd::utils::String> biorbd::rigidbody::Markers::markerNames() const{
455  // Extract the name of all the markers of a model
456  std::vector<biorbd::utils::String> names;
457  for (unsigned int i=0; i<nbMarkers(); ++i)
458  names.push_back(marker(i).biorbd::utils::Node::name());
459 
460  return names;
461 }
462 
463 std::vector<biorbd::utils::String> biorbd::rigidbody::Markers::technicalMarkerNames() const{
464  // Extract the name of all the technical markers of a model
465  std::vector<biorbd::utils::String> names;
466  for (unsigned int i=0; i<nbMarkers(); ++i)
467  if (marker(i).isTechnical())
468  names.push_back(marker(i).biorbd::utils::Node::name());
469 
470  return names;
471 }
472 
473 std::vector<biorbd::utils::String> biorbd::rigidbody::Markers::anatomicalMarkerNames() const{
474  // Extract the names of all the anatomical markers of a model
475  std::vector<biorbd::utils::String> names;
476  for (unsigned int i=0; i<nbMarkers(); ++i)
477  if (marker(i).isAnatomical())
478  names.push_back(marker(i).biorbd::utils::Node::name());
479 
480  return names;
481 }
biorbd::rigidbody::Markers::segmentMarkers
std::vector< biorbd::rigidbody::NodeSegment > segmentMarkers(const biorbd::rigidbody::GeneralizedCoordinates &Q, unsigned int idx, bool removeAxis=true, bool updateKin=true)
Return all the markers of the segment idx at a given Q in the global reference frame.
Definition: Markers.cpp:263
biorbd::rigidbody::Markers::inverseKinematics
bool inverseKinematics(const std::vector< biorbd::rigidbody::NodeSegment > &markers, const biorbd::rigidbody::GeneralizedCoordinates &Qinit, biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxes=true)
Performs an inverse kinematics.
Definition: Markers.cpp:347
biorbd::rigidbody::Markers::~Markers
virtual ~Markers()
Destroy class properly.
Definition: Markers.cpp:27
biorbd::utils::Vector3d
Wrapper around Eigen Vector3d and attach it to a parent.
Definition: Vector3d.h:24
biorbd::rigidbody::GeneralizedCoordinates
Class GeneralizedCoordinates.
Definition: GeneralizedCoordinates.h:15
biorbd::rigidbody::Markers::markersVelocity
std::vector< biorbd::rigidbody::NodeSegment > markersVelocity(const biorbd::rigidbody::GeneralizedCoordinates &Q, const biorbd::rigidbody::GeneralizedVelocity &Qdot, bool removeAxis=true, bool updateKin=true)
Return the velocity of all the markers.
Definition: Markers.cpp:182
biorbd::rigidbody::Markers::DeepCopy
biorbd::rigidbody::Markers DeepCopy() const
Deep copy of the markers.
Definition: Markers.cpp:32
biorbd::utils::Node::name
const biorbd::utils::String & name() const
Return the name of the node.
Definition: Node.cpp:59
biorbd::rigidbody::Markers::markersJacobian
std::vector< biorbd::utils::Matrix > markersJacobian(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxis=true, bool updateKin=true)
Return the jacobian of the markers.
Definition: Markers.cpp:308
biorbd::rigidbody::Markers
Holder for the marker set.
Definition: Markers.h:23
biorbd::rigidbody::Markers::addMarker
void addMarker(const biorbd::rigidbody::NodeSegment &pos, const biorbd::utils::String &name, const biorbd::utils::String &parentName, bool technical, bool anatomical, const biorbd::utils::String &axesToRemove, int id=-1)
Add a marker to the set.
Definition: Markers.cpp:47
biorbd::rigidbody::Markers::technicalMarkersJacobian
std::vector< biorbd::utils::Matrix > technicalMarkersJacobian(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxis=true, bool updateKin=true)
Return the jacobian of the technical markers.
Definition: Markers.cpp:316
biorbd::rigidbody::Markers::markerNames
std::vector< biorbd::utils::String > markerNames() const
Return the names of all the markers.
Definition: Markers.cpp:454
biorbd::rigidbody::Joints::nbQ
unsigned int nbQ() const
Return the number of generalized coordinates (Q)
Definition: Joints.cpp:108
biorbd::rigidbody::Markers::nbMarkers
unsigned int nbMarkers() const
Return the number of markers.
Definition: Markers.cpp:286
biorbd::rigidbody::Markers::anatomicalMarkers
std::vector< biorbd::rigidbody::NodeSegment > anatomicalMarkers(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxis=true, bool updateKin=true)
Return all the anatomical markers at a given Q in the global reference frame.
Definition: Markers.cpp:232
biorbd::rigidbody::Markers::markerVelocity
biorbd::rigidbody::NodeSegment markerVelocity(const biorbd::rigidbody::GeneralizedCoordinates &Q, const biorbd::rigidbody::GeneralizedVelocity &Qdot, unsigned int idx, bool removeAxis=true, bool updateKin=true)
Return the velocity of a marker.
Definition: Markers.cpp:157
biorbd::rigidbody::Joints
This is the core of the musculoskeletal model in biorbd.
Definition: Joints.h:40
biorbd::rigidbody::NodeSegment::isTechnical
bool isTechnical() const
Return if node is technical.
Definition: NodeSegment.cpp:106
biorbd::rigidbody::Markers::Markers
Markers()
Construct a marker set.
Definition: Markers.cpp:15
biorbd::rigidbody::Joints::segment
const biorbd::rigidbody::Segment & segment(unsigned int idx) const
Get a segment of index idx.
Definition: Joints.cpp:189
biorbd::utils::Matrix
A wrapper for the Eigen::MatrixXd.
Definition: Matrix.h:21
biorbd::rigidbody::NodeSegment::removeAxes
NodeSegment removeAxes() const
To remove axis.
Definition: NodeSegment.cpp:117
biorbd::rigidbody::Markers::technicalMarkerNames
std::vector< biorbd::utils::String > technicalMarkerNames() const
Return the names of all the technical markers.
Definition: Markers.cpp:463
biorbd::rigidbody::Markers::nbTechnicalMarkers
unsigned int nbTechnicalMarkers()
Return the number of technical markers.
Definition: Markers.cpp:413
biorbd::rigidbody::Markers::m_marks
std::shared_ptr< std::vector< biorbd::rigidbody::NodeSegment > > m_marks
The markers.
Definition: Markers.h:351
biorbd::rigidbody::Markers::anatomicalMarkerNames
std::vector< biorbd::utils::String > anatomicalMarkerNames() const
Return the names of all the anatomical markers.
Definition: Markers.cpp:473
biorbd::rigidbody::Markers::nbAnatomicalMarkers
unsigned int nbAnatomicalMarkers()
Return the number of anatomical markers.
Definition: Markers.cpp:443
biorbd::rigidbody::Markers::marker
const biorbd::rigidbody::NodeSegment & marker(unsigned int idx) const
Return the marker of index idx.
Definition: Markers.cpp:60
biorbd::utils::String
Wrapper around the std::string class with augmented functionality.
Definition: String.h:17
biorbd::rigidbody::NodeSegment
A point attached to a segment, generally speaking a skin marker.
Definition: NodeSegment.h:19
biorbd::rigidbody::Joints::UpdateKinematicsCustom
void UpdateKinematicsCustom(const biorbd::rigidbody::GeneralizedCoordinates *Q=nullptr, const biorbd::rigidbody::GeneralizedVelocity *Qdot=nullptr, const biorbd::rigidbody::GeneralizedAcceleration *Qddot=nullptr)
Update the kinematic variables such as body velocities and accelerations in the model to reflect the ...
Definition: Joints.cpp:1063
biorbd::utils::Node::parent
const biorbd::utils::String & parent() const
Return the parent name of the node.
Definition: Node.cpp:64
biorbd::rigidbody::GeneralizedVelocity
Class GeneralizedVelocity.
Definition: GeneralizedVelocity.h:15
biorbd::rigidbody::Markers::markers
std::vector< biorbd::rigidbody::NodeSegment > markers(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxis=true, bool updateKin=true)
Return all the markers at a given Q in the global reference frame.
Definition: Markers.cpp:129
biorbd::rigidbody::Markers::technicalMarkers
std::vector< biorbd::rigidbody::NodeSegment > technicalMarkers(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool removeAxis=true, bool updateKin=true)
Return all the technical markers at a given Q in the global reference frame.
Definition: Markers.cpp:202