SCL  1.0
Standard Control Library : Control, dynamics, physics, and simulation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Friends Macros Groups Pages
CDynamicsScl.hpp
1 /* This file is part of scl, a control and simulation library
2 for robots and biomechanical models.
3 
4 scl is free software; you can redistribute it and/or
5 modify it under the terms of the GNU Lesser General Public
6 License as published by the Free Software Foundation; either
7 version 3 of the License, or (at your option) any later version.
8 
9 Alternatively, you can redistribute it and/or
10 modify it under the terms of the GNU General Public License as
11 published by the Free Software Foundation; either version 2 of
12 the License, or (at your option) any later version.
13 
14 scl is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18 
19 You should have received a copy of the GNU Lesser General Public
20 License and a copy of the GNU General Public License along with
21 scl. If not, see <http://www.gnu.org/licenses/>.
22 */
23 /*
24  * CDynamicsScl.hpp
25  *
26  * Created on: Oct 7, 2013
27  *
28  * Copyright (C) 2013
29  *
30  * Author: Samir Menon <smenon@stanford.edu>
31  */
32 
33 #ifndef CDYNAMICSSCL_HPP_
34 #define CDYNAMICSSCL_HPP_
35 
36 #include <scl/dynamics/CDynamicsBase.hpp>
37 #include <scl/data_structs/SRigidBodyDyn.hpp>
38 
39 #include <sutil/CMappedTree.hpp>
40 #include <string>
41 
42 namespace scl
43 {
58  {
59 public:
60  /* *******************************************************************
61  * Overall Computational function.
62  * ******************************************************************* */
69  virtual sBool computeGCModel(
71  const SRobotSensors * arg_sensor_data,
74  SGcModel * arg_gc_model);
75 
76  /* *******************************************************************
77  * Coordinate Transformations
78  * ******************************************************************* */
82  inline virtual sBool computeTransformsForAllLinks(
86  const Eigen::VectorXd& arg_q);
87 
93  virtual sBool computeTransform(
95  SRigidBodyDyn& arg_link,
97  const Eigen::VectorXd& arg_q);
98 
112  Eigen::Affine3d& arg_T,
114  SRigidBodyDyn& arg_link,
117  const SRigidBodyDyn* arg_ancestor,
119  const Eigen::VectorXd& arg_q);
120 
121  /* *******************************************************************
122  * Compute Jacobians
123  * ******************************************************************* */
132  Eigen::MatrixXd& arg_J,
135  SRigidBodyDyn& arg_link,
137  const Eigen::VectorXd& arg_q,
139  const Eigen::Vector3d& arg_pos_local);
140 
148  virtual sBool computeJacobian(
150  Eigen::MatrixXd& arg_J,
152  const SRigidBodyDyn& arg_link,
154  const Eigen::VectorXd& arg_q,
156  const Eigen::Vector3d& arg_pos_local) const;
157 
166  Eigen::MatrixXd& arg_J,
168  SRigidBodyDyn& arg_link,
171  const SRigidBodyDyn* arg_ancestor,
173  const Eigen::VectorXd& arg_q,
175  const Eigen::Vector3d& arg_pos_local,
178  const bool arg_recompute_transforms=true)
179  { return false; }
180 
188  const Eigen::VectorXd& arg_q) const
189  {
190  if(false == has_been_init_){ return false; }
192  //Update the com Jacobians.
193  for(it = arg_tree.begin(), ite = arg_tree.end(); it!=ite; ++it)
194  { flag = flag && computeJacobian(it->J_com_,*it, arg_q, it->link_ds_->com_); }
195  return flag;
196  }
197 
198  /* *******************************************************************
199  * Dynamics Helper functions.
200  * ******************************************************************* */
204  inline sBool computeInertiaGC(Eigen::MatrixXd &ret_Mgc,
208  const Eigen::VectorXd& arg_q)
209  {
210  if(false == has_been_init_){ return false; }
212  int dof = robot_parsed_data_->dof_;
213  ret_Mgc.setZero(dof, dof);//Set the generalized inertia to zero.
214  for(it = arg_tree.begin(), ite = arg_tree.end(); it!=ite; ++it)
215  {//Compute each link's contribution to the overall Mgc
216  if(it->link_ds_->is_root_){ continue; }//Root doesn't move
217  ret_Mgc += it->link_ds_->mass_ * (it->J_com_.block(0,0,3,dof).transpose() * it->J_com_.block(0,0,3,dof));
218  ret_Mgc += it->J_com_.block(3,0,3,dof).transpose() * it->link_ds_->inertia_ * it->J_com_.block(3,0,3,dof);
219  }
220  return flag;
221  }
222 
226  inline sBool computeForceGravityGC(Eigen::VectorXd &ret_FgravGC,
230  const Eigen::VectorXd& arg_q)
231  {
232  if(false == has_been_init_){ return false; }
234  //Update the generalized gravity force vector.
235  int dof = robot_parsed_data_->dof_;
236  ret_FgravGC.setZero(dof);
237  for(it = arg_tree.begin(), ite = arg_tree.end(); it!=ite; ++it)
238  {
239  if(it->link_ds_->is_root_){ continue; }//Root doesn't experience gravity
240  ret_FgravGC += it->link_ds_->mass_ * (it->J_com_.block(0,0,3,dof).transpose() * robot_parsed_data_->gravity_);
241  }
242  return flag;
243  }
244 
245  /* *******************************************************************
246  * Dynamics State functions.
247  * ******************************************************************* */
253  const Eigen::VectorXd& arg_q,
255  const Eigen::VectorXd& arg_dq);
256 
262  const Eigen::VectorXd& arg_q);
263 
264  /* *******************************************************************
265  * Contact Resolution functions.
266  * ******************************************************************* */
269  virtual sBool computeImpulseForce(
271  const SRobotSensors &arg_sensor_data,
274  const SGcModel &arg_gc_model,
276  sFloat arg_coeff_rest,
278  sFloat arg_collision_time,
280  scl::SForceContact &ret_force_);
281 
282  /* *******************************************************************
283  * Initialization functions.
284  * ******************************************************************* */
287 
289  virtual ~CDynamicsScl(){}
290 
299  virtual sBool init(const SRobotParsed& arg_robot_data);
300 };
301 
302 } /* namespace scl */
303 #endif /* CDYNAMICSSCL_HPP_ */
virtual sBool computeJacobian(Eigen::MatrixXd &arg_J, const SRigidBodyDyn &arg_link, const Eigen::VectorXd &arg_q, const Eigen::Vector3d &arg_pos_local) const
Definition: CDynamicsScl.cpp:245
virtual sBool computeTransformToAncestor(Eigen::Affine3d &arg_T, SRigidBodyDyn &arg_link, const SRigidBodyDyn *arg_ancestor, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.cpp:135
Definition: SRobotParsed.hpp:51
Definition: SGcModel.hpp:53
virtual sBool computeJacobianWithTransforms(Eigen::MatrixXd &arg_J, SRigidBodyDyn &arg_link, const Eigen::VectorXd &arg_q, const Eigen::Vector3d &arg_pos_local)
Definition: CDynamicsScl.cpp:212
sBool computeJacobianComForAllLinks(sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q) const
Definition: CDynamicsScl.hpp:184
const SRobotParsed * robot_parsed_data_
Definition: CDynamicsBase.hpp:313
virtual sBool computeTransform(SRigidBodyDyn &arg_link, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.cpp:86
Definition: SForce.hpp:114
virtual sBool computeImpulseForce(const SRobotSensors &arg_sensor_data, const SGcModel &arg_gc_model, sFloat arg_coeff_rest, sFloat arg_collision_time, scl::SForceContact &ret_force_)
Definition: CDynamicsScl.cpp:404
Definition: CDynamicsBase.hpp:55
virtual sBool computeGCModel(const SRobotSensors *arg_sensor_data, SGcModel *arg_gc_model)
Definition: CDynamicsScl.cpp:46
sBool has_been_init_
Definition: CDynamicsBase.hpp:310
bool sBool
Definition: DataTypes.hpp:54
CDynamicsScl()
Definition: CDynamicsScl.hpp:286
sUInt dof_
Definition: SRobotParsed.hpp:83
virtual sBool init(const SRobotParsed &arg_robot_data)
Definition: CDynamicsScl.cpp:423
Eigen::Vector3d gravity_
Definition: SRobotParsed.hpp:86
virtual sFloat computeEnergyKinetic(sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q, const Eigen::VectorXd &arg_dq)
Definition: CDynamicsScl.cpp:333
iterator begin()
Definition: CMappedList.hpp:486
Definition: CDynamicsScl.hpp:57
sBool computeInertiaGC(Eigen::MatrixXd &ret_Mgc, sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.hpp:204
sBool computeForceGravityGC(Eigen::VectorXd &ret_FgravGC, sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.hpp:226
Definition: CMappedTree.hpp:66
virtual ~CDynamicsScl()
Definition: CDynamicsScl.hpp:289
virtual sFloat computeEnergyPotential(sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.cpp:373
double sFloat
Definition: DataTypes.hpp:72
virtual sBool computeTransformsForAllLinks(sutil::CMappedTree< std::string, SRigidBodyDyn > &arg_tree, const Eigen::VectorXd &arg_q)
Definition: CDynamicsScl.cpp:171
virtual sBool computeJacobian(Eigen::MatrixXd &arg_J, SRigidBodyDyn &arg_link, const SRigidBodyDyn *arg_ancestor, const Eigen::VectorXd &arg_q, const Eigen::Vector3d &arg_pos_local, const bool arg_recompute_transforms=true)
Definition: CDynamicsScl.hpp:164
Definition: SRigidBodyDyn.hpp:56
Definition: SRobotIO.hpp:49