aligator  0.9.0
A primal-dual augmented Lagrangian-type solver for nonlinear trajectory optimization.
Loading...
Searching...
No Matches
riccati-impl.hpp
Go to the documentation of this file.
1
3#pragma once
4
5#include "blk-matrix.hpp"
6
7#include <proxsuite-nlp/linalg/bunchkaufman.hpp>
8#include <Eigen/LU>
9#include <Eigen/Cholesky>
10
12
13#include <optional>
14
15namespace aligator {
16namespace gar {
17template <typename Scalar> struct LQRKnotTpl;
18template <typename Scalar> struct LQRProblemTpl;
19
21template <class T, class A>
22inline boost::span<T> make_span_from_indices(std::vector<T, A> &vec, size_t i0,
23 size_t i1) {
24 return boost::make_span(vec.data() + i0, i1 - i0);
25}
26
28template <class T, class A>
29inline boost::span<const T> make_span_from_indices(const std::vector<T, A> &vec,
30 size_t i0, size_t i1) {
31 return boost::make_span(vec.data() + i0, i1 - i0);
32}
33
35template <typename _Scalar> struct StageFactor {
36 using Scalar = _Scalar;
38 using RowMatrixXs = Eigen::Matrix<Scalar, -1, -1, Eigen::RowMajor>;
39
40 struct value_t {
41 MatrixXs Pmat; //< Riccati matrix
42 VectorXs pvec; //< Riccati bias
43 MatrixXs Vxx; //< "cost-to-go" matrix
44 VectorXs vx; //< "cost-to-go" gradient
45 MatrixXs Vxt;
46 MatrixXs Vtt;
47 VectorXs vt;
48
49 value_t(uint nx, uint nth)
50 : Pmat(nx, nx), pvec(nx), Vxx(nx, nx), vx(nx), Vxt(nx, nth),
51 Vtt(nth, nth), vt(nth) {
52 Vxt.setZero();
53 Vtt.setZero();
54 vt.setZero();
55 }
56 };
57
58 StageFactor(uint nx, uint nu, uint nc, uint nx2, uint nth);
59
60 MatrixXs Qhat;
61 MatrixXs Rhat;
62 MatrixXs Shat;
63 VectorXs qhat;
64 VectorXs rhat;
67
68 // Parametric
69 MatrixXs Gxhat;
70 MatrixXs Guhat;
71
72 BlkMatrix<VectorXs, 4, 1> ff; //< feedforward gains
73 BlkMatrix<RowMatrixXs, 4, 1> fb; //< feedback gains
74 BlkMatrix<RowMatrixXs, 4, 1> fth; //< parameter feedback gains
75 BlkMatrix<MatrixXs, 2, 2> kktMat; //< reduced KKT matrix buffer
76 Eigen::BunchKaufman<MatrixXs> kktChol; //< reduced KKT LDLT solver
77 Eigen::PartialPivLU<MatrixXs> Efact; //< LU decomp. of E matrix
78 VectorXs yff_pre;
79 MatrixXs A_pre;
80 MatrixXs Yth_pre;
81 MatrixXs Ptilde; //< product Et.inv P * E.inv
82 MatrixXs Einv; //< product P * E.inv
83 MatrixXs EinvP; //< product P * E.inv
84 MatrixXs schurMat; //< Dual-space Schur matrix
85 Eigen::LLT<MatrixXs> schurChol; //< Cholesky decomposition of Schur matrix
86 value_t vm; //< cost-to-go parameters
87};
88
91template <typename Scalar> struct ProximalRiccatiKernel {
93 using RowMatrixXs = Eigen::Matrix<Scalar, -1, -1, Eigen::RowMajor>;
94 using RowMatrixRef = Eigen::Ref<RowMatrixXs>;
95 using ConstRowMatrixRef = Eigen::Ref<const RowMatrixXs>;
99
100 struct kkt0_t {
104 Eigen::BunchKaufman<MatrixXs> chol{mat.rows()};
105 kkt0_t(uint nx, uint nc, uint nth)
106 : mat({nx, nc}, {nx, nc}), ff(mat.rowDims()),
107 fth(mat.rowDims(), {nth}) {}
108 };
109
110 inline static void terminalSolve(const KnotType &model, const Scalar mueq,
111 StageFactorType &d);
112
113 inline static bool backwardImpl(boost::span<const KnotType> stages,
114 const Scalar mudyn, const Scalar mueq,
116
118 inline static void
119 computeInitial(VectorRef x0, VectorRef lbd0, const kkt0_t &kkt0,
120 const std::optional<ConstVectorRef> &theta_);
121
122 inline static void stageKernelSolve(const KnotType &model, StageFactorType &d,
123 value_t &vn, const Scalar mudyn,
124 const Scalar mueq);
125
127 inline static bool
132 const std::optional<ConstVectorRef> &theta_ = std::nullopt);
133};
134
135} // namespace gar
136} // namespace aligator
137
138#include "./riccati-impl.hxx"
139
140#ifdef ALIGATOR_ENABLE_TEMPLATE_INSTANTIATION
141#include "./riccati-impl.txx"
142#endif
Block matrix class, with a fixed-size number of row and column blocks.
const row_dim_t & rowDims() const
unsigned int uint
Definition work.hpp:7
boost::span< T > make_span_from_indices(std::vector< T, A > &vec, size_t i0, size_t i1)
Create a boost::span object from a vector and two indices.
Main package namespace.
constexpr span< I > make_span(I *f, std::size_t c) noexcept
Definition make_span.hpp:17
Struct describing a stage of a constrained LQ problem.
Eigen::BunchKaufman< MatrixXs > chol
Kernel for use in Riccati-like algorithms for the proximal LQ subproblem.
static bool backwardImpl(boost::span< const KnotType > stages, const Scalar mudyn, const Scalar mueq, boost::span< StageFactorType > datas)
Eigen::Ref< const RowMatrixXs > ConstRowMatrixRef
static void terminalSolve(const KnotType &model, const Scalar mueq, StageFactorType &d)
static bool forwardImpl(boost::span< const KnotType > stages, boost::span< const StageFactorType > datas, boost::span< VectorXs > xs, boost::span< VectorXs > us, boost::span< VectorXs > vs, boost::span< VectorXs > lbdas, const std::optional< ConstVectorRef > &theta_=std::nullopt)
Forward sweep.
typename StageFactor< Scalar >::value_t value_t
static void stageKernelSolve(const KnotType &model, StageFactorType &d, value_t &vn, const Scalar mudyn, const Scalar mueq)
static void computeInitial(VectorRef x0, VectorRef lbd0, const kkt0_t &kkt0, const std::optional< ConstVectorRef > &theta_)
Solve initial stage.
Eigen::Matrix< Scalar, -1, -1, Eigen::RowMajor > RowMatrixXs
Eigen::Ref< RowMatrixXs > RowMatrixRef
Per-node struct for all computations in the factorization.
BlkMatrix< RowMatrixXs, 4, 1 > fth
Eigen::Matrix< Scalar, -1, -1, Eigen::RowMajor > RowMatrixXs
Eigen::PartialPivLU< MatrixXs > Efact
BlkMatrix< VectorXs, 4, 1 > ff
Eigen::LLT< MatrixXs > schurChol
BlkMatrix< RowMatrixXs, 4, 1 > fb
Eigen::BunchKaufman< MatrixXs > kktChol
BlkMatrix< MatrixXs, 2, 2 > kktMat
StageFactor(uint nx, uint nu, uint nc, uint nx2, uint nth)