mosek_quadprog.cpp 8.7 KB

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