active_set.h 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120
  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_ACTIVE_SET_H
  9. #define IGL_ACTIVE_SET_H
  10. #include "igl_inline.h"
  11. #include "SolverStatus.h"
  12. #include <Eigen/Core>
  13. #include <Eigen/Sparse>
  14. namespace igl
  15. {
  16. struct active_set_params;
  17. ///
  18. /// Minimize convex quadratic energy subject to linear inequality constraints
  19. ///
  20. /// min ½ Zᵀ A Z + Zᵀ B + constant
  21. /// Z
  22. /// subject to
  23. /// Aeq Z = Beq
  24. /// Aieq Z <= Bieq
  25. /// lx <= Z <= ux
  26. /// Z(known) = Y
  27. ///
  28. /// that Z(known) = Y, optionally also subject to the constraints Aeq*Z = Beq,
  29. /// and further optionally subject to the linear inequality constraints that
  30. /// Aieq*Z <= Bieq and constant inequality constraints lx <= x <= ux
  31. ///
  32. /// @param[in] A n by n matrix of quadratic coefficients
  33. /// @param[in] B n by 1 column of linear coefficients
  34. /// @param[in] known list of indices to known rows in Z
  35. /// @param[in] Y list of fixed values corresponding to known rows in Z
  36. /// @param[in] Aeq meq by n list of linear equality constraint coefficients
  37. /// @param[in] Beq meq by 1 list of linear equality constraint constant values
  38. /// @param[in] Aieq mieq by n list of linear inequality constraint coefficients
  39. /// @param[in] Bieq mieq by 1 list of linear inequality constraint constant values
  40. /// @param[in] lx n by 1 list of lower bounds [] implies -Inf
  41. /// @param[in] ux n by 1 list of upper bounds [] implies Inf
  42. /// @param[in] params struct of additional parameters (see below)
  43. /// @param[in,out] Z if not empty, is taken to be an n by 1 list of initial guess values. Set to solution on output.
  44. /// @return true on success, false on error
  45. ///
  46. /// \note Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  47. /// secs, igl/min_quad_with_fixed.h 7.1 secs
  48. ///
  49. /// \pre rows of [Aeq;Aieq] **must** be linearly independent. Should be
  50. /// using QR decomposition otherwise:
  51. /// https://v8doc.sas.com/sashtml/ormp/chap5/sect32.htm
  52. ///
  53. /// \warning This solver is fairly experimental. It works reasonably well for
  54. /// bbw problems but doesn't generalize well to other problems. NASOQ and
  55. /// OSQP are better general purpose solvers.
  56. template <
  57. typename AT,
  58. typename DerivedB,
  59. typename Derivedknown,
  60. typename DerivedY,
  61. typename AeqT,
  62. typename DerivedBeq,
  63. typename AieqT,
  64. typename DerivedBieq,
  65. typename Derivedlx,
  66. typename Derivedux,
  67. typename DerivedZ
  68. >
  69. IGL_INLINE igl::SolverStatus active_set(
  70. const Eigen::SparseMatrix<AT>& A,
  71. const Eigen::MatrixBase<DerivedB> & B,
  72. const Eigen::MatrixBase<Derivedknown> & known,
  73. const Eigen::MatrixBase<DerivedY> & Y,
  74. const Eigen::SparseMatrix<AeqT>& Aeq,
  75. const Eigen::MatrixBase<DerivedBeq> & Beq,
  76. const Eigen::SparseMatrix<AieqT>& Aieq,
  77. const Eigen::MatrixBase<DerivedBieq> & Bieq,
  78. const Eigen::MatrixBase<Derivedlx> & lx,
  79. const Eigen::MatrixBase<Derivedux> & ux,
  80. const igl::active_set_params & params,
  81. Eigen::PlainObjectBase<DerivedZ> & Z
  82. );
  83. };
  84. #include "EPS.h"
  85. /// Input parameters controling active_set
  86. ///
  87. /// \fileinfo
  88. struct igl::active_set_params
  89. {
  90. /// Auu_pd whether Auu is positive definite {false}
  91. bool Auu_pd;
  92. /// max_iter Maximum number of iterations (0 = Infinity, {100})
  93. int max_iter;
  94. /// inactive_threshold Threshold on Lagrange multiplier values to determine
  95. /// whether to keep constraints active {EPS}
  96. double inactive_threshold;
  97. /// constraint_threshold Threshold on whether constraints are violated (0
  98. /// is perfect) {EPS}
  99. double constraint_threshold;
  100. /// solution_diff_threshold Threshold on the squared norm of the difference
  101. /// between two consecutive solutions {EPS}
  102. double solution_diff_threshold;
  103. /// @private
  104. active_set_params():
  105. Auu_pd(false),
  106. max_iter(100),
  107. inactive_threshold(igl::DOUBLE_EPS),
  108. constraint_threshold(igl::DOUBLE_EPS),
  109. solution_diff_threshold(igl::DOUBLE_EPS)
  110. {};
  111. };
  112. #ifndef IGL_STATIC_LIBRARY
  113. # include "active_set.cpp"
  114. #endif
  115. #endif