tests/eiquadprog-fast.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2019 CNRS
3 //
4 // This file is part of eiquadprog.
5 //
6 // eiquadprog is free software: you can redistribute it and/or modify
7 // it under the terms of the GNU Lesser General Public License as published by
8 // the Free Software Foundation, either version 3 of the License, or
9 //(at your option) any later version.
10 
11 // eiquadprog is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 // GNU Lesser General Public License for more details.
15 
16 // You should have received a copy of the GNU Lesser General Public License
17 // along with eiquadprog. If not, see <https://www.gnu.org/licenses/>.
18 
20 
21 #include <Eigen/Core>
22 #include <boost/test/unit_test.hpp>
23 #include <iostream>
24 
25 using namespace eiquadprog::solvers;
26 
34 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
35 
36 // min ||x||^2
37 
38 BOOST_AUTO_TEST_CASE(test_unbiased) {
39  EiquadprogFast qp;
40  qp.reset(2, 0, 0);
41 
42  Eigen::MatrixXd Q(2, 2);
43  Q.setZero();
44  Q(0, 0) = 1.0;
45  Q(1, 1) = 1.0;
46 
47  Eigen::VectorXd C(2);
48  C.setZero();
49 
50  Eigen::MatrixXd Aeq(0, 2);
51 
52  Eigen::VectorXd Beq(0);
53 
54  Eigen::MatrixXd Aineq(0, 2);
55 
56  Eigen::VectorXd Bineq(0);
57 
58  Eigen::VectorXd x(2);
59 
60  Eigen::VectorXd solution(2);
61  solution.setZero();
62 
63  double val = 0.0;
64 
66 
67  EiquadprogFast_status status =
68  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
69 
70  BOOST_CHECK_EQUAL(status, expected);
71 
72  BOOST_CHECK_CLOSE(qp.getObjValue(), val, 1e-6);
73 
74  BOOST_CHECK(x.isApprox(solution));
75 }
76 
77 // min ||x-x_0||^2, x_0 = (1 1)^T
78 
79 BOOST_AUTO_TEST_CASE(test_biased) {
80  EiquadprogFast qp;
81  qp.reset(2, 0, 0);
82 
83  Eigen::MatrixXd Q(2, 2);
84  Q.setZero();
85  Q(0, 0) = 1.0;
86  Q(1, 1) = 1.0;
87 
88  Eigen::VectorXd C(2);
89  C(0) = -1.;
90  C(1) = -1.;
91 
92  Eigen::MatrixXd Aeq(0, 2);
93 
94  Eigen::VectorXd Beq(0);
95 
96  Eigen::MatrixXd Aineq(0, 2);
97 
98  Eigen::VectorXd Bineq(0);
99 
100  Eigen::VectorXd x(2);
101 
102  Eigen::VectorXd solution(2);
103  solution(0) = 1.;
104  solution(1) = 1.;
105 
106  double val = -1.;
107 
109 
110  EiquadprogFast_status status =
111  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
112 
113  BOOST_CHECK_EQUAL(status, expected);
114 
115  BOOST_CHECK_CLOSE(qp.getObjValue(), val, 1e-6);
116 
117  BOOST_CHECK(x.isApprox(solution));
118 }
119 
120 // min ||x||^2
121 // s.t.
122 // x[1] = 1 - x[0]
123 
124 BOOST_AUTO_TEST_CASE(test_equality_constraints) {
125  EiquadprogFast qp;
126  qp.reset(2, 1, 0);
127 
128  Eigen::MatrixXd Q(2, 2);
129  Q.setZero();
130  Q(0, 0) = 1.0;
131  Q(1, 1) = 1.0;
132 
133  Eigen::VectorXd C(2);
134  C.setZero();
135 
136  Eigen::MatrixXd Aeq(1, 2);
137  Aeq(0, 0) = 1.;
138  Aeq(0, 1) = 1.;
139 
140  Eigen::VectorXd Beq(1);
141  Beq(0) = -1.;
142 
143  Eigen::MatrixXd Aineq(0, 2);
144 
145  Eigen::VectorXd Bineq(0);
146 
147  Eigen::VectorXd x(2);
148 
149  Eigen::VectorXd solution(2);
150  solution(0) = 0.5;
151  solution(1) = 0.5;
152 
153  double val = 0.25;
154 
156 
157  EiquadprogFast_status status =
158  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
159 
160  BOOST_CHECK_EQUAL(status, expected);
161 
162  BOOST_CHECK_CLOSE(qp.getObjValue(), val, 1e-6);
163 
164  BOOST_CHECK(x.isApprox(solution));
165 }
166 
167 // min ||x||^2
168 // s.t.
169 // x[i] >= 1
170 
171 BOOST_AUTO_TEST_CASE(test_inequality_constraints) {
172  EiquadprogFast qp;
173  qp.reset(2, 0, 2);
174 
175  Eigen::MatrixXd Q(2, 2);
176  Q.setZero();
177  Q(0, 0) = 1.0;
178  Q(1, 1) = 1.0;
179 
180  Eigen::VectorXd C(2);
181  C.setZero();
182 
183  Eigen::MatrixXd Aeq(0, 2);
184 
185  Eigen::VectorXd Beq(0);
186 
187  Eigen::MatrixXd Aineq(2, 2);
188  Aineq.setZero();
189  Aineq(0, 0) = 1.;
190  Aineq(1, 1) = 1.;
191 
192  Eigen::VectorXd Bineq(2);
193  Bineq(0) = -1.;
194  Bineq(1) = -1.;
195 
196  Eigen::VectorXd x(2);
197 
198  Eigen::VectorXd solution(2);
199  solution(0) = 1.;
200  solution(1) = 1.;
201 
202  double val = 1.;
203 
205 
206  EiquadprogFast_status status =
207  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
208 
209  BOOST_CHECK_EQUAL(status, expected);
210 
211  BOOST_CHECK_CLOSE(qp.getObjValue(), val, 1e-6);
212 
213  BOOST_CHECK(x.isApprox(solution));
214 }
215 
216 // min ||x-x_0||^2, x_0 = (1 1)^T
217 // s.t.
218 // x[1] = 5 - x[0]
219 // x[1] >= 3
220 
222  EiquadprogFast qp;
223  qp.reset(2, 1, 1);
224 
225  Eigen::MatrixXd Q(2, 2);
226  Q.setZero();
227  Q(0, 0) = 1.0;
228  Q(1, 1) = 1.0;
229 
230  Eigen::VectorXd C(2);
231  C(0) = -1.;
232  C(1) = -1.;
233 
234  Eigen::MatrixXd Aeq(1, 2);
235  Aeq(0, 0) = 1.;
236  Aeq(0, 1) = 1.;
237 
238  Eigen::VectorXd Beq(1);
239  Beq(0) = -5.;
240 
241  Eigen::MatrixXd Aineq(1, 2);
242  Aineq.setZero();
243  Aineq(0, 1) = 1.;
244 
245  Eigen::VectorXd Bineq(1);
246  Bineq(0) = -3.;
247 
248  Eigen::VectorXd x(2);
249 
250  Eigen::VectorXd solution(2);
251  solution(0) = 2.;
252  solution(1) = 3.;
253 
254  double val = 1.5;
255 
257 
258  EiquadprogFast_status status =
259  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
260 
261  BOOST_CHECK_EQUAL(status, expected);
262 
263  BOOST_CHECK_CLOSE(qp.getObjValue(), val, 1e-6);
264 
265  BOOST_CHECK(x.isApprox(solution));
266 }
267 
268 // min ||x||^2
269 // s.t.
270 // x[0] = 1
271 // x[0] = -1
272 
273 BOOST_AUTO_TEST_CASE(test_unfeasible_equalities) {
274  EiquadprogFast qp;
275  qp.reset(2, 2, 0);
276 
277  Eigen::MatrixXd Q(2, 2);
278  Q.setZero();
279  Q(0, 0) = 1.0;
280  Q(1, 1) = 1.0;
281 
282  Eigen::VectorXd C(2);
283  C.setZero();
284 
285  Eigen::MatrixXd Aeq(2, 2);
286  Aeq.setZero();
287  Aeq(0, 0) = 1.;
288  Aeq(1, 0) = 1.;
289 
290  Eigen::VectorXd Beq(2);
291  Beq(0) = -1.;
292  Beq(1) = 1.;
293 
294  Eigen::MatrixXd Aineq(0, 2);
295 
296  Eigen::VectorXd Bineq(0);
297 
298  Eigen::VectorXd x(2);
299 
301 
302  EiquadprogFast_status status =
303  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
304 
305  BOOST_CHECK_EQUAL(status, expected);
306 }
307 
308 // min ||x||^2
309 // s.t.
310 // x[0] >= 1
311 // x[0] <= -1
312 //
313 // correctly fails, but returns wrong error code
314 
315 BOOST_AUTO_TEST_CASE(test_unfeasible_inequalities) {
316  EiquadprogFast qp;
317  qp.reset(2, 0, 2);
318 
319  Eigen::MatrixXd Q(2, 2);
320  Q.setZero();
321  Q(0, 0) = 1.0;
322  Q(1, 1) = 1.0;
323 
324  Eigen::VectorXd C(2);
325  C.setZero();
326 
327  Eigen::MatrixXd Aeq(0, 2);
328 
329  Eigen::VectorXd Beq(0);
330 
331  Eigen::MatrixXd Aineq(2, 2);
332  Aineq.setZero();
333  Aineq(0, 0) = 1.;
334  Aineq(1, 0) = -1.;
335 
336  Eigen::VectorXd Bineq(2);
337  Bineq(0) = -1;
338  Bineq(1) = -1;
339 
340  Eigen::VectorXd x(2);
341 
343 
344  EiquadprogFast_status status =
345  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
346 
347  BOOST_WARN_EQUAL(status, expected);
348  BOOST_CHECK(status != EIQUADPROG_FAST_OPTIMAL);
349 }
350 
351 // min ||x-x_0||^2, x_0 = (1 1)^T
352 // s.t.
353 // x[1] = 1 - x[0]
354 // x[0] <= 0
355 // x[1] <= 0
356 //
357 // correctly fails, but returns wrong error code
358 
359 BOOST_AUTO_TEST_CASE(test_unfeasible_constraints) {
360  EiquadprogFast qp;
361  qp.reset(2, 1, 2);
362 
363  Eigen::MatrixXd Q(2, 2);
364  Q.setZero();
365  Q(0, 0) = 1.0;
366  Q(1, 1) = 1.0;
367 
368  Eigen::VectorXd C(2);
369  C(0) = -1.;
370  C(1) = -1.;
371 
372  Eigen::MatrixXd Aeq(1, 2);
373  Aeq(0, 0) = 1.;
374  Aeq(0, 1) = 1.;
375 
376  Eigen::VectorXd Beq(1);
377  Beq(0) = -1.;
378 
379  Eigen::MatrixXd Aineq(2, 2);
380  Aineq.setZero();
381  Aineq(0, 0) = -1.;
382  Aineq(1, 1) = -1.;
383 
384  Eigen::VectorXd Bineq(2);
385  Bineq.setZero();
386 
387  Eigen::VectorXd x(2);
388 
390 
391  EiquadprogFast_status status =
392  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
393 
394  BOOST_WARN_EQUAL(status, expected);
395  BOOST_CHECK(status != EIQUADPROG_FAST_OPTIMAL);
396 }
397 
398 // min -||x||^2
399 // DOES NOT WORK!
400 
401 BOOST_AUTO_TEST_CASE(test_unbounded) {
402  EiquadprogFast qp;
403  qp.reset(2, 0, 0);
404 
405  Eigen::MatrixXd Q(2, 2);
406  Q.setZero();
407  Q(0, 0) = -1.0;
408  Q(1, 1) = -1.0;
409 
410  Eigen::VectorXd C(2);
411  C.setZero();
412 
413  Eigen::MatrixXd Aeq(0, 2);
414 
415  Eigen::VectorXd Beq(0);
416 
417  Eigen::MatrixXd Aineq(0, 2);
418 
419  Eigen::VectorXd Bineq(0);
420 
421  Eigen::VectorXd x(2);
422 
424 
425  EiquadprogFast_status status =
426  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
427 
428  BOOST_WARN_EQUAL(status, expected);
429  BOOST_WARN(status != EIQUADPROG_FAST_OPTIMAL); // SHOULD pass!
430 }
431 
432 // min -||x||^2
433 // s.t.
434 // 0<= x[0] <= 1
435 // 0<= x[1] <= 1
436 // DOES NOT WORK!
437 
438 BOOST_AUTO_TEST_CASE(test_nonconvex) {
439  EiquadprogFast qp;
440  qp.reset(2, 0, 4);
441 
442  Eigen::MatrixXd Q(2, 2);
443  Q.setZero();
444  Q(0, 0) = -1.0;
445  Q(1, 1) = -1.0;
446 
447  Eigen::VectorXd C(2);
448  C.setZero();
449 
450  Eigen::MatrixXd Aeq(0, 2);
451 
452  Eigen::VectorXd Beq(0);
453 
454  Eigen::MatrixXd Aineq(4, 2);
455  Aineq.setZero();
456  Aineq(0, 0) = 1.;
457  Aineq(1, 0) = -1.;
458  Aineq(2, 1) = 1.;
459  Aineq(3, 1) = -1.;
460 
461  Eigen::VectorXd Bineq(4);
462  Bineq(0) = 0.;
463  Bineq(1) = 1.;
464  Bineq(2) = 0.;
465  Bineq(3) = 1.;
466 
467  Eigen::VectorXd x(2);
468 
469  Eigen::VectorXd solution(2);
470  solution(0) = 1.;
471  solution(1) = 1.;
472 
473  double val = -1.;
474 
476 
477  EiquadprogFast_status status =
478  qp.solve_quadprog(Q, C, Aeq, Beq, Aineq, Bineq, x);
479 
480  BOOST_CHECK_EQUAL(status, expected);
481 
482  BOOST_WARN_CLOSE(qp.getObjValue(), val, 1e-6);
483 
484  BOOST_WARN(x.isApprox(solution));
485 }
486 
487 BOOST_AUTO_TEST_SUITE_END()
eiquadprog::solvers::EiquadprogFast::solve_quadprog
EiquadprogFast_status solve_quadprog(const MatrixXd &Hess, const VectorXd &g0, const MatrixXd &CE, const VectorXd &ce0, const MatrixXd &CI, const VectorXd &ci0, VectorXd &x)
Definition: src/eiquadprog-fast.cpp:189
eiquadprog::solvers::EIQUADPROG_FAST_REDUNDANT_EQUALITIES
@ EIQUADPROG_FAST_REDUNDANT_EQUALITIES
Definition: eiquadprog-fast.hpp:75
eiquadprog::solvers::EIQUADPROG_FAST_UNBOUNDED
@ EIQUADPROG_FAST_UNBOUNDED
Definition: eiquadprog-fast.hpp:73
BOOST_AUTO_TEST_CASE
BOOST_AUTO_TEST_CASE(test_unbiased)
Definition: tests/eiquadprog-fast.cpp:38
eiquadprog::solvers::EiquadprogFast
Definition: eiquadprog-fast.hpp:78
eiquadprog::solvers::EIQUADPROG_FAST_INFEASIBLE
@ EIQUADPROG_FAST_INFEASIBLE
Definition: eiquadprog-fast.hpp:72
eiquadprog::solvers::EiquadprogFast::getObjValue
double getObjValue() const
Definition: eiquadprog-fast.hpp:113
eiquadprog-fast.hpp
eiquadprog::solvers::EiquadprogFast_status
EiquadprogFast_status
Definition: eiquadprog-fast.hpp:70
eiquadprog::solvers
Definition: eiquadprog-fast.hpp:65
eiquadprog::solvers::EiquadprogFast::reset
void reset(size_t dim_qp, size_t num_eq, size_t num_ineq)
Definition: src/eiquadprog-fast.cpp:20
eiquadprog::solvers::EIQUADPROG_FAST_OPTIMAL
@ EIQUADPROG_FAST_OPTIMAL
Definition: eiquadprog-fast.hpp:71


eiquadprog
Author(s): Gabriele Buondonno, Andrea Del Prete, Luca Di Gaspero, Angelo Furfaro, Benjamin Stephens, Gael Guennebaud
autogenerated on Wed May 28 2025 02:55:57