Directory: | ./ |
---|---|
File: | include/crocoddyl/core/actions/lqr.hpp |
Date: | 2025-03-26 19:23:43 |
Exec | Total | Coverage | |
---|---|---|---|
Lines: | 18 | 35 | 51.4% |
Branches: | 32 | 66 | 48.5% |
Line | Branch | Exec | Source |
---|---|---|---|
1 | /////////////////////////////////////////////////////////////////////////////// | ||
2 | // BSD 3-Clause License | ||
3 | // | ||
4 | // Copyright (C) 2019-2025, LAAS-CNRS, University of Edinburgh | ||
5 | // Heriot-Watt University | ||
6 | // Copyright note valid unless otherwise stated in individual files. | ||
7 | // All rights reserved. | ||
8 | /////////////////////////////////////////////////////////////////////////////// | ||
9 | |||
10 | #ifndef CROCODDYL_CORE_ACTIONS_LQR_HPP_ | ||
11 | #define CROCODDYL_CORE_ACTIONS_LQR_HPP_ | ||
12 | |||
13 | #include <stdexcept> | ||
14 | |||
15 | #include "crocoddyl/core/action-base.hpp" | ||
16 | #include "crocoddyl/core/fwd.hpp" | ||
17 | #include "crocoddyl/core/states/euclidean.hpp" | ||
18 | |||
19 | namespace crocoddyl { | ||
20 | |||
21 | /** | ||
22 | * @brief Linear-quadratic regulator (LQR) action model | ||
23 | * | ||
24 | * A linear-quadratic regulator (LQR) action has a transition model of the form | ||
25 | * \f[ \begin{equation} | ||
26 | * \mathbf{x}^' = \mathbf{A x + B u + f}. | ||
27 | * \end{equation} \f] | ||
28 | * Its cost function is quadratic of the form: | ||
29 | * \f[ \begin{equation} | ||
30 | * \ell(\mathbf{x},\mathbf{u}) = \begin{bmatrix}1 | ||
31 | * \\ \mathbf{x} \\ \mathbf{u}\end{bmatrix}^T \begin{bmatrix}0 & | ||
32 | * \mathbf{q}^T & \mathbf{r}^T \\ \mathbf{q} & \mathbf{Q} | ||
33 | * & | ||
34 | * \mathbf{N}^T \\ | ||
35 | * \mathbf{r} & \mathbf{N} & \mathbf{R}\end{bmatrix} | ||
36 | * \begin{bmatrix}1 \\ \mathbf{x} \\ | ||
37 | * \mathbf{u}\end{bmatrix} | ||
38 | * \end{equation} \f] | ||
39 | * and the linear equality and inequality constraints has the form: | ||
40 | * \f[ \begin{aligned} | ||
41 | * \mathbf{g(x,u)} = \mathbf{G}\begin{bmatrix} \mathbf{x} \\ \mathbf{u} | ||
42 | * \end{bmatrix} [x,u] + \mathbf{g} \leq \mathbf{0} | ||
43 | * &\mathbf{h(x,u)} = \mathbf{H}\begin{bmatrix} \mathbf{x} \\ \mathbf{u} | ||
44 | * \end{bmatrix} [x,u] + \mathbf{h} \end{aligned} \f] | ||
45 | */ | ||
46 | template <typename _Scalar> | ||
47 | class ActionModelLQRTpl : public ActionModelAbstractTpl<_Scalar> { | ||
48 | public: | ||
49 | EIGEN_MAKE_ALIGNED_OPERATOR_NEW | ||
50 | ✗ | CROCODDYL_DERIVED_CAST(ActionModelBase, ActionModelLQRTpl) | |
51 | |||
52 | typedef _Scalar Scalar; | ||
53 | typedef ActionDataAbstractTpl<Scalar> ActionDataAbstract; | ||
54 | typedef ActionModelAbstractTpl<Scalar> Base; | ||
55 | typedef ActionDataLQRTpl<Scalar> Data; | ||
56 | typedef StateVectorTpl<Scalar> StateVector; | ||
57 | typedef MathBaseTpl<Scalar> MathBase; | ||
58 | typedef typename MathBase::VectorXs VectorXs; | ||
59 | typedef typename MathBase::MatrixXs MatrixXs; | ||
60 | |||
61 | /** | ||
62 | * @brief Initialize the LQR action model | ||
63 | * | ||
64 | * @param[in] A State matrix | ||
65 | * @param[in] B Input matrix | ||
66 | * @param[in] Q State weight matrix | ||
67 | * @param[in] R Input weight matrix | ||
68 | * @param[in] N State-input weight matrix | ||
69 | */ | ||
70 | ActionModelLQRTpl(const MatrixXs& A, const MatrixXs& B, const MatrixXs& Q, | ||
71 | const MatrixXs& R, const MatrixXs& N); | ||
72 | |||
73 | /** | ||
74 | * @brief Initialize the LQR action model | ||
75 | * | ||
76 | * @param[in] A State matrix | ||
77 | * @param[in] B Input matrix | ||
78 | * @param[in] Q State weight matrix | ||
79 | * @param[in] R Input weight matrix | ||
80 | * @param[in] N State-input weight matrix | ||
81 | * @param[in] f Dynamics drift | ||
82 | * @param[in] q State weight vector | ||
83 | * @param[in] r Input weight vector | ||
84 | */ | ||
85 | ActionModelLQRTpl(const MatrixXs& A, const MatrixXs& B, const MatrixXs& Q, | ||
86 | const MatrixXs& R, const MatrixXs& N, const VectorXs& f, | ||
87 | const VectorXs& q, const VectorXs& r); | ||
88 | |||
89 | /** | ||
90 | * @brief Initialize the LQR action model | ||
91 | * | ||
92 | * @param[in] A State matrix | ||
93 | * @param[in] B Input matrix | ||
94 | * @param[in] Q State weight matrix | ||
95 | * @param[in] R Input weight matrix | ||
96 | * @param[in] N State-input weight matrix | ||
97 | * @param[in] G State-input inequality constraint matrix | ||
98 | * @param[in] H State-input equality constraint matrix | ||
99 | * @param[in] f Dynamics drift | ||
100 | * @param[in] q State weight vector | ||
101 | * @param[in] r Input weight vector | ||
102 | * @param[in] g State-input inequality constraint bias | ||
103 | * @param[in] h State-input equality constraint bias | ||
104 | */ | ||
105 | ActionModelLQRTpl(const MatrixXs& A, const MatrixXs& B, const MatrixXs& Q, | ||
106 | const MatrixXs& R, const MatrixXs& N, const MatrixXs& G, | ||
107 | const MatrixXs& H, const VectorXs& f, const VectorXs& q, | ||
108 | const VectorXs& r, const VectorXs& g, const VectorXs& h); | ||
109 | |||
110 | /** | ||
111 | * @brief Initialize the LQR action model | ||
112 | * | ||
113 | * @param[in] nx Dimension of state vector | ||
114 | * @param[in] nu Dimension of control vector | ||
115 | * @param[in] drif_free Enable / disable the bias term of the linear dynamics | ||
116 | * (default true) | ||
117 | */ | ||
118 | ActionModelLQRTpl(const std::size_t nx, const std::size_t nu, | ||
119 | const bool drift_free = true); | ||
120 | |||
121 | /** @brief Copy constructor */ | ||
122 | ActionModelLQRTpl(const ActionModelLQRTpl& copy); | ||
123 | |||
124 | 420 | virtual ~ActionModelLQRTpl() = default; | |
125 | |||
126 | virtual void calc(const std::shared_ptr<ActionDataAbstract>& data, | ||
127 | const Eigen::Ref<const VectorXs>& x, | ||
128 | const Eigen::Ref<const VectorXs>& u) override; | ||
129 | virtual void calc(const std::shared_ptr<ActionDataAbstract>& data, | ||
130 | const Eigen::Ref<const VectorXs>& x) override; | ||
131 | virtual void calcDiff(const std::shared_ptr<ActionDataAbstract>& data, | ||
132 | const Eigen::Ref<const VectorXs>& x, | ||
133 | const Eigen::Ref<const VectorXs>& u) override; | ||
134 | virtual void calcDiff(const std::shared_ptr<ActionDataAbstract>& data, | ||
135 | const Eigen::Ref<const VectorXs>& x) override; | ||
136 | virtual std::shared_ptr<ActionDataAbstract> createData() override; | ||
137 | |||
138 | /** | ||
139 | * @brief Cast the LQR model to a different scalar type. | ||
140 | * | ||
141 | * It is useful for operations requiring different precision or scalar types. | ||
142 | * | ||
143 | * @tparam NewScalar The new scalar type to cast to. | ||
144 | * @return ActionModelLQRTpl<NewScalar> A action model with the | ||
145 | * new scalar type. | ||
146 | */ | ||
147 | template <typename NewScalar> | ||
148 | ActionModelLQRTpl<NewScalar> cast() const; | ||
149 | |||
150 | virtual bool checkData( | ||
151 | const std::shared_ptr<ActionDataAbstract>& data) override; | ||
152 | |||
153 | /** | ||
154 | * @brief Create a random LQR model | ||
155 | * | ||
156 | * @param[in] nx State dimension | ||
157 | * @param[in] nu Control dimension | ||
158 | * @param[in] ng Inequality constraint dimension (default 0) | ||
159 | * @param[in] nh Equality constraint dimension (defaul 0) | ||
160 | */ | ||
161 | static ActionModelLQRTpl Random(const std::size_t nx, const std::size_t nu, | ||
162 | const std::size_t ng = 0, | ||
163 | const std::size_t nh = 0); | ||
164 | |||
165 | /** @brief Return the state matrix */ | ||
166 | const MatrixXs& get_A() const; | ||
167 | |||
168 | /** @brief Return the input matrix */ | ||
169 | const MatrixXs& get_B() const; | ||
170 | |||
171 | /** @brief Return the dynamics drift */ | ||
172 | const VectorXs& get_f() const; | ||
173 | |||
174 | /** @brief Return the state weight matrix */ | ||
175 | const MatrixXs& get_Q() const; | ||
176 | |||
177 | /** @brief Return the input weight matrix */ | ||
178 | const MatrixXs& get_R() const; | ||
179 | |||
180 | /** @brief Return the state-input weight matrix */ | ||
181 | const MatrixXs& get_N() const; | ||
182 | |||
183 | /** @brief Return the state-input inequality constraint matrix */ | ||
184 | const MatrixXs& get_G() const; | ||
185 | |||
186 | /** @brief Return the state-input equality constraint matrix */ | ||
187 | const MatrixXs& get_H() const; | ||
188 | |||
189 | /** @brief Return the state weight vector */ | ||
190 | const VectorXs& get_q() const; | ||
191 | |||
192 | /** @brief Return the input weight vector */ | ||
193 | const VectorXs& get_r() const; | ||
194 | |||
195 | /** @brief Return the state-input inequality constraint bias */ | ||
196 | const VectorXs& get_g() const; | ||
197 | |||
198 | /** @brief Return the state-input equality constraint bias */ | ||
199 | const VectorXs& get_h() const; | ||
200 | |||
201 | /** | ||
202 | * @brief Modify the LQR action model | ||
203 | * | ||
204 | * @param[in] A State matrix | ||
205 | * @param[in] B Input matrix | ||
206 | * @param[in] Q State weight matrix | ||
207 | * @param[in] R Input weight matrix | ||
208 | * @param[in] N State-input weight matrix | ||
209 | * @param[in] G State-input inequality constraint matrix | ||
210 | * @param[in] H State-input equality constraint matrix | ||
211 | * @param[in] f Dynamics drift | ||
212 | * @param[in] q State weight vector | ||
213 | * @param[in] r Input weight vector | ||
214 | * @param[in] g State-input inequality constraint bias | ||
215 | * @param[in] h State-input equality constraint bias | ||
216 | */ | ||
217 | void set_LQR(const MatrixXs& A, const MatrixXs& B, const MatrixXs& Q, | ||
218 | const MatrixXs& R, const MatrixXs& N, const MatrixXs& G, | ||
219 | const MatrixXs& H, const VectorXs& f, const VectorXs& q, | ||
220 | const VectorXs& r, const VectorXs& g, const VectorXs& h); | ||
221 | |||
222 | ✗ | DEPRECATED("Use get_A", const MatrixXs& get_Fx() const { return get_A(); }) | |
223 | ✗ | DEPRECATED("Use get_B", const MatrixXs& get_Fu() const { return get_B(); }) | |
224 | ✗ | DEPRECATED("Use get_f", const VectorXs& get_f0() const { return get_f(); }) | |
225 | ✗ | DEPRECATED("Use get_q", const VectorXs& get_lx() const { return get_q(); }) | |
226 | ✗ | DEPRECATED("Use get_r", const VectorXs& get_lu() const { return get_r(); }) | |
227 | ✗ | DEPRECATED("Use get_Q", const MatrixXs& get_Lxx() const { return get_Q(); }) | |
228 | ✗ | DEPRECATED("Use get_R", const MatrixXs& get_Lxu() const { return get_R(); }) | |
229 | ✗ | DEPRECATED("Use get_N", const MatrixXs& get_Luu() const { return get_N(); }) | |
230 | ✗ | DEPRECATED( | |
231 | "Use set_LQR", void set_Fx(const MatrixXs& A) { | ||
232 | set_LQR(A, B_, Q_, R_, N_, G_, H_, f_, q_, r_, g_, h_); | ||
233 | }) | ||
234 | ✗ | DEPRECATED( | |
235 | "Use set_LQR", void set_Fu(const MatrixXs& B) { | ||
236 | set_LQR(A_, B, Q_, R_, N_, G_, H_, f_, q_, r_, g_, h_); | ||
237 | }) | ||
238 | ✗ | DEPRECATED( | |
239 | "Use set_LQR", void set_f0(const VectorXs& f) { | ||
240 | set_LQR(A_, B_, Q_, R_, N_, G_, H_, f, q_, r_, g_, h_); | ||
241 | }) | ||
242 | ✗ | DEPRECATED( | |
243 | "Use set_LQR", void set_lx(const VectorXs& q) { | ||
244 | set_LQR(A_, B_, Q_, R_, N_, G_, H_, f_, q, r_, g_, h_); | ||
245 | }) | ||
246 | ✗ | DEPRECATED( | |
247 | "Use set_LQR", void set_lu(const VectorXs& r) { | ||
248 | set_LQR(A_, B_, Q_, R_, N_, G_, H_, f_, q_, r, g_, h_); | ||
249 | }) | ||
250 | ✗ | DEPRECATED( | |
251 | "Use set_LQR", void set_Lxx(const MatrixXs& Q) { | ||
252 | set_LQR(A_, B_, Q, R_, N_, G_, H_, f_, q_, r_, g_, h_); | ||
253 | }) | ||
254 | ✗ | DEPRECATED( | |
255 | "Use set_LQR", void set_Luu(const MatrixXs& R) { | ||
256 | set_LQR(A_, B_, Q_, R, N_, G_, H_, f_, q_, r_, g_, h_); | ||
257 | }) | ||
258 | ✗ | DEPRECATED( | |
259 | "Use set_LQR", void set_Lxu(const MatrixXs& N) { | ||
260 | set_LQR(A_, B_, Q_, R_, N, G_, H_, f_, q_, r_, g_, h_); | ||
261 | }) | ||
262 | |||
263 | /** | ||
264 | * @brief Print relevant information of the LQR model | ||
265 | * | ||
266 | * @param[out] os Output stream object | ||
267 | */ | ||
268 | virtual void print(std::ostream& os) const override; | ||
269 | |||
270 | protected: | ||
271 | using Base::ng_; //!< Equality constraint dimension | ||
272 | using Base::nh_; //!< Inequality constraint dimension | ||
273 | using Base::nu_; //!< Control dimension | ||
274 | using Base::state_; //!< Model of the state | ||
275 | |||
276 | private: | ||
277 | MatrixXs A_; | ||
278 | MatrixXs B_; | ||
279 | MatrixXs Q_; | ||
280 | MatrixXs R_; | ||
281 | MatrixXs N_; | ||
282 | MatrixXs G_; | ||
283 | MatrixXs H_; | ||
284 | VectorXs f_; | ||
285 | VectorXs q_; | ||
286 | VectorXs r_; | ||
287 | VectorXs g_; | ||
288 | VectorXs h_; | ||
289 | MatrixXs L_; | ||
290 | bool drift_free_; | ||
291 | bool updated_lqr_; | ||
292 | }; | ||
293 | |||
294 | template <typename _Scalar> | ||
295 | struct ActionDataLQRTpl : public ActionDataAbstractTpl<_Scalar> { | ||
296 | typedef _Scalar Scalar; | ||
297 | typedef MathBaseTpl<Scalar> MathBase; | ||
298 | typedef ActionDataAbstractTpl<Scalar> Base; | ||
299 | typedef typename MathBase::VectorXs VectorXs; | ||
300 | |||
301 | template <template <typename Scalar> class Model> | ||
302 | 1648 | explicit ActionDataLQRTpl(Model<Scalar>* const model) | |
303 | : Base(model), | ||
304 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | R_u_tmp(VectorXs::Zero(static_cast<Eigen::Index>(model->get_nu()))), |
305 |
1/2✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
|
1648 | Q_x_tmp(VectorXs::Zero( |
306 |
4/8✓ Branch 2 taken 1617 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1617 times.
✗ Branch 6 not taken.
✓ Branch 9 taken 1617 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 1617 times.
✗ Branch 13 not taken.
|
3296 | static_cast<Eigen::Index>(model->get_state()->get_ndx()))) { |
307 | // Setting the linear model and quadratic cost as they are constant | ||
308 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 5 taken 1617 times.
✗ Branch 6 not taken.
|
1648 | const std::size_t nq = model->get_state()->get_nq(); |
309 |
1/2✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
|
1648 | const std::size_t nu = model->get_nu(); |
310 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | Fx = model->get_A(); |
311 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | Fu = model->get_B(); |
312 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | Lxx = model->get_Q(); |
313 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | Luu = model->get_R(); |
314 |
2/4✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
|
1648 | Lxu = model->get_N(); |
315 |
3/6✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1617 times.
✗ Branch 8 not taken.
|
1648 | Gx = model->get_G().leftCols(2 * nq); |
316 |
3/6✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1617 times.
✗ Branch 8 not taken.
|
1648 | Gu = model->get_G().rightCols(nu); |
317 |
3/6✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1617 times.
✗ Branch 8 not taken.
|
1648 | Hx = model->get_H().leftCols(2 * nq); |
318 |
3/6✓ Branch 1 taken 1617 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1617 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1617 times.
✗ Branch 8 not taken.
|
1648 | Hu = model->get_H().rightCols(nu); |
319 | 1648 | } | |
320 | 3205 | virtual ~ActionDataLQRTpl() = default; | |
321 | |||
322 | using Base::cost; | ||
323 | using Base::Fu; | ||
324 | using Base::Fx; | ||
325 | using Base::Gu; | ||
326 | using Base::Gx; | ||
327 | using Base::Hu; | ||
328 | using Base::Hx; | ||
329 | using Base::Lu; | ||
330 | using Base::Luu; | ||
331 | using Base::Lx; | ||
332 | using Base::Lxu; | ||
333 | using Base::Lxx; | ||
334 | using Base::r; | ||
335 | using Base::xnext; | ||
336 | |||
337 | VectorXs R_u_tmp; // Temporary variable for storing Hessian-vector product | ||
338 | // (size: nu) | ||
339 | VectorXs Q_x_tmp; // Temporary variable for storing Hessian-vector product | ||
340 | // (size: nx) | ||
341 | }; | ||
342 | |||
343 | } // namespace crocoddyl | ||
344 | |||
345 | /* --- Details -------------------------------------------------------------- */ | ||
346 | /* --- Details -------------------------------------------------------------- */ | ||
347 | /* --- Details -------------------------------------------------------------- */ | ||
348 | #include "crocoddyl/core/actions/lqr.hxx" | ||
349 | |||
350 | CROCODDYL_DECLARE_EXTERN_TEMPLATE_CLASS(crocoddyl::ActionModelLQRTpl) | ||
351 | CROCODDYL_DECLARE_EXTERN_TEMPLATE_STRUCT(crocoddyl::ActionDataLQRTpl) | ||
352 | |||
353 | #endif // CROCODDYL_CORE_ACTIONS_LQR_HPP_ | ||
354 |