ProxSDP is an open-source semidefinite programming (SDP) solver based on the paper "Exploiting Low-Rank Structure in Semidefinite Programming by Approximate Operator Splitting".
The main advantage of ProxSDP over other state-of-the-art solvers is the ability to exploit the low-rank structure inherent to several SDP problems.
- General conic convex optimization problems with the presence of the positive semidefinite cone, second-order cone and positive orthant;
- Semidefinite relaxation of nonconvex problems, for example, max-cut, binary MIMO, optimal power flow, sensor localization, sum-of-squares;
- Control theory problems with LMI constraints;
- Nuclear norm minimization problems, for example, matrix completion;
ProxSDP is licensed under the MIT License.
Install ProxSDP using Julia package manager:
import Pkg
Pkg.add("ProxSDP")
The online documentation is available at https://mariohsouto.github.io/ProxSDP.jl/latest/.
Consider the semidefinite programming relaxation of the max-cut problem, as in:
max 0.25 * W•X
s.t. diag(X) = 1,
X ≽ 0,
This problem can be solved by the following code using ProxSDP and JuMP.
# Load packages
using ProxSDP, JuMP, LinearAlgebra
# Number of vertices
n = 4
# Graph weights
W = [
18.0 -5.0 -7.0 -6.0
-5.0 6.0 0.0 -1.0
-7.0 0.0 8.0 -1.0
-6.0 -1.0 -1.0 8.0
]
# Build Max-Cut SDP relaxation via JuMP
model = Model(ProxSDP.Optimizer)
set_optimizer_attribute(model, "log_verbose", true)
set_optimizer_attribute(model, "tol_gap", 1e-4)
set_optimizer_attribute(model, "tol_feasibility", 1e-4)
@variable(model, X[1:n, 1:n], PSD)
@objective(model, Max, 0.25 * LinearAlgebra.dot(W, X))
@constraint(model, LinearAlgebra.diag(X) .== 1)
# Solve optimization problem with ProxSDP
optimize!(model)
# Retrieve solution
Xsol = value.(X)
Another alternative is to use ProxSDP via Convex.jl:
# Load packages
using Convex, ProxSDP
# Number of vertices
n = 4
# Graph weights
W = [
18.0 -5.0 -7.0 -6.0
-5.0 6.0 0.0 -1.0
-7.0 0.0 8.0 -1.0
-6.0 -1.0 -1.0 8.0
]
# Define optimization problem
X = Semidefinite(n)
problem = maximize(0.25 * dot(W, X), diag(X) == 1)
# Solve optimization problem with ProxSDP
solve!(problem, ProxSDP.Optimizer(log_verbose=true, tol_gap=1e-4, tol_feasibility=1e-4))
# Get the objective value
problem.optval
# Retrieve solution
evaluate(X)
The published version of the paper can be found here and the arXiv version here.
We kindly request that you cite the paper as:
@article{souto2020exploiting,
author = {Mario Souto and Joaquim D. Garcia and \'Alvaro Veiga},
title = {Exploiting low-rank structure in semidefinite programming by approximate operator splitting},
journal = {Optimization},
pages = {1-28},
year = {2020},
publisher = {Taylor & Francis},
doi = {10.1080/02331934.2020.1823387},
URL = {https://doi.org/10.1080/02331934.2020.1823387}
}
The preprint version of the paper can be found here.
- ProxSDP is a research software, therefore it should not be used in production.
- Please open an issue if you find any problems, developers will try to fix and find alternatives.
- There is no continuous development for 32-bit systems, the package should work, but might reach some issues.
- ProxSDP assumes primal and dual feasibility.
- Support for exponential and power cones
- Warm start