GroupsWithPropertyT/SL.jl

218 lines
5.2 KiB
Julia

using ArgParse
using GroupAlgebras
using PropertyT
using Mods
import Primes: isprime
import SCS.SCSSolver
function SL_generatingset(n::Int)
indexing = [(i,j) for i in 1:n for j in 1:n if i≠j]
S = [E(i,j,N=n) for (i,j) in indexing];
S = vcat(S, [convert(Array{Int,2},x') for x in S]);
S = vcat(S, [convert(Array{Int,2},inv(x)) for x in S]);
return unique(S)
end
function E(i::Int, j::Int; val=1, N::Int=3, mod=Inf)
@assert i≠j
m = eye(Int, N)
m[i,j] = val
if mod == Inf
return m
else
return [Mod(x,mod) for x in m]
end
end
function cofactor(i,j,M)
z1 = ones(Bool,size(M,1))
z1[i] = false
z2 = ones(Bool,size(M,2))
z2[j] = false
return M[z1,z2]
end
import Base.LinAlg.det
function det(M::Array{Mod,2})
if size(M,1) size(M,2)
d = Mod(0,M[1,1].mod)
elseif size(M,1) == 2
d = M[1,1]*M[2,2] - M[1,2]*M[2,1]
else
d = zero(eltype(M))
for i in 1:size(M,1)
d += (-1)^(i+1)*M[i,1]*det(cofactor(i,1,M))
end
end
# @show (M, d)
return d
end
function adjugate(M)
K = similar(M)
for i in 1:size(M,1), j in 1:size(M,2)
K[j,i] = (-1)^(i+j)*det(cofactor(i,j,M))
end
return K
end
import Base: inv, one, zero, *
one(::Type{Mod}) = 1
zero(::Type{Mod}) = 0
zero(x::Mod) = Mod(x.mod)
function inv(M::Array{Mod,2})
d = det(M)
d 0*d || thow(ArgumentError("Matrix is not invertible!"))
return inv(det(M))*adjugate(M)
return adjugate(M)
end
function SL_generatingset(n::Int, p::Int)
p == 0 && return SL_generatingset(n)
(p > 1 && n > 0) || throw(ArgumentError("Both n and p should be positive integers!"))
isprime(p) || throw(ArgumentError("p should be a prime number!"))
indexing = [(i,j) for i in 1:n for j in 1:n if i≠j]
S = [E(i,j, N=n, mod=p) for (i,j) in indexing]
S = vcat(S, [inv(s) for s in S])
S = vcat(S, [permutedims(x, [2,1]) for x in S]);
return unique(S)
end
function products{T}(U::AbstractVector{T}, V::AbstractVector{T})
result = Vector{T}()
for u in U
for v in V
push!(result, u*v)
end
end
return unique(result)
end
function ΔandSDPconstraints{T<:Number}(identity::Array{T,2}, S::Vector{Array{T,2}})
B₁ = vcat([identity], S)
B₂ = products(B₁, B₁);
B₃ = products(B₁, B₂);
B₄ = products(B₁, B₃);
@assert B₄[1:length(B₂)] == B₂
product_matrix = PropertyT.create_product_matrix(B₄,length(B₂));
sdp_constraints = PropertyT.constraints_from_pm(product_matrix, length(B₄))
L_coeff = PropertyT.splaplacian_coeff(S, B₂, length(B₄));
Δ = GroupAlgebraElement(L_coeff, product_matrix)
return Δ, sdp_constraints
end
ID(n::Int) = eye(Int, n)
function ID(n::Int, p::Int)
p==0 && return ID(n)
return [Mod(x,p) for x in eye(Int,N)]
end
#=
To use file property(T).jl (specifically: check_property_T function)
You need to define:
function ΔandSDPconstraints(identity, S):: (Δ, sdp_constraints)
=#
function cpuinfo_physicalcores()
maxcore = -1
for line in eachline("/proc/cpuinfo")
if startswith(line, "core id")
maxcore = max(maxcore, parse(Int, split(line, ':')[2]))
end
end
maxcore < 0 && error("failure to read core ids from /proc/cpuinfo")
return maxcore + 1
end
function parse_commandline()
s = ArgParseSettings()
@add_arg_table s begin
"--tol"
help = "set numerical tolerance for the SDP solver"
arg_type = Float64
default = 1e-9
"--iterations"
help = "set maximal number of iterations for the SDP solver"
arg_type = Int
default = 100000
"--upper-bound"
help = "Set an upper bound for the spectral gap"
arg_type = Float64
default = Inf
"--cpus"
help = "Set number of cpus used by solver"
arg_type = Int
required = false
"-N"
help = "Consider matrices of size N"
arg_type = Int
default = 3
"-p"
help = "Matrices over filed of p-elements (0 = over ZZ)"
arg_type = Int
default = 0
end
return parse_args(s)
end
function main()
parsed_args = parse_commandline()
println("Parsed args:")
# SL(3,Z)
# upper_bound = 0.28-1e-5
# tol = 1e-12
# iterations = 500000
# SL(4,Z)
# upper_bound = 1.31
# tol = 3e-11
# upper_bound=0.738 # (N,p) = (3,7)
tol = parsed_args["tol"]
iterations = parsed_args["iterations"]
# solver = MosekSolver(INTPNT_CO_TOL_REL_GAP=tol, QUIET=false)
solver = SCSSolver(eps=tol, max_iters=iterations, verbose=true)
N = parsed_args["N"]
upper_bound = parsed_args["upper-bound"]
p = parsed_args["p"]
if p == 0
name = "SL$(N)Z"
else
name = "SL$(N)_p"
end
S() = SL_generatingset(N, p)
if parsed_args["cpus"] nothing
if parsed_args["cpus"] > cpuinfo_physicalcores()
warn("Number of specified cores exceeds the physical core cound. Performance will suffer.")
end
Blas.set_num_threads(parsed_args["cpus"])
end
@time PropertyT.check_property_T(name, ID(N,p), S, solver, upper_bound, tol)
end
main()