Tries on AutF4 with semidirect product
This commit is contained in:
parent
e1949b5aaf
commit
5629a39ddd
120
AutF4.jl
Normal file
120
AutF4.jl
Normal file
@ -0,0 +1,120 @@
|
|||||||
|
using Combinatorics
|
||||||
|
|
||||||
|
using JuMP
|
||||||
|
import SCS: SCSSolver
|
||||||
|
import Mosek: MosekSolver
|
||||||
|
|
||||||
|
push!(LOAD_PATH, "./")
|
||||||
|
using SemiDirectProduct
|
||||||
|
using GroupAlgebras
|
||||||
|
include("property(T).jl")
|
||||||
|
|
||||||
|
const N = 4
|
||||||
|
|
||||||
|
const VERBOSE = true
|
||||||
|
|
||||||
|
function permutation_matrix(p::Vector{Int})
|
||||||
|
n = length(p)
|
||||||
|
sort(p) == collect(1:n) || throw(ArgumentError("Input array must be a permutation of 1:n"))
|
||||||
|
A = eye(n)
|
||||||
|
return A[p,:]
|
||||||
|
end
|
||||||
|
|
||||||
|
SymmetricGroup(n) = [nthperm(collect(1:n), k) for k in 1:factorial(n)]
|
||||||
|
|
||||||
|
# const SymmetricGroup = [permutation_matrix(x) for x in SymmetricGroup_perms]
|
||||||
|
|
||||||
|
function E(i, j; dim::Int=N)
|
||||||
|
@assert i≠j
|
||||||
|
k = eye(dim)
|
||||||
|
k[i,j] = 1
|
||||||
|
return k
|
||||||
|
end
|
||||||
|
|
||||||
|
function eltary_basis_vector(i; dim::Int=N)
|
||||||
|
result = zeros(dim)
|
||||||
|
if 0 < i ≤ dim
|
||||||
|
result[i] = 1
|
||||||
|
end
|
||||||
|
return result
|
||||||
|
end
|
||||||
|
|
||||||
|
v(i; dim=N) = eltary_basis_vector(i,dim=dim)
|
||||||
|
|
||||||
|
ϱ(i,j::Int,n=N) = SemiDirectProductElement(E(i,j,dim=n), v(j,dim=n))
|
||||||
|
λ(i,j::Int,n=N) = SemiDirectProductElement(E(i,j,dim=n), -v(j,dim=n))
|
||||||
|
|
||||||
|
function ɛ(i, n::Int=N)
|
||||||
|
result = eye(n)
|
||||||
|
result[i,i] = -1
|
||||||
|
return SemiDirectProductElement(result)
|
||||||
|
end
|
||||||
|
|
||||||
|
σ(permutation::Vector{Int}) =
|
||||||
|
SemiDirectProductElement(permutation_matrix(permutation))
|
||||||
|
|
||||||
|
# Standard generating set: 103 elements
|
||||||
|
|
||||||
|
function generatingset_ofAutF(n::Int=N)
|
||||||
|
indexing = [[i,j] for i in 1:n for j in 1:n if i≠j]
|
||||||
|
ϱs = [ϱ(ij...) for ij in indexing]
|
||||||
|
λs = [λ(ij...) for ij in indexing]
|
||||||
|
ɛs = [ɛ(i) for i in 1:N]
|
||||||
|
σs = [σ(perm) for perm in SymmetricGroup(n)]
|
||||||
|
S = vcat(ϱs, λs, ɛs, σs);
|
||||||
|
S = unique(vcat(S, [inv(x) for x in S]));
|
||||||
|
return S
|
||||||
|
end
|
||||||
|
|
||||||
|
#=
|
||||||
|
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!
|
||||||
|
=#
|
||||||
|
|
||||||
|
const ID = eye(N+1)
|
||||||
|
|
||||||
|
const S₁ = generatingset_ofAutF(N)
|
||||||
|
|
||||||
|
matrix_S₁ = [matrix_repr(x) for x in S₁]
|
||||||
|
|
||||||
|
const TOL=10.0^-7
|
||||||
|
|
||||||
|
matrix_S₁[1:10,:][:,1]
|
||||||
|
|
||||||
|
Δ, cm = prepare_Laplacian_and_constraints(matrix_S₁)
|
||||||
|
|
||||||
|
#solver = SCSSolver(eps=TOL, max_iters=ITERATIONS, verbose=true);
|
||||||
|
solver = MosekSolver(MSK_DPAR_INTPNT_CO_TOL_REL_GAP=TOL,
|
||||||
|
# MSK_DPAR_INTPNT_CO_TOL_PFEAS=1e-15,
|
||||||
|
# MSK_DPAR_INTPNT_CO_TOL_DFEAS=1e-15,
|
||||||
|
# MSK_IPAR_PRESOLVE_USE=0,
|
||||||
|
QUIET=!VERBOSE)
|
||||||
|
|
||||||
|
# κ, A = solve_for_property_T(S₁, solver, verbose=VERBOSE)
|
||||||
|
|
||||||
|
product_matrix = readdlm("SL₃Z.product_matrix", Int)
|
||||||
|
L = readdlm("SL₃Z.Δ.coefficients")[:, 1]
|
||||||
|
Δ = GroupAlgebraElement(L, product_matrix)
|
||||||
|
|
||||||
|
A = readdlm("matrix.A.Mosek")
|
||||||
|
κ = readdlm("kappa.Mosek")[1]
|
||||||
|
|
||||||
|
# @show eigvals(A)
|
||||||
|
@assert isapprox(eigvals(A), abs(eigvals(A)), atol=TOL)
|
||||||
|
@assert A == Symmetric(A)
|
||||||
|
|
||||||
|
|
||||||
|
const A_sqrt = real(sqrtm(A))
|
||||||
|
|
||||||
|
SOS_EOI_fp_L₁, Ω_fp_dist = check_solution(κ, A_sqrt, Δ)
|
||||||
|
|
||||||
|
κ_rational = rationalize(BigInt, κ;)
|
||||||
|
A_sqrt_rational = rationalize(BigInt, A_sqrt)
|
||||||
|
Δ_rational = rationalize(BigInt, Δ)
|
||||||
|
|
||||||
|
SOS_EOI_rat_L₁, Ω_rat_dist = check_solution(κ_rational, A_sqrt_rational, Δ_rational)
|
88
SemiDirectProduct.jl
Normal file
88
SemiDirectProduct.jl
Normal file
@ -0,0 +1,88 @@
|
|||||||
|
module SemiDirectProduct
|
||||||
|
|
||||||
|
import Base: convert, show, isequal, ==, size, inv
|
||||||
|
import Base: +, -, *, //
|
||||||
|
|
||||||
|
export SemiDirectProductElement, matrix_repr
|
||||||
|
|
||||||
|
"""
|
||||||
|
Implements elements of a semidirect product of groups H and N, where N is normal in the product. Usually written as H ⋉ N.
|
||||||
|
The multiplication inside semidirect product is defined as
|
||||||
|
(h₁, n₁) ⋅ (h₂, n₂) = (h₁h₂, n₁φ(h₁)(n₂)),
|
||||||
|
where φ:H → Aut(N) is a homomorphism.
|
||||||
|
|
||||||
|
In the case below we implement H = GL(n,K) and N = Kⁿ, the Affine Group (i.e. GL(n,K) ⋉ Kⁿ) where elements of GL(n,K) act on vectors in Kⁿ via matrix multiplication.
|
||||||
|
# Arguments:
|
||||||
|
* `h::Array{T,2}` : square invertible matrix (element of GL(n,K))
|
||||||
|
* `n::Vector{T,1}` : vector in Kⁿ
|
||||||
|
* `φ = φ(h,n) = φ(h)(n)` :2-argument function which defines the action of GL(n,K) on Kⁿ; matrix-vector multiplication by default.
|
||||||
|
"""
|
||||||
|
immutable SemiDirectProductElement{T<:Number}
|
||||||
|
h::Array{T,2}
|
||||||
|
n::Vector{T}
|
||||||
|
φ::Function
|
||||||
|
|
||||||
|
function SemiDirectProductElement(h::Array{T,2},n::Vector{T},φ::Function)
|
||||||
|
# size(h,1) == size(h,2)|| throw(ArgumentError("h has to be square matrix"))
|
||||||
|
det(h) ≠ 0 || throw(ArgumentError("h has to be invertible!"))
|
||||||
|
new(h,n,φ)
|
||||||
|
end
|
||||||
|
end
|
||||||
|
|
||||||
|
SemiDirectProductElement{T}(h::Array{T,2}, n::Vector{T}, φ) =
|
||||||
|
SemiDirectProductElement{T}(h,n,φ)
|
||||||
|
|
||||||
|
SemiDirectProductElement{T}(h::Array{T,2}, n::Vector{T}) =
|
||||||
|
SemiDirectProductElement(h,n,*)
|
||||||
|
|
||||||
|
SemiDirectProductElement{T}(h::Array{T,2}) =
|
||||||
|
SemiDirectProductElement(h,zeros(h[:,1]))
|
||||||
|
|
||||||
|
SemiDirectProductElement{T}(n::Vector{T}) =
|
||||||
|
SemiDirectProductElement(eye(eltype(n), n))
|
||||||
|
|
||||||
|
convert{T<:Number}(::Type{T}, X::SemiDirectProductElement) =
|
||||||
|
SemiDirectProductElement(convert(Array{T,2},X.h),
|
||||||
|
convert(Vector{T},X.n),
|
||||||
|
X.φ)
|
||||||
|
|
||||||
|
size(X::SemiDirectProductElement) = (size(X.h), size(X.n))
|
||||||
|
|
||||||
|
matrix_repr{T}(X::SemiDirectProductElement{T}) =
|
||||||
|
[X.h X.n; zeros(T, 1, size(X.h,2)) [1]]
|
||||||
|
|
||||||
|
show{T}(io::IO, X::SemiDirectProductElement{T}) = print(io,
|
||||||
|
"Element of SemiDirectProduct over $T of size $(size(X)):\n",
|
||||||
|
matrix_repr(X))
|
||||||
|
|
||||||
|
function isequal{T}(X::SemiDirectProductElement{T}, Y::SemiDirectProductElement{T})
|
||||||
|
X.h == Y.h || return false
|
||||||
|
X.n == Y.n || return false
|
||||||
|
X.φ == Y.φ || return false
|
||||||
|
return true
|
||||||
|
end
|
||||||
|
|
||||||
|
function isequal{T,S}(X::SemiDirectProductElement{T}, Y::SemiDirectProductElement{S})
|
||||||
|
W = promote_type(T,S)
|
||||||
|
warn("Comparing elements with different coefficients! trying to promoting to $W")
|
||||||
|
X = convert(W, X)
|
||||||
|
Y = convert(W, Y)
|
||||||
|
return isequal(X,Y)
|
||||||
|
end
|
||||||
|
|
||||||
|
(==)(X::SemiDirectProductElement, Y::SemiDirectProductElement) = isequal(X, Y)
|
||||||
|
|
||||||
|
function semidirect_multiplication{T}(X::SemiDirectProductElement{T},
|
||||||
|
Y::SemiDirectProductElement{T})
|
||||||
|
size(X) == size(Y) || throw(ArgumentError("trying to multiply elements from different groups!"))
|
||||||
|
return SemiDirectProductElement(X.h*Y.h, X.n + X.φ(X.h, Y.n))
|
||||||
|
end
|
||||||
|
|
||||||
|
(*){T}(X::SemiDirectProductElement{T}, Y::SemiDirectProductElement{T}) =
|
||||||
|
semidirect_multiplication(X,Y)
|
||||||
|
|
||||||
|
inv{T}(X::SemiDirectProductElement{T}) =
|
||||||
|
SemiDirectProductElement(inv(X.h), X.φ(inv(X.h), -X.n))
|
||||||
|
|
||||||
|
|
||||||
|
end
|
Loading…
Reference in New Issue
Block a user