diff --git a/Project.toml b/Project.toml index 5db198c..4b36bb7 100644 --- a/Project.toml +++ b/Project.toml @@ -1,7 +1,7 @@ name = "PeriodicGraphs" uuid = "18c5b727-b240-4874-878a-f2e242435bab" authors = ["Lionel Zoubritzky lionel.zoubritzky@gmail.com"] -version = "0.5.0" +version = "0.5.1" [deps] Graphs = "86223c79-3864-5bf0-83f7-82e725a168b6" diff --git a/src/PeriodicGraphs.jl b/src/PeriodicGraphs.jl index bdd6d64..5b25e20 100644 --- a/src/PeriodicGraphs.jl +++ b/src/PeriodicGraphs.jl @@ -223,8 +223,8 @@ const PeriodicGraph1D = PeriodicGraph{1} const PeriodicGraph2D = PeriodicGraph{2} const PeriodicGraph3D = PeriodicGraph{3} -function PeriodicGraph{N}(nv::Integer, t, s) where N - return PeriodicGraph{N}(Ref(nv), t, s, Ref(-1//1)) +function PeriodicGraph{N}(ne::Integer, t, s) where N + return PeriodicGraph{N}(Ref(ne), t, s, Ref(-1//1)) end """ @@ -233,9 +233,81 @@ end Construct a `PeriodicGraph{N}` with `nv` vertices and 0 edge. """ function PeriodicGraph{N}(n::Integer = 0) where N - @assert n >= 0 return PeriodicGraph{N}(0, [PeriodicVertex{N}[] for _ in 1:n], [1 for _ in 1:n]) end + +function from_edges(nv, t::AbstractVector{PeriodicEdge{N}}) where {N} + isempty(t) && return PeriodicGraph{N}(nv) + + directedgestart = ones(Int, nv) + counterdirect = zeros(Int, nv) + numdifferentindirect = zeros(Int, nv) + last_e = first(t) + numdifferentindirect[last_e.dst.v] = 1 + for e in t + s = e.src + d = e.dst.v + directedgestart[d] += 1 + counterdirect[s] += 1 + numdifferentindirect[d] += ((s != last_e.src) | (d != last_e.dst.v)) + last_e = e + end + + nlist = Vector{Vector{PeriodicVertex{N}}}(undef, nv) + counterindirect = Vector{Vector{Int}}(undef, nv) + for i in 1:nv + nlist[i] = Vector{PeriodicVertex{N}}(undef, counterdirect[i] + directedgestart[i] - 1) + counterdirect[i] = 0 + counterindirect[i] = zeros(Int, numdifferentindirect[i] + 1) + end + + # First pass: set the direct edges + idx_indirect = ones(Int, nv) + last_e = first(t) + for e in t + s = e.src + d = e.dst.v + nlist[s][directedgestart[s] + counterdirect[s]] = e.dst + counterdirect[s] += 1 + if ((s != last_e.src) | (d != last_e.dst.v)) + last_d = last_e.dst.v + idx_indirect[last_d] += 1 + val = idx_indirect[last_d] + counterindirect[last_d][val] += counterindirect[last_d][val-1] + end + counterindirect[d][idx_indirect[d]] += 1 + last_e = e + end + + @simd for i in 1:nv; idx_indirect[i] = 1; end + + # Second pass: set the indirect edges + last_e = first(t) + current_idx = idx_indirect[last_e.dst.v] + for e in t + s = e.src + d = e.dst.v + if ((s != last_e.src) | (d != last_e.dst.v)) + idx_indirect[last_e.dst.v] += 1 + current_idx = idx_indirect[d] + end + nlist[d][counterindirect[d][current_idx]] = PeriodicVertex(s, .- e.dst.ofs) + counterindirect[d][current_idx] -= 1 + last_e = e + end + + return PeriodicGraph{N}(length(t), nlist, directedgestart) +end + +function from_edges(t::AbstractVector{PeriodicEdge{N}}) where N + @static if VERSION < v"1.6-" + isempty(t) && return from_edges(0, t) + return from_edges(maximum(max(e.src, e.dst.v) for e in t), t) + else + return from_edges(maximum(max(e.src, e.dst.v) for e in t; init=0), t) + end +end + """ PeriodicGraph([nv::Integer, ]edge_list::AbstractVector{PeriodicEdge{N}}) PeriodicGraph{N}([nv::Integer, ]edge_list::AbstractVector{PeriodicEdge{N}}) @@ -260,21 +332,22 @@ julia> ne(g) ``` """ function PeriodicGraph{N}(nv::Integer, t::AbstractVector{PeriodicEdge{N}}) where N - sort!(t); unique!(t) - ne = length(t) - g = PeriodicGraph{N}(nv) - for e in t - add_edge!(g, e) + for (i, e) in enumerate(t) + if isindirectedge(e) + t[i] = reverse(e) + end end - return g + sort!(t); unique!(t) + return from_edges(nv, t) end + PeriodicGraph(nv::Integer, t::AbstractVector{PeriodicEdge{N}}) where {N} = PeriodicGraph{N}(nv, t) function PeriodicGraph{N}(t::AbstractVector{PeriodicEdge{N}}) where N @static if VERSION < v"1.6-" - isempty(t) && return PeriodicGraph(0, t) - return PeriodicGraph(maximum(max(e.src, e.dst.v) for e in t), t) + isempty(t) && return PeriodicGraph{N}(0, t) + return PeriodicGraph{N}(maximum(max(e.src, e.dst.v) for e in t), t) else - return PeriodicGraph(maximum(max(e.src, e.dst.v) for e in t; init=0), t) + return PeriodicGraph{N}(maximum(max(e.src, e.dst.v) for e in t; init=0), t) end end PeriodicGraph(t::AbstractVector{PeriodicEdge{N}}) where {N} = PeriodicGraph{N}(t) @@ -325,6 +398,19 @@ end Base.isempty(k::KeyString) = isempty(SubString(k.x, k.start[])) Base.IteratorSize(::Type{<:KeyString}) = Base.SizeUnknown() +function edges_from_string(key::KeyString{Int}, ::Val{N}) where N + edgs = PeriodicEdge{N}[] + while !isempty(key) + src = popfirst!(key) + dst = popfirst!(key) + ofs = SVector{N,Int}([popfirst!(key) for _ in 1:N]) + push!(edgs, PeriodicEdge{N}(src, dst, ofs)) + end + return edgs +end + + + """ PeriodicGraph(key::AbstractString) PeriodicGraph{N}(key::AbstractString) @@ -341,6 +427,10 @@ and for all `i` between `1` and `m`, the number of edges, the `i`-th edge is des This compact representation of a graph can be obtained simply by `print`ing the graph or with `string`. +!!! note + Use `parse(PeriodicGraph, key)` or `parse(PeriodicGraph{N}, key)` for a faster + implementation if `key` was obtained from `string(g)` with `g` a `PeriodicGraph{N}`. + ## Examples ```jldoctest julia> PeriodicGraph("2 1 2 0 0 2 1 1 0 1 1 0 -1") @@ -352,7 +442,7 @@ PeriodicGraph3D(1, PeriodicEdge3D[(1, 1, (0,0,1)), (1, 1, (0,1,0)), (1, 1, (1,0, julia> string(ans) "3 1 1 0 0 1 1 1 0 1 0 1 1 1 0 0" -julia> string(PeriodicGraph(ans)) == ans +julia> string(parse(PeriodicGraph3D, ans)) == ans true ``` """ @@ -360,24 +450,24 @@ function PeriodicGraph{N}(s::AbstractString) where N key = KeyString{Int}(s) M = popfirst!(key) M != N && throw(DimensionMismatch("Cannot construct a $N-dimensional graph from a $M-dimensional key")) - edgs = PeriodicEdge{N}[] - while !isempty(key) - src = popfirst!(key) - dst = popfirst!(key) - ofs = SVector{N,Int}([popfirst!(key) for _ in 1:N]) - push!(edgs, PeriodicEdge{N}(src, dst, ofs)) - end - return PeriodicGraph{N}(edgs) + return PeriodicGraph{N}(edges_from_string(key, Val(N))) end function PeriodicGraph(s::AbstractString) - N = first(KeyString{Int}(s)) - return PeriodicGraph{N}(s) -end -function PeriodicGraph{N}(g::PeriodicGraph{N}) where N - PeriodicGraph{N}(g.ne[], [copy(x) for x in g.nlist], copy(g.directedgestart)) + key = KeyString{Int}(s) + N = popfirst!(key) + return PeriodicGraph{N}(edges_from_string(key, Val(N))) end -PeriodicGraph(g::PeriodicGraph{N}) where {N} = PeriodicGraph{N}(g) +function Base.parse(::Type{PeriodicGraph{N}}, s::AbstractString) where N + key = KeyString{Int}(s) + popfirst!(key) # No verification is done for this function + return from_edges(edges_from_string(key, Val(N))) +end +function Base.parse(::Type{PeriodicGraph}, s::AbstractString) + key = KeyString{Int}(s) + N = popfirst!(key) + return from_edges(edges_from_string(key, Val(N))) +end function show(io::IO, g::PeriodicGraph{N}) where N if get(io, :typeinfo, Any) != PeriodicGraph{N} @@ -392,6 +482,13 @@ function print(io::IO, g::PeriodicGraph{N}) where N end end + +function PeriodicGraph{N}(g::PeriodicGraph{N}) where N + PeriodicGraph{N}(g.ne[], [copy(x) for x in g.nlist], copy(g.directedgestart)) +end +PeriodicGraph(g::PeriodicGraph{N}) where {N} = PeriodicGraph{N}(g) +Base.copy(g::PeriodicGraph{N}) where {N} = PeriodicGraph{N}(g) + """ ==(g1::PeriodicGraph, g2::PeriodicGraph) @@ -925,8 +1022,7 @@ function graph_width!(g::PeriodicGraph{N}) where N g.width[] = width == nv(g)+1 ? Rational(maxa) : width end -function Graphs._neighborhood(g::Union{PeriodicGraph{0},PeriodicGraph{1},PeriodicGraph{2},PeriodicGraph{3}}, v::Integer, d::Real, distmx::AbstractMatrix{U}, neighborfn::Function) where U <: Real - @assert typeof(neighborfn) === typeof(outneighbors) +function Graphs._neighborhood(g::Union{PeriodicGraph{0},PeriodicGraph{1},PeriodicGraph{2},PeriodicGraph{3}}, v::Integer, d::Real, distmx::AbstractMatrix{U}, ::typeof(outneighbors)) where U <: Real N = ndims(g) Q = Tuple{PeriodicVertex{N}, U}[] d < zero(U) && return Q @@ -957,8 +1053,7 @@ function Graphs._neighborhood(g::Union{PeriodicGraph{0},PeriodicGraph{1},Periodi return Q end -function Graphs._neighborhood(g::PeriodicGraph{N}, v::Integer, d::Real, distmx::AbstractMatrix{U}, neighborfn::Function) where {N,U <: Real} - @assert typeof(neighborfn) === typeof(outneighbors) +function Graphs._neighborhood(g::PeriodicGraph{N}, v::Integer, d::Real, distmx::AbstractMatrix{U}, ::typeof(outneighbors)) where {N,U <: Real} Q = Tuple{PeriodicVertex, U}[] d < zero(U) && return Q start_vertex = PeriodicVertex{N}(v) diff --git a/src/precompile.jl b/src/precompile.jl index 69c4d30..f8d76e2 100644 --- a/src/precompile.jl +++ b/src/precompile.jl @@ -21,7 +21,11 @@ function _precompile_() @enforce Base.precompile(Tuple{Type{PeriodicGraph{i}}}) @enforce Base.precompile(Tuple{Type{PeriodicGraph{i}},Vector{PeriodicEdge{i}}}) @enforce Base.precompile(Tuple{Type{PeriodicGraph{i}},String}) + @enforce Base.precompile(Tuple{Type{PeriodicGraph{i}},SubString{String}}) @enforce Base.precompile(Tuple{Type{PeriodicGraph},PeriodicGraph{i}}) + @enforce Base.precompile(Tuple{typeof(copy),PeriodicGraph{i}}) + @enforce Base.precompile(Tuple{typeof(Base.parse),Type{PeriodicGraph{i}},String}) + @enforce Base.precompile(Tuple{typeof(Base.parse),Type{PeriodicGraph{i}},SubString{String}}) @enforce Base.precompile(Tuple{typeof(==),PeriodicVertex{i},PeriodicVertex{i}}) @enforce Base.precompile(Tuple{typeof(<),PeriodicVertex{i},PeriodicVertex{i}}) @@ -33,7 +37,7 @@ function _precompile_() @enforce Base.precompile(Tuple{typeof(==),PeriodicGraph{i},PeriodicGraph{i}}) @enforce Base.precompile(Tuple{typeof(hash),PeriodicGraph{i},UInt}) - @enforce Base.precompile(Tuple{typeof(Graphs._neighborhood),PeriodicGraph{i},Int,Int,Graphs.DefaultDistance,Function}) + @enforce Base.precompile(Tuple{typeof(Graphs._neighborhood),PeriodicGraph{i},Int,Int,Graphs.DefaultDistance,typeof(outneighbors)}) @enforce Base.precompile(Tuple{typeof(Graphs.add_edge!),PeriodicGraph{i},PeriodicEdge{i}}) @enforce Base.precompile(Tuple{typeof(add_vertex!),PeriodicGraph{i}}) @enforce Base.precompile(Tuple{typeof(Graphs.add_vertices!),PeriodicGraph{i},Int}) diff --git a/test/runtests.jl b/test/runtests.jl index 849d602..b7de3b3 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -7,12 +7,6 @@ import StaticArrays: SVector # using Aqua # Aqua.test_all(PeriodicGraphs) -@static if VERSION < v"1.4-" - function only(t) - length(t) == 1 || throw(ArgumentError("Not exactly 1 element")) - return t[1] - end -end @testset "PeriodicVertex" begin @test_throws MethodError PeriodicVertex2D() # Always require a vertex identifier @@ -22,6 +16,9 @@ end @test PeriodicVertex{4}(5) == PeriodicVertex{4}(5) @test PeriodicVertex1D(1, SVector{1,Int}(3)) == PeriodicVertex1D(1, [3]) @test PeriodicVertex(3, (1,2)) == PeriodicVertex2D(3, (1,2)) == PeriodicVertex(3, SVector{2,Int}(1,2)) + h = hash(PeriodicVertex(3, (1,2))) + @test h == hash(PeriodicVertex(3, SVector{2,Int}(1,2))) + @test hash(PeriodicVertex2D(3)) != h != hash(PeriodicVertex3D(3, (1,2,0))) != hash(PeriodicVertex(2, (1,2))) @test PeriodicVertex[(1, (-1,0))] == [PeriodicVertex2D(1, (-1,0))] == PeriodicVertex2D[(1, (-1,0))] @test_throws MethodError PeriodicVertex(2, [4]) # Performance trap if the dimension is unknown @test_throws DimensionMismatch PeriodicVertex2D(1, (2,)) @@ -63,6 +60,8 @@ end @test PeriodicEdge(1, 1, SVector{1,Int}(2)) == PeriodicEdge1D(1, 1, SVector{1,Int}(2)) @test PeriodicEdge(2, 1, ()) == PeriodicEdge{0}((2, 1, SVector{0,Int}())) @test PeriodicEdge((1, 2, (1,0,0))) == PeriodicEdge3D((1, 2, (1,0,0))) + @test hash(PeriodicEdge((1, 2, (1,0,0)))) == hash(PeriodicEdge3D((1, 2, (1,0,0)))) + @test hash(PeriodicEdge((1, 2, (1,)))) != hash(PeriodicEdge((2, 1, (-1,)))) @test convert(PeriodicEdge, (3, PeriodicVertex{1}(2, (1,)))) == convert(PeriodicEdge1D, (3, PeriodicVertex{1}(2, (1,)))) @test PeriodicEdge1D[(1,2,SVector{1,Int}(3))] == PeriodicEdge[(1,2,(3,))] == [PeriodicEdge(1, 2, (3,))] @static if VERSION < v"1.6.0-DEV" @@ -122,12 +121,53 @@ end @test !has_vertex(g, 1) end + +function naive_periodicgraph(nv::Integer, t::AbstractVector{PeriodicEdge{N}}) where N + sort!(t); unique!(t) + g = PeriodicGraph{N}(nv) + for e in t + add_edge!(g, e) + end + return g +end + +function is_well_formed(g::PeriodicGraph{N}) where N + isempty(g.nlist) || sum(length, g.nlist) == 2*g.ne[] || return false + for (i, l) in enumerate(g.nlist) + if isempty(l) + g.directedgestart[i] == 1 || return false + continue + end + issorted(l) || return false + allunique(l) || return false + any(==(PeriodicVertex{N}(i)), l) && return false + dedgestart = g.directedgestart[i] + if dedgestart > length(l) + dedgestart == length(l) + 1 || return false + PeriodicGraphs.isindirectedge(PeriodicEdge{N}(i, l[end])) || return false + continue + end + PeriodicGraphs.isindirectedge(PeriodicEdge{N}(i, l[dedgestart])) && return false + if dedgestart != 1 && !PeriodicGraphs.isindirectedge(PeriodicEdge{N}(i, l[dedgestart-1])) + return false + end + for j in dedgestart:length(l) + rev = PeriodicVertex{N}(i, .- l[j].ofs) + any(==(rev), g.nlist[l[j].v]) || return false + end + end + g.width[] == -1 || g.width[] > 0 || return false + return true +end + @testset "Basic graph construction and modification" begin @test_throws MethodError PeriodicGraph(2) # missing dimension g = PeriodicGraph3D() + @test is_well_formed(g) @testset "add_vertex!, rem_vertex!" begin @test add_vertex!(g) + @test is_well_formed(g) @test ne(g) == 0 && nv(g) == 1 @test g == PeriodicGraph3D(1) @test add_vertex!(g) @@ -136,6 +176,7 @@ end @test_throws ArgumentError rem_vertex!(g, 0) @test_throws ArgumentError rem_vertex!(g, 3) @test rem_vertex!(g, 1) + @test is_well_formed(g) @test ne(g) == 0 && nv(g) == 1 @test g == PeriodicGraph3D(1) @test rem_vertex!(g, 1) @@ -155,17 +196,20 @@ end @test only(vertices(g)) == 1 @test rem_vertices!(g, [1]) == Int[] @test add_vertices!(g, 5) == 5 + @test is_well_formed(g) @test ne(g) == 0 && nv(g) == 5 && g == PeriodicGraph3D(5) @test vertices(g) == 1:5 @test issetequal(rem_vertices!(g, [2,3]), [1,4,5]) @test ne(g) == 0 && nv(g) == 3 && g == PeriodicGraph3D(3) @test vertices(g) == 1:3 @test rem_vertex!(g, 2) + @test is_well_formed(g) @test ne(g) == 0 && nv(g) == 2 && g == PeriodicGraph3D(2) @test vertices(g) == [1,2] @test rem_vertices!(g, [1,2]) == Int[] @test add_vertices!(g, 5) == 5 @test rem_vertices!(g, [2,4], true) == [1,3,5] + @test is_well_formed(g) @test g == PeriodicGraph3D(3) @test rem_vertices!(g, [1,2,3], true) == Int[] @test g == PeriodicGraph3D() @@ -184,6 +228,7 @@ end @test !has_edge(g, 1, 1) @test !has_edge(g, PeriodicEdge3D(1, 1, (0,1,1))) @test add_edge!(g, 1, PeriodicVertex(1, (-1,0,0))) + @test is_well_formed(g) @test has_edge(g, 1, 1) @test !add_edge!(g, PeriodicEdge(1, 1, (-1,0,0))) @test has_edge(g, 1, 1) @@ -192,6 +237,7 @@ end @test g == PeriodicGraph3D(PeriodicEdge3D[(1, 1, (1,0,0))]) @test g != PeriodicGraph3D(1) @test rem_edge!(g, PeriodicEdge(1, 1, (-1,0,0))) + @test is_well_formed(g) @test g == PeriodicGraph3D(1) @test !has_edge(g, 1, 1) @test !has_edge(g, 1, PeriodicVertex(1, (1,0,0))) @@ -238,6 +284,52 @@ end @test ne(g) == 0 @test nv(g) == 1 end + + @testset "Multiple edges" begin + edgs1 = PeriodicEdge2D[(1, 1, (0,1)), (1, 1, (-1,0)), (1, 1, (1,0))] + g1 = PeriodicGraph(edgs1) + @test nv(g1) == 1 + @test ne(g1) == 2 + @test is_well_formed(g1) + @test g1 == naive_periodicgraph(1, edgs1) == PeriodicGraph(1, edgs1) + _edgs1 = collect(edges(g1)) + @test g1 == PeriodicGraph(_edgs1) + + edgs2 = PeriodicEdge2D[(3, 1, (0,0)), (2, 3, (0,0)), (1, 1, (2,0)), (1, 3, (0,1))] + g2 = PeriodicGraph(edgs2) + @test nv(g2) == 3 + @test ne(g2) == 4 + @test is_well_formed(g2) + @test g2 == naive_periodicgraph(3, edgs2) == PeriodicGraph(3, edgs2) + _edgs2 = collect(edges(g2)) + @test g2 == PeriodicGraph(_edgs2) + sort!(_edgs2) + @test g2 == PeriodicGraphs.from_edges(_edgs2) == PeriodicGraphs.from_edges(3, _edgs2) + + edgs3 = PeriodicEdge2D[(1, 1, (0,1)), (1, 1, (1,0)), (1, 1, (2,0)), (1, 2, (0,1)), + (1, 3, (-1,0)), (1, 3, (0,1)), (2, 3, (0,0)), (3, 4, (0,0))] + g3 = PeriodicGraph(edgs3) + @test nv(g3) == 4 + @test ne(g3) == 8 + @test is_well_formed(g3) + @test g3 == naive_periodicgraph(4, edgs3) == PeriodicGraph(4, edgs3) + _edgs3 = collect(edges(g3)) + @test g3 == PeriodicGraph(_edgs3) + sort!(_edgs3) + @test g3 == PeriodicGraphs.from_edges(_edgs3) == PeriodicGraphs.from_edges(4, _edgs3) + + edgs4 = PeriodicEdge3D[(2, 3, (-1,0,0)), (2, 3, (0,1,0)), (4, 4, (1,-1,0)), + (1, 3, (0,0,1)), (3, 2, (0,0,0))] + g4 = PeriodicGraph(edgs4) + @test nv(g4) == 4 + @test ne(g4) == 5 + @test is_well_formed(g4) + @test g4 == naive_periodicgraph(4, edgs4) == PeriodicGraph(4, edgs4) + _edgs4 = collect(edges(g4)) + @test g4 == PeriodicGraph(_edgs4) + sort!(_edgs4) + @test g4 == PeriodicGraphs.from_edges(_edgs4) == PeriodicGraphs.from_edges(4, _edgs4) + end end @testset "String graph construction" begin @@ -257,7 +349,12 @@ end @test_throws ArgumentError PeriodicGraph{4}("4 1") @test_throws ArgumentError PeriodicGraph("2 1 1 0 1 0") @test_throws LoopException PeriodicGraph("3 1 1 0 0 0") - @test string(PeriodicGraph("2 1 2 0 0 2 3 1 0")) == "2 1 2 0 0 2 3 1 0" + g = PeriodicGraph("2 1 2 0 0 2 3 1 0") + sg = string(g) + @test sg == "2 1 2 0 0 2 3 1 0" + @test parse(PeriodicGraph, sg) == g + @test parse(PeriodicGraph2D, sg) == g + @test PeriodicGraph(sg) == g @test string(PeriodicGraph("1 23 4 5 10 9 43")) == "1 4 23 -5 9 10 -43" @test string(PeriodicGraph{5}()) == "5" end @@ -272,8 +369,10 @@ end g2 = PeriodicGraph(g1) g3 = deepcopy(g1) @test g1 == g2 == g3 + @test hash(g1) == hash(g2) == hash(g3) rem_vertex!(g1, 1) @test g1 != g2 == g3 + @test hash(g1) != hash(g2) == hash(g3) @test add_edge!(g2, PeriodicEdge(1, 2, (-1,0))) @test add_edge!(g2, 1, PeriodicVertex(1, (1,1))) @test add_edge!(g2, PeriodicEdge(1, 3, (0,0))) @@ -310,6 +409,7 @@ end @test_throws ArgumentError offset_representatives!(gg, [2]) gcopy = deepcopy(g) @test offset_representatives!(g, [0,0,0,0]) == gcopy == g + @test is_well_formed(g) ggg = offset_representatives!(gg, [(1,-1,1),0]) @test ggg == gg @test collect(edges(gg)) == PeriodicEdge3D[(1, 2, (1,-1,0)), (2, 2, (1,0,0))] @@ -318,6 +418,7 @@ end @test gg == gcopy @test_throws DimensionMismatch swap_axes!(g, [1,3]) ggg = swap_axes!(gg, (1,2,3)) + @test is_well_formed(ggg) @test gcopy == gg == ggg ggg = swap_axes!(gg, (2,1,3)) @test gcopy != gg == ggg @@ -330,6 +431,7 @@ end g = PeriodicGraph3D(1) @test add_vertices!(g, 2) == 2 @test add_edge!(g, PeriodicEdge3D(2, 3, (0, 0, 1))) + @test is_well_formed(g) @test only(neighbors(g, 2)) == PeriodicVertex3D(3, (0, 0, 1)) @test only(neighbors(g, 3)) == PeriodicVertex3D(2, (0, 0, -1)) @test !has_edge(g, 2, 2) @@ -467,8 +569,10 @@ end PeriodicEdge3D(4, 4, (1, 1, 0))]) gg = deepcopy(g) @test rem_vertices!(g, Int[]) == [1, 2, 3, 4] + @test is_well_formed(g) @test g == gg @test rem_vertices!(g, [2, 4], false) == [1, 3] + @test is_well_formed(g) @test nv(g) == 2 @test ne(g) == 1 @test neighbors(g, 1) == PeriodicVertex3D[(2, (0, 0, 1))] @@ -488,6 +592,7 @@ end 2 5 0 0 4 5 1 1") g2 = deepcopy(g1) + @test g2 == parse(PeriodicGraph, string(g1)) @test rem_vertices!(g1, [2,4], false) == [1, 5, 3] @test g1 == PeriodicGraph("2 1 1 1 0 1 3 0 0") @test rem_vertices!(g2, [2,4], true) == [1, 3, 5]