Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

lattice(lat; kw...) #95

Merged
merged 3 commits into from
Sep 21, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 6 additions & 3 deletions src/convert.jl
Original file line number Diff line number Diff line change
Expand Up @@ -14,13 +14,16 @@ Base.convert(::Type{T}, l::Bravais) where T<:Bravais = T(l)

# Constructors for conversion

Sublat{E,T}(s::Sublat, name = s.name) where {E,T} =
Sublat{E,T,V}(s::Sublat, name = s.name) where {E,T,V<:Vector} =
Sublat([padright(site, zero(T), Val(E)) for site in s.sites], name)

# We need this to promote different sublats into common dimensionality and type to combine
# into a lattice, while neglecting orbital dimension
Base.promote(ss::Sublat{E,T}...) where {E,T} = ss
Base.promote_rule(::Type{Sublat{E1,T1}}, ::Type{Sublat{E2,T2}}) where {E1,E2,T1,T2} =
Sublat{max(E1, E2), promote_type(T1, T2)}
function Base.promote_rule(::Type{Sublat{E1,T1,Vector{SVector{E1,T1}}}}, ::Type{Sublat{E2,T2,Vector{SVector{E2,T2}}}}) where {E1,E2,T1,T2}
E´ = max(E1, E2)
T´ = promote_type(T1, T2)
return Sublat{E´, T´, Vector{SVector{E´,T´}}}
end

Bravais{E,L,T}(b::Bravais) where {E,L,T} = Bravais(padtotype(b.matrix, SMatrix{E,L,T}))
2 changes: 1 addition & 1 deletion src/hamiltonian.jl
Original file line number Diff line number Diff line change
Expand Up @@ -263,7 +263,7 @@ julia> h[(3,3)][[1,2],[1,2]] .= Ref(@SMatrix[1 2; 2 1])
# See also:
`onsite`, `hopping`, `bloch`, `bloch!`
"""
hamiltonian(lat, ts...; orbitals = missing, kw...) =
hamiltonian(lat::AbstractLattice, ts...; orbitals = missing, kw...) =
_hamiltonian(lat, sanitize_orbs(orbitals, lat.unitcell.names), ts...; kw...)
_hamiltonian(lat::AbstractLattice, orbs; kw...) = _hamiltonian(lat, orbs, TightbindingModel(); kw...)
_hamiltonian(lat::AbstractLattice, orbs, m::TightbindingModel; type::Type = Complex{numbertype(lat)}, kw...) =
Expand Down
66 changes: 43 additions & 23 deletions src/lattice.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,8 +3,8 @@ abstract type AbstractLattice{E,L,T} end
#######################################################################
# Sublat (sublattice)
#######################################################################
struct Sublat{E,T}
sites::Vector{SVector{E,T}}
struct Sublat{E,T,V<:AbstractVector{SVector{E,T}}}
sites::V
name::NameType
end

Expand Down Expand Up @@ -39,12 +39,18 @@ Sublat{2,Float64} : sublattice of Float64-typed sites in 2D space
Name : :A
```
"""
sublat(sites::Vector{<:SVector}; name = :_, kw...) =
sublat(sites::AbstractVector{<:SVector}; name = :_, kw...) =
Sublat(sites, nametype(name))
sublat(vs::Union{Tuple,AbstractVector{<:Number}}...; kw...) = sublat(toSVectors(vs...); kw...)

toSVectors(vs...) = [promote(toSVector.(vs)...)...]

dims(::NTuple{N,Sublat{E,T}}) where {N,E,T} = E

sublatnames(ss::NTuple{N,Sublat{E,T}}) where {N,E,T} = (s -> s.name).(ss)

Base.eltype(::NTuple{N,Sublat{E,T}}) where {N,E,T} = T

#######################################################################
# Unitcell
#######################################################################
Expand All @@ -54,28 +60,33 @@ struct Unitcell{E,T,N}
offsets::Vector{Int} # Linear site number offsets for each sublat
end # so that diff(offset) == sublatsites


Unitcell(sublats::Sublat...; kw...) = Unitcell(promote(sublats...); kw...)
Unitcell(sublats::NTuple{N,Sublat{E,T}};
dim = Val(E), type = float(T), names = [s.name for s in sublats], kw...) where {N,E,T} =
_Unitcell(sublats, dim, type, names)

Unitcell(s; dim = Val(dims(s)), type = float(eltype(s)), names = sublatnames(s)) =
_unitcell(s, dim, type, names)

# Dynamic dispatch
_Unitcell(sublats, dim::Integer, type, names) = _Unitcell(sublats, Val(dim), type, names)
_unitcell(sublats, dim::Integer, type, names) = _unitcell(sublats, Val(dim), type, names)

function _Unitcell(sublats::NTuple{N,Sublat}, dim::Val{E}, type::Type{T}, names) where {N,E,T}
function _unitcell(sublats::NTuple{N,Sublat}, dim::Val{E}, type::Type{T}, names) where {N,E,T}
sites = SVector{E,T}[]
offsets = [0] # length(offsets) == length(sublats) + 1
for s in eachindex(sublats)
for site in sublats[s].sites
push!(sites, padright(site, Val(E)))
push!(sites, padtotype(site, SVector{E,T}))
end
push!(offsets, length(sites))
end
names´ = uniquenames(sanitize_names(names, Val(N)))
return Unitcell(sites, names´, offsets)
end

_unitcell(u::Unitcell{E,T,N}, dim::Val{E}, type::Type{T}, names) where {E,T,N} =
Unitcell(u.sites, uniquenames(sanitize_names(names, Val(N))), u.offsets)

_unitcell(u::Unitcell{E,T,N}, dim::Val{E2}, type::Type{T2}, names) where {E,T,E2,T2,N} =
Unitcell(padtotype.(u.sites, SVector{E2,T2}), uniquenames(sanitize_names(names, Val(N))), u.offsets)

sanitize_names(name::Union{NameType,Int}, ::Val{N}) where {N} = ntuple(_ -> NameType(name), Val(N))
sanitize_names(names::AbstractVector, ::Val{N}) where {N} = ntuple(i -> NameType(names[i]), Val(N))
sanitize_names(names::NTuple{N,Union{NameType,Int}}, ::Val{N}) where {N} = NameType.(names)
Expand Down Expand Up @@ -116,15 +127,6 @@ function sublat(u::Unitcell, siteidx)
return l
end

# function boundingbox(u::Unitcell{E,T}) where {E,T}
# min´ = max´ = first(u.sites)
# for r in u.sites
# min´ = min.(min´, r)
# max´ = max.(max´, r)
# end
# return min´, max´
# end

sublatsites(u::Unitcell) = diff(u.offsets)

nsublats(u::Unitcell) = length(u.names)
Expand All @@ -133,8 +135,14 @@ sublats(u::Unitcell) = 1:nsublats(u)

sublatname(u::Unitcell, s) = u.names[s]

sublatnames(u::Unitcell) = u.names

transform!(f::Function, u::Unitcell) = (u.sites .= f.(u.sites); u)

dims(::Unitcell{E}) where {E} = E

Base.eltype(::Unitcell{E,T}) where {E,T} = T

Base.copy(u::Unitcell) = Unitcell(copy(u.sites), u.names, copy(u.offsets))

Base.isequal(u1::Unitcell, u2::Unitcell) =
Expand Down Expand Up @@ -224,6 +232,12 @@ corresponds to a bounded lattice with no Bravais vectors.
A keyword `names` can be used to rename `sublats`. Given names can be replaced to ensure
that all sublattice names are unique.

lattice(lat::AbstractLattice; bravais = missing, dim = missing, type = missing, names = missing)

Create a new lattice by applying any non-missing `kw` to `lat`. For performance, allocations
will be avoided if possible (depends on `kw`), so the result can share memory of `lat`. To
avoid that, do `lattice(copy(lat); kw...)`.

See also `LatticePresets` for built-in lattices.

# Examples
Expand Down Expand Up @@ -255,17 +269,23 @@ Lattice{3,2,Float64} : 2D lattice in 3D space
`LatticePresets`, `bravais`, `sublat`, `supercell`, `intracell`
"""
function lattice(ss::Sublat...; bravais = (), kw...)
ucell = Unitcell(ss...; kw...)
br = Bravais(bravais, ucell)
return lattice(br, ucell)
u = Unitcell(ss...; kw...)
b = Bravais(bravais, u)
return lattice(u, b)
end

function lattice(bravais::B, unitcell::U) where {E2,L2,E,T,B<:Bravais{E2,L2},U<:Unitcell{E,T}}
function lattice(unitcell::U, bravais::B) where {E2,L2,E,T,B<:Bravais{E2,L2},U<:Unitcell{E,T}}
L = min(E,L2) # L should not exceed E
bravais´ = convert(Bravais{E,L,T}, bravais)
return Lattice(bravais´, unitcell)
end

function lattice(lat::Lattice; bravais = bravais(lat), kw...)
u = Unitcell(lat.unitcell; kw...)
b = Bravais(bravais, u)
return lattice(u, b)
end

#######################################################################
# Supercell
#######################################################################
Expand Down
17 changes: 16 additions & 1 deletion test/test_lattice.jl
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
using Quantica: nsites, Sublat, Bravais, Lattice, Superlattice
using Quantica: nsites, Sublat, Bravais, Lattice, Superlattice, allsitepositions
using Random
using LinearAlgebra: I

Expand All @@ -19,6 +19,21 @@ end
@test lattice(s; bravais = br, type = t, dim = Val(e)) isa Lattice{e,min(l,e),t}
@test lattice(s; bravais = br, type = t, dim = e) isa Lattice{e,min(l,e),t}
end
lat = lattice(sublat((0,0,0)), sublat((1,1,1f0)); bravais = SMatrix{3,3}(I))
lat2 = lattice(lat, bravais = ())
@test lat2 isa Lattice{3,0}
@test allsitepositions(lat) === allsitepositions(lat2)
lat2 = lattice(lat, bravais = (), names = :A)
@test lat2 isa Lattice{3,0}
@test allsitepositions(lat) === allsitepositions(lat2)
lat2 = lattice(lat, dim = Val(2))
@test lat2 isa Lattice{2,2} # must be L <= E
@test allsitepositions(lat) !== allsitepositions(lat2)
lat2 = lattice(lat, type = Float64)
@test lat2 isa Lattice{3,3}
@test allsitepositions(lat) !== allsitepositions(lat2)
lat2 = lattice(lat, dim = Val(2), bravais = SA[1 2; 3 4])
@test bravais(lat2) == SA[1 2; 3 4]
end

@testset "lattice presets" begin
Expand Down