mosek_quadprog.cpp 9.9 KB

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