mosek_quadprog.cpp 10.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "mosek_quadprog.h"
  9. #include "mosek_guarded.h"
  10. #include <cstdio>
  11. #include "../find.h"
  12. #include "../verbose.h"
  13. #include "../speye.h"
  14. #include "../matrix_to_list.h"
  15. #include "../list_to_matrix.h"
  16. #include "../harwell_boeing.h"
  17. #include "../EPS.h"
  18. igl::mosek::MosekData::MosekData()
  19. {
  20. // These are the default settings that worked well for BBW. Your miles may
  21. // very well be kilometers.
  22. // >1e0 NONSOLUTION
  23. // 1e-1 artifacts in deformation
  24. // 1e-3 artifacts in isolines
  25. // 1e-4 seems safe
  26. // 1e-8 MOSEK DEFAULT SOLUTION
  27. douparam[MSK_DPAR_INTPNT_TOL_REL_GAP]=1e-8;
  28. #if MSK_VERSION_MAJOR >= 8
  29. douparam[MSK_DPAR_INTPNT_QO_TOL_REL_GAP]=1e-12;
  30. #endif
  31. // Force using multiple threads, not sure if MOSEK is properly destroying
  32. //extra threads...
  33. #if MSK_VERSION_MAJOR >= 7
  34. intparam[MSK_IPAR_NUM_THREADS] = 6;
  35. #elif MSK_VERSION_MAJOR == 6
  36. intparam[MSK_IPAR_INTPNT_NUM_THREADS] = 6;
  37. #endif
  38. #if MSK_VERSION_MAJOR == 6
  39. // Force turn off data check
  40. intparam[MSK_IPAR_DATA_CHECK]=MSK_OFF;
  41. #endif
  42. // Turn off presolving
  43. // intparam[MSK_IPAR_PRESOLVE_USE] = MSK_PRESOLVE_MODE_OFF;
  44. // Force particular matrix reordering method
  45. // MSK_ORDER_METHOD_NONE cuts time in half roughly, since half the time is
  46. // usually spent reordering the matrix
  47. // !! WARNING Setting this parameter to anything but MSK_ORDER_METHOD_FREE
  48. // seems to have the effect of setting it to MSK_ORDER_METHOD_NONE
  49. // *Or maybe Mosek is spending a bunch of time analyzing the matrix to
  50. // choose the right ordering method when really any of them are
  51. // instantaneous
  52. intparam[MSK_IPAR_INTPNT_ORDER_METHOD] = MSK_ORDER_METHOD_NONE;
  53. // 1.0 means optimizer is very lenient about declaring model infeasible
  54. douparam[MSK_DPAR_INTPNT_TOL_INFEAS] = 1e-8;
  55. // Hard to say if this is doing anything, probably nothing dramatic
  56. douparam[MSK_DPAR_INTPNT_TOL_PSAFE]= 1e2;
  57. // Turn off convexity check
  58. intparam[MSK_IPAR_CHECK_CONVEXITY] = MSK_CHECK_CONVEXITY_NONE;
  59. }
  60. template <typename Index, typename Scalar>
  61. IGL_INLINE bool igl::mosek::mosek_quadprog(
  62. const Index n,
  63. std::vector<Index> & Qi,
  64. std::vector<Index> & Qj,
  65. std::vector<Scalar> & Qv,
  66. const std::vector<Scalar> & c,
  67. const Scalar cf,
  68. const Index m,
  69. std::vector<Scalar> & Av,
  70. std::vector<Index> & Ari,
  71. const std::vector<Index> & Acp,
  72. const std::vector<Scalar> & lc,
  73. const std::vector<Scalar> & uc,
  74. const std::vector<Scalar> & lx,
  75. const std::vector<Scalar> & ux,
  76. MosekData & mosek_data,
  77. std::vector<Scalar> & x)
  78. {
  79. // I J V vectors of Q should all be same length
  80. assert(Qv.size() == Qi.size());
  81. assert(Qv.size() == Qj.size());
  82. // number of columns in linear constraint matrix must be ≤ number of
  83. // variables
  84. assert( (int)Acp.size() == (n+1));
  85. // linear bound vectors must be size of number of constraints or empty
  86. assert( ((int)lc.size() == m) || ((int)lc.size() == 0));
  87. assert( ((int)uc.size() == m) || ((int)uc.size() == 0));
  88. // constant bound vectors must be size of number of variables or empty
  89. assert( ((int)lx.size() == n) || ((int)lx.size() == 0));
  90. assert( ((int)ux.size() == n) || ((int)ux.size() == 0));
  91. // allocate space for solution in x
  92. x.resize(n);
  93. // variables for mosek task, env and result code
  94. MSKenv_t env;
  95. MSKtask_t task;
  96. // Create the MOSEK environment
  97. #if MSK_VERSION_MAJOR >= 7
  98. mosek_guarded(MSK_makeenv(&env,NULL));
  99. #elif MSK_VERSION_MAJOR == 6
  100. mosek_guarded(MSK_makeenv(&env,NULL,NULL,NULL,NULL));
  101. #endif
  102. ///* Directs the log stream to the 'printstr' function. */
  103. //// Little function mosek needs in order to know how to print to std out
  104. //const auto & printstr = [](void *handle, char str[])
  105. //{
  106. // printf("%s",str);
  107. //}
  108. //mosek_guarded(MSK_linkfunctoenvstream(env,MSK_STREAM_LOG,NULL,printstr));
  109. // initialize mosek environment
  110. #if MSK_VERSION_MAJOR <= 7
  111. mosek_guarded(MSK_initenv(env));
  112. #endif
  113. // Create the optimization task
  114. mosek_guarded(MSK_maketask(env,m,n,&task));
  115. verbose("Creating task with %ld linear constraints and %ld variables...\n",m,n);
  116. //// Tell mosek how to print to std out
  117. //mosek_guarded(MSK_linkfunctotaskstream(task,MSK_STREAM_LOG,NULL,printstr));
  118. // Give estimate of number of variables
  119. mosek_guarded(MSK_putmaxnumvar(task,n));
  120. if(m>0)
  121. {
  122. // Give estimate of number of constraints
  123. mosek_guarded(MSK_putmaxnumcon(task,m));
  124. // Give estimate of number of non zeros in A
  125. mosek_guarded(MSK_putmaxnumanz(task,Av.size()));
  126. }
  127. // Give estimate of number of non zeros in Q
  128. mosek_guarded(MSK_putmaxnumqnz(task,Qv.size()));
  129. if(m>0)
  130. {
  131. // Append 'm' empty constraints, the constrainst will initially have no
  132. // bounds
  133. #if MSK_VERSION_MAJOR >= 7
  134. mosek_guarded(MSK_appendcons(task,m));
  135. #elif MSK_VERSION_MAJOR == 6
  136. mosek_guarded(MSK_append(task,MSK_ACC_CON,m));
  137. #endif
  138. }
  139. // Append 'n' variables
  140. #if MSK_VERSION_MAJOR >= 7
  141. mosek_guarded(MSK_appendvars(task,n));
  142. #elif MSK_VERSION_MAJOR == 6
  143. mosek_guarded(MSK_append(task,MSK_ACC_VAR,n));
  144. #endif
  145. // add a contant term to the objective
  146. mosek_guarded(MSK_putcfix(task,cf));
  147. // loop over variables
  148. for(int j = 0;j<n;j++)
  149. {
  150. if(c.size() > 0)
  151. {
  152. // Set linear term c_j in the objective
  153. mosek_guarded(MSK_putcj(task,j,c[j]));
  154. }
  155. // Set constant bounds on variable j
  156. if(lx[j] == ux[j])
  157. {
  158. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_FX,lx[j],ux[j]));
  159. }else
  160. {
  161. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_RA,lx[j],ux[j]));
  162. }
  163. if(m>0)
  164. {
  165. // Input column j of A
  166. #if MSK_VERSION_MAJOR >= 7
  167. mosek_guarded(
  168. MSK_putacol(
  169. task,
  170. j,
  171. Acp[j+1]-Acp[j],
  172. &Ari[Acp[j]],
  173. &Av[Acp[j]])
  174. );
  175. #elif MSK_VERSION_MAJOR == 6
  176. mosek_guarded(
  177. MSK_putavec(
  178. task,
  179. MSK_ACC_VAR,
  180. j,
  181. Acp[j+1]-Acp[j],
  182. &Ari[Acp[j]],
  183. &Av[Acp[j]])
  184. );
  185. #endif
  186. }
  187. }
  188. // loop over constraints
  189. for(int i = 0;i<m;i++)
  190. {
  191. // put bounds on constraints
  192. mosek_guarded(MSK_putbound(task,MSK_ACC_CON,i,MSK_BK_RA,lc[i],uc[i]));
  193. }
  194. // Input Q for the objective (REMEMBER Q SHOULD ONLY BE LOWER TRIANGLE)
  195. mosek_guarded(MSK_putqobj(task,Qv.size(),&Qi[0],&Qj[0],&Qv[0]));
  196. // Set up task parameters
  197. for(
  198. std::map<MSKiparame,int>::iterator pit = mosek_data.intparam.begin();
  199. pit != mosek_data.intparam.end();
  200. pit++)
  201. {
  202. mosek_guarded(MSK_putintparam(task,pit->first,pit->second));
  203. }
  204. for(
  205. std::map<MSKdparame,double>::iterator pit = mosek_data.douparam.begin();
  206. pit != mosek_data.douparam.end();
  207. pit++)
  208. {
  209. mosek_guarded(MSK_putdouparam(task,pit->first,pit->second));
  210. }
  211. // Now the optimizer has been prepared
  212. MSKrescodee trmcode;
  213. // run the optimizer
  214. mosek_guarded(MSK_optimizetrm(task,&trmcode));
  215. //// Print a summary containing information about the solution for debugging
  216. //// purposes
  217. //MSK_solutionsummary(task,MSK_STREAM_LOG);
  218. // Get status of solution
  219. MSKsolstae solsta;
  220. #if MSK_VERSION_MAJOR >= 7
  221. MSK_getsolsta (task,MSK_SOL_ITR,&solsta);
  222. #elif MSK_VERSION_MAJOR == 6
  223. MSK_getsolutionstatus(task,MSK_SOL_ITR,NULL,&solsta);
  224. #endif
  225. bool success = false;
  226. switch(solsta)
  227. {
  228. case MSK_SOL_STA_OPTIMAL:
  229. case MSK_SOL_STA_NEAR_OPTIMAL:
  230. MSK_getsolutionslice(task,MSK_SOL_ITR,MSK_SOL_ITEM_XX,0,n,&x[0]);
  231. //printf("Optimal primal solution\n");
  232. //for(size_t j=0; j<n; ++j)
  233. //{
  234. // printf("x[%ld]: %g\n",j,x[j]);
  235. //}
  236. success = true;
  237. break;
  238. case MSK_SOL_STA_DUAL_INFEAS_CER:
  239. case MSK_SOL_STA_PRIM_INFEAS_CER:
  240. case MSK_SOL_STA_NEAR_DUAL_INFEAS_CER:
  241. case MSK_SOL_STA_NEAR_PRIM_INFEAS_CER:
  242. //printf("Primal or dual infeasibility certificate found.\n");
  243. break;
  244. case MSK_SOL_STA_UNKNOWN:
  245. //printf("The status of the solution could not be determined.\n");
  246. break;
  247. default:
  248. //printf("Other solution status.");
  249. break;
  250. }
  251. MSK_deletetask(&task);
  252. MSK_deleteenv(&env);
  253. return success;
  254. }
  255. IGL_INLINE bool igl::mosek::mosek_quadprog(
  256. const Eigen::SparseMatrix<double> & Q,
  257. const Eigen::VectorXd & c,
  258. const double cf,
  259. const Eigen::SparseMatrix<double> & A,
  260. const Eigen::VectorXd & lc,
  261. const Eigen::VectorXd & uc,
  262. const Eigen::VectorXd & lx,
  263. const Eigen::VectorXd & ux,
  264. MosekData & mosek_data,
  265. Eigen::VectorXd & x)
  266. {
  267. using namespace Eigen;
  268. using namespace std;
  269. typedef int Index;
  270. typedef double Scalar;
  271. // Q should be square
  272. assert(Q.rows() == Q.cols());
  273. // Q should be symmetric
  274. #ifdef EIGEN_HAS_A_BUG_AND_FAILS_TO_LET_ME_COMPUTE_Q_MINUS_Q_TRANSPOSE
  275. assert( (Q-Q.transpose()).sum() < FLOAT_EPS);
  276. #endif
  277. // Only keep lower triangular part of Q
  278. SparseMatrix<Scalar> QL;
  279. //QL = Q.template triangularView<Lower>();
  280. QL = Q.triangularView<Lower>();
  281. VectorXi Qi,Qj;
  282. VectorXd Qv;
  283. find(QL,Qi,Qj,Qv);
  284. vector<Index> vQi = matrix_to_list(Qi);
  285. vector<Index> vQj = matrix_to_list(Qj);
  286. vector<Scalar> vQv = matrix_to_list(Qv);
  287. // Convert linear term
  288. vector<Scalar> vc = matrix_to_list(c);
  289. assert(lc.size() == A.rows());
  290. assert(uc.size() == A.rows());
  291. // Convert A to harwell boeing format
  292. vector<Scalar> vAv;
  293. vector<Index> vAr,vAc;
  294. Index nr;
  295. harwell_boeing(A,nr,vAv,vAr,vAc);
  296. assert(lx.size() == Q.rows());
  297. assert(ux.size() == Q.rows());
  298. vector<Scalar> vlc = matrix_to_list(lc);
  299. vector<Scalar> vuc = matrix_to_list(uc);
  300. vector<Scalar> vlx = matrix_to_list(lx);
  301. vector<Scalar> vux = matrix_to_list(ux);
  302. vector<Scalar> vx;
  303. bool ret = mosek_quadprog<Index,Scalar>(
  304. Q.rows(),vQi,vQj,vQv,
  305. vc,
  306. cf,
  307. nr,
  308. vAv, vAr, vAc,
  309. vlc,vuc,
  310. vlx,vux,
  311. mosek_data,
  312. vx);
  313. list_to_matrix(vx,x);
  314. return ret;
  315. }
  316. #ifdef IGL_STATIC_LIBRARY
  317. // Explicit template declarations
  318. #endif