11// VolEsti (volume computation and sampling library)
22
3- // Copyright (c) 2012-2020 Vissarion Fisikopoulos
4- // Copyright (c) 2018-2020 Apostolos Chalkis
3+ // Copyright (c) 2012-2024 Vissarion Fisikopoulos
4+ // Copyright (c) 2018-2024 Apostolos Chalkis
55// Copyright (c) 2021 Vaibhav Thakkar
66
77// Contributed and/or modified by Alexandros Manochis, as part of Google Summer of Code 2020 program.
1515
1616#include < utility>
1717#include < Eigen/Eigen>
18+ #include " Spectra/include/Spectra/SymEigsSolver.h"
19+ #include " Spectra/include/Spectra/Util/SelectionRule.h"
20+ #include " Spectra/include/Spectra/MatOp/DenseSymMatProd.h"
1821
1922
2023/*
3134 tolerance parameters tol, reg
3235
3336 Output: center of the ellipsoid y
34- matrix V = E_transpose * E
37+ matrix E2^{-1} = E_transpose * E
3538*/
3639
3740// using Custom_MT as to deal with both dense and sparse matrices, MT will be the type of result matrix
41+ // TODO: Change the return data structure to std::tuple
3842template <typename MT, typename Custom_MT, typename VT, typename NT>
3943std::pair<std::pair<MT, VT>, bool > max_inscribed_ellipsoid (Custom_MT A, VT b, VT const & x0,
40- unsigned int const & maxiter,
41- NT const & tol, NT const & reg)
44+ unsigned int const & maxiter,
45+ NT const & tol, NT const & reg)
4246{
4347 typedef Eigen::DiagonalMatrix<NT, Eigen::Dynamic> Diagonal_MT;
4448
@@ -49,8 +53,8 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
4953 last_r1 = std::numeric_limits<NT>::lowest (),
5054 last_r2 = std::numeric_limits<NT>::lowest (),
5155 prev_obj = std::numeric_limits<NT>::lowest (),
52- gap, rmu, res, objval, r1, r2 ,r3, rel, Rel ,
53- astep, ax, ay, az, tau;
56+ gap, rmu, res, objval, r1, r2 ,r3, astep, ax ,
57+ ay, az, tau, logdetE2 ;
5458
5559 NT const reg_lim = std::pow (10.0 , -10.0 ), tau0 = 0.75 , minmu = std::pow (10.0 , -8.0 );
5660
@@ -74,9 +78,10 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
7478
7579 Y = y.asDiagonal ();
7680
77- E2 .noalias () = MT (A_trans * Y * A).inverse ();
81+ E2 .noalias () = MT (A_trans * Y * A);
82+ Eigen::LLT<MT> llt (E2 );
7883
79- Q.noalias () = A * E2 * A_trans;
84+ Q.noalias () = A * llt. solve ( A_trans) ;
8085 h = Q.diagonal ();
8186 h = h.cwiseSqrt ();
8287
@@ -115,19 +120,38 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
115120
116121 res = std::max (r1, r2);
117122 res = std::max (res, r3);
118- objval = std::log (E2 .determinant ()) / 2.0 ;
119-
120- Eigen::SelfAdjointEigenSolver <MT> eigensolver (E2 ); // E2 is positive definite matrix
121- // computing eigenvalues of E2
122- rel = eigensolver.eigenvalues ().minCoeff ();
123- Rel = eigensolver.eigenvalues ().maxCoeff ();
123+ logdetE2 = llt.matrixL ().toDenseMatrix ().diagonal ().array ().log ().sum ();
124+ objval = logdetE2; // logdet of E2 is already divided by 2
124125
125126 if (i % 10 == 0 ) {
127+
128+ NT rel, Rel;
129+
130+ // computing eigenvalues of E2
131+ Spectra::DenseSymMatProd<NT> op (E2 );
132+ // The value of ncv is chosen empirically
133+ Spectra::SymEigsSolver<NT, Spectra::SELECT_EIGENVALUE::BOTH_ENDS,
134+ Spectra::DenseSymMatProd<NT>> eigs (&op, 2 , std::min (std::max (10 , n/5 ), n));
135+ eigs.init ();
136+ int nconv = eigs.compute ();
137+ if (eigs.info () == Spectra::COMPUTATION_INFO::SUCCESSFUL) {
138+ Rel = 1.0 / eigs.eigenvalues ().coeff (1 );
139+ rel = 1.0 / eigs.eigenvalues ().coeff (0 );
140+ } else {
141+ Eigen::SelfAdjointEigenSolver<MT> eigensolver (E2 ); // E2 is positive definite matrix
142+ if (eigensolver.info () == Eigen::ComputationInfo::Success) {
143+ Rel = 1.0 / eigensolver.eigenvalues ().coeff (0 );
144+ rel = 1.0 / eigensolver.eigenvalues ().template tail <1 >().value ();
145+ } else {
146+ std::runtime_error (" Computations failed." );
147+ }
148+ }
126149
127150 if (std::abs ((last_r1 - r1) / std::min (NT (std::abs (last_r1)), NT (std::abs (r1)))) < 0.01 &&
128151 std::abs ((last_r2 - r2) / std::min (NT (abs (last_r2)), NT (std::abs (r2)))) < 0.01 &&
129152 Rel / rel > 100.0 &&
130153 reg > reg_lim) {
154+
131155 converged = false ;
132156 // Stopped making progress
133157 break ;
@@ -138,9 +162,10 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
138162
139163 // stopping criterion
140164 if ((res < tol * (1.0 + bnrm) && rmu <= minmu) ||
141- (i > 4 && prev_obj != std::numeric_limits<NT>::lowest () &&
142- ((std::abs (objval - prev_obj) <= tol * objval && std::abs (objval - prev_obj) <= tol * prev_obj) ||
143- (prev_obj >= (1.0 - tol) * objval || objval <= (1.0 - tol) * prev_obj) ) ) ) {
165+ (i > 1 && prev_obj != std::numeric_limits<NT>::lowest () &&
166+ (std::abs (objval - prev_obj) <= tol * std::min (std::abs (objval), std::abs (prev_obj)) ||
167+ std::abs (prev_obj - objval) <= tol) ) ) {
168+
144169 // converged
145170 x += x0;
146171 converged = true ;
@@ -151,7 +176,7 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
151176 YQ.noalias () = Y * Q;
152177 G = YQ.cwiseProduct (YQ.transpose ());
153178 y2h = 2.0 * yh;
154- YA = Y * A;
179+ YA. noalias () = Y * A;
155180
156181 vec_iter1 = y2h.data ();
157182 vec_iter2 = z.data ();
@@ -165,10 +190,9 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
165190
166191 G.diagonal () += y2h_z;
167192 h_z = h + z;
193+ Eigen::PartialPivLU<MT> luG (G);
194+ T.noalias () = luG.solve (h_z.asDiagonal ()*YA);
168195
169- for (int j = 0 ; j < n; ++j) {
170- T.col (j) = G.colPivHouseholderQr ().solve ( VT (YA.col (j).cwiseProduct (h_z)) );
171- }
172196 ATP.noalias () = MT (y2h.asDiagonal ()*T - YA).transpose ();
173197
174198 vec_iter1 = R3.data ();
@@ -184,11 +208,11 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
184208 R23 = R2 - R3Dy;
185209 ATP_A.noalias () = ATP * A;
186210 ATP_A.diagonal () += ones_m * reg;
187- dx = ATP_A.colPivHouseholderQr ().solve (R1 + ATP * R23); // predictor step
211+ dx = ATP_A.lu ().solve (R1 + ATP * R23); // predictor step
188212
189213 // corrector and combined step & length
190214 Adx.noalias () = A * dx;
191- dyDy = G. colPivHouseholderQr () .solve (y2h.cwiseProduct (Adx-R23));
215+ dyDy = luG .solve (y2h.cwiseProduct (Adx-R23));
192216
193217 dy = y.cwiseProduct (dyDy);
194218 dz = R3Dy - z.cwiseProduct (dyDy);
@@ -228,9 +252,13 @@ std::pair<std::pair<MT, VT>, bool> max_inscribed_ellipsoid(Custom_MT A, VT b, VT
228252 i++;
229253 }
230254
231- return std::pair<std::pair<MT, VT>, bool >(std::pair<MT, VT>(E2 , x), converged);
255+ if (!converged) {
256+ x += x0;
257+ }
232258
259+ return std::pair<std::pair<MT, VT>, bool >(std::pair<MT, VT>(E2 , x), converged);
233260}
234261
235262
236263#endif
264+
0 commit comments