GroupsWithPropertyT/AutFN.jl

242 lines
6.8 KiB
Julia
Raw Normal View History

using ArgParse
2017-03-13 16:18:42 +01:00
using Groups
using GroupAlgebras
using PropertyT
import SCS.SCSSolver
2017-03-13 16:18:42 +01:00
#=
Note that the element
α(i,j,k) = ϱ(i,j)*ϱ(i,k)*inv(ϱ(i,j))*inv(ϱ(i,k)),
which surely belongs to ball of radius 4 in Aut(F₄) becomes trivial under the representation
Aut(F₄) GL₄()ℤ⁴ GL₅().
Moreover, due to work of Potapchik and Rapinchuk [1] every real representation of Aut(Fₙ) into GLₘ() (for m 2n-2) factors through GLₙ()ℤⁿ, so will have the same problem.
We need a different approach: Here we actually compute in Aut(𝔽₄)
=#
import Combinatorics.nthperm
SymmetricGroup(n) = [nthperm(collect(1:n), k) for k in 1:factorial(n)]
function generating_set_of_AutF(N::Int)
indexing = [[i,j] for i in 1:N for j in 1:N if i≠j]
σs = [symmetric_AutSymbol(perm) for perm in SymmetricGroup(N)[2:end]];
ϱs = [rmul_AutSymbol(i,j) for (i,j) in indexing]
λs = [lmul_AutSymbol(i,j) for (i,j) in indexing]
ɛs = [flip_AutSymbol(i) for i in 1:N];
S = vcat(ϱs,λs)
S = vcat(S..., σs..., ɛs)
S = vcat(S..., [inv(g) for g in S])
return Vector{AutWord}(unique(S)), one(AutWord)
2017-03-13 16:18:42 +01:00
end
function generating_set_of_OutF(N::Int)
indexing = [[i,j] for i in 1:N for j in 1:N if i≠j]
ϱs = [rmul_AutSymbol(i,j) for (i,j) in indexing]
λs = [lmul_AutSymbol(i,j) for (i,j) in indexing]
ɛs = [flip_AutSymbol(i) for i in 1:N];
S = ϱs
push!(S, λs..., ɛs...)
push!(S,[inv(g) for g in S]...)
return Vector{AutWord}(unique(S)), one(AutWord)
2017-03-13 16:18:42 +01:00
end
2017-04-10 20:02:17 +02:00
function generating_set_of_SOutF(N::Int)
indexing = [[i,j] for i in 1:N for j in 1:N if i≠j]
ϱs = [rmul_AutSymbol(i,j) for (i,j) in indexing]
λs = [lmul_AutSymbol(i,j) for (i,j) in indexing]
S = ϱs
push!(S, λs...)
push!(S,[inv(g) for g in S]...)
return Vector{AutWord}(unique(S)), one(AutWord)
end
2017-03-13 16:18:42 +01:00
function generating_set_of_Sym(N::Int)
σs = [symmetric_AutSymbol(perm) for perm in SymmetricGroup(N)[2:end]];
S = σs
push!(S, [inv(s) for s in S]...)
return Vector{AutWord}(unique(S)), one(AutWord)
2017-03-13 16:18:42 +01:00
end
function products(S1::Vector{AutWord}, S2::Vector{AutWord})
result = Vector{AutWord}()
seen = Set{Vector{FGWord}}()
n = length(S1)
for (i,x) in enumerate(S1)
for y in S2
z::AutWord = x*y
v::Vector{FGWord} = z(domain)
if !in(v, seen)
push!(seen, v)
push!(result, z)
end
end
end
return result
end
function products_images(S1::Vector{AutWord}, S2::Vector{AutWord})
result = Vector{Vector{FGWord}}()
seen = Set{Vector{FGWord}}()
n = length(S1)
for (i,x) in enumerate(S1)
z = x(domain)
for y in S2
v = y(z)
if !in(v, seen)
push!(seen, v)
push!(result, v)
end
end
end
return result
end
function hashed_product{T}(image::T, B, images_dict::Dict{T, Int})
n = size(B,1)
column = zeros(Int,n)
Threads.@threads for j in 1:n
w = (B[j])(image)
k = images_dict[w]
k 0 || throw(ArgumentError(
"($i,$j): $(x^-1)*$y don't seem to be supported on basis!"))
column[j] = k
end
return column
end
2017-04-10 21:43:56 +02:00
function create_product_matrix(images, basis::Vector{AutWord})
2017-03-13 16:18:42 +01:00
n = length(basis)
product_matrix = zeros(Int, (n, n));
print("Creating hashtable of images...")
@time images_dict = Dict{Vector{FGWord}, Int}(x => i
for (i,x) in enumerate(images))
for i in 1:n
z = (inv(basis[i]))(domain)
product_matrix[i,:] = hashed_product(z, basis, images_dict)
end
return product_matrix
end
2017-04-10 21:43:56 +02:00
function generate_balls{T}(S::Vector{T}, Id::T; radius=4)
sizes = Vector{Int}()
S = vcat([Id], S)
B = [Id]
for i in 1:radius
B = products(B, S);
push!(sizes, length(B))
2017-03-13 16:18:42 +01:00
end
2017-04-10 21:43:56 +02:00
return B, sizes
end
2017-03-13 16:18:42 +01:00
2017-04-10 21:43:56 +02:00
function ΔandSDPconstraints(Id::AutWord, S::Vector{AutWord}, r::Int=2)
B, sizes = generate_balls(S, Id, radius=2*r)
basis = B[1:sizes[r]]
B_images = unique([f(domain) for f in B])
println("Generated balls of sizes $sizes")
2017-03-13 16:18:42 +01:00
println("Creating product matrix...")
2017-04-10 21:43:56 +02:00
@time pm = create_product_matrix(B_images, basis)
2017-03-13 16:18:42 +01:00
println("Creating sdp_constratints...")
2017-03-14 16:13:05 +01:00
@time sdp_constraints = PropertyT.constraints_from_pm(pm)
2017-03-13 16:18:42 +01:00
2017-04-10 21:43:56 +02:00
L_coeff = PropertyT.splaplacian_coeff(S, basis, length(B_images))
2017-03-14 16:13:05 +01:00
Δ = PropertyT.GroupAlgebraElement(L_coeff, Array{Int,2}(pm))
2017-03-13 16:18:42 +01:00
return Δ, sdp_constraints
end
const symbols = [FGSymbol("x₁",1), FGSymbol("x₂",1), FGSymbol("x₃",1), FGSymbol("x₄",1), FGSymbol("x₅",1), FGSymbol("x₆",1)]
const TOL=1e-8
const N = 4
const domain = Vector{FGWord}(symbols[1:N])
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 (default: 1e-5)"
arg_type = Float64
default = 1e-5
"--iterations"
help = "set maximal number of iterations for the SDP solver (default: 20000)"
arg_type = Int
default = 20000
"--upper-bound"
help = "Set an upper bound for the spectral gap (default: Inf)"
arg_type = Float64
default = Inf
"--cpus"
help = "Set number of cpus used by solver (default: auto)"
arg_type = Int
required = false
"-N"
help = "Consider automorphisms of free group on N generators (default: N=3)"
arg_type = Int
default = 3
end
return parse_args(s)
end
2017-03-13 16:18:42 +01:00
# const name = "SYM$N"
# const upper_bound=factorial(N)-TOL^(1/5)
# S() = generating_set_of_Sym(N)
# name = "AutF$N"
# S() = generating_set_of_AutF(N)
function main()
parsed_args = parse_commandline()
tol = parsed_args["tol"]
iterations = parsed_args["iterations"]
solver = SCSSolver(eps=tol, max_iters=iterations, verbose=true, linearsolver=SCS.Indirect)
N = parsed_args["N"]
upper_bound = parsed_args["upper-bound"]
2017-04-10 20:02:17 +02:00
name = "SOutF$N"
name = name*"-$(string(upper_bound))"
2017-04-10 20:02:17 +02:00
S() = generating_set_of_SOutF(N)
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
2017-04-10 20:01:56 +02:00
@time PropertyT.check_property_T(name, S, solver, upper_bound, tol, 2)
return 0
end
2017-03-13 16:18:42 +01:00
main()