robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
unconstr_state_equation.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_UNCONSTR_STATE_EQUATION_HPP_
2#define ROBOTOC_UNCONSTR_STATE_EQUATION_HPP_
3
4#include "Eigen/Core"
5
9
10
11namespace robotoc {
12
21void linearizeUnconstrForwardEuler(const double dt, const SplitSolution& s,
22 const SplitSolution& s_next,
23 SplitKKTMatrix& kkt_matrix,
24 SplitKKTResidual& kkt_residual);
25
32 SplitKKTResidual& kkt_residual);
33
44void linearizeUnconstrBackwardEuler(const double dt,
45 const Eigen::VectorXd& q_prev,
46 const Eigen::VectorXd& v_prev,
47 const SplitSolution& s,
48 const SplitSolution& s_next,
49 SplitKKTMatrix& kkt_matrix,
50 SplitKKTResidual& kkt_residual);
51
62 const Eigen::VectorXd& q_prev,
63 const Eigen::VectorXd& v_prev,
64 const SplitSolution& s,
65 SplitKKTMatrix& kkt_matrix,
66 SplitKKTResidual& kkt_residual);
67
75void evalUnconstrForwardEuler(const double dt, const SplitSolution& s,
76 const SplitSolution& s_next,
77 SplitKKTResidual& kkt_residual);
78
87void evalUnconstrBackwardEuler(const double dt, const Eigen::VectorXd& q_prev,
88 const Eigen::VectorXd& v_prev,
89 const SplitSolution& s,
90 SplitKKTResidual& kkt_residual);
91
92} // namespace robotoc
93
94#endif // ROBOTOC_UNCONSTR_STATE_EQUATION_HPP_
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
Definition: constraint_component_base.hpp:17
void evalUnconstrBackwardEuler(const double dt, const Eigen::VectorXd &q_prev, const Eigen::VectorXd &v_prev, const SplitSolution &s, SplitKKTResidual &kkt_residual)
Computes the residual in the state equation of backward Euler.
void linearizeUnconstrBackwardEulerTerminal(const double dt, const Eigen::VectorXd &q_prev, const Eigen::VectorXd &v_prev, const SplitSolution &s, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual)
Linearizes the state equation of backward Euler at the terminal stage.
void linearizeUnconstrForwardEuler(const double dt, const SplitSolution &s, const SplitSolution &s_next, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual)
Linearizes the state equation of forward Euler.
void linearizeUnconstrForwardEulerTerminal(const SplitSolution &s, SplitKKTResidual &kkt_residual)
Linearizes the state equation of forward Euler.
void linearizeUnconstrBackwardEuler(const double dt, const Eigen::VectorXd &q_prev, const Eigen::VectorXd &v_prev, const SplitSolution &s, const SplitSolution &s_next, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual)
Linearizes the state equation of backward Euler.
void evalUnconstrForwardEuler(const double dt, const SplitSolution &s, const SplitSolution &s_next, SplitKKTResidual &kkt_residual)
Computes the residual in the state equation of forward Euler.