Basic types and definitions
Introduction
A PeriodicGraph{N}
g
is the representation of an N
-periodic graph.
Each vertex has a unique representative, indexed from 1 to n = nv(g)
. Each vertex x
of the graph is represented by a PeriodicVertex{N}
containing both the representative v
and the offset o
between the unit cell containing the vertex and a reference unit cell, accessible through the syntax v, o = x
, or with v = first(x)
and o = last(x)
. The offset is a N
-uplet of integers, stored as a SVector{N,Int}
.
Each edge e
is represented by a PeriodicEdge{N}
defined by its source vertex inside the reference unit cell and of representative s
, and its destination vertex d
, a PeriodicVertex{N}
, which can be accessed through the syntax s, d = e
, or with s = first(e)
and d = last(e)
(or d = Graphs.dst(e)
).
For convenience, aliases are exported for 1D, 2D and 3D (N = 1
, N = 2
and N = 3
) under the names PeriodicGraph1D
, PeriodicEdge2D
, PeriodicVertex3D
, etc.
PeriodicVertex
A PeriodicVertex
can be built like so:
julia> PeriodicVertex1D(5, (-1,))
PeriodicVertex1D(5, (-1,))
julia> PeriodicVertex{4}(2) # shorthand for the vertices of the reference cell
PeriodicVertex{4}(2, (0,0,0,0))
PeriodicEdge
A PeriodicEdge
, can be defined from src
and dst
or, equivalently, be the identifiers of both source and destination vertices, and the cell offset between source and destination. For example, in 3D, the edge between vertex 1
and vertex 4
in cell (0, 1, 0)
is PeriodicEdge3D(1, PeriodicVertex3D(4, (0,1,0)))
, or, equivalently, PeriodicEdge3D(1, 4, (0,1,0))
. Since PeriodicEdge(u, v, ofs)
and PeriodicEdge(v, u, .-ofs)
represent the same edge, one of them is called the direct edge, when it has either u < v
or u == v && ofs > zero(ofs)
, and the other is the indirect edge. Functions isdirectedge
and directedge
are used for this.
More examples:
julia> e = PeriodicEdge(3, PeriodicVertex(2, (0,1,0,0))) # dimensions are inferred from the input
PeriodicEdge{4}(3, 2, (0,1,0,0))
julia> isdirectedge(e)
false
julia> directedge(e)
PeriodicEdge{4}(2, 3, (0,-1,0,0))
julia> src, (dstv, ofs) = PeriodicEdge3D(5, 6, (1,0,2));
julia> src
5
julia> dstv
6
julia> ofs
3-element StaticArrays.SVector{3, Int64} with indices SOneTo(3):
1
0
2
julia> PeriodicEdge2D(5, 5, (0,0))
ERROR: LoopException: a loop from vertex 5 to itself in the same unit cell is a forbidden edges. Maybe the offset is wrong?
Note that loops (that is, edges of the form (u, u, (0,0,0,...,0))
) are forbidden and will throw a PeriodicGraphs.LoopException
if created. To bypass this check, use the unexported PeriodicGraphs.unsafe_edge
function.
PeriodicGraph
Finally, N
-periodic graphs, represented by the type PeriodicGraph{N}
, are defined by the number of vertices in the reference cell and the set of edges starting from the reference cell. When the number of vertices is simply the highest number appearing in the list of edges, it can be omitted. Periodic graphs can be built through several methods:
julia> PeriodicGraph{5}() # create the empty graph
PeriodicGraph{5}(0, PeriodicEdge{5}[])
julia> PeriodicGraph1D(4) # create a graph with 4 vertices but no edge
PeriodicGraph1D(4, PeriodicEdge1D[])
julia> PeriodicGraph(2, PeriodicEdge{4}[(1, 1, (0,0,1,1)), (1, 1, (0,1,0,-1))]) # the dimension can be inferred
PeriodicGraph{4}(2, PeriodicEdge{4}[(1, 1, (0,0,1,1)), (1, 1, (0,1,0,-1))])
julia> PeriodicGraph3D(PeriodicEdge3D[(1, 3, (0,1,0)), (2, 2, (0,0,-1)), (1, 2, (1,0,0)), (2, 3, (0,1,1))])
PeriodicGraph3D(3, PeriodicEdge3D[(1, 2, (1,0,0)), (1, 3, (0,1,0)), (2, 2, (0,0,1)), (2, 3, (0,1,1))])
julia> parse(PeriodicGraph3D, "3 1 2 1 0 0 1 3 0 1 0 2 2 0 0 1 2 3 0 1 1") # compact representation of the previous graph
PeriodicGraph3D(3, PeriodicEdge3D[(1, 2, (1,0,0)), (1, 3, (0,1,0)), (2, 2, (0,0,1)), (2, 3, (0,1,1))])
julia> string(ans) # to obtain the compact representation from the graph
"3 1 2 1 0 0 1 3 0 1 0 2 2 0 0 1 2 3 0 1 1"
julia> string(PeriodicGraph("2 1 2 0 0 2 3 0 0 3 1 0 0 3 1 1 0 2 1 0 1 2 3 0 1"))
"2 1 2 0 -1 1 2 0 0 1 3 -1 0 1 3 0 0 2 3 0 0 2 3 0 1"
API
Type definitions
PeriodicGraphs.PeriodicVertex
— TypePeriodicVertex{N}
Vertex type for an N
-periodic graph.
A vertex is uniquely determined by the identifier of its representative in the a fixed initial cell, and the offset of the cell containing the vertex compared to the the initial cell. Vertex identifiers start at 1.
PeriodicGraphs.PeriodicEdge
— TypePeriodicEdge{N} <: Graphs.SimpleGraphs.AbstractSimpleEdge{Int}
Edge type for an N
-periodic graph.
An edge is uniquely determined by the vertex identifiers of its source and destination, and the cell offset between the source vertex and the destination vertex.
PeriodicGraphs.PeriodicGraph
— TypePeriodicGraph{N} <: AbstractGraph{Int}
Type representing an undirected N
-periodic graph.
PeriodicGraphs.PeriodicGraph
— MethodPeriodicGraph{N}(nv::Integer=0)
Construct a PeriodicGraph{N}
with nv
vertices and 0 edge.
PeriodicGraphs.PeriodicGraph
— MethodPeriodicGraph([nv::Integer, ]edge_list::AbstractVector{PeriodicEdge{N}})
PeriodicGraph{N}([nv::Integer, ]edge_list::AbstractVector{PeriodicEdge{N}})
Construct a PeriodicGraph{N}
from a vector of edges. If nv
is unspecified, the number of vertices is the highest that is used in an edge in edge_list
.
Implementation Notes
This constructor works the fastest when edge_list
is sorted by the lexical ordering and does not contain any duplicates.
Examples
julia> el = PeriodicEdge2D[(1, 1, (1,0)), (1, 3, (0,1)), (3, 1, (0,-1))];
julia> g = PeriodicGraph(el)
PeriodicGraph2D(3, PeriodicEdge2D[(1, 1, (1,0)), (1, 3, (0,1))])
julia> ne(g)
2
PeriodicGraphs.PeriodicGraph
— MethodPeriodicGraph(key::AbstractString)
PeriodicGraph{N}(key::AbstractString)
Construct a PeriodicGraph{N}
from a key
, which is a string of whitespace-separated values of the form "N src1 dst1 ofs1_1 ofs1_2 ... ofs1_N src2 dst2 ofs2_1 ofs2_2 ... ofs2_N ... srcm dstm ofsm_1 ofsm_2 ... ofsm_N"
where N
is the number of repeating dimensions of the graph, m
is the number of edges and for all i
between 1
and m
, the number of edges, the i
-th edge is described as
srci
, the vertex identifier of the source vertex,dsti
, the vertex identifier of the destination vertex and(ofsi_1, ofsi_2, ..., ofsi_N)
the offset of the edge.
This compact representation of a graph can be obtained simply by print
ing the graph or with string
.
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
julia> PeriodicGraph("2 1 2 0 0 2 1 1 0 1 1 0 -1")
PeriodicGraph2D(2, PeriodicEdge2D[(1, 1, (0,1)), (1, 2, (-1,0)), (1, 2, (0,0))])
julia> PeriodicGraph3D("3 1 1 0 0 1 1 1 0 1 0 1 1 1 0 0")
PeriodicGraph3D(1, PeriodicEdge3D[(1, 1, (0,0,1)), (1, 1, (0,1,0)), (1, 1, (1,0,0))])
julia> string(ans)
"3 1 1 0 0 1 1 1 0 1 0 1 1 1 0 0"
julia> string(parse(PeriodicGraph3D, ans)) == ans
true
Base.parse
— Methodparse(::Type{PeriodicGraph}, key::AbstractString)
parse(::Type{PeriodicGraph{N}}, key::AbstractString)
Parse a string representation of a PeriodicGraph
back to a PeriodicGraph
.
See PeriodicGraph(key::AbstractString)
or PeriodicGraph{N}(key::AbstractString)
for details on the string representations of a PeriodicGraph
.
This function assumes that the string is a valid representation of a PeriodicGraph
with all its edges direct, sorted in lexicographical order and unique, as obtained from string(g)
or print(g)
where g
is a PeriodicGraph
. No check is performed to ensure this condition.
If the input string may not obey these conditions, use PeriodicGraph(key)
or PeriodicGraph{N}(key)
instead.
Other basic functions
PeriodicGraphs.isdirectedge
— Functionisdirectedge(e::PeriodicEdge)
Return true
if e
is direct, in the sense being of the form (u, v, ofs)
with either u < v
or u == v && ofs > zero(ofs)
.
An edge e
is indirect iff reverse(e)
is not.
Examples
julia> isdirectedge(PeriodicEdge1D(3, 4, (0,)))
true
julia> isdirectedge(PeriodicEdge2D(5, 2, (0,0)))
false
julia> isdirectedge(PeriodicEdge3D(3, 3, (0,-1,2)))
false
See also directedge
PeriodicGraphs.directedge
— Functiondirectedge(e::PeriodicEdge{D}) where D
directedge(src::PeriodicVertex{D}, dst::PeriodicVertex{D}) where D
directedge(src, dst, ofs::Union{SVector{D,T},NTuple{D,T}}) where {D,T}
Return the direct edge corresponding to e = PeriodicEdge{D}(src, dst)
, i.e. e
itself if e
is direct, or reverse(e)
otherwise.
Examples
julia> directedge(PeriodicEdge1D(3, 4, (0,)))
PeriodicEdge1D(3, 4, (0,))
julia> directedge(PeriodicEdge2D(5, 2, (0,0)))
PeriodicEdge2D(2, 5, (0,0))
julia> directedge(PeriodicEdge3D(3, 3, (0,-1,2)))
PeriodicEdge3D(3, 3, (0,1,-2))
See also isdirectedge
PeriodicGraphs.unsafe_edge
— Typeunsafe_edge{N}
Internal constructor for PeriodicEdge{N}
that bypasses the loop check.
PeriodicGraphs.LoopException
— TypeLoopException <: Exception
Error type for constructing an invalid PeriodicEdge{N}
of the form (u, u, zeros(Int,N))
. Loops are not expected in the algorithms implemented in PeriodicGraphs.jl. If you still want to construct them, use the unsafe_edge{N}
constructor instead of PeriodicEdge{N}
.