mosek_linprog.cpp 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 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_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. using namespace std;
  47. // number of constraints
  48. const int m = A.rows();
  49. // number of variables
  50. const int n = A.cols();
  51. vector<double> vAv;
  52. vector<int> vAri,vAcp;
  53. int nr;
  54. harwell_boeing(A,nr,vAv,vAri,vAcp);
  55. MSKtask_t task;
  56. // Create the optimization task
  57. mosek_guarded(MSK_maketask(env,m,n,&task));
  58. // no threads
  59. mosek_guarded(MSK_putintparam(task,MSK_IPAR_NUM_THREADS,1));
  60. if(m>0)
  61. {
  62. // Append 'm' empty constraints, the constrainst will initially have no
  63. // bounds
  64. mosek_guarded(MSK_appendcons(task,m));
  65. }
  66. mosek_guarded(MSK_appendvars(task,n));
  67. const auto & key = [](const double lxj, const double uxj) ->
  68. MSKboundkeye
  69. {
  70. MSKboundkeye k = MSK_BK_FR;
  71. if(isfinite(lxj) && isfinite(uxj))
  72. {
  73. if(lxj == uxj)
  74. {
  75. k = MSK_BK_FX;
  76. }else{
  77. k = MSK_BK_RA;
  78. }
  79. }else if(isfinite(lxj))
  80. {
  81. k = MSK_BK_LO;
  82. }else if(isfinite(uxj))
  83. {
  84. k = MSK_BK_UP;
  85. }
  86. return k;
  87. };
  88. // loop over variables
  89. for(int j = 0;j<n;j++)
  90. {
  91. if(c.size() > 0)
  92. {
  93. // Set linear term c_j in the objective
  94. mosek_guarded(MSK_putcj(task,j,c(j)));
  95. }
  96. // Set constant bounds on variable j
  97. const double lxj = lx.size()>0?lx[j]:-numeric_limits<double>::infinity();
  98. const double uxj = ux.size()>0?ux[j]: numeric_limits<double>::infinity();
  99. mosek_guarded(MSK_putvarbound(task,j,key(lxj,uxj),lxj,uxj));
  100. if(m>0)
  101. {
  102. // Input column j of A
  103. mosek_guarded(
  104. MSK_putacol(
  105. task,
  106. j,
  107. vAcp[j+1]-vAcp[j],
  108. &vAri[vAcp[j]],
  109. &vAv[vAcp[j]])
  110. );
  111. }
  112. }
  113. // loop over constraints
  114. for(int i = 0;i<m;i++)
  115. {
  116. // Set constraint bounds for row i
  117. const double lci = lc.size()>0?lc[i]:-numeric_limits<double>::infinity();
  118. const double uci = uc.size()>0?uc[i]: numeric_limits<double>::infinity();
  119. mosek_guarded(MSK_putconbound(task,i,key(lci,uci),lci,uci));
  120. }
  121. // Now the optimizer has been prepared
  122. MSKrescodee trmcode;
  123. // run the optimizer
  124. mosek_guarded(MSK_optimizetrm(task,&trmcode));
  125. // Get status
  126. MSKsolstae solsta;
  127. MSK_getsolsta (task,MSK_SOL_ITR,&solsta);
  128. bool success = false;
  129. switch(solsta)
  130. {
  131. case MSK_SOL_STA_OPTIMAL:
  132. case MSK_SOL_STA_NEAR_OPTIMAL:
  133. x.resize(n);
  134. /* Request the basic solution. */
  135. MSK_getxx(task,MSK_SOL_BAS,x.data());
  136. success = true;
  137. break;
  138. case MSK_SOL_STA_DUAL_INFEAS_CER:
  139. case MSK_SOL_STA_PRIM_INFEAS_CER:
  140. case MSK_SOL_STA_NEAR_DUAL_INFEAS_CER:
  141. case MSK_SOL_STA_NEAR_PRIM_INFEAS_CER:
  142. //printf("Primal or dual infeasibility certificate found.\n");
  143. break;
  144. case MSK_SOL_STA_UNKNOWN:
  145. //printf("The status of the solution could not be determined.\n");
  146. break;
  147. default:
  148. //printf("Other solution status.");
  149. break;
  150. }
  151. MSK_deletetask(&task);
  152. return success;
  153. }