mosek_quadprog.h 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 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. #ifndef IGL_MOSEK_MOSEK_QUADPROG_H
  9. #define IGL_MOSEK_MOSEK_QUADPROG_H
  10. #include "../igl_inline.h"
  11. #include <vector>
  12. #include <map>
  13. #include <mosek.h>
  14. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  15. #include <Eigen/Dense>
  16. #include <Eigen/Sparse>
  17. namespace igl
  18. {
  19. namespace mosek
  20. {
  21. struct MosekData
  22. {
  23. // Integer parameters
  24. std::map<MSKiparame,int> intparam;
  25. // Double parameters
  26. std::map<MSKdparame,double> douparam;
  27. // Default values
  28. IGL_INLINE MosekData();
  29. };
  30. // Solve a convex quadratic optimization problem with linear and constant
  31. // bounds, that is:
  32. //
  33. // Minimize: ½ * xT * Q⁰ * x + cT * x + cf
  34. //
  35. // Subject to: lc ≤ Ax ≤ uc
  36. // lx ≤ x ≤ ux
  37. //
  38. // where we are trying to find the optimal vector of values x.
  39. //
  40. // Note: Q⁰ must be symmetric and the ½ is a convention of MOSEK
  41. //
  42. // Note: Because of how MOSEK accepts different parts of the system, Q
  43. // should be stored in IJV (aka Coordinate) format and should only include
  44. // entries in the lower triangle. A should be stored in Column compressed
  45. // (aka Harwell Boeing) format. As described:
  46. // http://netlib.org/linalg/html_templates/node92.html
  47. // or
  48. // http://en.wikipedia.org/wiki/Sparse_matrix
  49. // #Compressed_sparse_column_.28CSC_or_CCS.29
  50. //
  51. //
  52. // Templates:
  53. // Index type for index variables
  54. // Scalar type for floating point variables (gets cast to double?)
  55. // Input:
  56. // n number of variables, i.e. size of x
  57. // Qi vector of qnnz row indices of non-zeros in LOWER TRIANGLE ONLY of
  58. // Q⁰
  59. // Qj vector of qnnz column indices of non-zeros in LOWER TRIANGLE ONLY
  60. // of Q⁰
  61. // Qv vector of qnnz values of non-zeros in LOWER TRIANGLE ONLY of Q⁰,
  62. // such that:
  63. //
  64. // Q⁰(Qi[k],Qj[k]) = Qv[k] for k ∈ [0,Qnnz-1], where Qnnz is the
  65. //
  66. // number of non-zeros in Q⁰
  67. // c (optional) vector of n values of c, transpose of coefficient row
  68. // vector of linear terms, EMPTY means c == 0
  69. // cf (ignored) value of constant term in objective, 0 means cf == 0, so
  70. // optional only in the sense that it is mandatory
  71. // m number of constraints, therefore also number of rows in linear
  72. // constraint coefficient matrix A, and in linear constraint bound
  73. // vectors lc and uc
  74. // Av vector of non-zero values of A, in column compressed order
  75. // Ari vector of row indices corresponding to non-zero values of A,
  76. // Acp vector of indices into Ari and Av of the first entry for each
  77. // column of A, size(Acp) = (# columns of A) + 1 = n + 1
  78. // lc vector of m linear constraint lower bounds
  79. // uc vector of m linear constraint upper bounds
  80. // lx vector of n constant lower bounds
  81. // ux vector of n constant upper bounds
  82. // Output:
  83. // x vector of size n to hold output of optimization
  84. // Return:
  85. // true only if optimization was successful with no errors
  86. //
  87. // Note: All indices are 0-based
  88. //
  89. template <typename Index, typename Scalar>
  90. IGL_INLINE bool mosek_quadprog(
  91. const Index n,
  92. /* mosek won't allow this to be const*/ std::vector<Index> & Qi,
  93. /* mosek won't allow this to be const*/ std::vector<Index> & Qj,
  94. /* mosek won't allow this to be const*/ std::vector<Scalar> & Qv,
  95. const std::vector<Scalar> & c,
  96. const Scalar cf,
  97. const Index m,
  98. /* mosek won't allow this to be const*/ std::vector<Scalar> & Av,
  99. /* mosek won't allow this to be const*/ std::vector<Index> & Ari,
  100. const std::vector<Index> & Acp,
  101. const std::vector<Scalar> & lc,
  102. const std::vector<Scalar> & uc,
  103. const std::vector<Scalar> & lx,
  104. const std::vector<Scalar> & ux,
  105. MosekData & mosek_data,
  106. std::vector<Scalar> & x);
  107. // Wrapper with Eigen elements
  108. //
  109. // Inputs:
  110. // Q n by n square quadratic coefficients matrix **only lower triangle
  111. // is used**.
  112. // c n-long vector of linear coefficients
  113. // cf constant coefficient
  114. // A m by n square linear coefficienst matrix of inequality constraints
  115. // lc m-long vector of lower bounds for linear inequality constraints
  116. // uc m-long vector of upper bounds for linear inequality constraints
  117. // lx n-long vector of lower bounds
  118. // ux n-long vector of upper bounds
  119. // mosek_data parameters struct
  120. // Outputs:
  121. // x n-long solution vector
  122. // Returns true only if optimization finishes without error
  123. //
  124. IGL_INLINE bool mosek_quadprog(
  125. const Eigen::SparseMatrix<double> & Q,
  126. const Eigen::VectorXd & c,
  127. const double cf,
  128. const Eigen::SparseMatrix<double> & A,
  129. const Eigen::VectorXd & lc,
  130. const Eigen::VectorXd & uc,
  131. const Eigen::VectorXd & lx,
  132. const Eigen::VectorXd & ux,
  133. MosekData & mosek_data,
  134. Eigen::VectorXd & x);
  135. }
  136. }
  137. #ifndef IGL_STATIC_LIBRARY
  138. # include "mosek_quadprog.cpp"
  139. #endif
  140. #endif