@@ -12,6 +12,152 @@ using T = double;
1212using namespace proxsuite ;
1313using namespace proxsuite ::proxqp;
1414
15+ DOCTEST_TEST_CASE (" sparse random strongly convex qp with inequality constraints"
16+ " and empty equality constraints" )
17+ {
18+ std::cout << " ---testing sparse random strongly convex qp with inequality "
19+ " constraints "
20+ " and empty equality constraints---"
21+ << std::endl;
22+ double sparsity_factor = 0.15 ;
23+ T eps_abs = T (1e-9 );
24+ utils::rand::set_seed (1 );
25+ dense::isize dim = 10 ;
26+
27+ dense::isize n_eq (0 );
28+ dense::isize n_in (dim / 4 );
29+ T strong_convexity_factor (1 .e -2 );
30+ proxqp::dense::Model<T> qp_random = proxqp::utils::dense_strongly_convex_qp (
31+ dim, n_eq, n_in, sparsity_factor, strong_convexity_factor);
32+
33+ dense::QP<T> qp{ dim, n_eq, n_in }; // creating QP object
34+ qp.settings .eps_abs = eps_abs;
35+ qp.settings .eps_rel = 0 ;
36+
37+ // Testing with empty but properly sized matrix A of size (0, 10)
38+ std::cout << " Solving QP with" << std::endl;
39+ std::cout << " dim: " << dim << std::endl;
40+ std::cout << " n_eq: " << n_eq << std::endl;
41+ std::cout << " n_in: " << n_in << std::endl;
42+ std::cout << " H : " << qp_random.H << std::endl;
43+ std::cout << " g : " << qp_random.g << std::endl;
44+ std::cout << " A : " << qp_random.A << std::endl;
45+ std::cout << " A.cols() : " << qp_random.A .cols () << std::endl;
46+ std::cout << " A.rows() : " << qp_random.A .rows () << std::endl;
47+ std::cout << " b : " << qp_random.b << std::endl;
48+ std::cout << " C : " << qp_random.C << std::endl;
49+ std::cout << " u : " << qp_random.u << std::endl;
50+ std::cout << " l : " << qp_random.l << std::endl;
51+
52+ qp.init (qp_random.H ,
53+ qp_random.g ,
54+ std::nullopt ,
55+ qp_random.b ,
56+ qp_random.C ,
57+ qp_random.l ,
58+ qp_random.u );
59+ qp.solve ();
60+
61+ T pri_res = (dense::positive_part (qp_random.C * qp.results .x - qp_random.u ) +
62+ dense::negative_part (qp_random.C * qp.results .x - qp_random.l ))
63+ .lpNorm <Eigen::Infinity>();
64+ T dua_res = (qp_random.H * qp.results .x + qp_random.g +
65+ qp_random.C .transpose () * qp.results .z )
66+ .lpNorm <Eigen::Infinity>();
67+ DOCTEST_CHECK (pri_res <= eps_abs);
68+ DOCTEST_CHECK (dua_res <= eps_abs);
69+
70+ std::cout << " ------using API solving qp with dim: " << dim
71+ << " neq: " << n_eq << " nin: " << n_in << std::endl;
72+ std::cout << " primal residual: " << pri_res << std::endl;
73+ std::cout << " dual residual: " << dua_res << std::endl;
74+ std::cout << " total number of iteration: " << qp.results .info .iter
75+ << std::endl;
76+ std::cout << " setup timing " << qp.results .info .setup_time << " solve time "
77+ << qp.results .info .solve_time << std::endl;
78+
79+ // Testing with empty matrix A of size (0, 0)
80+ qp_random.A = Eigen::MatrixXd ();
81+ qp_random.b = Eigen::VectorXd ();
82+
83+ dense::QP<T> qp2{ dim, n_eq, n_in }; // creating QP object
84+ qp2.settings .eps_abs = eps_abs;
85+ qp2.settings .eps_rel = 0 ;
86+
87+ std::cout << " Solving QP with" << std::endl;
88+ std::cout << " dim: " << dim << std::endl;
89+ std::cout << " n_eq: " << n_eq << std::endl;
90+ std::cout << " n_in: " << n_in << std::endl;
91+ std::cout << " H : " << qp_random.H << std::endl;
92+ std::cout << " g : " << qp_random.g << std::endl;
93+ std::cout << " A : " << qp_random.A << std::endl;
94+ std::cout << " A.cols() : " << qp_random.A .cols () << std::endl;
95+ std::cout << " A.rows() : " << qp_random.A .rows () << std::endl;
96+ std::cout << " b : " << qp_random.b << std::endl;
97+ std::cout << " C : " << qp_random.C << std::endl;
98+ std::cout << " u : " << qp_random.u << std::endl;
99+ std::cout << " l : " << qp_random.l << std::endl;
100+
101+ qp2.init (qp_random.H ,
102+ qp_random.g ,
103+ qp_random.A ,
104+ qp_random.b ,
105+ qp_random.C ,
106+ qp_random.l ,
107+ qp_random.u );
108+ qp2.solve ();
109+
110+ pri_res = (dense::positive_part (qp_random.C * qp.results .x - qp_random.u ) +
111+ dense::negative_part (qp_random.C * qp.results .x - qp_random.l ))
112+ .lpNorm <Eigen::Infinity>();
113+ dua_res = (qp_random.H * qp.results .x + qp_random.g +
114+ qp_random.C .transpose () * qp.results .z )
115+ .lpNorm <Eigen::Infinity>();
116+ DOCTEST_CHECK (pri_res <= eps_abs);
117+ DOCTEST_CHECK (dua_res <= eps_abs);
118+
119+ std::cout << " ------using API solving qp with dim: " << dim
120+ << " neq: " << n_eq << " nin: " << n_in << std::endl;
121+ std::cout << " primal residual: " << pri_res << std::endl;
122+ std::cout << " dual residual: " << dua_res << std::endl;
123+ std::cout << " total number of iteration: " << qp.results .info .iter
124+ << std::endl;
125+ std::cout << " setup timing " << qp.results .info .setup_time << " solve time "
126+ << qp.results .info .solve_time << std::endl;
127+
128+ // Testing with nullopt
129+ dense::QP<T> qp3{ dim, n_eq, n_in }; // creating QP object
130+ qp3.settings .eps_abs = eps_abs;
131+ qp3.settings .eps_rel = 0 ;
132+
133+ qp3.init (qp_random.H ,
134+ qp_random.g ,
135+ std::nullopt ,
136+ std::nullopt ,
137+ qp_random.C ,
138+ qp_random.l ,
139+ qp_random.u );
140+ qp3.solve ();
141+
142+ pri_res = (dense::positive_part (qp_random.C * qp.results .x - qp_random.u ) +
143+ dense::negative_part (qp_random.C * qp.results .x - qp_random.l ))
144+ .lpNorm <Eigen::Infinity>();
145+ dua_res = (qp_random.H * qp.results .x + qp_random.g +
146+ qp_random.C .transpose () * qp.results .z )
147+ .lpNorm <Eigen::Infinity>();
148+ DOCTEST_CHECK (pri_res <= eps_abs);
149+ DOCTEST_CHECK (dua_res <= eps_abs);
150+
151+ std::cout << " ------using API solving qp with dim: " << dim
152+ << " neq: " << n_eq << " nin: " << n_in << std::endl;
153+ std::cout << " primal residual: " << pri_res << std::endl;
154+ std::cout << " dual residual: " << dua_res << std::endl;
155+ std::cout << " total number of iteration: " << qp.results .info .iter
156+ << std::endl;
157+ std::cout << " setup timing " << qp.results .info .setup_time << " solve time "
158+ << qp.results .info .solve_time << std::endl;
159+ }
160+
15161DOCTEST_TEST_CASE (" sparse random strongly convex qp with equality and "
16162 " inequality constraints: test update H" )
17163{
0 commit comments