Skip to content

Commit 525867b

Browse files
authored
Merge pull request #125 from jcarpent/devel
Fix dimension check issues
2 parents c35a6d3 + 7f59eee commit 525867b

File tree

5 files changed

+91
-22
lines changed

5 files changed

+91
-22
lines changed

include/proxsuite/proxqp/dense/helpers.hpp

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -183,25 +183,25 @@ update(optional<MatRef<T>> H,
183183
{
184184
// check the model is valid
185185
if (g != nullopt) {
186-
PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().rows(),
186+
PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().size(),
187187
model.dim,
188188
"the dimension wrt the primal variable x "
189189
"variable for updating g is not valid.");
190190
}
191191
if (b != nullopt) {
192-
PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().rows(),
192+
PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().size(),
193193
model.n_eq,
194194
"the dimension wrt equality constrained "
195195
"variables for updating b is not valid.");
196196
}
197197
if (u != nullopt) {
198-
PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().rows(),
198+
PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().size(),
199199
model.n_in,
200200
"the dimension wrt inequality constrained "
201201
"variables for updating u is not valid.");
202202
}
203203
if (l != nullopt) {
204-
PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().rows(),
204+
PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().size(),
205205
model.n_in,
206206
"the dimension wrt inequality constrained "
207207
"variables for updating l is not valid.");

include/proxsuite/proxqp/dense/preconditioner/ruiz.hpp

Lines changed: 6 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -89,8 +89,8 @@ ruiz_scale_qp_in_place( //
8989
T aux = sqrt(std::max({
9090
infty_norm(H.col(k).head(k)),
9191
infty_norm(H.row(k).tail(n - k)),
92-
infty_norm(A.col(k)),
93-
infty_norm(C.col(k)),
92+
n_eq > 0 ? infty_norm(A.col(k)) : T(0),
93+
n_in > 0 ? infty_norm(C.col(k)) : T(0),
9494
}));
9595
if (aux == T(0)) {
9696
delta(k) = T(1);
@@ -104,8 +104,8 @@ ruiz_scale_qp_in_place( //
104104
T aux = sqrt(std::max({
105105
infty_norm(H.col(k).head(k)),
106106
infty_norm(H.col(k).tail(n - k)),
107-
infty_norm(A.col(k)),
108-
infty_norm(C.col(k)),
107+
n_eq > 0 ? infty_norm(A.col(k)) : T(0),
108+
n_in > 0 ? infty_norm(C.col(k)) : T(0),
109109
}));
110110
if (aux == T(0)) {
111111
delta(k) = T(1);
@@ -118,8 +118,8 @@ ruiz_scale_qp_in_place( //
118118

119119
T aux = sqrt(std::max({
120120
infty_norm(H.col(k)),
121-
infty_norm(A.col(k)),
122-
infty_norm(C.col(k)),
121+
n_eq > 0 ? infty_norm(A.col(k)) : T(0),
122+
n_in > 0 ? infty_norm(C.col(k)) : T(0),
123123
}));
124124
if (aux == T(0)) {
125125
delta(k) = T(1);

include/proxsuite/proxqp/dense/wrapper.hpp

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -136,7 +136,7 @@ struct QP
136136
// check the model is valid
137137
if (g != nullopt && g.value().size() != 0) {
138138
PROXSUITE_CHECK_ARGUMENT_SIZE(
139-
g.value().rows(),
139+
g.value().size(),
140140
model.dim,
141141
"the dimension wrt the primal variable x variable for initializing g "
142142
"is not valid.");
@@ -145,7 +145,7 @@ struct QP
145145
}
146146
if (b != nullopt && b.value().size() != 0) {
147147
PROXSUITE_CHECK_ARGUMENT_SIZE(
148-
b.value().rows(),
148+
b.value().size(),
149149
model.n_eq,
150150
"the dimension wrt equality constrained variables for initializing b "
151151
"is not valid.");
@@ -154,7 +154,7 @@ struct QP
154154
}
155155
if (u != nullopt && u.value().size() != 0) {
156156
PROXSUITE_CHECK_ARGUMENT_SIZE(
157-
u.value().rows(),
157+
u.value().size(),
158158
model.n_in,
159159
"the dimension wrt inequality constrained variables for initializing u "
160160
"is not valid.");
@@ -163,7 +163,7 @@ struct QP
163163
}
164164
if (l != nullopt && l.value().size() != 0) {
165165
PROXSUITE_CHECK_ARGUMENT_SIZE(
166-
l.value().rows(),
166+
l.value().size(),
167167
model.n_in,
168168
"the dimension wrt inequality constrained variables for initializing l "
169169
"is not valid.");

include/proxsuite/proxqp/sparse/wrapper.hpp

Lines changed: 8 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -179,7 +179,7 @@ struct QP
179179
}
180180
if (g != nullopt && g.value().size() != 0) {
181181
PROXSUITE_CHECK_ARGUMENT_SIZE(
182-
g.value().rows(),
182+
g.value().size(),
183183
model.dim,
184184
"the dimension wrt the primal variable x variable for initializing g "
185185
"is not valid.");
@@ -188,7 +188,7 @@ struct QP
188188
}
189189
if (b != nullopt && b.value().size() != 0) {
190190
PROXSUITE_CHECK_ARGUMENT_SIZE(
191-
b.value().rows(),
191+
b.value().size(),
192192
model.n_eq,
193193
"the dimension wrt equality constrained variables for initializing b "
194194
"is not valid.");
@@ -197,7 +197,7 @@ struct QP
197197
}
198198
if (u != nullopt && u.value().size() != 0) {
199199
PROXSUITE_CHECK_ARGUMENT_SIZE(
200-
u.value().rows(),
200+
u.value().size(),
201201
model.n_in,
202202
"the dimension wrt inequality constrained variables for initializing u "
203203
"is not valid.");
@@ -206,7 +206,7 @@ struct QP
206206
}
207207
if (l != nullopt && l.value().size() != 0) {
208208
PROXSUITE_CHECK_ARGUMENT_SIZE(
209-
l.value().rows(),
209+
l.value().size(),
210210
model.n_in,
211211
"the dimension wrt inequality constrained variables for initializing l "
212212
"is not valid.");
@@ -384,25 +384,25 @@ struct QP
384384

385385
// check the model is valid
386386
if (g != nullopt) {
387-
PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().rows(),
387+
PROXSUITE_CHECK_ARGUMENT_SIZE(g.value().size(),
388388
model.dim,
389389
"the dimension wrt the primal variable x "
390390
"variable for updating g is not valid.");
391391
}
392392
if (b != nullopt) {
393-
PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().rows(),
393+
PROXSUITE_CHECK_ARGUMENT_SIZE(b.value().size(),
394394
model.n_eq,
395395
"the dimension wrt equality constrained "
396396
"variables for updating b is not valid.");
397397
}
398398
if (u != nullopt) {
399-
PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().rows(),
399+
PROXSUITE_CHECK_ARGUMENT_SIZE(u.value().size(),
400400
model.n_in,
401401
"the dimension wrt inequality constrained "
402402
"variables for updating u is not valid.");
403403
}
404404
if (l != nullopt) {
405-
PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().rows(),
405+
PROXSUITE_CHECK_ARGUMENT_SIZE(l.value().size(),
406406
model.n_in,
407407
"the dimension wrt inequality constrained "
408408
"variables for updating l is not valid.");

test/src/dense_qp_solve.cpp

Lines changed: 69 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -12,6 +12,75 @@ using T = double;
1212
using namespace proxsuite;
1313
using namespace proxsuite::proxqp;
1414

15+
DOCTEST_TEST_CASE("proxqp::dense: test init with fixed sizes matrices")
16+
{
17+
double sparsity_factor = 0.15;
18+
T eps_abs = T(1e-9);
19+
utils::rand::set_seed(1);
20+
dense::isize dim = 10;
21+
22+
dense::isize n_eq(5), n_in(2);
23+
T strong_convexity_factor(1.e-2);
24+
proxqp::dense::Model<T> qp = proxqp::utils::dense_strongly_convex_qp(
25+
dim, n_eq, n_in, sparsity_factor, strong_convexity_factor);
26+
27+
Eigen::Matrix<T, 10, 10> H = qp.H;
28+
Eigen::Matrix<T, 10, 1> g = qp.g;
29+
Eigen::Matrix<T, 5, 10> A = qp.A;
30+
Eigen::Matrix<T, 5, 1> b = qp.b;
31+
Eigen::Matrix<T, 2, 10> C = qp.C;
32+
Eigen::Matrix<T, 2, 1> l = qp.l;
33+
Eigen::Matrix<T, 2, 1> u = qp.u;
34+
35+
{
36+
Results<T> results = dense::solve<T>(
37+
H, g, A, b, C, l, u, nullopt, nullopt, nullopt, eps_abs, 0);
38+
39+
T pri_res = std::max((qp.A * results.x - qp.b).lpNorm<Eigen::Infinity>(),
40+
(helpers::positive_part(qp.C * results.x - qp.u) +
41+
helpers::negative_part(qp.C * results.x - qp.l))
42+
.lpNorm<Eigen::Infinity>());
43+
T dua_res = (qp.H * results.x + qp.g + qp.A.transpose() * results.y +
44+
qp.C.transpose() * results.z)
45+
.lpNorm<Eigen::Infinity>();
46+
DOCTEST_CHECK(pri_res <= eps_abs);
47+
DOCTEST_CHECK(dua_res <= eps_abs);
48+
49+
std::cout << "------using API solving qp with dim: " << dim
50+
<< " neq: " << n_eq << " nin: " << n_in << std::endl;
51+
std::cout << "primal residual: " << pri_res << std::endl;
52+
std::cout << "dual residual: " << dua_res << std::endl;
53+
std::cout << "total number of iteration: " << results.info.iter
54+
<< std::endl;
55+
std::cout << "setup timing " << results.info.setup_time << " solve time "
56+
<< results.info.solve_time << std::endl;
57+
}
58+
59+
{
60+
dense::QP<T> qp_problem(dim, n_eq, 0);
61+
qp_problem.init(H, g, A, b, nullopt, nullopt, nullopt);
62+
qp_problem.settings.eps_abs = eps_abs;
63+
qp_problem.solve();
64+
65+
const Results<T>& results = qp_problem.results;
66+
67+
T pri_res = (qp.A * results.x - qp.b).lpNorm<Eigen::Infinity>();
68+
T dua_res = (qp.H * results.x + qp.g + qp.A.transpose() * results.y)
69+
.lpNorm<Eigen::Infinity>();
70+
DOCTEST_CHECK(pri_res <= eps_abs);
71+
DOCTEST_CHECK(dua_res <= eps_abs);
72+
73+
std::cout << "------using API solving qp with dim: " << dim
74+
<< " neq: " << n_eq << " nin: " << n_in << std::endl;
75+
std::cout << "primal residual: " << pri_res << std::endl;
76+
std::cout << "dual residual: " << dua_res << std::endl;
77+
std::cout << "total number of iteration: " << results.info.iter
78+
<< std::endl;
79+
std::cout << "setup timing " << results.info.setup_time << " solve time "
80+
<< results.info.solve_time << std::endl;
81+
}
82+
}
83+
1584
DOCTEST_TEST_CASE("sparse random strongly convex qp with equality and "
1685
"inequality constraints: test solve function")
1786
{

0 commit comments

Comments
 (0)