- 3.0.2 optimal control module.
LqrTest.h
Go to the documentation of this file.
1 /**********************************************************************************************************************
2 This file is part of the Control Toolbox (https://github.com/ethz-adrl/control-toolbox), copyright by ETH Zurich.
3 Licensed under the BSD-2 license (see LICENSE file in main directory)
4 **********************************************************************************************************************/
5 
12 #pragma once
13 
14 #include <chrono>
15 
16 #ifdef MATLAB
17 #include <matlabCppInterface/Engine.hpp>
18 #endif
19 
20 // Bring in gtest
21 #include <gtest/gtest.h>
22 
23 
24 namespace ct {
25 namespace optcon {
26 namespace example {
27 
28 TEST(LQRTest, DARETest)
29 {
30  const size_t stateDim = 2;
31  const size_t controlDim = 1;
32 
33  Eigen::Matrix<double, stateDim, stateDim> A;
34  Eigen::Matrix<double, stateDim, controlDim> B;
35  Eigen::Matrix<double, stateDim, stateDim> Q;
36  Eigen::Matrix<double, controlDim, controlDim> R;
37  Eigen::Matrix<double, controlDim, stateDim> K;
38 
39  A << 1, 1, 1, 0;
40  B << 0, 1;
41  Q << 1, 0, 0, 1;
42  R << 1;
43 
45  Eigen::Matrix<double, stateDim, stateDim> P = dare.computeSteadyStateRiccatiMatrix(Q, R, A, B, K, true);
46  Eigen::Matrix<double, stateDim, stateDim> P_test;
47  P_test << 6.932484752255643, 4.332273119899151, 4.332273119899151, 4.55195134961773;
48  ASSERT_LT((P - P_test).array().abs().maxCoeff(), 1e-12);
49 }
50 
51 
52 TEST(LQRTest, quadTest)
53 {
54  // std::cout << "QUADROTOR TEST"<<std::endl;
55  // std::cout << "==================================="<<std::endl;
56  // std::cout << "==================================="<<std::endl << std::endl << std::endl;
57 
58  const size_t stateDim = 12;
59  const size_t controlDim = 4;
60 
61  Eigen::Matrix<double, stateDim, stateDim> A;
62  Eigen::Matrix<double, stateDim, controlDim> B;
63  Eigen::Matrix<double, stateDim, stateDim> Q;
64  Eigen::Matrix<double, controlDim, controlDim> R;
65  Eigen::Matrix<double, controlDim, stateDim> K;
66  Eigen::Matrix<double, controlDim, stateDim> Kiterative;
67 
69 
70  Q << 10, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 10, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2000, 0, 0, 0, 0, 0, 0, 0, 0,
71  0, 0, 0, 0, 0.1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.5, 0, 0, 0, 0,
72  0, 0, 0, 0, 0, 0, 0, 0, 0.1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0,
73  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
74  0, 0.02;
75 
76  R << 100, 0, 0, 0, 0, 1000, 0, 0, 0, 0, 1000, 0, 0, 0, 0, 100;
77 
78 
79  A << 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0,
80  0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0,
81  0, 0, 9.81, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -9.81, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -0, -0, 0, 0, 0, 0, 0, 0, 0,
82  0, 0, 0, 0, 0, -0, 0, 0, 0, 0, -0, -0, 0, 0, 0, 0, -0, -0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -0, 0,
83  0;
84 
85 
86  B << 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -0, 0, 0, 0, 1.39665, 0, 0,
87  0, 0, 142.857, -0, 0, 0, 0, 142.857, 0, 0, -0, 0, 83.3333;
88 
89  Kiterative.setZero();
90  K.setZero();
91 
92 
93  bool foundSolutionIterative = lqr.compute(Q, R, A, B, Kiterative, false, true);
94  ASSERT_TRUE(foundSolutionIterative);
95 
96 #ifdef CT_USE_LAPACK
97  bool foundSolutionDirect = lqr.compute(Q, R, A, B, K, false);
98  ASSERT_TRUE(foundSolutionDirect);
99  ASSERT_LT((K - Kiterative).array().abs().maxCoeff(), 1e-4);
100 #endif
101 
102  int nTests = 100;
103 #ifdef CT_USE_LAPACK
104  auto start1 = std::chrono::system_clock::now();
105  for (int i = 0; i < nTests; i++)
106  {
107  lqr.compute(Q, R, A, B, K, false);
108  }
109  auto end1 = std::chrono::system_clock::now();
110  auto elapsed1 = std::chrono::duration_cast<std::chrono::milliseconds>(end1 - start1);
111  std::cout << "solved " << nTests << " lqr problems with state dimension " << stateDim << " in " << elapsed1.count()
112  << " ms (average: " << elapsed1.count() / static_cast<double>(nTests) << " ms / lqr)" << std::endl;
113 #endif
114 
115  auto start2 = std::chrono::system_clock::now();
116  for (int i = 0; i < nTests; i++)
117  {
118  lqr.compute(Q, R, A, B, Kiterative, false, true);
119  }
120  auto end2 = std::chrono::system_clock::now();
121  auto elapsed2 = std::chrono::duration_cast<std::chrono::milliseconds>(end2 - start2);
122  std::cout << "solved " << nTests << " lqr problems iteratively with state dimension " << stateDim << " in "
123  << elapsed2.count() << " ms (average: " << elapsed2.count() / static_cast<double>(nTests) << " ms / lqr)"
124  << std::endl;
125 }
126 
127 #ifdef MATLAB
128 TEST(LQRTest, matlabTest)
129 {
130  matlab::Engine engine(true);
131  ASSERT_TRUE(engine.good());
132 
133  const size_t stateDim = 5;
134  std::string stateDimString = std::to_string(stateDim);
135 
136  Eigen::MatrixXd Ad;
137  Eigen::MatrixXd Bd;
138  Eigen::MatrixXd Qd;
139  Eigen::MatrixXd Rd;
140  Eigen::MatrixXd K_Matlab;
141 
142  Eigen::Matrix<double, stateDim, stateDim> A;
143  Eigen::Matrix<double, stateDim, stateDim> B;
144  Eigen::Matrix<double, stateDim, stateDim> Q;
145  Eigen::Matrix<double, stateDim, stateDim> R;
146  Eigen::Matrix<double, stateDim, stateDim> K_Cpp;
147  Eigen::Matrix<double, stateDim, stateDim> K_Cpp_iteratively;
148 
150 
151  // std::cout << "ARTIFICIAL TEST "<<std::endl;
152  // std::cout << "==================================="<<std::endl;
153  // std::cout << "==================================="<<std::endl << std::endl << std::endl;
154 
155  for (int i = 0; i < 10; i++)
156  {
157  std::cout << "Test " << std::to_string(i) << std::endl;
158  std::cout << "===================================" << std::endl;
159 
160  std::cout << "1. Generating problem in Matlab" << std::endl;
161  engine.executeCommand("A = magic(" + stateDimString + ");");
162  engine.executeCommand("B = magic(" + stateDimString + ");");
163  engine.executeCommand("Q = diag(100*rand(" + stateDimString + ",1));");
164  engine.executeCommand("R = diag(100*rand(" + stateDimString + ",1));");
165  engine.executeCommand("N = zeros(" + stateDimString + ");");
166 
167  std::cout << "2. Computing LQR in Matlab" << std::endl;
168  std::cout << engine.executeCommand("[K,S,E] = lqr(A,B,Q,R,N);");
169 
170  std::cout << "3. Obtaining problem from Matlab" << std::endl;
171  engine.get("A", Ad);
172  A = Ad;
173  engine.get("B", Bd);
174  B = Bd;
175  engine.get("Q", Qd);
176  Q = Qd;
177  engine.get("R", Rd);
178  R = Rd;
179 
180  std::cout << "4. Obtaining LQR solution from Matlab" << std::endl;
181  engine.get("K", K_Matlab);
182 
183  std::cout << "5. Computing LQR solution in C++" << std::endl;
184 
185  bool foundSolutionDirect = lqr.compute(Q, R, A, B, K_Cpp, false);
186  ASSERT_TRUE(foundSolutionDirect);
187 
188  bool foundSolutionIterative = lqr.compute(Q, R, A, B, K_Cpp_iteratively, false, true);
189  ASSERT_TRUE(foundSolutionIterative);
190 
191 
192  std::cout << "7. Comparing both solutions" << std::endl;
193  ASSERT_LT((K_Matlab - K_Cpp).array().abs().maxCoeff(), 1e-4);
194  ASSERT_LT((K_Matlab - K_Cpp_iteratively).array().abs().maxCoeff(), 1e-4);
195 
196  std::cout << std::endl << std::endl << std::endl;
197  }
198 }
199 #endif //MATLAB
200 
201 } // namespace example
202 } // namespace optcon
203 } // namespace ct
continuous-time infinite-horizon LQR
Definition: LQR.hpp:32
TEST(ConstraintComparison, comparisonAnalyticAD)
Definition: ConstraintComparison.h:158
for i
Definition: mpc_unittest_plotting.m:14
Discrete-Time Algebraic Riccati Equation.
Definition: DARE.hpp:24
bool compute(const state_matrix_t &Q, const control_matrix_t &R, const state_matrix_t &A, const control_gain_matrix_t &B, control_feedback_t &K, bool RisDiagonal=false, bool solveRiccatiIteratively=false)
design the infinite-horizon LQR controller.
Definition: LQR-impl.hpp:16
state_matrix_t computeSteadyStateRiccatiMatrix(const state_matrix_t &Q, const control_matrix_t &R, const state_matrix_t &A, const control_gain_matrix_t &B, control_feedback_t &K, bool verbose=false, const SCALAR eps=1e-6, size_t maxIter=1000)
Definition: DARE-impl.hpp:18