robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
state_equation.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_STATE_EQUATION_HPP_
2#define ROBOTOC_STATE_EQUATION_HPP_
3
4#include "Eigen/Core"
5
13
14
15namespace robotoc {
16
26void evalStateEquation(const Robot& robot, const double dt,
27 const SplitSolution& s,
28 const Eigen::VectorXd& q_next,
29 const Eigen::VectorXd& v_next,
30 SplitKKTResidual& kkt_residual);
31
40void evalStateEquation(const Robot& robot, const double dt,
41 const SplitSolution& s,
42 const SplitSolution& s_next,
43 SplitKKTResidual& kkt_residual);
44
56void linearizeStateEquation(const Robot& robot, const double dt,
57 const Eigen::VectorXd& q_prev,
58 const SplitSolution& s, const SplitSolution& s_next,
59 StateEquationData& data, SplitKKTMatrix& kkt_matrix,
60 SplitKKTResidual& kkt_residual);
61
73void correctLinearizeStateEquation(const Robot& robot, const double dt,
74 const SplitSolution& s,
75 const SplitSolution& s_next,
76 StateEquationData& data,
77 SplitKKTMatrix& kkt_matrix,
78 SplitKKTResidual& kkt_residual);
79
86
98 const Eigen::VectorXd& q0,
99 const Eigen::VectorXd& v0,
100 const SplitSolution& s0,
101 const StateEquationData& data,
102 SplitDirection& d0);
103
104} // namespace robotoc
105
106#endif // ROBOTOC_STATE_EQUATION_HPP_
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
Newton direction of the solution to the optimal control problem split into a time stage.
Definition: split_direction.hpp:20
The KKT matrix split into a time stage.
Definition: split_kkt_matrix.hpp:18
KKT residual split into each time stage.
Definition: split_kkt_residual.hpp:18
Solution to the optimal control problem split into a time stage.
Definition: split_solution.hpp:20
Data for the state equations.
Definition: state_equation_data.hpp:16
Definition: constraint_component_base.hpp:17
void linearizeStateEquation(const Robot &robot, const double dt, const Eigen::VectorXd &q_prev, const SplitSolution &s, const SplitSolution &s_next, StateEquationData &data, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual)
Linearizes the state equation.
void evalStateEquation(const Robot &robot, const double dt, const SplitSolution &s, const Eigen::VectorXd &q_next, const Eigen::VectorXd &v_next, SplitKKTResidual &kkt_residual)
Computes the residual in the state equation.
void correctLinearizeStateEquation(const Robot &robot, const double dt, const SplitSolution &s, const SplitSolution &s_next, StateEquationData &data, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual)
Corrects the linearized state equation using the Jacobian of the Lie group.
void correctCostateDirection(StateEquationData &data, SplitDirection &d)
Corrects the costate direction using the Jacobian of the Lie group.
void computeInitialStateDirection(const Robot &robot, const Eigen::VectorXd &q0, const Eigen::VectorXd &v0, const SplitSolution &s0, const StateEquationData &data, SplitDirection &d0)
Computes the initial state direction using the result of linearizeStateEquation() and correctLinear...