Skip to content

Commit 7ec10c4

Browse files
author
Antoine Bambade
committed
add examples on how initializing QP model in the documentation
1 parent 7918a9c commit 7ec10c4

File tree

6 files changed

+117
-2
lines changed

6 files changed

+117
-2
lines changed

doc/2-PROXQP_API/2-ProxQP_api.md

Lines changed: 17 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -97,7 +97,7 @@ For loading ProxQP with dense backend it is as simple as the following code belo
9797
</tr>
9898
</table>
9999

100-
The dimensions of the problem (i.e., n is the dimension of primal variable x, n_eq the number of equality constraints, and n_in the number of inequality constraints) are used for allocating the space needed for the Qp object. The dense Qp object is templated by the floatting precision of the QP model (in the example above in C++ a double precision).
100+
The dimensions of the problem (i.e., n is the dimension of primal variable x, n_eq the number of equality constraints, and n_in the number of inequality constraints) are used for allocating the space needed for the Qp object. The dense Qp object is templated by the floatting precision of the QP model (in the example above in C++ a double precision). Note that for model to be valid, the primal dimension (i.e., n) must be strictly positive. If it is not the case an assertion will be raised precising this issue.
101101

102102
For loading ProxQP with sparse backend they are two possibilities:
103103
* one can use as before the dimensions of the QP problem (i.e., n, n_eq and n_in)
@@ -139,7 +139,22 @@ Once you have defined a Qp object, the init method enables you setting up the QP
139139
</tr>
140140
</table>
141141

142-
Note that with its dense backend, ProxQP solver can manipulate both matrices in dense and sparse representations (in the example above the matrices are in sparse format).
142+
Note that with its dense backend, ProxQP solver can manipulate both matrices in dense and sparse representations (in the example above the matrices are in sparse format). Note that if some elements of your QP model are not defined (for example a QP without linear cost or inequality constraints), you can either pass a None argument, or a matrix with zero shape for specifying it. We provide an example below in cpp and python (for the dense case, it is similar with sparse backend).
143+
144+
<table class="manual">
145+
<tr>
146+
<th>examples/cpp/initializing_with_none.cpp</th>
147+
<th>examples/python/initializing_with_none.py</th>
148+
</tr>
149+
<tr>
150+
<td valign="top">
151+
\include initializing_with_none.cpp
152+
</td>
153+
<td valign="top">
154+
\include initializing_with_none.py
155+
</td>
156+
</tr>
157+
</table>
143158

144159
With the init method, you can also setting-up on the same time some other parameters in the following order:
145160
* compute_preconditioner: a boolean parameter for executing or not the preconditioner. The preconditioner is an algorithm used (for the moment we use [Ruiz equilibrator](https://cds.cern.ch/record/585592/files/CM-P00040415.pdf)) for reducing the ill-conditioning of the QP problem, and hence speeding-up the solver and increasing its accuracy. It consists mostly of an heuristic involving linear scalings. Note that for very ill-conditioned QP problem, when one asks for a very accurate solution, the unscaling procedure can become less precise (we provide some remarks about this subject in section 6.D of the [following paper](https://hal.inria.fr/hal-03683733/file/Yet_another_QP_solver_for_robotics_and_beyond.pdf)). By default its value is set to true.

doc/3-ProxQP_solve.md

Lines changed: 17 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -130,3 +130,20 @@ Finally, note that in C++, if you want to change one option, the order described
130130
</td>
131131
</tr>
132132
</table>
133+
134+
Note that if some elements of your QP model are not defined (for example a QP without linear cost or inequality constraints), you can either pass a None argument, or a matrix with zero shape for specifying it. We provide an example below in cpp and python (for the dense case, it is similar with sparse backend).
135+
136+
<table class="manual">
137+
<tr>
138+
<th>examples/cpp/initializing_with_none_without_api.cpp</th>
139+
<th>examples/python/initializing_with_none_without_api.py</th>
140+
</tr>
141+
<tr>
142+
<td valign="top">
143+
\include initializing_with_none_without_api.cpp
144+
</td>
145+
<td valign="top">
146+
\include initializing_with_none_without_api.py
147+
</td>
148+
</tr>
149+
</table>
Lines changed: 28 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,28 @@
1+
#include "proxsuite/proxqp/dense/dense.hpp"
2+
#include <proxsuite/proxqp/utils/random_qp_problems.hpp> // used for generating a random convex Qp
3+
4+
using namespace proxsuite::proxqp;
5+
using T = double;
6+
int
7+
main()
8+
{
9+
dense::isize dim = 10;
10+
dense::isize n_eq(0);
11+
dense::isize n_in(0);
12+
dense::QP<T> Qp(dim, n_eq, n_in);
13+
T strong_convexity_factor(0.1);
14+
T sparsity_factor(0.15);
15+
// we generate a qp, so the function used from helpers.hpp is
16+
// in proxqp namespace. The qp is in dense eigen format and
17+
// you can control its sparsity ratio and strong convexity factor.
18+
dense::Model<T> qp = utils::dense_strongly_convex_qp(
19+
dim, n_eq, n_in, sparsity_factor, strong_convexity_factor);
20+
21+
Qp.init(qp.H, qp.g,qp.A,qp.b,qp.C,qp.u,qp.l); // initialization with zero shape matrices
22+
//it is equivalent to do Qp.init(qp.H, qp.g, std::nullopt,std::nullopt,std::nullopt,std::nullopt,std::nullopt);
23+
Qp.solve();
24+
// print an optimal solution x,y and z
25+
std::cout << "optimal x: " << Qp.results.x << std::endl;
26+
std::cout << "optimal y: " << Qp.results.y << std::endl;
27+
std::cout << "optimal z: " << Qp.results.z << std::endl;
28+
}
Lines changed: 27 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,27 @@
1+
#include "proxsuite/proxqp/dense/dense.hpp"
2+
#include <proxsuite/proxqp/utils/random_qp_problems.hpp> // used for generating a random convex Qp
3+
4+
using namespace proxsuite::proxqp;
5+
using T = double;
6+
int
7+
main()
8+
{
9+
dense::isize dim = 10;
10+
dense::isize n_eq(0);
11+
dense::isize n_in(0);
12+
T strong_convexity_factor(0.1);
13+
T sparsity_factor(0.15);
14+
// we generate a qp, so the function used from helpers.hpp is
15+
// in proxqp namespace. The qp is in dense eigen format and
16+
// you can control its sparsity ratio and strong convexity factor.
17+
dense::Model<T> qp = utils::dense_strongly_convex_qp(
18+
dim, n_eq, n_in, sparsity_factor, strong_convexity_factor);
19+
20+
Results<T> results =
21+
dense::solve<T>(qp.H, qp.g,qp.A,qp.b,qp.C,qp.u,qp.l); // initialization with zero shape matrices
22+
//it is equivalent to do dense::solve<T>(qp.H, qp.g, std::nullopt,std::nullopt,std::nullopt,std::nullopt,std::nullopt);
23+
// print an optimal solution x,y and z
24+
std::cout << "optimal x: " << results.x << std::endl;
25+
std::cout << "optimal y: " << results.y << std::endl;
26+
std::cout << "optimal z: " << results.z << std::endl;
27+
}
Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,15 @@
1+
import numpy as np
2+
import proxsuite
3+
4+
H = np.array([[65.0, -22.0, -16.0], [-22.0, 14.0, 7.0], [-16.0, 7.0, 5.0]])
5+
g = np.array([-13.0, 15.0, 7.0])
6+
A = None
7+
b = None
8+
C = None
9+
u = None
10+
l = None
11+
12+
Qp = proxsuite.proxqp.dense.QP(3,0,0)
13+
Qp.init(H, g, A, b, C, u, l) # it is equivalent to do as well Qp.init(H, g)
14+
Qp.solve()
15+
print("optimal x: {}".format(Qp.results.x))
Lines changed: 13 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,13 @@
1+
import numpy as np
2+
import proxsuite
3+
4+
H = np.array([[65.0, -22.0, -16.0], [-22.0, 14.0, 7.0], [-16.0, 7.0, 5.0]])
5+
g = np.array([-13.0, 15.0, 7.0])
6+
A = None
7+
b = None
8+
C = None
9+
u = None
10+
l = None
11+
12+
results = proxsuite.proxqp.dense.solve(H, g, A, b, C, u, l) # it is equivalent to do as well proxsuite.proxqp.dense.solve(H, g)
13+
print("optimal x: {}".format(results.x))

0 commit comments

Comments
 (0)