mosek_quadprog.h 4.6 KB

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