Biorbd
RotoTransNodes.cpp
1 #define BIORBD_API_EXPORTS
2 #include "RigidBody/RotoTransNodes.h"
3 
4 #include <rbdl/Model.h>
5 #include <rbdl/Kinematics.h>
6 #include "Utils/String.h"
7 #include "Utils/Matrix.h"
8 #include "Utils/Rotation.h"
9 #include "Utils/RotoTransNode.h"
10 #include "RigidBody/GeneralizedCoordinates.h"
11 #include "RigidBody/Joints.h"
12 #include "RigidBody/Segment.h"
13 
15  m_RTs(std::make_shared<std::vector<biorbd::utils::RotoTransNode>>())
16 {
17  //ctor
18 }
19 
21 {
22  m_RTs = other.m_RTs;
23 }
24 
26 {
27 
28 }
29 
31 {
33  copy.DeepCopy(*this);
34  return copy;
35 }
36 
38 {
39  m_RTs->resize(other.m_RTs->size());
40  for (unsigned int i=0; i<other.m_RTs->size(); ++i)
41  (*m_RTs)[i] = (*other.m_RTs)[i].DeepCopy();
42 }
43 
45 {
46  m_RTs->push_back(biorbd::utils::RotoTransNode());
47 }
48 
49 // Add a new marker to the existing pool of markers
51  const biorbd::utils::RotoTransNode &RotoTrans)
52 {
53  m_RTs->push_back(biorbd::utils::RotoTransNode(RotoTrans));
54 }
55 
57 {
58  return static_cast<unsigned int>(m_RTs->size());
59 }
60 
61 
62 // Get the markers in the global reference
63 const std::vector<biorbd::utils::RotoTransNode>& biorbd::rigidbody::RotoTransNodes::RTs() const
64 {
65  return *m_RTs;
66 }
67 
68 std::vector<biorbd::utils::RotoTransNode> biorbd::rigidbody::RotoTransNodes::RTs(
69  const biorbd::utils::String& segmentName)
70 {
71  std::vector<biorbd::utils::RotoTransNode> pos;
72  for (unsigned int i=0; i<nbRTs(); ++i) // Scan through all the markers and select the good ones
73  if (!RT(i).parent().compare(segmentName))
74  pos.push_back(RT(i));
75  return pos;
76 }
77 
79  unsigned int idx)
80 {
81  return (*m_RTs)[idx];
82 }
83 
84 // Get the RotoTransNodes at the position given by Q
85 std::vector<biorbd::utils::RotoTransNode> biorbd::rigidbody::RotoTransNodes::RTs(
87  bool updateKin)
88 {
89  // Assuming that this is also a Joints type (via BiorbdModel)
90  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
91  if (updateKin) {
92  model.UpdateKinematicsCustom (&Q);
93  }
94 
95  std::vector<biorbd::utils::RotoTransNode> pos;
96  for (unsigned int i=0; i<nbRTs(); ++i)
97  pos.push_back(RT(Q, i, false));
98 
99  return pos;
100 }
101 
102 // Get an IMU at the position given by Q
105  unsigned int idx,
106  bool updateKin)
107 {
108  // Assuming that this is also a Joints type (via BiorbdModel)
109  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
110  if (updateKin) {
111  model.UpdateKinematicsCustom (&Q);
112  }
113 
114  biorbd::utils::RotoTransNode node = RT(idx);
115  unsigned int id = static_cast<unsigned int>(model.GetBodyBiorbdId(node.parent()));
116 
117  return model.globalJCS(id) * node;
118 }
119 
120 std::vector<biorbd::utils::RotoTransNode> biorbd::rigidbody::RotoTransNodes::segmentRTs(
122  unsigned int idx,
123  bool updateKin)
124 {
125  // Assuming that this is also a Joints type (via BiorbdModel)
126  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
127  if (updateKin) {
128  model.UpdateKinematicsCustom (&Q);
129  }
130 
131  // Segment name to find
132  biorbd::utils::String name(model.segment(idx).name());
133 
134  std::vector<biorbd::utils::RotoTransNode> pos;
135  for (unsigned int i=0; i<nbRTs(); ++i) // scan all the markers and select the right ones
136  if (!((*m_RTs)[i]).parent().compare(name))
137  pos.push_back(RT(Q,i,false));
138 
139  return pos;
140 }
141 
142 // Get the Jacobian of the markers
143 std::vector<biorbd::utils::Matrix> biorbd::rigidbody::RotoTransNodes::RTsJacobian(
145  bool updateKin)
146 {
147  // Assuming that this is also a Joints type (via BiorbdModel)
148  biorbd::rigidbody::Joints &model = dynamic_cast<biorbd::rigidbody::Joints &>(*this);
149  if (updateKin) {
150  model.UpdateKinematicsCustom (&Q);
151  }
152 
153  std::vector<biorbd::utils::Matrix> G;
154 
155  for (unsigned int idx=0; idx<nbRTs(); ++idx){
156  // Actual marker
157  biorbd::utils::RotoTransNode node = RT(idx);
158 
159  unsigned int id = model.GetBodyId(node.parent().c_str());
160  biorbd::utils::Matrix G_tp(biorbd::utils::Matrix::Zero(9,model.dof_count));
161 
162  // Calculate the Jacobian of this Tag
163  model.CalcMatRotJacobian(Q, id, node.rot(), G_tp, false); // False for speed
164 
165  G.push_back(G_tp);
166  }
167 
168  return G;
169 }
170 
171 std::vector<biorbd::utils::String> biorbd::rigidbody::RotoTransNodes::RTsNames()
172 {
173  // Extract the name of all the markers of a model
174  std::vector<biorbd::utils::String> names;
175  for (unsigned int i=0; i<nbRTs(); ++i)
176  names.push_back(RT(i).biorbd::utils::Node::name());
177  return names;
178 }
biorbd::rigidbody::RotoTransNodes::RotoTransNodes
RotoTransNodes()
Construct RT set.
Definition: RotoTransNodes.cpp:14
biorbd::rigidbody::RotoTransNodes::m_RTs
std::shared_ptr< std::vector< biorbd::utils::RotoTransNode > > m_RTs
All the RTs.
Definition: RotoTransNodes.h:146
biorbd::rigidbody::GeneralizedCoordinates
Class GeneralizedCoordinates.
Definition: GeneralizedCoordinates.h:15
biorbd::utils::RotoTransNode
A RotoTrans which is attached to a segment.
Definition: RotoTransNode.h:16
biorbd::utils::Node::name
const biorbd::utils::String & name() const
Return the name of the node.
Definition: Node.cpp:59
biorbd::rigidbody::Joints::globalJCS
biorbd::utils::RotoTrans globalJCS(const biorbd::rigidbody::GeneralizedCoordinates &Q, const biorbd::utils::String &name)
Return the joint coordinate system (JCS) for the segment in global reference frame at a given Q.
Definition: Joints.cpp:267
biorbd::rigidbody::Joints::GetBodyBiorbdId
int GetBodyBiorbdId(const biorbd::utils::String &segmentName) const
Return the biorbd body identification.
Definition: Joints.cpp:245
biorbd::rigidbody::RotoTransNodes::RTs
const std::vector< biorbd::utils::RotoTransNode > & RTs() const
Return all the RTs in the local reference of the segment.
Definition: RotoTransNodes.cpp:63
biorbd::rigidbody::Joints
This is the core of the musculoskeletal model in biorbd.
Definition: Joints.h:40
biorbd::rigidbody::RotoTransNodes::~RotoTransNodes
virtual ~RotoTransNodes()
Destroy the class properly.
Definition: RotoTransNodes.cpp:25
biorbd::rigidbody::Joints::segment
const biorbd::rigidbody::Segment & segment(unsigned int idx) const
Get a segment of index idx.
Definition: Joints.cpp:189
biorbd::rigidbody::RotoTransNodes
Hold a set of RotoTransNodes.
Definition: RotoTransNodes.h:23
biorbd::utils::Matrix
A wrapper for the Eigen::MatrixXd.
Definition: Matrix.h:21
biorbd::rigidbody::Joints::CalcMatRotJacobian
void CalcMatRotJacobian(const biorbd::rigidbody::GeneralizedCoordinates &Q, unsigned int segmentIdx, const RigidBodyDynamics::Math::Matrix3d &rotation, RigidBodyDynamics::Math::MatrixNd &G, bool updateKin)
Calculate the jacobian matrix of a rotation matrix.
Definition: Joints.cpp:1072
biorbd::rigidbody::RotoTransNodes::addRT
void addRT()
Add a new RT to the set.
Definition: RotoTransNodes.cpp:44
biorbd::rigidbody::RotoTransNodes::nbRTs
unsigned int nbRTs() const
Return the number of RTs in the set.
Definition: RotoTransNodes.cpp:56
biorbd::rigidbody::RotoTransNodes::RTsNames
std::vector< biorbd::utils::String > RTsNames()
Return the names of the RTs.
Definition: RotoTransNodes.cpp:171
biorbd::utils::String
Wrapper around the std::string class with augmented functionality.
Definition: String.h:17
biorbd::rigidbody::RotoTransNodes::RT
const biorbd::utils::RotoTransNode & RT(unsigned int idx)
Return the RTs of a specified index.
Definition: RotoTransNodes.cpp:78
biorbd::rigidbody::RotoTransNodes::DeepCopy
biorbd::rigidbody::RotoTransNodes DeepCopy() const
Deep copy of the RTs data.
Definition: RotoTransNodes.cpp:30
biorbd::rigidbody::RotoTransNodes::segmentRTs
std::vector< biorbd::utils::RotoTransNode > segmentRTs(const biorbd::rigidbody::GeneralizedCoordinates &Q, unsigned int idx, bool updateKin=true)
Return all the RTs on a specified segment.
Definition: RotoTransNodes.cpp:120
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::RotoTransNodes::RTsJacobian
std::vector< biorbd::utils::Matrix > RTsJacobian(const biorbd::rigidbody::GeneralizedCoordinates &Q, bool updateKin=true)
Return the jacobian of the RTs.
Definition: RotoTransNodes.cpp:143
biorbd::utils::RotoTrans::rot
biorbd::utils::Rotation rot() const
Return the rotation matrix.
Definition: RotoTrans.cpp:97