libigl v2.5.0
Loading...
Searching...
No Matches
igl::mosek Namespace Reference

Classes

struct  MosekData
 Structure for holding MOSEK data for solving a quadratic program. More...
 

Functions

template<typename DerivedV , typename DerivedEle , typename Derivedb , typename Derivedbc , typename DerivedW >
bool bbw (const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedEle > &Ele, const Eigen::PlainObjectBase< Derivedb > &b, const Eigen::PlainObjectBase< Derivedbc > &bc, igl::BBWData &data, igl::mosek::MosekData &mosek_data, Eigen::PlainObjectBase< DerivedW > &W)
 Compute Bounded Biharmonic Weights on a given domain (V,Ele) with a given set of boundary conditions.
 
MSKrescodee mosek_guarded (const MSKrescodee r)
 Little function to wrap around mosek call to handle errors.
 
bool mosek_linprog (const Eigen::VectorXd &c, const Eigen::SparseMatrix< double > &A, const Eigen::VectorXd &lc, const Eigen::VectorXd &uc, const Eigen::VectorXd &lx, const Eigen::VectorXd &ux, Eigen::VectorXd &x)
 Solve a linear program using mosek.
 
bool mosek_linprog (const Eigen::VectorXd &c, const Eigen::SparseMatrix< double > &A, const Eigen::VectorXd &lc, const Eigen::VectorXd &uc, const Eigen::VectorXd &lx, const Eigen::VectorXd &ux, const MSKenv_t &env, Eigen::VectorXd &x)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename Index , typename Scalar >
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)
 
bool mosek_quadprog (const Eigen::SparseMatrix< double > &Q, const Eigen::VectorXd &c, const double cf, const Eigen::SparseMatrix< double > &A, const Eigen::VectorXd &lc, const Eigen::VectorXd &uc, const Eigen::VectorXd &lx, const Eigen::VectorXd &ux, MosekData &mosek_data, Eigen::VectorXd &x)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 

Function Documentation

◆ bbw()

template<typename DerivedV , typename DerivedEle , typename Derivedb , typename Derivedbc , typename DerivedW >
bool igl::mosek::bbw ( const Eigen::PlainObjectBase< DerivedV > &  V,
const Eigen::PlainObjectBase< DerivedEle > &  Ele,
const Eigen::PlainObjectBase< Derivedb > &  b,
const Eigen::PlainObjectBase< Derivedbc > &  bc,
igl::BBWData data,
igl::mosek::MosekData mosek_data,
Eigen::PlainObjectBase< DerivedW > &  W 
)

Compute Bounded Biharmonic Weights on a given domain (V,Ele) with a given set of boundary conditions.

Template Parameters
DerivedVderived type of eigen matrix for V (e.g. MatrixXd)
DerivedFderived type of eigen matrix for F (e.g. MatrixXi)
Derivedbderived type of eigen matrix for b (e.g. VectorXi)
Derivedbcderived type of eigen matrix for bc (e.g. MatrixXd)
DerivedWderived type of eigen matrix for W (e.g. MatrixXd)
Parameters
[in]V#V by dim vertex positions
[in]Ele#Elements by simplex-size list of element indices
[in]b#b boundary indices into V
[in]bc#b by #W list of boundary values
[in]dataobject containing options, initial guess --> solution and results
[in]mosek_dataobject containing mosek options
[out]W#V by #W list of unnormalized weights to normalize use igl::normalize_row_sums(W,W);
Returns
true on success, false on failure

◆ mosek_guarded()

MSKrescodee igl::mosek::mosek_guarded ( const MSKrescodee  r)

Little function to wrap around mosek call to handle errors.

Parameters
[in]rmosek error code returned from mosek call
Returns
r untouched

◆ mosek_linprog() [1/2]

bool igl::mosek::mosek_linprog ( const Eigen::VectorXd &  c,
const Eigen::SparseMatrix< double > &  A,
const Eigen::VectorXd &  lc,
const Eigen::VectorXd &  uc,
const Eigen::VectorXd &  lx,
const Eigen::VectorXd &  ux,
Eigen::VectorXd &  x 
)

Solve a linear program using mosek.

Given in the form:

min c'x
s.t. lc <= A x <= uc
     lx <= x <= ux
Parameters
[in]c#x list of linear objective coefficients
[in]A#A by #x matrix of linear inequality constraint coefficients
[in]lc#A list of lower constraint bounds
[in]uc#A list of upper constraint bounds
[in]lx#x list of lower variable bounds
[in]ux#x list of upper variable bounds
[out]x#x list of solution values
Returns
true iff success.

◆ mosek_linprog() [2/2]

bool igl::mosek::mosek_linprog ( const Eigen::VectorXd &  c,
const Eigen::SparseMatrix< double > &  A,
const Eigen::VectorXd &  lc,
const Eigen::VectorXd &  uc,
const Eigen::VectorXd &  lx,
const Eigen::VectorXd &  ux,
const MSKenv_t &  env,
Eigen::VectorXd &  x 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Wrapper that keeps mosek environment alive (if licence checking is becoming a bottleneck)

Parameters
[in]envmosek environment

◆ mosek_quadprog() [1/2]

template<typename Index , typename Scalar >
bool igl::mosek::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 
)

◆ mosek_quadprog() [2/2]

bool igl::mosek::mosek_quadprog ( const Eigen::SparseMatrix< double > &  Q,
const Eigen::VectorXd &  c,
const double  cf,
const Eigen::SparseMatrix< double > &  A,
const Eigen::VectorXd &  lc,
const Eigen::VectorXd &  uc,
const Eigen::VectorXd &  lx,
const Eigen::VectorXd &  ux,
MosekData mosek_data,
Eigen::VectorXd &  x 
)

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Parameters
[in]Qn by n square quadratic coefficients matrix only lower triangle is used.
[in]cn-long vector of linear coefficients
[in]cfconstant coefficient
[in]Am by n square linear coefficienst matrix of inequality constraints
[in]lcm-long vector of lower bounds for linear inequality constraints
[in]ucm-long vector of upper bounds for linear inequality constraints
[in]lxn-long vector of lower bounds
[in]uxn-long vector of upper bounds
[in]mosek_dataparameters struct
[out]xn-long solution vector
Returns
true only if optimization finishes without error