mosek_quadprog.cpp 9.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336
  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. mosek_guarded(MSK_initenv(env));
  108. // Create the optimization task
  109. mosek_guarded(MSK_maketask(env,m,n,&task));
  110. verbose("Creating task with %ld linear constraints and %ld variables...\n",m,n);
  111. //// Tell mosek how to print to std out
  112. //mosek_guarded(MSK_linkfunctotaskstream(task,MSK_STREAM_LOG,NULL,printstr));
  113. // Give estimate of number of variables
  114. mosek_guarded(MSK_putmaxnumvar(task,n));
  115. if(m>0)
  116. {
  117. // Give estimate of number of constraints
  118. mosek_guarded(MSK_putmaxnumcon(task,m));
  119. // Give estimate of number of non zeros in A
  120. mosek_guarded(MSK_putmaxnumanz(task,Av.size()));
  121. }
  122. // Give estimate of number of non zeros in Q
  123. mosek_guarded(MSK_putmaxnumqnz(task,Qv.size()));
  124. if(m>0)
  125. {
  126. // Append 'm' empty constraints, the constrainst will initially have no
  127. // bounds
  128. #if MSK_VERSION_MAJOR >= 7
  129. mosek_guarded(MSK_appendcons(task,m));
  130. #elif MSK_VERSION_MAJOR == 6
  131. mosek_guarded(MSK_append(task,MSK_ACC_CON,m));
  132. #endif
  133. }
  134. // Append 'n' variables
  135. #if MSK_VERSION_MAJOR >= 7
  136. mosek_guarded(MSK_appendvars(task,n));
  137. #elif MSK_VERSION_MAJOR == 6
  138. mosek_guarded(MSK_append(task,MSK_ACC_VAR,n));
  139. #endif
  140. // add a contant term to the objective
  141. mosek_guarded(MSK_putcfix(task,cf));
  142. // loop over variables
  143. for(int j = 0;j<n;j++)
  144. {
  145. if(c.size() > 0)
  146. {
  147. // Set linear term c_j in the objective
  148. mosek_guarded(MSK_putcj(task,j,c[j]));
  149. }
  150. // Set constant bounds on variable j
  151. if(lx[j] == ux[j])
  152. {
  153. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_FX,lx[j],ux[j]));
  154. }else
  155. {
  156. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_RA,lx[j],ux[j]));
  157. }
  158. if(m>0)
  159. {
  160. // Input column j of A
  161. #if MSK_VERSION_MAJOR >= 7
  162. mosek_guarded(
  163. MSK_putacol(
  164. task,
  165. j,
  166. Acp[j+1]-Acp[j],
  167. &Ari[Acp[j]],
  168. &Av[Acp[j]])
  169. );
  170. #elif MSK_VERSION_MAJOR == 6
  171. mosek_guarded(
  172. MSK_putavec(
  173. task,
  174. MSK_ACC_VAR,
  175. j,
  176. Acp[j+1]-Acp[j],
  177. &Ari[Acp[j]],
  178. &Av[Acp[j]])
  179. );
  180. #endif
  181. }
  182. }
  183. // loop over constraints
  184. for(int i = 0;i<m;i++)
  185. {
  186. // put bounds on constraints
  187. mosek_guarded(MSK_putbound(task,MSK_ACC_CON,i,MSK_BK_RA,lc[i],uc[i]));
  188. }
  189. // Input Q for the objective (REMEMBER Q SHOULD ONLY BE LOWER TRIANGLE)
  190. mosek_guarded(MSK_putqobj(task,Qv.size(),&Qi[0],&Qj[0],&Qv[0]));
  191. // Set up task parameters
  192. for(
  193. std::map<MSKiparame,int>::iterator pit = mosek_data.intparam.begin();
  194. pit != mosek_data.intparam.end();
  195. pit++)
  196. {
  197. mosek_guarded(MSK_putintparam(task,pit->first,pit->second));
  198. }
  199. for(
  200. std::map<MSKdparame,double>::iterator pit = mosek_data.douparam.begin();
  201. pit != mosek_data.douparam.end();
  202. pit++)
  203. {
  204. mosek_guarded(MSK_putdouparam(task,pit->first,pit->second));
  205. }
  206. // Now the optimizer has been prepared
  207. MSKrescodee trmcode;
  208. // run the optimizer
  209. mosek_guarded(MSK_optimizetrm(task,&trmcode));
  210. //// Print a summary containing information about the solution for debugging
  211. //// purposes
  212. //MSK_solutionsummary(task,MSK_STREAM_LOG);
  213. // Get status of solution
  214. MSKsolstae solsta;
  215. #if MSK_VERSION_MAJOR >= 7
  216. MSK_getsolsta (task,MSK_SOL_ITR,&solsta);
  217. #elif MSK_VERSION_MAJOR == 6
  218. MSK_getsolutionstatus(task,MSK_SOL_ITR,NULL,&solsta);
  219. #endif
  220. bool success = false;
  221. switch(solsta)
  222. {
  223. case MSK_SOL_STA_OPTIMAL:
  224. case MSK_SOL_STA_NEAR_OPTIMAL:
  225. MSK_getsolutionslice(task,MSK_SOL_ITR,MSK_SOL_ITEM_XX,0,n,&x[0]);
  226. //printf("Optimal primal solution\n");
  227. //for(size_t j=0; j<n; ++j)
  228. //{
  229. // printf("x[%ld]: %g\n",j,x[j]);
  230. //}
  231. success = true;
  232. break;
  233. case MSK_SOL_STA_DUAL_INFEAS_CER:
  234. case MSK_SOL_STA_PRIM_INFEAS_CER:
  235. case MSK_SOL_STA_NEAR_DUAL_INFEAS_CER:
  236. case MSK_SOL_STA_NEAR_PRIM_INFEAS_CER:
  237. //printf("Primal or dual infeasibility certificate found.\n");
  238. break;
  239. case MSK_SOL_STA_UNKNOWN:
  240. //printf("The status of the solution could not be determined.\n");
  241. break;
  242. default:
  243. //printf("Other solution status.");
  244. break;
  245. }
  246. MSK_deletetask(&task);
  247. MSK_deleteenv(&env);
  248. return success;
  249. }
  250. IGL_INLINE bool igl::mosek::mosek_quadprog(
  251. const Eigen::SparseMatrix<double> & Q,
  252. const Eigen::VectorXd & c,
  253. const double cf,
  254. const Eigen::SparseMatrix<double> & A,
  255. const Eigen::VectorXd & lc,
  256. const Eigen::VectorXd & uc,
  257. const Eigen::VectorXd & lx,
  258. const Eigen::VectorXd & ux,
  259. MosekData & mosek_data,
  260. Eigen::VectorXd & x)
  261. {
  262. using namespace Eigen;
  263. using namespace std;
  264. // Q should be square
  265. assert(Q.rows() == Q.cols());
  266. // Q should be symmetric
  267. #ifdef EIGEN_HAS_A_BUG_AND_FAILS_TO_LET_ME_COMPUTE_Q_MINUS_Q_TRANSPOSE
  268. assert( (Q-Q.transpose()).sum() < FLOAT_EPS);
  269. #endif
  270. // Only keep lower triangular part of Q
  271. SparseMatrix<double> QL;
  272. //QL = Q.template triangularView<Lower>();
  273. QL = Q.triangularView<Lower>();
  274. VectorXi Qi,Qj;
  275. VectorXd Qv;
  276. find(QL,Qi,Qj,Qv);
  277. vector<int> vQi = matrix_to_list(Qi);
  278. vector<int> vQj = matrix_to_list(Qj);
  279. vector<double> vQv = matrix_to_list(Qv);
  280. // Convert linear term
  281. vector<double> vc = matrix_to_list(c);
  282. assert(lc.size() == A.rows());
  283. assert(uc.size() == A.rows());
  284. // Convert A to harwell boeing format
  285. vector<double> vAv;
  286. vector<int> vAr,vAc;
  287. int nr;
  288. harwell_boeing(A,nr,vAv,vAr,vAc);
  289. assert(lx.size() == Q.rows());
  290. assert(ux.size() == Q.rows());
  291. vector<double> vlc = matrix_to_list(lc);
  292. vector<double> vuc = matrix_to_list(uc);
  293. vector<double> vlx = matrix_to_list(lx);
  294. vector<double> vux = matrix_to_list(ux);
  295. vector<double> vx;
  296. bool ret = mosek_quadprog(
  297. Q.rows(),vQi,vQj,vQv,
  298. vc,
  299. cf,
  300. nr,
  301. vAv, vAr, vAc,
  302. vlc,vuc,
  303. vlx,vux,
  304. mosek_data,
  305. vx);
  306. list_to_matrix(vx,x);
  307. return ret;
  308. }
  309. #ifdef IGL_STATIC_LIBRARY
  310. // Explicit template declarations
  311. #endif