robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
unconstr_split_backward_correction.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_UNCONSTR_SPLIT_BACKWARD_CORRECTION_HPP_
2#define ROBOTOC_UNCONSTR_SPLIT_BACKWARD_CORRECTION_HPP_
3
4#include <vector>
5
6#include "Eigen/Core"
7
14
15
16namespace robotoc {
17
23public:
29
34
39
44 const UnconstrSplitBackwardCorrection&) = default;
45
50 const UnconstrSplitBackwardCorrection&) = default;
51
56 UnconstrSplitBackwardCorrection&&) noexcept = default;
57
62 UnconstrSplitBackwardCorrection&&) noexcept = default;
63
73 void coarseUpdate(const Eigen::MatrixXd& aux_mat_next,
74 const double dt, const SplitKKTMatrix& kkt_matrix,
75 const SplitKKTResidual& kkt_residual,
76 const SplitSolution& s, SplitSolution& s_new);
77
87 void coarseUpdate(const Eigen::Block<Eigen::MatrixXd>& aux_mat_next,
88 const double dt, const SplitKKTMatrix& kkt_matrix,
89 const SplitKKTResidual& kkt_residual,
90 const SplitSolution& s, SplitSolution& s_new);
91
101 void coarseUpdate(const double dt, const SplitKKTMatrix& kkt_matrix,
102 const SplitKKTResidual& kkt_residual,
103 const SplitSolution& s, SplitSolution& s_new);
104
109 const Eigen::Block<const Eigen::MatrixXd> auxMat() const;
110
118 const SplitSolution& s_new_next,
119 SplitSolution& s_new);
120
126
135 const SplitSolution& s_new_prev,
136 SplitSolution& s_new);
137
143
150 static void computeDirection(const SplitSolution& s,
151 const SplitSolution& s_new,
152 SplitDirection& d);
153
154private:
155 int dimv_, dimx_, dimkkt_;
156 UnconstrKKTMatrixInverter kkt_mat_inverter_;
157 Eigen::MatrixXd H_, kkt_mat_inv_;
158 Eigen::VectorXd kkt_res_, d_, x_res_, dx_;
159
160
161 template <typename MatrixType>
162 void coarseUpdate_impl(const Eigen::MatrixBase<MatrixType>& aux_mat_next,
163 const double dt, const SplitKKTMatrix& kkt_matrix,
164 const SplitKKTResidual& kkt_residual,
165 const SplitSolution& s, SplitSolution& s_new) {
166 assert(aux_mat_next.rows() == dimx_);
167 assert(aux_mat_next.cols() == dimx_);
168 H_.topLeftCorner(dimv_, dimv_) = kkt_matrix.Qaa;
169 H_.bottomLeftCorner(dimx_, dimv_) = kkt_matrix.Qxu; // This is actually Qxa
170 H_.bottomRightCorner(dimx_, dimx_) = kkt_matrix.Qxx.transpose();
171 H_.bottomRightCorner(dimx_, dimx_).noalias() += aux_mat_next;
172 kkt_mat_inverter_.invert(dt, H_, kkt_mat_inv_);
173 kkt_res_.head(dimx_) = kkt_residual.Fx;
174 kkt_res_.segment(dimx_, dimv_) = kkt_residual.la;
175 kkt_res_.tail(dimx_) = kkt_residual.lx;
176 d_.noalias() = kkt_mat_inv_ * kkt_res_;
177 s_new.lmd = s.lmd - d_.head(dimv_);
178 s_new.gmm = s.gmm - d_.segment(dimv_, dimv_);
179 s_new.a = s.a - d_.segment(2*dimv_, dimv_);
180 s_new.q = s.q - d_.segment(3*dimv_, dimv_);
181 s_new.v = s.v - d_.tail(dimv_);
182 }
183
184};
185
186} // namespace robotoc
187
188#endif // ROBOTOC_UNCONSTR_SPLIT_BACKWARD_CORRECTION_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
Schur complement for SplitKKTMatrix for UnconstrParNMPC.
Definition: unconstr_kkt_matrix_inverter.hpp:18
void invert(const double dt, const Eigen::MatrixBase< MatrixType1 > &H, const Eigen::MatrixBase< MatrixType2 > &KKT_mat_inv)
Computes the inverse of the split KKT matrix of the time stage.
Definition: unconstr_kkt_matrix_inverter.hxx:40
Split backward correction.
Definition: unconstr_split_backward_correction.hpp:22
UnconstrSplitBackwardCorrection(const Robot &robot)
Construct split backward correction.
static void computeDirection(const SplitSolution &s, const SplitSolution &s_new, SplitDirection &d)
Computes the direction from the results of the forward correction.
void forwardCorrectionParallel(SplitSolution &s_new)
Performs the parallel part of the forward correction.
UnconstrSplitBackwardCorrection(const UnconstrSplitBackwardCorrection &)=default
Default copy constructor.
UnconstrSplitBackwardCorrection()
Default constructor.
void backwardCorrectionParallel(SplitSolution &s_new)
Performs the parallel part of the backward correction.
void backwardCorrectionSerial(const SplitSolution &s_next, const SplitSolution &s_new_next, SplitSolution &s_new)
Performs the serial part of the backward correction.
void forwardCorrectionSerial(const SplitSolution &s_prev, const SplitSolution &s_new_prev, SplitSolution &s_new)
Performs the serial part of the forward correction.
void coarseUpdate(const Eigen::MatrixXd &aux_mat_next, const double dt, const SplitKKTMatrix &kkt_matrix, const SplitKKTResidual &kkt_residual, const SplitSolution &s, SplitSolution &s_new)
Coarse updates the split solution of this time stage.
UnconstrSplitBackwardCorrection(UnconstrSplitBackwardCorrection &&) noexcept=default
Default move constructor.
UnconstrSplitBackwardCorrection & operator=(const UnconstrSplitBackwardCorrection &)=default
Default copy operator.
const Eigen::Block< const Eigen::MatrixXd > auxMat() const
Auxiliary matrix of this time stage.
Definition: constraint_component_base.hpp:17