mosek_linprog.cpp 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <[email protected]>
  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_linprog.h"
  9. #include "../mosek/mosek_guarded.h"
  10. #include "../harwell_boeing.h"
  11. #include <limits>
  12. #include <cmath>
  13. #include <vector>
  14. IGL_INLINE bool igl::mosek::mosek_linprog(
  15. const Eigen::VectorXd & c,
  16. const Eigen::SparseMatrix<double> & A,
  17. const Eigen::VectorXd & lc,
  18. const Eigen::VectorXd & uc,
  19. const Eigen::VectorXd & lx,
  20. const Eigen::VectorXd & ux,
  21. Eigen::VectorXd & x)
  22. {
  23. // variables for mosek task, env and result code
  24. MSKenv_t env;
  25. // Create the MOSEK environment
  26. mosek_guarded(MSK_makeenv(&env,NULL));
  27. // initialize mosek environment
  28. #if MSK_VERSION_MAJOR <= 7
  29. mosek_guarded(MSK_initenv(env));
  30. #endif
  31. const bool ret = mosek_linprog(c,A,lc,uc,lx,ux,env,x);
  32. MSK_deleteenv(&env);
  33. return ret;
  34. }
  35. IGL_INLINE bool igl::mosek::mosek_linprog(
  36. const Eigen::VectorXd & c,
  37. const Eigen::SparseMatrix<double> & A,
  38. const Eigen::VectorXd & lc,
  39. const Eigen::VectorXd & uc,
  40. const Eigen::VectorXd & lx,
  41. const Eigen::VectorXd & ux,
  42. const MSKenv_t & env,
  43. Eigen::VectorXd & x)
  44. {
  45. // following http://docs.mosek.com/7.1/capi/Linear_optimization.html
  46. // number of constraints
  47. const int m = A.rows();
  48. // number of variables
  49. const int n = A.cols();
  50. vector<double> vAv;
  51. vector<int> vAri,vAcp;
  52. int nr;
  53. harwell_boeing(A,nr,vAv,vAri,vAcp);
  54. MSKtask_t task;
  55. // Create the optimization task
  56. mosek_guarded(MSK_maketask(env,m,n,&task));
  57. // no threads
  58. mosek_guarded(MSK_putintparam(task,MSK_IPAR_NUM_THREADS,1));
  59. if(m>0)
  60. {
  61. // Append 'm' empty constraints, the constrainst will initially have no
  62. // bounds
  63. mosek_guarded(MSK_appendcons(task,m));
  64. }
  65. mosek_guarded(MSK_appendvars(task,n));
  66. const auto & key = [](const double lxj, const double uxj) ->
  67. MSKboundkeye
  68. {
  69. MSKboundkeye k = MSK_BK_FR;
  70. if(isfinite(lxj) && isfinite(uxj))
  71. {
  72. if(lxj == uxj)
  73. {
  74. k = MSK_BK_FX;
  75. }else{
  76. k = MSK_BK_RA;
  77. }
  78. }else if(isfinite(lxj))
  79. {
  80. k = MSK_BK_LO;
  81. }else if(isfinite(uxj))
  82. {
  83. k = MSK_BK_UP;
  84. }
  85. return k;
  86. };
  87. // loop over variables
  88. for(int j = 0;j<n;j++)
  89. {
  90. if(c.size() > 0)
  91. {
  92. // Set linear term c_j in the objective
  93. mosek_guarded(MSK_putcj(task,j,c(j)));
  94. }
  95. // Set constant bounds on variable j
  96. const double lxj = lx.size()>0?lx[j]:-numeric_limits<double>::infinity();
  97. const double uxj = ux.size()>0?ux[j]: numeric_limits<double>::infinity();
  98. mosek_guarded(MSK_putvarbound(task,j,key(lxj,uxj),lxj,uxj));
  99. if(m>0)
  100. {
  101. // Input column j of A
  102. mosek_guarded(
  103. MSK_putacol(
  104. task,
  105. j,
  106. vAcp[j+1]-vAcp[j],
  107. &vAri[vAcp[j]],
  108. &vAv[vAcp[j]])
  109. );
  110. }
  111. }
  112. // loop over constraints
  113. for(int i = 0;i<m;i++)
  114. {
  115. // Set constraint bounds for row i
  116. const double lci = lc.size()>0?lc[i]:-numeric_limits<double>::infinity();
  117. const double uci = uc.size()>0?uc[i]: numeric_limits<double>::infinity();
  118. mosek_guarded(MSK_putconbound(task,i,key(lci,uci),lci,uci));
  119. }
  120. // Now the optimizer has been prepared
  121. MSKrescodee trmcode;
  122. // run the optimizer
  123. mosek_guarded(MSK_optimizetrm(task,&trmcode));
  124. // Get status
  125. MSKsolstae solsta;
  126. MSK_getsolsta (task,MSK_SOL_ITR,&solsta);
  127. bool success = false;
  128. switch(solsta)
  129. {
  130. case MSK_SOL_STA_OPTIMAL:
  131. #if MSK_VERSION_MAJOR <= 8
  132. case MSK_SOL_STA_NEAR_OPTIMAL:
  133. #endif
  134. x.resize(n);
  135. /* Request the basic solution. */
  136. MSK_getxx(task,MSK_SOL_BAS,x.data());
  137. success = true;
  138. break;
  139. case MSK_SOL_STA_DUAL_INFEAS_CER:
  140. case MSK_SOL_STA_PRIM_INFEAS_CER:
  141. #if MSK_VERSION_MAJOR <= 8
  142. case MSK_SOL_STA_NEAR_DUAL_INFEAS_CER:
  143. case MSK_SOL_STA_NEAR_PRIM_INFEAS_CER:
  144. #endif
  145. //printf("Primal or dual infeasibility certificate found.\n");
  146. break;
  147. case MSK_SOL_STA_UNKNOWN:
  148. //printf("The status of the solution could not be determined.\n");
  149. break;
  150. default:
  151. //printf("Other solution status.");
  152. break;
  153. }
  154. MSK_deletetask(&task);
  155. return success;
  156. }