aligator  0.9.0
A primal-dual augmented Lagrangian-type solver for nonlinear trajectory optimization.
Loading...
Searching...
No Matches
dense-riccati.hpp
Go to the documentation of this file.
1
3#pragma once
4
5#include <proxsuite-nlp/linalg/bunchkaufman.hpp>
6
7#include "blk-matrix.hpp"
8#include "riccati-base.hpp"
9
10namespace aligator::gar {
11
16template <typename _Scalar>
17class RiccatiSolverDense : public RiccatiSolverBase<_Scalar> {
18public:
19 using Scalar = _Scalar;
26
27 std::vector<BlkMat44> kkts;
28 std::vector<BlkVec4> ffs;
29 std::vector<BlkRowMat41> fbs;
30 std::vector<BlkRowMat41> fts;
31 std::vector<Eigen::BunchKaufman<MatrixXs>> ldls;
32 std::vector<MatrixXs> Pxx;
33 std::vector<MatrixXs> Pxt;
34 std::vector<MatrixXs> Ptt;
35 std::vector<VectorXs> px;
36 std::vector<VectorXs> pt;
37 struct {
40 BlkMatrix<MatrixXs, 2, 1> fth; // parametric rhs
41 Eigen::BunchKaufman<MatrixXs> ldl;
43 VectorXs thGrad;
44 MatrixXs thHess;
45
46 explicit RiccatiSolverDense(const LQRProblemTpl<Scalar> &problem);
47
48 bool backward(const Scalar mudyn, const Scalar mueq);
49
50 bool forward(std::vector<VectorXs> &xs, std::vector<VectorXs> &us,
51 std::vector<VectorXs> &vs, std::vector<VectorXs> &lbdas,
52 const std::optional<ConstVectorRef> &theta = std::nullopt) const;
53
54 void cycleAppend(const KnotType &knot);
55 VectorRef getFeedforward(size_t i) { return ffs[i].matrix(); }
56 RowMatrixRef getFeedback(size_t i) { return fbs[i].matrix(); }
57
58protected:
60 void initialize();
62};
63
64} // namespace aligator::gar
65
66#ifdef ALIGATOR_ENABLE_TEMPLATE_INSTANTIATION
67#include "dense-riccati.txx"
68#endif
Block matrix class, with a fixed-size number of row and column blocks.
A stagewise-dense Riccati solver. This algorithm uses a dense Bunch-Kaufman factorization at every st...
void cycleAppend(const KnotType &knot)
BlkMatrix< VectorXs, 2, 1 > ff
std::vector< BlkRowMat41 > fts
std::vector< Eigen::BunchKaufman< MatrixXs > > ldls
void init_factor(const LQRKnotTpl< Scalar > &knot)
bool forward(std::vector< VectorXs > &xs, std::vector< VectorXs > &us, std::vector< VectorXs > &vs, std::vector< VectorXs > &lbdas, const std::optional< ConstVectorRef > &theta=std::nullopt) const
BlkMatrix< MatrixXs, 2, 2 > mat
RiccatiSolverDense(const LQRProblemTpl< Scalar > &problem)
BlkMatrix< MatrixXs, 2, 1 > fth
struct aligator::gar::RiccatiSolverDense::@1 kkt0
std::vector< BlkRowMat41 > fbs
Eigen::BunchKaufman< MatrixXs > ldl
const LQRProblemTpl< Scalar > * problem_
RowMatrixRef getFeedback(size_t i)
bool backward(const Scalar mudyn, const Scalar mueq)
Struct describing a stage of a constrained LQ problem.