mosek_quadprog.h 5.6 KB

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