6#include <Eigen/SparseCore>
12template <
typename InType,
typename OutScalar>
14 const Eigen::DenseBase<InType> &input,
15 Eigen::SparseMatrix<OutScalar> &out,
bool update) {
16 assert(i0 + input.rows() <= out.rows() &&
"Inconsistent rows");
17 assert(j0 + input.cols() <= out.cols() &&
"Inconsistent cols");
19 for (Index i = 0; i < input.rows(); i++) {
20 for (Index j = 0; j < input.cols(); j++) {
22 out.coeffRef(i0 + i, j0 + j) = input(i, j);
24 out.insert(i0 + i, j0 + j) = input(i, j);
29template <
typename Scalar>
31 Eigen::SparseMatrix<Scalar> &out,
bool update) {
33 assert(i0 <= i1 &&
"i0 should be lesser than i1. Can't assign empty range.");
34 assert(i1 <= out.rows() &&
"Inconsistent rows");
35 assert(i1 <= out.cols() &&
"Inconsistent cols");
36 for (Index kk = i0; kk < i1; kk++) {
38 out.coeffRef(kk, kk) = value;
40 out.insert(kk, kk) = value;
45template <
typename Scalar>
47 const Scalar mudyn,
const Scalar mueq,
48 Eigen::SparseMatrix<Scalar> &mat,
49 Eigen::Matrix<Scalar, -1, 1> &rhs,
bool update);
51template <
typename Scalar>
54 boost::span<
const typename math_types<Scalar>::VectorXs> xs,
55 boost::span<
const typename math_types<Scalar>::VectorXs> us,
56 boost::span<
const typename math_types<Scalar>::VectorXs> vs,
57 boost::span<
const typename math_types<Scalar>::VectorXs> lbdas,
58 const Scalar mudyn,
const Scalar mueq,
59 const std::optional<
typename math_types<Scalar>::ConstVectorRef> &theta_,
60 bool verbose =
false);
65template <
typename Scalar>
67 Scalar mueq,
typename math_types<Scalar>::MatrixXs &mat,
68 typename math_types<Scalar>::VectorXs &rhs);
71template <
typename Scalar>
73 const auto &knots = problem.
stages;
75 const size_t N = knots.size() - 1UL;
77 for (
size_t t = 0; t <= N; t++) {
78 const auto &model = knots[t];
79 nrows += model.nx + model.nu + model.nc;
87template <
typename Scalar>
91 decltype(
auto) knots = problem.
stages;
92 using MatrixXs =
typename math_types<Scalar>::MatrixXs;
93 using VectorXs =
typename math_types<Scalar>::VectorXs;
96 MatrixXs mat(nrows, nrows);
100 fmt::print(
"{:s} WARNING! Problem was not initialized.", __FUNCTION__);
102 return std::make_pair(mat, rhs);
106template <
typename Scalar>
109 const typename math_types<Scalar>::ConstVectorRef solution,
110 std::vector<
typename math_types<Scalar>::VectorXs> &xs,
111 std::vector<
typename math_types<Scalar>::VectorXs> &us,
112 std::vector<
typename math_types<Scalar>::VectorXs> &vs,
113 std::vector<
typename math_types<Scalar>::VectorXs> &lbdas) {
119 const uint nc0 = problem.
nc0();
121 lbdas[0] = solution.head(nc0);
124 for (
size_t t = 0; t <= N; t++) {
127 auto seg = solution.segment(idx, n);
128 xs[t] = seg.head(knot.
nx);
129 us[t] = seg.segment(knot.
nx, knot.
nu);
130 vs[t] = seg.segment(knot.
nx + knot.
nu, knot.
nc);
133 lbdas[t + 1] = solution.segment(idx, knot.
nx2);
139template <
typename Scalar>
141 using VectorXs =
typename math_types<Scalar>::VectorXs;
143 std::vector<VectorXs> xs;
144 std::vector<VectorXs> us;
145 std::vector<VectorXs> vs;
146 std::vector<VectorXs> lbdas;
154 lbdas[0].setZero(problem.
nc0());
155 for (
uint i = 0; i <= N; i++) {
156 const knot_t &kn = problem.
stages[i];
157 xs[i].setZero(kn.nx);
158 us[i].setZero(kn.nu);
159 vs[i].setZero(kn.nc);
162 lbdas[i + 1].setZero(kn.nx2);
164 if (problem.
stages.back().nu == 0) {
167 return std::make_tuple(std::move(xs), std::move(us), std::move(vs),
171#ifdef ALIGATOR_ENABLE_TEMPLATE_INSTANTIATION
173 const LQRProblemTpl<context::Scalar> &problem,
const context::Scalar mudyn,
175 context::VectorXs &rhs,
bool update);
176extern template std::array<context::Scalar, 3>
178 const LQRProblemTpl<context::Scalar> &,
182 const std::optional<context::ConstVectorRef> &,
bool);
void sparseAssignDiagonal(Eigen::Index i0, Eigen::Index i1, Scalar value, Eigen::SparseMatrix< Scalar > &out, bool update)
void sparseAssignDenseBlock(Eigen::Index i0, Eigen::Index j0, const Eigen::DenseBase< InType > &input, Eigen::SparseMatrix< OutScalar > &out, bool update)
Helper to assign a dense matrix into a range of coefficients of a sparse matrix.
auto lqrInitializeSolution(const LQRProblemTpl< Scalar > &problem)
bool lqrDenseMatrix(const LQRProblemTpl< Scalar > &problem, Scalar mudyn, Scalar mueq, typename math_types< Scalar >::MatrixXs &mat, typename math_types< Scalar >::VectorXs &rhs)
Fill in a KKT constraint matrix and vector for the given LQ problem with the given dual-regularizatio...
template auto lqrDenseMatrix< context::Scalar >(const LQRProblemTpl< context::Scalar > &, context::Scalar, context::Scalar)
void lqrDenseSolutionToTraj(const LQRProblemTpl< Scalar > &problem, const typename math_types< Scalar >::ConstVectorRef solution, std::vector< typename math_types< Scalar >::VectorXs > &xs, std::vector< typename math_types< Scalar >::VectorXs > &us, std::vector< typename math_types< Scalar >::VectorXs > &vs, std::vector< typename math_types< Scalar >::VectorXs > &lbdas)
Convert dense RHS solution to its trajectory [x,u,v,lambda] solution.
void lqrCreateSparseMatrix(const LQRProblemTpl< Scalar > &problem, const Scalar mudyn, const Scalar mueq, Eigen::SparseMatrix< Scalar > &mat, Eigen::Matrix< Scalar, -1, 1 > &rhs, bool update)
uint lqrNumRows(const LQRProblemTpl< Scalar > &problem)
Compute the number of rows in the problem matrix.
template std::array< context::Scalar, 3 > lqrComputeKktError< context::Scalar >(const LQRProblemTpl< context::Scalar > &, boost::span< const context::VectorXs >, boost::span< const context::VectorXs >, boost::span< const context::VectorXs >, boost::span< const context::VectorXs >, const context::Scalar, const context::Scalar, const std::optional< context::ConstVectorRef > &, bool)
std::array< Scalar, 3 > lqrComputeKktError(const LQRProblemTpl< Scalar > &problem, boost::span< const typename math_types< Scalar >::VectorXs > xs, boost::span< const typename math_types< Scalar >::VectorXs > us, boost::span< const typename math_types< Scalar >::VectorXs > vs, boost::span< const typename math_types< Scalar >::VectorXs > lbdas, const Scalar mudyn, const Scalar mueq, const std::optional< typename math_types< Scalar >::ConstVectorRef > &theta_, bool verbose=false)
template void lqrCreateSparseMatrix< context::Scalar >(const LQRProblemTpl< context::Scalar > &problem, const context::Scalar mudyn, const context::Scalar mueq, Eigen::SparseMatrix< context::Scalar > &mat, context::VectorXs &rhs, bool update)
Struct describing a stage of a constrained LQ problem.
uint nc0() const noexcept
int horizon() const noexcept