API


Types

Brillouin.WignerSeitz.CellType
struct Cell{D} <: AbstractArray{Array{StaticArraysCore.SArray{Tuple{D}, Float64, 1, D}, 1}, 1}
  • verts::Array{StaticArraysCore.SVector{D, Float64}, 1} where D

  • faces::Vector{Vector{Int64}}

  • basis::StaticArraysCore.SArray{Tuple{D}, StaticArraysCore.SVector{D, Float64}, 1, D} where D

  • setting::Base.RefValue{Brillouin.BasisEnum}

source
Brillouin.KPaths.KPathType
struct KPath{D} <: Brillouin.KPaths.AbstractPath{Pair{Symbol, StaticArraysCore.SArray{Tuple{D}, Float64, 1, D}}}
  • points::Dict{Symbol, StaticArraysCore.SVector{D, Float64}} where D

  • paths::Vector{Vector{Symbol}}

  • basis::Bravais.ReciprocalBasis

  • setting::Base.RefValue{Brillouin.BasisEnum}

source
Brillouin.KPaths.KPathInterpolantType
struct KPathInterpolant{D} <: Brillouin.KPaths.AbstractPath{StaticArraysCore.SArray{Tuple{D}, Float64, 1, D}}
  • kpaths::Array{Array{StaticArraysCore.SVector{D, Float64}, 1}, 1} where D

  • labels::Vector{Dict{Int64, Symbol}}

  • basis::Bravais.ReciprocalBasis

  • setting::Base.RefValue{Brillouin.BasisEnum}

source

Exported methods

Brillouin.basisMethod
basis(x::Union{KPath, KPathInterpolant, Cell})

Return the (reciprocal or direct) lattice basis associated with x, in Cartesian coordinates.

Methods in Brillouin by default return points in the lattice basis, i.e., points are referred to basis(x). This corresponds to the setting(x) == LATTICE. Coordinates may instead be referred to a Cartesian basis, corresponding to setting(x) == CARTESIAN by using cartesianize. The result of basis(x), however, is invariant to this and always refers to the lattice basis in Cartesian coordinates.

source
Brillouin.settingMethod
setting(x::Union{KPath, KPathInterpolant, Cell})

Return the basis setting of coordinates in x. The returned value is a member of the BasisEnum enum with member values LATTICE (i.e. coordinates in the basis of the lattice vectors) or CARTESIAN (i.e. coordinates in the Cartesian basis). By default, methods in Brillouin return coordinates in the LATTICE setting.

source
Brillouin.WignerSeitz.reduce_to_wignerseitzMethod
reduce_to_wignerseitz(v::StaticVector, Vs::BasisLike)  -->  v′

Return the periodic image v′ of the point v in the basis Vs.

v is assumed to be provided in the lattice basis (i.e., relative to Vs) and v′ is returned in similar fashion.

The returned point v′ lies in the Wigner-Seitz cell (or its boundary) defined by Vs, has the least possible norm among all equivalent images of v, and differs from v at most by integer lattice translations such that mod(v, 1) ≈ mod(v′, 1).

source
Brillouin.WignerSeitz.wignerseitzMethod
wignerseitz(basis::AbstractVector{<:SVector{D}}; merge::Bool = true, Nmax = 3)
wignerseitz(basis::AbstractVector{<:AbstractVector}; merge::Bool = true, Nmax = 3)
                                                            --> Cell{D}

Given a provided basis, return a Cell{D} containing the vertices and associated (outward oriented) faces of the Wigner-Seitz cell defined by basis in D dimensions. The returned vertices are given in the the basis of basis (see cartesianize! for conversion).

Keyword arguments

  • merge (default, true): if true, co-planar faces are merged to form polygonal planar faces (e.g., triangles, quadrilaterals, and ngons generally). If false, raw "unprocessed" triangles (D=3) and segments (D=2) are returned instead. merge has no impact for D=1.
  • Nmax (default, 3): includes -Nmax:Nmax points in the initial lattice used to generate the underlying Voronoi tesselation. It is unwise to set this to anything lower than 3 without explicitly testing convergence; and probably unnecessary to increase it beyond 3 as well.
source
Bravais.cartesianize!Method
cartesianize!(kpi::KPathInterpolant)

In-place transform an interpolated k-path kpi in a lattice basis to a Cartesian basis with (primitive) reciprocal lattice vectors basis(kpi).

source
Bravais.cartesianize!Method
cartesianize!(kp::KPath{D})

Transform a k-path kp in a lattice basis to a Cartesian basis with (primitive) reciprocal basis vectors basis(kp). Modifies kp in-place.

source
Bravais.cartesianizeMethod
cartesianize(kp::KPath{D})

Transform a k-path kp in a lattice basis to a Cartesian basis with (primitive) reciprocal basis vectors basis(kp).

source
Bravais.latticize!Method
latticize!(kpi::KPathInterpolant{D})

In-place transform an interpolated k-path kpi in a Cartesian basis to a lattice basis with (primitive) reciprocal lattice vectors basis(kpi).

source
Bravais.latticize!Method
latticize!(kp::KPath)

Transform a k-path kp in a Cartesian basis to a lattice basis with (primitive) reciprocal lattice vectors basis(kp). Modifies kp in-place.

source
Bravais.latticizeMethod
latticize(kp::KPath)

Transform a k-path kp in a Cartesian basis to a lattice basis with (primitive) reciprocal lattice vectors basis(kp).

source
Bravais.latticizeMethod
latticize(kpi::KPathInterpolant{D}, basis::AbstractVector{<:AbstractVector{<:Real}})

Transform an interpolated k-path kpi in a Cartesian basis to a lattice basis with (primitive) reciprocal lattice vectors basis.

If kpi is not in a Cartesian basis (i.e., if setting(kpi) == LATTICE), kpi is returned as-is.

source
Bravais.latticizeMethod
latticize(kp::KPath{D}, basis::AbstractVector{<:AbstractVector{<:Real}})

Transform a k-path kp in a Cartesian basis to a lattice basis with (primitive) reciprocal lattice vectors basis.

If kp is not in a Cartesian basis (i.e., if setting(kp) == LATTICE), kp is returned as-is.

source
Brillouin.KPaths.interpolateMethod
interpolate(kvs::AbstractVector{<:AbstractVector{<:Real}}, N::Integer)
                                                --> Vector{Vector{<:Real}}

Return an interpolated k-path between discrete k-points in kvs, with approximately N interpolation points in total (typically fewer).

Note that, in general, it is not possible to do this so that all interpolated k-points are equidistant; samples are however exactly equidistant across each linear segment defined by points in kvs and approximately equidistant across all segments.

See also interpolate(::KPath, ::Integer) and splice.

Future deprecation

This method signature is likely to be deprecated in future versions of Brillouin.jl.

source
Brillouin.KPaths.interpolateMethod
interpolate(kp::KPath, N::Integer)
interpolate(kp::KPath; N::Integer, density::Real) --> KPathInterpolant

Return an interpolant of kp with N points distributed approximately equidistantly across the full k-path (equidistance is measured in a Cartesian metric).

Note that the interpolant may contain slightly fewer or more points than N (typically fewer) in order to improve equidistance. N can also be provided as a keyword argument.

As an alternative to specifying the desired total number of interpolate points via N, a desired density per unit (reciprocal) length can be specified via the keyword argument density.

source
Brillouin.KPaths.irrfbz_pathMethod
irrfbz_path(sgnum::Integer, Rs, [::Union{Val(D), Integer},]=Val(3))  -->  ::KPath{D}

Returns a k-path (::KPath) in the (primitive) irreducible Brillouin zone for a space group with number sgnum, (conventional) direct lattice vectors Rs, and dimension D. The path includes all distinct high-symmetry lines and points as well as relevant parts of the Brillouin zone boundary.

The dimension D (1, 2, or 3) is specified as the third input argument, preferably as a static Val{D} type parameter (or, type-unstably, as an <:Integer). Defaults to Val(3).

Rs refers to the direct basis of the conventional unit cell, i.e., not the primitive direct basis vectors. The setting of Rs must agree with the conventional setting choices in the International Tables of Crystallography, Volume A (the "ITA conventional setting"). If Rs is a subtype of a StaticVector or NTuple, the dimension can be inferred from its (static) size; in this case, this dimension will take precedence (i.e. override, if different) over any dimension specified in the third input argument.

Notes

  • The returned k-points are given in the basis of the primitive reciprocal basis in the CDML setting. To obtain the associated transformation matrices between the ITA conventional setting and the CDML primitive setting, see primitivebasismatrix of Bravais.jl](https://thchr.github.io/Crystalline.jl/stable/bravais/) (or, equivalently, the relations defined Table 2 of [1]). To transform to a Cartesian basis, see cartesianize!.
  • To interpolate a KPath, see interpolate(::KPath, ::Integer) and splice(::KPath, ::Integer).
  • All paths currently assume time-reversal symmetry (or, equivalently, inversion symmetry). If neither are present, include the "inverted" -k paths manually.

Data and referencing

3D paths are sourced from the SeeK-path publication: please cite the original work [2].

References

source
Brillouin.KPaths.pathsMethod
paths(kp::KPath) -> Vector{Vector{Symbol}}

Return a vector of vectors, with each vector describing a connected path between between k-points referenced in kp (see also points(::KPath)).

source
Brillouin.KPaths.pointsMethod
points(kp::KPath{D}) -> Dict{Symbol, SVector{D,Float64}}

Return a dictionary of the k-points (values) and associated k-labels (keys) referenced in kp.

source
Brillouin.KPaths.spliceMethod
splice(kvs::AbstractVector{<:AbstractVector{<:Real}}, N::Integer)
                                                --> Vector{Vector{<:Real}}

Return an interpolated k-path between the discrete k-points in kvs, with N interpolation points inserted in each segment defined by pairs of adjacent k-points.

See also splice(::KPath, ::Integer) and interpolate.

Future deprecation

This method signature is likely to be deprecated in future versions of Brillouin.jl.

source
Brillouin.KPaths.spliceMethod
splice(kp::KPath, N::Integer) --> KPathInterpolant

Return an interpolant of kp with N points inserted into each k-path segment of kp.

source