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