2016-12-19 15:44:52 +01:00
|
|
|
|
using JuMP
|
2017-01-09 00:59:40 +01:00
|
|
|
|
import Base: rationalize
|
2017-01-14 15:24:16 +01:00
|
|
|
|
using GroupAlgebras
|
2016-12-19 15:44:52 +01:00
|
|
|
|
|
2017-02-11 13:31:01 +01:00
|
|
|
|
function products{T}(U::AbstractVector{T}, V::AbstractVector{T})
|
|
|
|
|
result = Vector{T}()
|
|
|
|
|
for u in U
|
|
|
|
|
for v in V
|
|
|
|
|
push!(result, u*v)
|
2016-12-21 10:00:22 +01:00
|
|
|
|
end
|
|
|
|
|
end
|
2017-01-18 17:51:58 +01:00
|
|
|
|
return unique(result)
|
2016-12-21 10:00:22 +01:00
|
|
|
|
end
|
2016-12-19 15:44:52 +01:00
|
|
|
|
|
2017-02-11 13:33:35 +01:00
|
|
|
|
function create_product_matrix(basis, limit)
|
2017-01-13 18:04:20 +01:00
|
|
|
|
product_matrix = zeros(Int, (limit,limit))
|
2016-12-21 10:00:22 +01:00
|
|
|
|
for i in 1:limit
|
2017-02-11 13:33:35 +01:00
|
|
|
|
x_inv::eltype(basis) = inv(basis[i])
|
2017-03-06 11:57:42 +01:00
|
|
|
|
Threads.@threads for j in 1:limit
|
2017-01-13 18:02:34 +01:00
|
|
|
|
w = x_inv*basis[j]
|
|
|
|
|
index = findfirst(basis, w)
|
2017-02-11 13:33:35 +01:00
|
|
|
|
index ≠ 0 || throw(ArgumentError("Product is not supported on basis: $w"))
|
|
|
|
|
product_matrix[i,j] = index
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
|
|
|
|
end
|
2017-02-11 13:33:35 +01:00
|
|
|
|
return product_matrix
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
|
|
|
|
|
2017-02-11 13:30:17 +01:00
|
|
|
|
function constraints_from_pm(pm, total_length=maximum(pm))
|
|
|
|
|
n = size(pm,1)
|
|
|
|
|
constraints = constraints = [Array{Int,1}[] for x in 1:total_length]
|
|
|
|
|
for j in 1:n
|
|
|
|
|
Threads.@threads for i in 1:n
|
|
|
|
|
idx = pm[i,j]
|
|
|
|
|
push!(constraints[idx], [i,j])
|
|
|
|
|
end
|
|
|
|
|
end
|
|
|
|
|
return constraints
|
|
|
|
|
end
|
2016-12-21 16:02:03 +01:00
|
|
|
|
|
2017-02-11 13:34:28 +01:00
|
|
|
|
function splaplacian_coeff(S, basis, n=length(basis))
|
|
|
|
|
result = spzeros(n)
|
2016-12-21 16:02:03 +01:00
|
|
|
|
result[1] = length(S)
|
|
|
|
|
for s in S
|
2017-02-11 13:34:28 +01:00
|
|
|
|
ind = findfirst(basis, s)
|
2016-12-21 16:02:03 +01:00
|
|
|
|
result[ind] += -1
|
|
|
|
|
end
|
|
|
|
|
return result
|
|
|
|
|
end
|
|
|
|
|
|
2017-02-11 13:34:28 +01:00
|
|
|
|
function laplacian_coeff(S, basis)
|
|
|
|
|
return full(splaplacian_coeff(S,basis))
|
2016-12-21 16:02:03 +01:00
|
|
|
|
end
|
|
|
|
|
|
2017-01-09 01:01:31 +01:00
|
|
|
|
|
|
|
|
|
function create_SDP_problem(matrix_constraints, Δ::GroupAlgebraElement)
|
2016-12-21 16:03:19 +01:00
|
|
|
|
N = size(Δ.product_matrix,1)
|
2017-01-09 01:01:31 +01:00
|
|
|
|
const Δ² = Δ*Δ
|
2016-12-19 15:44:52 +01:00
|
|
|
|
@assert length(Δ) == length(matrix_constraints)
|
2017-02-11 13:38:02 +01:00
|
|
|
|
m = JuMP.Model();
|
|
|
|
|
JuMP.@variable(m, A[1:N, 1:N], SDP)
|
|
|
|
|
JuMP.@SDconstraint(m, A >= zeros(size(A)))
|
|
|
|
|
JuMP.@variable(m, κ >= 0.0)
|
2017-03-06 11:53:25 +01:00
|
|
|
|
JuMP.@constraint(m, κ <= 0.26)
|
2017-02-11 13:38:02 +01:00
|
|
|
|
JuMP.@objective(m, Max, κ)
|
2017-03-06 11:53:25 +01:00
|
|
|
|
JuMP.@constraint(m, sum(A[i] for i in eachindex(A)) == 0)
|
2016-12-19 15:44:52 +01:00
|
|
|
|
|
2016-12-23 00:51:06 +01:00
|
|
|
|
for (pairs, δ², δ) in zip(matrix_constraints, Δ².coefficients, Δ.coefficients)
|
2017-02-11 13:38:02 +01:00
|
|
|
|
JuMP.@constraint(m, sum(A[i,j] for (i,j) in pairs) == δ² - κ*δ)
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
|
|
|
|
return m
|
|
|
|
|
end
|
|
|
|
|
|
2017-02-11 13:41:03 +01:00
|
|
|
|
function solve_SDP(sdp_constraints, Δ, solver; verbose=true)
|
|
|
|
|
SDP_problem = create_SDP_problem(sdp_constraints, Δ);
|
|
|
|
|
verbose && @show solver
|
|
|
|
|
|
|
|
|
|
JuMP.setsolver(SDP_problem, solver);
|
|
|
|
|
verbose && @show SDP_problem
|
|
|
|
|
# @time MathProgBase.writeproblem(SDP_problem, "/tmp/SDP_problem")
|
|
|
|
|
solution_status = JuMP.solve(SDP_problem);
|
2017-01-09 01:01:31 +01:00
|
|
|
|
verbose && @show solution_status
|
|
|
|
|
|
|
|
|
|
if solution_status != :Optimal
|
2017-02-26 13:48:31 +01:00
|
|
|
|
warn("The solver did not solve the problem successfully!")
|
2017-01-09 01:01:31 +01:00
|
|
|
|
end
|
2017-02-26 13:48:31 +01:00
|
|
|
|
|
2017-03-06 11:53:25 +01:00
|
|
|
|
κ = JuMP.getvalue(JuMP.getvariable(SDP_problem, :κ))
|
|
|
|
|
A = JuMP.getvalue(JuMP.getvariable(SDP_problem, :A))
|
|
|
|
|
@show sum(A)
|
2017-01-09 01:01:31 +01:00
|
|
|
|
return κ, A
|
|
|
|
|
end
|
|
|
|
|
|
|
|
|
|
function EOI{T<:Number}(Δ::GroupAlgebraElement{T}, κ::T)
|
|
|
|
|
return Δ*Δ - κ*Δ
|
|
|
|
|
end
|
|
|
|
|
|
2017-02-11 13:43:18 +01:00
|
|
|
|
function square_as_elt(vector, elt)
|
2016-12-23 00:51:06 +01:00
|
|
|
|
zzz = zeros(elt.coefficients)
|
2017-01-14 15:24:16 +01:00
|
|
|
|
zzz[1:length(vector)] = vector
|
|
|
|
|
# new_base_elt = GroupAlgebraElement(zzz, elt.product_matrix)
|
|
|
|
|
# return (new_base_elt*new_base_elt).coefficients
|
|
|
|
|
return GroupAlgebras.algebra_multiplication(zzz, zzz, elt.product_matrix)
|
|
|
|
|
end
|
|
|
|
|
|
|
|
|
|
function compute_SOS{T<:Number}(sqrt_matrix::Array{T,2},
|
|
|
|
|
elt::GroupAlgebraElement{T})
|
2017-02-11 13:44:51 +01:00
|
|
|
|
n = size(sqrt_matrix,2)
|
|
|
|
|
# result = zeros(T, length(elt.coefficients))
|
|
|
|
|
result = @parallel (+) for i in 1:n
|
|
|
|
|
square_as_elt(sqrt_matrix[:,i], elt)
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
2016-12-22 22:12:52 +01:00
|
|
|
|
return GroupAlgebraElement{T}(result, elt.product_matrix)
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
|
|
|
|
|
|
|
|
|
function correct_to_augmentation_ideal{T<:Rational}(sqrt_matrix::Array{T,2})
|
|
|
|
|
sqrt_corrected = similar(sqrt_matrix)
|
|
|
|
|
l = size(sqrt_matrix,2)
|
|
|
|
|
for i in 1:l
|
|
|
|
|
col = view(sqrt_matrix,:,i)
|
|
|
|
|
sqrt_corrected[:,i] = col - sum(col)//l
|
2017-01-09 01:01:31 +01:00
|
|
|
|
# @assert sum(sqrt_corrected[:,i]) == 0
|
2016-12-19 15:44:52 +01:00
|
|
|
|
end
|
|
|
|
|
return sqrt_corrected
|
|
|
|
|
end
|
2017-01-09 01:01:31 +01:00
|
|
|
|
|
2017-02-26 13:51:20 +01:00
|
|
|
|
function check_solution{T<:Number}(κ::T, sqrt_matrix::Array{T,2}, Δ::GroupAlgebraElement{T}; verbose=true, augmented=false)
|
2017-01-14 15:24:16 +01:00
|
|
|
|
result = compute_SOS(sqrt_matrix, Δ)
|
2017-02-11 13:45:56 +01:00
|
|
|
|
if augmented
|
|
|
|
|
@assert GroupAlgebras.ɛ(result) == 0//1
|
|
|
|
|
end
|
|
|
|
|
SOS_diff = EOI(Δ, κ) - result
|
|
|
|
|
|
|
|
|
|
eoi_SOS_L₁_dist = norm(SOS_diff,1)
|
|
|
|
|
|
|
|
|
|
if verbose
|
2017-02-11 13:52:32 +01:00
|
|
|
|
@show κ
|
2017-02-11 13:45:56 +01:00
|
|
|
|
if augmented
|
|
|
|
|
println("ɛ(Δ² - κΔ - ∑ξᵢ*ξᵢ) = ", GroupAlgebras.ɛ(SOS_diff))
|
|
|
|
|
else
|
|
|
|
|
ɛ_dist = Float64(round(GroupAlgebras.ɛ(SOS_diff),12))
|
|
|
|
|
println("ɛ(Δ² - κΔ - ∑ξᵢ*ξᵢ) ≈ $ɛ_dist")
|
|
|
|
|
end
|
|
|
|
|
L₁_dist = Float64(round(eoi_SOS_L₁_dist, 12))
|
|
|
|
|
println("‖Δ² - κΔ - ∑ξᵢ*ξᵢ‖₁ ≈ $L₁_dist")
|
|
|
|
|
end
|
|
|
|
|
|
2017-03-06 11:57:10 +01:00
|
|
|
|
distance_to_cone = κ - 2^2*eoi_SOS_L₁_dist
|
2017-02-11 13:45:56 +01:00
|
|
|
|
return distance_to_cone
|
2017-01-09 01:01:31 +01:00
|
|
|
|
end
|
|
|
|
|
|
2017-01-09 00:59:40 +01:00
|
|
|
|
function rationalize{T<:Integer, S<:Real}(::Type{T},
|
|
|
|
|
X::AbstractArray{S}; tol::Real=eps(eltype(X)))
|
|
|
|
|
r(x) = rationalize(T, x, tol=tol)
|
|
|
|
|
return r.(X)
|
2017-01-13 18:07:41 +01:00
|
|
|
|
end;
|
2017-02-11 13:46:22 +01:00
|
|
|
|
|
|
|
|
|
ℚ(x, tol::Real) = rationalize(BigInt, x, tol=tol)
|
|
|
|
|
|
2017-02-11 13:48:10 +01:00
|
|
|
|
|
2017-02-26 13:51:20 +01:00
|
|
|
|
function ℚ_distance_to_positive_cone(Δ::GroupAlgebraElement, κ, A;
|
|
|
|
|
tol=10.0^-7, verbose=true)
|
2017-02-11 13:48:10 +01:00
|
|
|
|
|
2017-02-26 13:51:20 +01:00
|
|
|
|
isapprox(eigvals(A), abs(eigvals(A)), atol=tol) ||
|
|
|
|
|
warn("The solution matrix doesn't seem to be positive definite!")
|
2017-02-11 13:48:10 +01:00
|
|
|
|
@assert A == Symmetric(A)
|
|
|
|
|
A_sqrt = real(sqrtm(A))
|
2017-02-11 13:52:32 +01:00
|
|
|
|
|
2017-02-26 13:51:20 +01:00
|
|
|
|
println("")
|
2017-02-11 13:48:10 +01:00
|
|
|
|
println("Checking in floating-point arithmetic...")
|
2017-02-26 13:51:20 +01:00
|
|
|
|
@time fp_distance = check_solution(κ, A_sqrt, Δ, verbose=verbose)
|
2017-03-06 11:55:40 +01:00
|
|
|
|
println("Floating point distance (to positive cone) ≈ $(Float64(trunc(fp_distance,8)))")
|
2017-02-11 13:48:10 +01:00
|
|
|
|
println("-------------------------------------------------------------")
|
|
|
|
|
println("")
|
|
|
|
|
|
2017-03-06 11:55:40 +01:00
|
|
|
|
if fp_distance ≤ 0
|
|
|
|
|
return fp_distance
|
|
|
|
|
end
|
|
|
|
|
|
2017-02-11 13:52:32 +01:00
|
|
|
|
println("Checking in rational arithmetic...")
|
2017-02-26 13:51:20 +01:00
|
|
|
|
κ_ℚ = ℚ(trunc(κ,Int(abs(log10(tol)))), tol)
|
|
|
|
|
A_sqrt_ℚ, Δ_ℚ = ℚ(A_sqrt, tol), ℚ(Δ, tol)
|
|
|
|
|
@time ℚ_distance = check_solution(κ_ℚ, A_sqrt_ℚ, Δ_ℚ, verbose=verbose)
|
2017-02-11 13:48:10 +01:00
|
|
|
|
@assert isa(ℚ_distance, Rational)
|
2017-03-06 11:55:40 +01:00
|
|
|
|
println("Rational distance (to positive cone) ≈ $(Float64(trunc(ℚ_distance,8)))")
|
2017-02-11 13:48:10 +01:00
|
|
|
|
println("-------------------------------------------------------------")
|
|
|
|
|
println("")
|
2017-03-06 11:55:40 +01:00
|
|
|
|
if ℚ_distance ≤ 0
|
|
|
|
|
return ℚ_distance
|
|
|
|
|
end
|
2017-02-11 13:52:32 +01:00
|
|
|
|
|
2017-02-11 13:48:10 +01:00
|
|
|
|
println("Projecting columns of A_sqrt to the augmentation ideal...")
|
|
|
|
|
A_sqrt_ℚ_aug = correct_to_augmentation_ideal(A_sqrt_ℚ)
|
2017-02-26 13:51:20 +01:00
|
|
|
|
@time ℚ_dist_to_Σ² = check_solution(κ_ℚ, A_sqrt_ℚ_aug, Δ_ℚ, verbose=verbose, augmented=true)
|
2017-02-11 13:48:10 +01:00
|
|
|
|
@assert isa(ℚ_dist_to_Σ², Rational)
|
2017-03-06 11:55:40 +01:00
|
|
|
|
println("Augmentation-projected rational distance (to positive cone)")
|
|
|
|
|
println("$(Float64(trunc(ℚ_dist_to_Σ²,8))) ≤ κ(G,S)")
|
|
|
|
|
println("-------------------------------------------------------------")
|
2017-02-11 13:48:10 +01:00
|
|
|
|
return ℚ_dist_to_Σ²
|
|
|
|
|
end
|