mosek_quadprog.cpp 9.4 KB

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