mosek_quadprog.cpp 7.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226
  1. #include "mosek_quadprog.h"
  2. #include "mosek_guarded.h"
  3. #include <cstdio>
  4. #include "../verbose.h"
  5. // Little function mosek needs in order to know how to print to std out
  6. static void MSKAPI printstr(void *handle, char str[])
  7. {
  8. printf("%s",str);
  9. }
  10. template <typename Index, typename Scalar>
  11. IGL_INLINE bool igl::mosek_quadprog(
  12. const Index n,
  13. const std::vector<Index> & Qi,
  14. const std::vector<Index> & Qj,
  15. const std::vector<Scalar> & Qv,
  16. const std::vector<Scalar> & c,
  17. const Scalar cf,
  18. const Index m,
  19. const std::vector<Index> & Ari,
  20. const std::vector<Index> & Acp,
  21. const std::vector<Scalar> & Av,
  22. const std::vector<Scalar> & lc,
  23. const std::vector<Scalar> & uc,
  24. const std::vector<Scalar> & lx,
  25. const std::vector<Scalar> & ux,
  26. std::vector<Scalar> & x)
  27. {
  28. // I J V vectors of Q should all be same length
  29. assert(Qv.size() == Qi.size());
  30. assert(Qv.size() == Qj.size());
  31. // number of columns in linear constraint matrix must be ≤ number of
  32. // variables
  33. assert(Acp.size() == (n+1));
  34. // linear bound vectors must be size of number of constraints or empty
  35. assert( (lc.size() == m) || (lc.size() == 0));
  36. assert( (uc.size() == m) || (uc.size() == 0));
  37. // constant bound vectors must be size of number of variables or empty
  38. assert( (lx.size() == n) || (lx.size() == 0));
  39. assert( (ux.size() == n) || (ux.size() == 0));
  40. // allocate space for solution in x
  41. x.resize(n);
  42. // variables for mosek task, env and result code
  43. MSKenv_t env;
  44. MSKtask_t task;
  45. // Create the MOSEK environment
  46. mosek_guarded(MSK_makeenv(&env,NULL,NULL,NULL,NULL));
  47. /* Directs the log stream to the 'printstr' function. */
  48. mosek_guarded(MSK_linkfunctoenvstream(env,MSK_STREAM_LOG,NULL,printstr));
  49. // initialize mosek environment
  50. mosek_guarded(MSK_initenv(env));
  51. // Create the optimization task
  52. mosek_guarded(MSK_maketask(env,m,n,&task));
  53. verbose("Creating task with %ld linear constraints and %ld variables...\n",m,n);
  54. // Tell mosek how to print to std out
  55. mosek_guarded(MSK_linkfunctotaskstream(task,MSK_STREAM_LOG,NULL,printstr));
  56. // Give estimate of number of variables
  57. mosek_guarded(MSK_putmaxnumvar(task,n));
  58. if(m>0)
  59. {
  60. // Give estimate of number of constraints
  61. mosek_guarded(MSK_putmaxnumcon(task,m));
  62. // Give estimate of number of non zeros in A
  63. mosek_guarded(MSK_putmaxnumanz(task,Av.size()));
  64. }
  65. // Give estimate of number of non zeros in Q
  66. mosek_guarded(MSK_putmaxnumqnz(task,Qv.size()));
  67. if(m>0)
  68. {
  69. // Append 'm' empty constraints, the constrainst will initially have no
  70. // bounds
  71. mosek_guarded(MSK_append(task,MSK_ACC_CON,m));
  72. }
  73. // Append 'n' variables
  74. mosek_guarded(MSK_append(task,MSK_ACC_VAR,n));
  75. // add a contant term to the objective
  76. mosek_guarded(MSK_putcfix(task,cf));
  77. // loop over variables
  78. for(int j = 0;j<n;j++)
  79. {
  80. if(c.size() > 0)
  81. {
  82. // Set linear term c_j in the objective
  83. mosek_guarded(MSK_putcj(task,j,c[j]));
  84. }
  85. // Set constant bounds on variable j
  86. if(lx[j] == ux[j])
  87. {
  88. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_FX,lx[j],ux[j]));
  89. }else
  90. {
  91. mosek_guarded(MSK_putbound(task,MSK_ACC_VAR,j,MSK_BK_RA,lx[j],ux[j]));
  92. }
  93. if(m>0)
  94. {
  95. // Input column j of A
  96. mosek_guarded(
  97. MSK_putavec(
  98. task,
  99. MSK_ACC_VAR,
  100. j,
  101. Acp[j+1]-Acp[j],
  102. &Ari[Acp[j]],
  103. &Av[Acp[j]])
  104. );
  105. }
  106. }
  107. // loop over constraints
  108. for(size_t i = 0;i<m;i++)
  109. {
  110. // put bounds on constraints
  111. mosek_guarded(MSK_putbound(task,MSK_ACC_CON,i,MSK_BK_RA,lc[i],uc[i]));
  112. }
  113. // Input Q for the objective (REMEMBER Q SHOULD ONLY BE LOWER TRIANGLE)
  114. mosek_guarded(MSK_putqobj(task,Qv.size(),&Qi[0],&Qj[0],&Qv[0]));
  115. //// Set up task parameters
  116. // set tolerance
  117. //mosek_guarded(
  118. // MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_DFEAS,1e-8));
  119. //mosek_guarded(
  120. // MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_DSAFE,1.0));
  121. // 1.0 means optimizer is very leniant about declaring model infeasible
  122. //mosek_guarded(
  123. // MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_INFEAS,1e-8));
  124. //mosek_guarded(
  125. // MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_PATH,1e-8));
  126. //mosek_guarded(
  127. // MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_PFEAS,1e-8));
  128. // Hard to say if this is doing anything, probably nothing dramatic
  129. mosek_guarded(MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_PSAFE,1e2));
  130. // >1e0 NONSOLUTION
  131. // 1e-1 artifacts in deformation
  132. // 1e-3 artifacts in isolines
  133. // 1e-4 seems safe
  134. // 1e-8 MOSEK DEFAULT SOLUTION
  135. mosek_guarded(MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_REL_GAP,1e-8));
  136. //mosek_guarded(MSK_putdouparam(task,MSK_DPAR_INTPNT_TOL_REL_STEP,0.9999));
  137. // Turn off presolving
  138. //mosek_guarded(
  139. // MSK_putintparam(task,MSK_IPAR_PRESOLVE_USE,MSK_PRESOLVE_MODE_OFF));
  140. // Force particular matrix reordering method
  141. // MSK_ORDER_METHOD_NONE cuts time in half roughly, since half the time is
  142. // usually spent reordering the matrix
  143. // !! WARNING Setting this parameter to anything but MSK_ORDER_METHOD_FREE
  144. // seems to have the effect of setting it to MSK_ORDER_METHOD_NONE
  145. // *Or maybe Mosek is spending a bunch of time analyzing the matrix to
  146. // choose the right ordering method when really any of them are
  147. // instantaneous
  148. mosek_guarded(
  149. MSK_putintparam(task,MSK_IPAR_INTPNT_ORDER_METHOD,MSK_ORDER_METHOD_NONE));
  150. // Turn off convexity check
  151. mosek_guarded(
  152. MSK_putintparam(task,MSK_IPAR_CHECK_CONVEXITY,MSK_CHECK_CONVEXITY_NONE));
  153. // Force using multiple threads, not sure if MOSEK is properly destorying
  154. //extra threads...
  155. mosek_guarded(MSK_putintparam(task,MSK_IPAR_INTPNT_NUM_THREADS,6));
  156. // Force turn off data check
  157. mosek_guarded(MSK_putintparam(task,MSK_IPAR_DATA_CHECK,MSK_OFF));
  158. // Now the optimizer has been prepared
  159. MSKrescodee trmcode;
  160. // run the optimizer
  161. mosek_guarded(MSK_optimizetrm(task,&trmcode));
  162. // Print a summary containing information about the solution for debugging
  163. // purposes
  164. MSK_solutionsummary(task,MSK_STREAM_LOG);
  165. // Get status of solution
  166. MSKsolstae solsta;
  167. MSK_getsolutionstatus(task,MSK_SOL_ITR,NULL,&solsta);
  168. bool success = false;
  169. switch(solsta)
  170. {
  171. case MSK_SOL_STA_OPTIMAL:
  172. case MSK_SOL_STA_NEAR_OPTIMAL:
  173. MSK_getsolutionslice(task,MSK_SOL_ITR,MSK_SOL_ITEM_XX,0,n,&x[0]);
  174. printf("Optimal primal solution\n");
  175. //for(size_t j=0; j<n; ++j)
  176. //{
  177. // printf("x[%ld]: %g\n",j,x[j]);
  178. //}
  179. success = true;
  180. break;
  181. case MSK_SOL_STA_DUAL_INFEAS_CER:
  182. case MSK_SOL_STA_PRIM_INFEAS_CER:
  183. case MSK_SOL_STA_NEAR_DUAL_INFEAS_CER:
  184. case MSK_SOL_STA_NEAR_PRIM_INFEAS_CER:
  185. printf("Primal or dual infeasibility certificate found.\n");
  186. break;
  187. case MSK_SOL_STA_UNKNOWN:
  188. printf("The status of the solution could not be determined.\n");
  189. break;
  190. default:
  191. printf("Other solution status.");
  192. break;
  193. }
  194. MSK_deletetask(&task);
  195. MSK_deleteenv(&env);
  196. return success;
  197. }
  198. #ifndef IGL_HEADER_ONLY
  199. // Explicit template specialization
  200. #endif