aligator  0.6.1
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#pragma once
2
3#include <proxsuite-nlp/linalg/bunchkaufman.hpp>
4
5#include "blk-matrix.hpp"
6#include "lqr-problem.hpp"
7#include "riccati-base.hpp"
8#include <tracy/Tracy.hpp>
9
10namespace aligator::gar {
11
13template <typename _Scalar>
14class RiccatiSolverDense : public RiccatiSolverBase<_Scalar> {
15public:
16 using Scalar = _Scalar;
18 using Base = RiccatiSolverBase<Scalar>;
19 using BlkMat44 = BlkMatrix<MatrixXs, 4, 4>;
20 using BlkRowMat41 = BlkMatrix<RowMatrixXs, 4, 1>;
21 using BlkVec4 = BlkMatrix<VectorXs, 4, 1>;
22 using KnotType = LQRKnotTpl<Scalar>;
23
24 struct FactorData {
29 Eigen::BunchKaufman<MatrixXs> ldl;
30 };
31
32 static FactorData init_factor(const LQRKnotTpl<Scalar> &knot) {
33 std::array<long, 4> dims = {knot.nu, knot.nc, knot.nx2, knot.nx2};
34 using ldl_t = decltype(FactorData::ldl);
35 long ntot = std::accumulate(dims.begin(), dims.end(), 0);
36 uint nth = knot.nth;
37 return FactorData{BlkMat44(dims, dims), BlkVec4(dims, {1}),
38 BlkRowMat41(dims, {knot.nx}), BlkRowMat41(dims, {nth}),
39 ldl_t{ntot}};
40 }
41
42 std::vector<FactorData> datas;
43 std::vector<MatrixXs> Pxx;
44 std::vector<MatrixXs> Pxt;
45 std::vector<MatrixXs> Ptt;
46 std::vector<VectorXs> px;
47 std::vector<VectorXs> pt;
48 struct {
49 BlkMatrix<MatrixXs, 2, 2> mat;
50 BlkMatrix<VectorXs, 2, 1> ff;
51 BlkMatrix<MatrixXs, 2, 1> fth; // parametric rhs
52 Eigen::BunchKaufman<MatrixXs> ldl;
54 VectorXs thGrad;
55 MatrixXs thHess;
56
57 explicit RiccatiSolverDense(const LQRProblemTpl<Scalar> &problem);
58
59 bool backward(const Scalar mudyn, const Scalar mueq);
60
61 bool forward(std::vector<VectorXs> &xs, std::vector<VectorXs> &us,
62 std::vector<VectorXs> &vs, std::vector<VectorXs> &lbdas,
63 const std::optional<ConstVectorRef> &theta = std::nullopt) const;
64
65 VectorRef getFeedforward(size_t i) { return datas[i].ff.matrix(); }
66 RowMatrixRef getFeedback(size_t i) { return datas[i].fb.matrix(); }
67
68protected:
69 void initialize();
70 const LQRProblemTpl<Scalar> *problem_;
71};
72
73} // namespace aligator::gar
74
75#include "dense-riccati.hxx"
76
77#ifdef ALIGATOR_ENABLE_TEMPLATE_INSTANTIATION
78#include "dense-riccati.txx"
79#endif
A stagewise-dense Riccati solver.
BlkMatrix< VectorXs, 4, 1 > BlkVec4
struct aligator::gar::RiccatiSolverDense::@2 kkt0
BlkMatrix< VectorXs, 2, 1 > ff
BlkMatrix< RowMatrixXs, 4, 1 > BlkRowMat41
static FactorData 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, 4, 4 > BlkMat44
BlkMatrix< MatrixXs, 2, 2 > mat
RiccatiSolverDense(const LQRProblemTpl< Scalar > &problem)
BlkMatrix< MatrixXs, 2, 1 > fth
std::vector< FactorData > datas
Eigen::BunchKaufman< MatrixXs > ldl
const LQRProblemTpl< Scalar > * problem_
RowMatrixRef getFeedback(size_t i)
bool backward(const Scalar mudyn, const Scalar mueq)
unsigned int uint
Definition logger.hpp:11