robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
riccati_factorizer.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_RICCATI_FACTORIZER_HPP_
2#define ROBOTOC_RICCATI_FACTORIZER_HPP_
3
4#include "Eigen/Core"
5#include "Eigen/LU"
6
17
18#include <limits>
19#include <cmath>
20
21
22namespace robotoc {
23
29public:
37 RiccatiFactorizer(const Robot& robot, const double max_dts0=0.1);
38
43
48
53
58
62 RiccatiFactorizer(RiccatiFactorizer&&) noexcept = default;
63
67 RiccatiFactorizer& operator=(RiccatiFactorizer&&) noexcept = default;
68
75 void setRegularization(const double max_dts0);
76
86 SplitKKTMatrix& kkt_matrix,
87 SplitKKTResidual& kkt_residual,
89 LQRPolicy& lqr_policy);
90
103 SplitKKTMatrix& kkt_matrix,
104 SplitKKTResidual& kkt_residual,
106 LQRPolicy& lqr_policy, const bool sto,
107 const bool has_next_sto_phase);
108
118 const SplitRiccatiFactorization& riccati,
119 SplitRiccatiFactorization& riccati_m, STOPolicy& sto_policy,
120 const bool has_next_sto_phase) const;
121
130 SplitKKTMatrix& kkt_matrix,
131 SplitKKTResidual& kkt_residual,
133
143 SplitKKTMatrix& kkt_matrix,
144 SplitKKTResidual& kkt_residual,
146 const bool sto);
147
148private:
149 bool has_floating_base_;
150 int dimv_, dimu_;
151 double max_dts0_, eps_;
152 Eigen::LLT<Eigen::MatrixXd> llt_, llt_s_;
153 LQRPolicy lqr_policy_;
154 BackwardRiccatiRecursionFactorizer backward_recursion_;
156
157};
158
159
173 const SplitKKTResidual& kkt_residual,
174 const LQRPolicy& lqr_policy, SplitDirection& d,
175 SplitDirection& d_next, const bool sto,
176 const bool has_next_sto_phase);
177
187 const SplitKKTResidual& kkt_residual,
188 const SplitDirection& d,
189 SplitDirection& d_next);
190
200 const bool has_prev_sto_phase);
201
211 SplitDirection& d, const bool sto,
212 const bool has_next_sto_phase);
213
221 SplitDirection& d, const bool sto);
222
233 const SplitRiccatiFactorization& riccati, SplitDirection& d,
234 const bool sto, const bool has_next_sto_phase);
235
236} // namespace robotoc
237
238#endif // ROBOTOC_RICCATI_FACTORIZER_HPP_
Factorizer of the backward Riccati recursion.
Definition: backward_riccati_recursion_factorizer.hpp:19
The state feedback and feedforward policy of LQR subproblem at a time stage.
Definition: lqr_policy.hpp:16
Riccati factorizer.
Definition: riccati_factorizer.hpp:28
RiccatiFactorizer(const RiccatiFactorizer &)=default
Default copy constructor.
RiccatiFactorizer()
Default constructor.
void backwardRiccatiRecursion(const SplitRiccatiFactorization &riccati_next, SplitKKTMatrix &kkt_matrix, SplitKKTResidual &kkt_residual, SplitRiccatiFactorization &riccati, LQRPolicy &lqr_policy)
Performs the backward Riccati recursion.
void setRegularization(const double max_dts0)
Sets the regularization on the STO.
RiccatiFactorizer & operator=(const RiccatiFactorizer &)=default
Default copy operator.
RiccatiFactorizer(RiccatiFactorizer &&) noexcept=default
Default move constructor.
void backwardRiccatiRecursionPhaseTransition(const SplitRiccatiFactorization &riccati, SplitRiccatiFactorization &riccati_m, STOPolicy &sto_policy, const bool has_next_sto_phase) const
Performs the backward Riccati recursion for the phase transition.
RiccatiFactorizer(const Robot &robot, const double max_dts0=0.1)
Constructs a factorizer.
~RiccatiFactorizer()
Destructor.
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
The state feedback and feedforward policy of the switching time optimization (STO).
Definition: sto_policy.hpp:16
Riccati factorization matrix and vector for the switching constraint.
Definition: split_constrained_riccati_factorization.hpp:17
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
Riccati factorization matrix and vector for a time stage.
Definition: split_riccati_factorization.hpp:15
Definition: constraint_component_base.hpp:17
void computeSwitchingTimeDirection(const STOPolicy &sto_policy, SplitDirection &d, const bool has_prev_sto_phase)
Computes the switching time direction.
void computeLagrangeMultiplierDirection(const SplitRiccatiFactorization &riccati, SplitDirection &d, const bool sto, const bool has_next_sto_phase)
Computes the Newton direction of the Lagrange multiplier with respect to the switching constraint.
void computeCostateDirection(const SplitRiccatiFactorization &riccati, SplitDirection &d, const bool sto, const bool has_next_sto_phase)
Computes the Newton direction of the costate.
void forwardRiccatiRecursion(const SplitKKTMatrix &kkt_matrix, const SplitKKTResidual &kkt_residual, const LQRPolicy &lqr_policy, SplitDirection &d, SplitDirection &d_next, const bool sto, const bool has_next_sto_phase)
Performs the forward Riccati recursion and computes the state direction.