Public API


Exported types

Crystalline.BandRepType
struct BandRep <: AbstractVector{Int64}
  • wyckpos::String

  • sitesym::String

  • label::String

  • dim::Int64

  • decomposable::Bool

  • spinful::Bool

  • irvec::Vector{Int64}

  • irlabs::Vector{String}

source
Crystalline.BandRepSetType
struct BandRepSet <: AbstractVector{BandRep}
  • sgnum::Int64

  • bandreps::Vector{BandRep}

  • kvs::Vector{<:KVec}

  • klabs::Vector{String}

  • irlabs::Vector{String}

  • allpaths::Bool

  • spinful::Bool

  • timereversal::Bool

source
Crystalline.CharacterTableType
struct CharacterTable{D} <: Crystalline.AbstractCharacterTable
  • ops::Array{SymOperation{D}, 1} where D

  • irlabs::Vector{String}

  • table::Matrix{ComplexF64}

  • tag::String

source
Crystalline.KVecMethod
KVec{D}(str::AbstractString) --> KVec{D}
KVec(str::AbstractString)    --> KVec
KVec(::AbstractVector, ::AbstractMatrix) --> KVec

Return a KVec by parsing the string representations str, supplied in one of the following formats:

"($1,$2,$3)"
"[$1,$2,$3]"
"$1,$2,$3"

where the coordinates $1,$2, and $3 are strings that may contain fractions, decimal numbers, and "free" parameters {'α','β','γ'} (or, alternatively and equivalently, {'u','v','w'} or {'x','y','z'}).

Fractions such as 1/2 and decimal numbers can be parsed: but use of any other special operator besides / will produce undefined behavior (e.g. do not use *).

Example

julia> KVec("0.25,α,0")
[1/4, α, 0]
source
Crystalline.LGIrrepType
struct LGIrrep{D} <: Crystalline.AbstractIrrep{D}
  • cdml::String

  • g::LittleGroup

  • matrices::Vector{Matrix{ComplexF64}}

  • translations::Vector{Vector{Float64}}

  • reality::Reality

  • iscorep::Bool

source
Crystalline.LittleGroupType
struct LittleGroup{D} <: Crystalline.AbstractGroup{D, SymOperation{D}}
  • num::Int64

  • kv::KVec

  • klab::String

  • operations::Array{SymOperation{D}, 1} where D

source
Crystalline.MSpaceGroupType
struct MSpaceGroup{D} <: Crystalline.AbstractGroup{D, MSymOperation{D}}
  • num::Tuple{Int64, Int64}

  • operations::Array{MSymOperation{D}, 1} where D

source
Crystalline.ModulatedFourierLatticeType
ModulatedFourierLattice{D} <: AbstractFourierLattice{D}

A D-dimensional concrete Fourier (plane wave) lattice, derived from a UnityFourierLattice{D} by scaling and modulating its orbit coefficients by complex numbers; in general, the coefficients do not have unit norm.

source
Crystalline.MultTableMethod
MultTable(ops::AbstractVector, modτ=true)

Compute the multiplication (or Cayley) table of ops, an iterable of SymOperations.

A MultTable is returned, which contains symmetry operations resulting from composition of row and col operators; the table of indices give the symmetry operators relative to the ordering of ops.

Keyword arguments

  • modτ (default: true): whether composition of operations is taken modulo lattice

vectors (true) or not (false).

source
Crystalline.PGIrrepType
struct PGIrrep{D} <: Crystalline.AbstractIrrep{D}
  • cdml::String

  • g::PointGroup

  • matrices::Vector{Matrix{ComplexF64}}

  • reality::Reality

  • iscorep::Bool

source
Crystalline.PointGroupType
struct PointGroup{D} <: Crystalline.AbstractGroup{D, SymOperation{D}}
  • num::Int64

  • label::String

  • operations::Array{SymOperation{D}, 1} where D

source
Crystalline.RVecMethod
RVec{D}(str::AbstractString) --> RVec{D}
RVec(str::AbstractString)    --> RVec
RVec(::AbstractVector, ::AbstractMatrix) --> RVec

Return a RVec by parsing the string representations str, supplied in one of the following formats:

"($1,$2,$3)"
"[$1,$2,$3]"
"$1,$2,$3"

where the coordinates $1,$2, and $3 are strings that may contain fractions, decimal numbers, and "free" parameters {'α','β','γ'} (or, alternatively and equivalently, {'u','v','w'} or {'x','y','z'}).

Fractions such as 1/2 and decimal numbers can be parsed: but use of any other special operator besides / will produce undefined behavior (e.g. do not use *).

Example

julia> RVec("0.25,α,0")
[1/4, α, 0]
source
Crystalline.SiteGroupType
struct SiteGroup{D} <: Crystalline.AbstractGroup{D, SymOperation{D}}
  • num::Int64

  • wp::WyckoffPosition

  • operations::Array{SymOperation{D}, 1} where D

  • cosets::Array{SymOperation{D}, 1} where D

source
Crystalline.SpaceGroupType
struct SpaceGroup{D} <: Crystalline.AbstractGroup{D, SymOperation{D}}
  • num::Int64

  • operations::Array{SymOperation{D}, 1} where D

source
Crystalline.SubperiodicGroupType
struct SubperiodicGroup{D, P} <: Crystalline.AbstractGroup{D, SymOperation{D}}
  • num::Int64

  • operations::Array{SymOperation{D}, 1} where D

A subperiodic group of embedding dimension D and periodicity dimension P.

Fields:

  • operations: the SymOperations of the finite factor group $G/T$, where $G$ is the

subperiodic group and $T$ is the translation group of the associated lattice.

  • num: the canonical number of the group, following the International Tables for

Crystallography, Volume E.

source
Crystalline.SymOperationType
struct SymOperation{D} <: Crystalline.AbstractOperation{D}
  • rotation::Crystalline.SquareStaticMatrices.SqSMatrix{D, Float64} where D

  • translation::StaticArraysCore.SVector{D, Float64} where D

source

Exported methods

Base.invMethod
inv(op::SymOperation{D}) --> SymOperation{D}

Compute the inverse {W|w}⁻¹ ≡ {W⁻¹|-W⁻¹w} of an operator op ≡ {W|w}.

source
Base.positionMethod
position(x::Union{AbstractGroup, AbstractIrrep})

If a position is associated with x, return it; if no position is associated, return nothing.

Applicable cases include LittleGroup (return the associated k-vector) and SiteGroup (returns the associated Wyckoff position), as well as their associated irrep types (LGIrrep and SiteIrrep).

source
Bravais.cartesianizeMethod
cartesianize(op::SymOperation{D}, Rs::DirectBasis{D}) --> SymOperation{D}

Converts opˡ from a lattice basis to a Cartesian basis, by computing the transformed operators opᶜ = 𝐑*opˡ*𝐑⁻¹ via the Cartesian basis matrix 𝐑 (whose columns are the DirectBasis vectors Rs[i]).

Note 1

The matrix 𝐑 maps vectors coefficients in a lattice basis 𝐯ˡ to coefficients in a Cartesian basis 𝐯ᶜ as 𝐯ˡ = 𝐑⁻¹𝐯ᶜ and vice versa as 𝐯ᶜ = 𝐑𝐯ˡ. Since a general transformation P transforms an "original" vectors with coefficients 𝐯 to new coefficients 𝐯′ via 𝐯′ = P⁻¹𝐯 and since we here here consider the lattice basis as the "original" basis we have P = 𝐑⁻¹. As such, the transformation of the operator op transforms as opᶜ = P⁻¹*opˡ*P, i.e. opᶜ = transform(opˡ,P) = transform(opˡ,𝐑⁻¹).

Note 2

The display (e.g. Seitz and xyzt notation) of SymOperations e.g. in the REPL implicitly assumes integer coefficients for its point-group matrix: as a consequence, displaying SymOperations in a Cartesian basis may produce undefined behavior. The matrix representation remains valid, however.

source
Bravais.centeringMethod
centering(g::AbstractGroup) --> Char

Return the conventional centering type of a group.

For groups without lattice structure (e.g., point groups), return nothing.

source
Bravais.conventionalizeMethod
conventionalize(flat::AbstractFourierLattice, cntr::Char) --> ::typeof(flat′)

Given flat referred to a primitive basis with centering cntr, compute the derived (but physically equivalent) lattice flat′ referred to the associated conventional basis.

See also the complementary method primitivize(::AbstractFourierLattice, ::Char) for additional details.

source
Bravais.conventionalizeMethod
conventionalize(v′::AbstractVec, cntr::Char)  -->  v::typeof(v′)

Transforms a primitive coordinate vector v′ back to a standard conventional basis (specified by the centering type cntr), returning the conventional coordinate vector v.

See also primitivize and transform.

source
Bravais.conventionalizeMethod
conventionalize(op′::SymOperation, cntr::Char, modw::Bool=true) --> typeof(op)

Return a symmetry operation op $= \{W|w\}$ in a conventional setting, transformed from an input symmetry operation op′ $≡ \{W'|w'\}$ in a primitive setting.

See primitivize(::SymOperation, ::Char, ::Bool) for description of the centering argument cntr and optional argument modw.

source
Bravais.primitivizeMethod
primitivize(flat::AbstractFourierLattice, cntr::Char) --> ::typeof(flat)

Given flat referred to a conventional basis with centering cntr, compute the derived (but physically equivalent) lattice flat′ referred to the associated primitive basis.

Specifically, if flat refers to a direct conventional basis Rs $≡ (\mathbf{a} \mathbf{b} \mathbf{c})$ [with coordinate vectors $\mathbf{r} ≡ (r_1, r_2, r_3)^{\mathrm{T}}$] then flat′ refers to a direct primitive basis Rs′ $≡ (\mathbf{a}' \mathbf{b}' \mathbf{c}') ≡ (\mathbf{a} \mathbf{b} \mathbf{c})\mathbf{P}$ [with coordinate vectors $\mathbf{r}' ≡ (r_1', r_2', r_3')^{\mathrm{T}} = \mathbf{P}^{-1}\mathbf{r}$], where $\mathbf{P}$ denotes the basis-change matrix obtained from primitivebasismatrix(...).

To compute the associated primitive basis vectors, see primitivize(::DirectBasis, ::Char) [specifically, Rs′ = primitivize(Rs, cntr)].

Examples

A centered ('c') lattice from plane group 5 in 2D, plotted in its conventional and primitive basis (requires using PyPlot):

julia> sgnum = 5; D = 2; cntr = centering(sgnum, D)  # 'c' (body-centered)

julia> Rs   = directbasis(sgnum, Val(D))     # conventional basis (rectangular)
julia> flat = levelsetlattice(sgnum, Val(D)) # Fourier lattice in basis of Rs
julia> flat = modulate(flat)                 # modulate the lattice coefficients
julia> plot(flat, Rs)

julia> Rs′   = primitivize(Rs, cntr)    # primitive basis (oblique)
julia> flat′ = primitivize(flat, cntr)  # Fourier lattice in basis of Rs′

julia> using PyPlot
julia> plot(flat′, Rs′)
source
Bravais.primitivizeMethod
primitivize(v::AbstractVec, cntr::Char)  -->  v′::typeof(v)

Transforms a conventional coordinate vector v to a standard primitive basis (specified by the centering type cntr), returning the primitive coordinate vector v′.

Note that a basis change matrix $\mathbf{P}$ (as returned e.g. by Bravais.primitivebasismatrix) transforms direct coordinate vectors (RVec) as $\mathbf{r}' = \mathbf{P}^{-1}\mathbf{r}$ but transforms reciprocal coordinates (KVec) as $\mathbf{k}' = \mathbf{P}^{\text{T}}\mathbf{k}$ [ITA6]. Recall also the distinction between transforming a basis and the coordinates of a vector.

source
Bravais.primitivizeMethod
primitivize(op::SymOperation, cntr::Char, modw::Bool=true) --> typeof(op)

Return a symmetry operation op′ $≡ \{W'|w'\}$ in a primitive setting, transformed from an input symmetry operation op $= \{W|w\}$ in a conventional setting. The operations $\{W'|w'\}$ and $\{W|w\}$ are related by a transformation $\{P|p\}$ via (cf. Section 1.5.2.3 of [ITA6]):

\[\{W'|w'\} = \{P|p\}⁻¹\{W|w\}\{P|p\}.\]

where $P$ and $p$ are the basis change matrix and origin shifts, respectively. The relevant transformation $\{P|p\}$ is inferred from the centering type, as provided by cntr (see also Bravais.centering).

By default, translation parts of op′, i.e. $w'$ are reduced modulo 1 (modw = true); to disable this, set modw = false.

source
Bravais.transformMethod
transform(v::AbstractVec, P::AbstractMatrix)  -->  v′::typeof(v)

Return a transformed coordinate vector v′ from an original coordinate vector v using a basis change matrix P.

Note that a basis change matrix $\mathbf{P}$ transforms direct coordinate vectors (RVec) as $\mathbf{r}' = \mathbf{P}^{-1}\mathbf{r}$ but transforms reciprocal coordinates (KVec) as $\mathbf{k}' = \mathbf{P}^{\mathrm{T}}\mathbf{k}$ [ITA6]

source
Bravais.transformMethod
transform(op::SymOperation, P::AbstractMatrix{<:Real}, 
          p::Union{AbstractVector{<:Real}, Nothing}=nothing,
          modw::Bool=true)                           --> SymOperation

Transforms a op $= \{\mathbf{W}|\mathbf{w}\}$ by a rotation matrix P and a translation vector p (can be nothing for zero-translations), producing a new symmetry operation op′ $= \{\mathbf{W}'|\mathbf{w}'\}$ (cf. Section 1.5.2.3 of [ITA6])

\[\{\mathbf{W}'|\mathbf{w}'\} = \{\mathbf{P}|\mathbf{p}\}^{-1}\{\mathbf{W}|\mathbf{w}\} \{\mathbf{P}|\mathbf{p}\}\]

with

\[\mathbf{W}' = \mathbf{P}^{-1}\mathbf{W}\mathbf{P} \text{ and } \mathbf{w}' = \mathbf{P}^{-1}(\mathbf{w}+\mathbf{W}\mathbf{p}-\mathbf{p})\]

By default, the translation part of op′, i.e. $\mathbf{w}'$, is reduced to the range $[0,1)$, i.e. computed modulo 1. This can be disabled by setting modw = false (default, modw = true).

See also Bravais.primitivize(::SymOperation, ::Char, ::Bool) and Bravais.conventionalize(::SymOperation, ::Char, ::Bool).

source
Crystalline.:⊕Method
⊕(ir1::T, ir2::T, ir3::T...) where T<:AbstractIrrep --> T

Compute the representation obtained from direct sum of the irreps ir1, ir2, ir3, etc. The resulting representation is reducible and has dimension irdim(ir1) + irdim(ir2) + irdim(ir3) + ….

The groups of the provided irreps must be identical. If T isa LGIrrep, the irrep translation factors must also be identical (due to an implementation detail of the LGIrrep type).

Also provided via Base.:+.

source
Crystalline.bandrepsFunction
bandreps(sgnum::Integer, D::Integer=3; 
         allpaths::Bool=false, spinful::Bool=false, timereversal::Bool=true)

Returns the elementary band representations (EBRs) as a BandRepSet for space group sgnum and dimension D.

Keyword arguments

  • allpaths: include a minimal sufficient set (false, default) or all (true) k-vectors.
  • spinful: single- (false, default) or double-valued (true) irreps, as appropriate for spinless and spinful particles, respectively. Only available for D=3.
  • timereversal: assume presence (true, default) or absence (false) of time-reversal symmetry.

References

3D EBRs are obtained from the Bilbao Crystallographic Server's BANDREP program; please reference the original research papers noted there if used in published work.

source
Crystalline.basisdimMethod
basisdim(BRS::BandRepSet) --> Int

Return the dimension of the (linearly independent parts) of a band representation set. This is $d^{\text{bs}} = d^{\text{ai}}$ in the notation of Po, Watanabe, & Vishwanath, Nature Commun. 8, 50 (2017), or equivalently, the rank of matrix(BRS) over the ring of integers. This is the number of linearly independent basis vectors that span the expansions of a band structure viewed as symmetry data.

source
Crystalline.calc_realityMethod
calc_reality(lgir::LGIrrep, 
             sgops::AbstractVector{SymOperation{D}},
             αβγ::Union{Vector{<:Real},Nothing}=nothing) --> ::(Enum Reality)

Compute and return the reality of a lgir::LGIrrep using the Herring criterion.

The computed value is one of three integers in ${1,-1,0}$. In practice, this value is returned via a member of the Enum Reality, which has instances REAL = 1, PSEUDOREAL = -1, and COMPLEX = 0.

Optional arguments

As a sanity check, a value of αβγ can be provided to check for invariance along a symmetry symmetry line/plane/general point in k-space. The reality must be invariant to this choice.

Note

The provided space group operations sgops must be the set reduced by primitive translation vectors; i.e. using spacegroup(...) directly is not allowable in general (since the irreps we reference only include these "reduced" operations). This reduced set of operations can be obtained e.g. from the Γ point irreps of ISOTROPY's dataset, or alternatively, from reduce_ops(spacegroup(...), true).

Implementation

The Herring criterion evaluates the following sum

$[∑ χ({β|b}²)]/[g_0/M(k)]$

over symmetry operations ${β|b}$ that take $k → -k$. Here $g_0$ is the order of the point group of the space group and $M(k)$ is the order of star($k$) [both in a primitive basis].

See e.g. Cornwell, p. 150-152 & 187-188 (which we mainly followed), Inui Eq. (13.48), Dresselhaus, p. 618, or Herring's original paper.

source
Crystalline.charactersMethod
characters(irs::AbstractVector{<:AbstractIrrep}, αβγ=nothing)

Compute the character table associated with vector of AbstractIrreps irs, returning a CharacterTable.

Optional arguments

Optionally, an αβγ::AbstractVector{<:Real} variable can be passed to evaluate the irrep (and associated characters) with concrete free parameters (e.g., for LGIrreps, a concrete k-vector sampled from a "line-irrep"). Defaults to nothing, indicating it being either irrelevant (e.g., for PGIrreps) or all free parameters implicitly set to zero.

source
Crystalline.classcharactersMethod
classcharacters(irs::AbstractVector{<:AbstractIrrep}, αβγ=nothing)

Compute the character table associated with the conjugacy classes of a vector of AbstractIrreps irs, returning a ClassCharacterTable.

Since characters depend only on the conjugacy class (this is not true for ray, or projective, irreps), the class-specific characters often more succintly communicate the same information as the characters for each operation (as returned by characters).

See also: classes.

Optional arguments

Optionally, an αβγ::AbstractVector{<:Real} variable can be passed to evaluate the irrep (and associated characters) with concrete free parameters (e.g., for LGIrreps, a concrete k-vector sampled from a "line-irrep"). Defaults to nothing, indicating it being either irrelevant (e.g., for PGIrreps) or all free parameters implicitly set to zero.

source
Crystalline.classesMethod
classes(ops::AbstractVector{SymOperation{D}}, [cntr::Union{Char, Nothing}])
                                                -->  Vector{Vector{SymOperation{D}}}

Return the conjugacy classes of a group $G$ defined by symmetry operations ops.

Definitions

Two elements $a$ and $b$ in $G$ are considered conjugate if there exists a $g ∈ G$ such that $gag^{-1} = b$. This defines an equivalence relation $\sim$, i.e., we say that $a \sim b$ if $a$ and $b$ are conjugate. The conjugacy classes of $G$ are the distinct equivalence classes that can be identified under this equivalence relation, i.e. the grouping of $G$ into subsets that are equivalent under conjugacy.

Extended help

If ops describe operations in a crystal system that is not primitive (i.e., if its centering type is not p or P) but is presented in a conventional setting, the centering symbol cntr must be given. If ops is not in a centered crystal system, or if ops is already reduced to a primitive setting, cntr should be given as nothing (default behavior) or, alternatively, as P or p (depending on dimensionality).

A single-argument calls to classes with SpaceGroup or LittleGroup types will assume that ops is provided in a conventional setting, i.e., will forward the method call to classes(ops, centering(ops, dim(ops))). To avoid this behavior (if ops was already reduced to a primitive setting prior to calling classes), cntr should be provided explicitly as nothing.

source
Crystalline.classificationMethod
classification(BRS_or_F::Union{BandRepSet, Smith}) --> String

Return the symmetry indicator group $X^{\text{BS}}$ of an EBR basis F_or_BRS, provided as a BandRepSet or Smith decomposition.

Technically, the calculation answers the question "what direct product of $\mathbb{Z}_n$ groups is the the quotient group $X^{\text{BS}} = \{\text{BS}\}/\{\text{AI}\}$ isomorphic to?" (see Po, Watanabe, & Vishwanath, Nature Commun. 8, 50 (2017) for more information).

source
Crystalline.composeMethod
compose(op::SymOperation, kv::KVec[, checkabc::Bool=true])  -->  KVec

Return the composition op $= \{\mathbf{W}|\mathbf{w}\}$ and a reciprocal-space vector kv $= \mathbf{k}$.

The operation is taken to act directly, returning

\[ \mathbf{k}' = \{\mathbf{W}|\mathbf{w}\}\mathbf{k} = \mathbf{W}^{\text{T}}\mathbf{k}.\]

Note the transposition of $\mathbf{W}$, arising as a result of the implicit real-space basis of $\{\mathbf{W}|\mathbf{w}\}$ versus the reciprocal-space basis specification of $\mathbf{k}$. Note also that the composition of $\{\mathbf{W}|\mathbf{w}\}$ with $\mathbf{k}$ is invariant under $\mathbf{w}$, i.e., translations do not act in reciprocal space.

Extended help

If checkabc = false, the free part of KVec is not transformed (can be improve performance in situations when kabc is zero, and several transformations are requested).

source
Crystalline.composeMethod
compose(op::SymOperation, rv::RVec)  -->  RVec

Return the composition of op $= \{\mathbf{W}|\mathbf{w}\}$ and a real-space vector rv $= \mathbf{r}$.

The operation is taken to act directly, returning

\[ \mathbf{r}' = \{\mathbf{W}|\mathbf{w}\}\mathbf{r} = \mathbf{W}\mathbf{r} + \mathbf{w}.\]

The corresponding inverse action $\{\mathbf{W}|\mathbf{w}\}^{-1}\mathbf{r} = \mathbf{W}^{-1}\mathbf{r} - \mathbf{W}^{-1}\mathbf{w}$ can be obtained via compose(inv(op), rv).

source
Crystalline.composeMethod
compose(op1::T, op2::T, modτ::Bool=true) where T<:SymOperation

Compose two symmetry operations op1 $= \{W₁|w₁\}$ and op2 $= \{W₂|w₂\}$ using the composition rule (in Seitz notation)

$\{W₁|w₁\}\{W₂|w₂\} = \{W₁W₂|w₁+W₁w₂\}$

By default, the translation part of the $\{W₁W₂|w₁+W₁w₂\}$ is reduced to the range $[0,1[$, i.e. computed modulo 1. This can be toggled off (or on) by the Boolean flag modτ (enabled, i.e. true, by default). Returns another SymOperation.

The multiplication operator * is overloaded for SymOperations to call compose, in the manner op1 * op2 = compose(op1, op2, modτ=true).

source
Crystalline.cosetsMethod
cosets(g::SiteGroup) -> Array{SymOperation{D}, 1} where D

Return the cosets of a SiteGroup g.

The cosets generate the orbit of the Wyckoff position position(g) (see also orbit(::SiteGroup)) and furnish a left-coset decomposition of the underlying space group, jointly with the operations in g itself.

source
Crystalline.dimMethod
dim(BR::BandRep) --> Int

Return the number of bands included in the provided BandRep.

If the bands are "nondetachable" (i.e. if BR.decomposable = false), this is equal to a band connectivity μ.

source
Crystalline.find_isomorphic_parent_pointgroupMethod
find_isomorphic_parent_pointgroup(g::AbstractVector{SymOperation{D}}) 
                                                --> PointGroup{D}, Vector{Int}, Bool

Given a group g (or a collection of operators, defining a group), identifies a "parent" point group that is isomorphic to g.

Three variables are returned:

  • pg: the identified "parent" point group, with operators sorted to match the sorting of g's operators.
  • Iᵖ²ᵍ: a permutation vector which transforms the standard sorting of point group operations (as returned by pointgroup(::String)) to the operator sorting of g.
  • equal: a boolean, identifying whether the point group parts of g operations are identical (true) or merely isomorphic to the point group operations in g. In practice, this indicates whether pg and g are in the same setting or not.

Implementation

The identification is made partly on the basis of comparison of operators (this is is sufficient for the equal = true case) and partly on the basis of comparison of multiplication tables (equal = false case); the latter can be combinatorially slow if the sorting of operators is unlucky (i.e., if permutation between sortings in g and pg differ by many pairwise permutations).

Beyond mere isomorphisms of multiplication tables, the search also guarantees that all rotation orders are shared between pg and g. This disambiguates point groups that are intrinsically isomorphic to eachother, e.g. "m" and "-1", but which still differ in their spatial interpretation.

Properties

The following properties hold for g, pg, and Iᵖ²ᵍ:

pg, Iᵖ²ᵍ, equal = find_isomorphic_parent_pointgroup(g)
@assert MultTable(pg) == MultTable(pointgroup(g))
pg′ = pointgroup(label(pg), dim(pg)) # "standard" sorting
@assert pg′[Iᵖ²ᵍ] == pg

If equal = true, the following also holds:

pointgroup(g) == operations(pg) == operations(pg′)[Iᵖ²ᵍ]

Example

sgnum = 141
wp    = wyckoffs(sgnum, Val(3))[end] # 4a Wyckoff position
sg    = spacegroup(sgnum, Val(3))
siteg = sitegroup(sg, wp)
pg, Iᵖ²ᵍ, equal = find_isomorphic_parent_pointgroup(siteg)
source
Crystalline.find_representationFunction
find_representation(symvals::AbstractVector{Number}, 
                    lgirs::AbstractVector{<:AbstractIrrep},
                    αβγ::Union{AbstractVector{<:Real},Nothing}=nothing,
                    assert_return_T::Type{<:Union{Integer, AbstractFloat}}=Int),
                    atol::Real=DEFAULT_ATOL,
                    maxresnorm::Real=1e-3,
                    verbose::Bool=false)

                    --> Vector{assert_return_T}

From a vector (or vector of vectors) of symmetry eigenvalues symvals sampled along all the operations of a group gᵢ, whose irreps are contained in irs (evaluated with optional free parameters αβγ), return the multiplicities of each irrep.

Optionally, the multiciplities' element type can be specified via the assert_return_T argument (performing checked conversion; returns nothing if representation in assert_return_T is impossible). This can be useful if one suspects a particular band to transform like a fraction of an irrep (i.e., the specified symmetry data is incomplete).

If no valid set of multiplicities exist (i.e., is solvable, and has real-valued and assert_return_T representible type), the sentinel value nothing is returned. Optional debugging information can in this case be shown by setting verbose=true.

Extended help

Effectively, this applies the projection operator P⁽ʲ⁾ of each irrep's character set χ⁽ʲ⁾(gᵢ) (j = 1, ... , Nⁱʳʳ) to the symmetry data sᵢ ≡ symvals:

P⁽ʲ⁾  ≡ (dⱼ/|g|) ∑ᵢ χ⁽ʲ⁾(gᵢ)*gᵢ         [characters χ⁽ʲ⁾(gᵢ), irrep dimension dⱼ]
P⁽ʲ⁾s = (dⱼ/|g|) ∑ᵢ χ⁽ʲ⁾(gᵢ)*sᵢ = nⱼ,   [number of bands that transform like jth irrep]

returning the irrep multiplicities mⱼ ≡ nⱼ/dⱼ.

source
Crystalline.findmaximalMethod
findmaximal(sitegs::AbstractVector{<:SiteGroup})

Given a vector of SiteGroups associated with the Wyckoff positions of a space group, return those SiteGroups that are associated with a maximal Wyckoff positions.

Results are returned as a view into the input vector (i.e. as an AbstractVector{<:SiteGroup}). The associated Wyckoff positions can be retrieved via position.

Definition

A Wyckoff position is maximal if its site symmetry group has higher order than the site symmetry groups of any "nearby" Wyckoff positions (i.e. Wyckoff positions that can be connected, i.e. made equivalent, through parameter variation to the considered Wyckoff position).

Example

julia> sgnum = 5;

julia> D = 2;

julia> wps = wyckoffs(sgnum, Val(D));

julia> sg  = spacegroup(sgnum, Val(D));


julia> sitegs = sitegroup.(Ref(sg), wps)
2-element Vector{SiteGroup{2}}:
 [1] (4b: [α, β])
 [1, m₁₀] (2a: [0, β])

julia> only(findmaximal(sitegs))
SiteGroup{2} ⋕5 (c1m1) at 2a = [0, β] with 2 operations:
 1
 m₁₀
source
Crystalline.generateMethod
generate(
    gens::AbstractArray{SymOperation{D}, 1};
    cntr,
    modτ,
    Nmax
) -> Crystalline.GenericGroup

Return the group generated from a finite set of generators gens.

Keyword arguments

  • cntr (default, nothing): check equivalence of operations modulo primitive lattice vectors (see also isapprox(::SymOperation, ::SymOperation, cntr::Union{Nothing, Char}); only nonequivalent operations are included in the returned group.
  • modτ (default, true): the group composition operation can either be taken modulo lattice vectors (true) or not (false, useful e.g. for site symmetry groups). In this case, the provided generators will also be taken modulo integer lattice translations.
  • Nmax (default, 256): the maximum size of the generated group. This is essentially a cutoff set to ensure halting of execution in case the provided set of generators do not define a finite group (especially relevant if modτ=false). If more operations than Nmax are generated, the method throws an overflow error.
source
Crystalline.generatorsMethod
generators(num::Integer, T::Type{AbstractGroup{D}}[, optargs])
generators(pgiuc::String, T::PointGroup{D}})              -->  Vector{SymOperation{D}}

Return the generators of the group type T which may be a SpaceGroup{D} or a PointGroup{D} parameterized by its dimensionality D. Depending on T, the group is determined by inputting as the first argument:

  • SpaceGroup{D}: the space group number num::Integer.
  • PointGroup{D}: the point group IUC label pgiuc::String (see also [pointgroup(::String)) or the canonical point group number num::Integer, which can optionally be supplemented by an integer-valued setting choice setting::Integer (see also pointgroup(::Integer, ::Integer, ::Integer)]).
  • SubperiodicGroup{D}: the subperiodic group number num::Integer.

Setting choices match those in spacegroup, pointgroup, and subperiodicgroup.

Iterated composition of the returned symmetry operations will generate all operations of the associated space or point group (see generate). As an example, generate(generators(num,SpaceGroup{D}))andspacegroup(num, D)` return identical operations (with different sorting typically); and similarly so for point and subperiodic groups.

Example

Generators of space group 200:

julia> generators(200, SpaceGroup{3})
4-element Vector{SymOperation{3}}:
 2₀₀₁
 2₀₁₀
 3₁₁₁⁺
 -1

Generators of point group m-3m:

julia> generators("2/m", PointGroup{3})
2-element Vector{SymOperation{3}}:
 2₀₁₀
 -1

Generators of the Frieze group 𝓅2mg:

julia> generators(7, SubperiodicGroup{2, 1})
2-element Vector{SymOperation{2}}:
 2
 {m₁₀|½,0}

Citing

Please cite the original data sources if used in published work:

Extended help

Note that the returned generators are not guaranteed to be the smallest possible set of generators; i.e., there may exist other generators with fewer elements (an example is space group 168 (P6), for which the returned generators are [2₀₀₁, 3₀₀₁⁺] even though the group could be generated by just [6₀₀₁⁺]). The returned generators, additionally, are not guaranteed to be minimal, i.e., they may include proper subsets that generate the group themselves (e.g., in space group 75 (P4), the returned generators are [2₀₀₁, 4₀₀₁⁺] although the subset [4₀₀₁⁺] is sufficient to generate the group). The motivation for this is to expose as similar generators as possible for similar crystal systems (see e.g. Section 8.3.5 of the International Tables of Crystallography, Vol. A, Ed. 5 (ITA) for further background).

Note also that, contrary to conventions in ITA, the identity operation is excluded among the returned generators (except in space group 1) since it composes trivially and adds no additional context.

source
Crystalline.generatorsMethod
generators(num::Integer, ::Type{SubperiodicGroup{D,P}})  -->  ::Vector{SymOperation{D}}

Return a canonical set of generators for the subperiodic group num of embedding dimension D and periodicity dimension P. See also subperiodicgroup.

See also generators(::Integer, ::Type{SpaceGroup}) and information therein.

Example

julia> generators(7, SubperiodicGroup{2, 1})
2-element Vector{SymOperation{2}}:
 2
 {m₁₀|½,0}

Data sources

The generators returned by this function were originally retrieved from the Bilbao Crystallographic Database, SUBPERIODIC GENPOS.

source
Crystalline.is_abelianMethod
is_abelian(ops::AbstractVector{SymOperation}, [cntr::Union{Char, Nothing}])  -->  Bool

Return the whether the group composed of the elements ops is Abelian.

A group $G$ is Abelian if all its elements commute mutually, i.e., if $g = hgh^{-1}$ for all $g,h ∈ G$.

See discussion of the setting argument cntr in classes.

source
Crystalline.isnormalMethod
isnormal(opsᴳ::AbstractVector{<:SymOperation},
         opsᴴ::AbstractVector{<:SymOperation};
         verbose::Bool=false)                    --> Bool

Determine whether the operations in group $H$ are normal in the group $G$ (each with operations opsᴳ and opsᴴ), in the sense that

\[ghg⁻¹ ∈ H, ∀ g∈G, ∀ h∈H\]

Returns a Boolean answer (true if normal, false if not).

Note

This compares space groups rather than space group types, i.e. the comparison assumes a matching setting choice between $H$ and $G$. To compare space group types with different conventional settings, they must first be transformed to a shared setting.

source
Crystalline.israyrepFunction
israyrep(lgir::LGIrrep, αβγ=nothing) -> (::Bool, ::Matrix)

Computes whether a given little group irrep ir is a ray representation by computing the coefficients αᵢⱼ in DᵢDⱼ=αᵢⱼDₖ; if any αᵢⱼ differ from unity, we consider the little group irrep a ray representation (as opposed to the simpler "vector" representations where DᵢDⱼ=Dₖ). The function returns a boolean (true => ray representation) and the coefficient matrix αᵢⱼ.

source
Crystalline.issubgroupMethod
issubgroup(opsᴳ::T, opsᴴ::T′) where T⁽′⁾<:AbstractVector{SymOperation} --> Bool

Determine whether the operations in group $H$ are a subgroup of the group $G$ (each with operations opsᴳ and opsᴴ, respectively), i.e. whether $H<G$. Specifically, this requires that $G$ and $H$ are both groups and that for every $h∈H$ there exists an element $g∈G$ such that $h=g$.

Returns a Boolean answer (true if normal, false if not).

Note

This compares space groups rather than space group types, i.e. the comparison assumes a matching setting choice between $H$ and $G$. To compare space group types with different conventional settings, they must first be transformed to a shared setting.

source
Crystalline.issymmorphFunction
issymmorph(sgnum::Integer, D::Integer=3) --> Bool

Return whether the space group with number sgnum and dimensionality D is symmorphic (true) or nonsymmorphic (false).

Equivalent to issymmorph(spacegroup(sgnum, D)) but uses memoization for performance.

source
Crystalline.issymmorphMethod
issymmorph(sg::Union{SpaceGroup, LittleGroup}) --> Bool

Return whether a given space group sg is symmorphic (true) or nonsymmorphic (false).

source
Crystalline.issymmorphMethod
issymmorph(op::SymOperation, cntr::Char) --> Bool

Return whether a given symmetry operation op is symmorphic (true) or nonsymmorphic (false).

The operation is assumed provided in conventional basis with centering type cntr: checking symmorphism is then equivalent to checking whether the operation's translation part is zero or a lattice vector in the associated primitive basis.

source
Crystalline.iucFunction
iuc(sgnum::Integer, D::Integer=3) --> String

Return the IUC (International Union of Crystallography) notation for space group number sgnum in dimension D (1, 2, or 3), as used in the International Tables of Crystallography.

The notation is sometimes also known as the Hermann-Mauguin notation.

source
Crystalline.levelsetlatticeMethod
levelsetlattice(sgnum::Integer, D::Integer=2, idxmax::NTuple=ntuple(i->2,D))
    --> UnityFourierLattice{D}

Compute a "neutral"/uninitialized Fourier lattice basis, a UnityFourierLattice, consistent with the symmetries of the space group sgnum in dimension D. The resulting lattice flat is expanded in a Fourier basis split into symmetry-derived orbits, with intra-orbit coefficients constrained by the symmetries of the space-group. The inter-orbit coefficients are, however, free and unconstrained.

The Fourier resolution along each reciprocal lattice vector is controlled by idxmax: e.g., if D = 2 and idxmax = (2, 3), the resulting Fourier lattice may contain reciprocal lattice vectors (k₁, k₂) with k₁∈[0,±1,±2] and k₂∈[0,±1,±2,±3], referred to a 𝐆-basis.

This "neutral" lattice can, and usually should, be subsequently modulated by modulate (which modulates the inter-orbit coefficients, which may eliminate "synthetic symmetries" that can exist in the "neutral" configuration, due to all inter-orbit coefficients being set to unity).

Examples

Compute a UnityFourierLattice, modulate it with random inter-orbit coefficients via modulate, and finally plot it (via PyPlot.jl):

julia> uflat = levelsetlattice(16, Val(2))
julia> flat  = modulate(uflat)
julia> Rs    = directbasis(16, Val(2)) 
julia> using PyPlot
julia> plot(flat, Rs)
source
Crystalline.lgirrepsMethod
lgirreps(sgnum::Integer, D::Union{Val{Int}, Integer}=Val(3))
                                                -> Dict{String, Vector{LGIrrep{D}}}

For given space group number sgnum and dimension D, return the associated little group (or "small") irreps (LGIrrep{D}s) at high-symmetry k-points, lines, and planes.

Returns a Dict with little group k-point labels as keys and vectors of LGIrrep{D}s as values.

Notes

  • The returned irreps are complex in general. Real irreps (as needed in time-reversal invariant settings) can subsequently be obtained with the realify method.
  • Returned irreps are spinless.
  • The irrep labelling follows CDML conventions.
  • Irreps along lines or planes may depend on free parameters αβγ that parametrize the k point. To evaluate the irreps at a particular value of αβγ and return the associated matrices, use (lgir::LGIrrep)(αβγ). If αβγ is an empty tuple in this call, the matrices associated with lgir will be evaluated assuming αβγ = [0,0,...].

References

The underlying data is sourced from the ISOTROPY ISO-IR dataset. Please cite the original reference material associated with ISO-IR:

  1. Stokes, Hatch, & Campbell, ISO-IR, ISOTROPY Software Suite.
  2. Stokes, Campbell, & Cordes, Acta Cryst. A. 69, 388-395 (2013).

The ISO-IR dataset is occasionally missing some k-points that lie outside the basic domain but still resides in the representation domain (i.e. k-points with postscripted 'A', 'B', etc. labels, such as 'ZA'). In such cases, the missing irreps may instead have been manually sourced from the Bilbao Crystallographic Database.

source
Crystalline.littlegroupFunction
littlegroup(sg::SpaceGroup, kv::KVec) -> LittleGroup
littlegroup(
    sg::SpaceGroup,
    kv::KVec,
    klab::String
) -> LittleGroup

Return the little group associated with space group sg at the k-vector kv.

Optionally, an associated k-vector label klab can be provided; if not provided, the empty string is used as label.

source
Crystalline.littlegroupsMethod
littlegroups(sgnum::Integer, D::Union{Val{Int}, Integer}=Val(3)) 
                                                    -> Dict{String, LittleGroup{D}}

For given space group number sgnum and dimension D, return the associated little groups (LittleGroups{D}s) at high-symmetry k-points, lines, and planes (see also lgirreps).

Returns a Dict with little group k-point labels as keys and vectors of LittleGroup{D}s as values.

Notes

A conventional crystallographic setting is assumed (as in spacegroup).

Unlike spacegroup, "centering"-copies of symmetry operations are not included in the returned LittleGroups; as an example, space group 110 (body-centered, with centering symbol 'I') has a centering translation [1/2,1/2,1/2] in the conventional setting: the symmetry operations returned by spacegroup thus includes e.g. both {1|0} and {1|½,½,½} while the symmetry operations returned by littlegroups only include {1|0} (and so on).

Currently, only D = 3 is supported.

References

The underlying data is sourced from the ISOTROPY dataset: see also lgirreps.

source
Crystalline.matrixMethod
matrix(BRS::BandRepSet; includedim::Bool=true)

Return a matrix representation of BRS::BandRepSet, with band representations as columns and irreps over rows.

By default, the last row will give the "filling" of each BandRep (or, more precisely, number of included bands per BandRep, i.e. dim.(BRS). To toggle this off, set the keyword argument includedim to false (default is includedim = true).

source
Crystalline.matrixMethod
matrix(op::AbstractOperation{D}) --> SMatrix{D, D+1, Float64}

Return the D×D+1 matrix representation of op.

source
Crystalline.maximal_subgroupsMethod
maximal_subgroups(num::Integer, AG::Type{<:AbstractGroup}=SpaceGrop{3}; kind)

Returns the graph structure of the maximal subgroups as a GroupRelationGraph for a group of type AG and number num.

Visualization

The resulting group structure can be plotted using Makie.jl (e.g., GLMakie.jl) using plot(::GroupRelationGraph):

julia> using Crystalline
julia> gr = maximal_subgroups(112, SpaceGroup{3})
julia> using GraphMakie, GLMakie
julia> plot(gr)

Keyword arguments

  • kind (default, Crystalline.TRANSLATIONENGLEICHE): to return klassengleiche relations, set kind = Crystalline.KLASSENGLEICHE). For klassengleiche relationships, only a selection of reasonably low-index relationships are returned.

Data sources

The group relationships returned by this function were retrieved from the Bilbao Crystallographic Server's MAXSUB program. Please cite the original reference work associated with MAXSUB:

source
Crystalline.minimal_supergroupsMethod
minimal_supergroups(num::Integer, AG::Type{<:AbstractGroup}=SpaceGrop{3}; kind)

Returns the graph structure of the minimal supergroups as a GroupRelationGraph for a group of type AG and number num.

Visualization

The resulting group structure can be plotted using Makie.jl (e.g., GLMakie.jl) using plot(::GroupRelationGraph):

julia> using Crystalline
julia> gr = minimal_supergroups(112, SpaceGroup{3})
julia> using GraphMakie, GLMakie
julia> plot(gr)

Keyword arguments

  • kind (default, Crystalline.TRANSLATIONENGLEICHE): to return klassengleiche relations, set kind = Crystalline.KLASSENGLEICHE). For klassengleiche relationships, only a selection of reasonably low-index relationships are returned.

Data sources

The group relationships returned by this function were retrieved from the Bilbao Crystallographic Server's MAXSUB program. Please cite the original reference work associated with MAXSUB:

source
Crystalline.modulateMethod
modulate(flat::UnityFourierLattice{D},
modulation::AbstractVector{ComplexF64}=rand(ComplexF64, length(getcoefs(flat))),
expon::Union{Nothing, Real}=nothing, Gs::Union{ReciprocalBasis{D}, Nothing}=nothing)
                        --> ModulatedFourierLattice{D}

Derive a concrete, modulated Fourier lattice from a UnityFourierLattice flat (containing the interrelations between orbit coefficients), by multiplying the "normalized" orbit coefficients by a modulation, a complex modulating vector (in general, should be complex; otherwise restores unintended symmetry to the lattice). Distinct modulation vectors produce distinct realizations of the same lattice described by the original flat. By default, a random complex vector is used.

An exponent expon can be provided, which introduces a penalty term to short- wavelength features (i.e. high-|G| orbits) by dividing the orbit coefficients by |G|^expon; producing a "simpler" and "smoother" lattice boundary when expon > 0 (reverse for expon < 0). This basically amounts to a continuous "simplifying" operation on the lattice (it is not necessarily a smoothing operation; it simply suppresses "high-frequency" components). If expon = nothing, no rescaling is performed. If Gs is provided as nothing, the orbit norm is computed in the reciprocal lattice basis (and, so, may not strictly speaking be a norm if the lattice basis is not cartesian); to account for the basis explicitly, Gs must be provided as a ReciprocalBasis, see also normscale.

source
Crystalline.mspacegroupMethod
mspacegroup(BNS₁, BNS₂)
mspacegroup(OG₃)        --> MSpaceGroup{3}

Return the magnetic space group with BNS numbers (BNS₁, BNS₂) or the sequential OG number OG₃ (from 1 to 1651).

Data sources

The data underlying this function was retrieved from ISOTROPY's compilation of Daniel Litvin's magnetic space group tables (http://www.bk.psu.edu/faculty/litvin/Download.html).

source
Crystalline.mullikenMethod
mulliken(pgir::PGIrrep{D}) -> String

Return the Mulliken label of a point group irrep pgir.

Notes

This functionality is a simple mapping between the tabulated CDML point group irrep labels and associated Mulliken labels [1], using the listings from the Bilbao Crystallographic Database [2].

Ignoring subscript, the rough rules associated with assignment of Mulliken labels are:

  1. Irrep dimensionality:
    • 1D irreps: if a real irrep, assign A or B (B if antisymmetric under a principal rotation); if a complex irrep, assigned label ¹E or ²E.
    • 2D irreps: assign label E.
    • 3D irreps: assign label T.
  2. u and g subscripts: if the group contains inversion, indicate whether irrep is symmetric (g ~ gerade) or antisymmetric (u ~ ungerade) under inversion.
  3. Prime superscripts: if the group contains a mirror m aligned with a principal rotation axis, but does not contain inversion, indicate whether irrep is symmetric (′) or antisymmetric (′′) under this mirror.
  4. Numeral subscripts: the rules for assignment of numeral subscripts are too complicated in general - and indeed, we are unaware of a general coherent rule – to describe here.

References

source
Crystalline.nontrivial_factorsMethod
nontrivial_factors(F::Smith) -> Any

Return the nontrivial (i.e., ≠ {0,1}) elementary factors of an EBR basis, provided as a BandRepSet or Smith decomposition.

source
Crystalline.normscale!Method
normscale!(flat::ModulatedFourierLattice, expon::Real,
           Gs::Union{ReciprocalBasis, Nothing} = nothing) --> ModulatedFourierLattice

In-place equivalent of normscale: mutates flat.

source
Crystalline.normscaleMethod
normscale(flat::ModulatedFourierLattice, expon::Real, 
          Gs::Union{ReciprocalBasis, Nothing} = nothing)  --> ModulatedFourierLattice

Applies inverse-orbit norm rescaling of expansion coefficients with a norm exponent expon. If Gs is nothing, the orbit norm is computed in the lattice basis (and, so, is not strictly a norm); by providing Gs as ReciprocalBasis, the norm is evaluated correctly in cartesian setting. See further discussion in modulate.

An in-place equivalent is provided in normscale!.

source
Crystalline.orbitMethod
orbit(g::SiteGroup)  -->  Vector{WyckoffPosition}

Compute the orbit of the Wyckoff position associated with the site symmetry group g.

Extended help

The orbit of a Wyckoff position $\mathbf{r}$ in a space group $G$ is defined as the set of inequivalent points in the unit cell that can be obtained by applying the elements of $G$ to $\mathbf{r}$. Equivalently, every element of the orbit of $\mathbf{r}$ can be written as the composition of a coset representative of the Wyckoff position's site group in $G$ with $\mathbf{r}$.

source
Crystalline.orbitMethod
orbit(g::AbstractVector{<:SymOperation}, kv::KVec, cntr::Char)  -->  Vector{KVec{D}}
orbit(lg::LittleGroup)
orbit(lgir::LGIrrep)

Return the orbit of of the reciprocal-space vector kv under the action of the group g, also known as the star of k.

The orbit of kv in g is the set of inequivalent k-points obtained by composition of all the symmetry operations of g with kv. Two reciprocal vectors $\mathbf{k}$ and $\mathbf{k}'$ are equivalent if they differ by a primitive reciprocal lattice vector.

If kv and g are specified in a conventional basis but refer to a non-primitive lattice, the centering type cntr must be provided to ensure that only equivalence by primitive (not conventional) reciprocal lattice vectors are considered. If the centering type of the group g can be inferred from g (e.g., if g is a SpaceGroup), orbit will assume a conventional setting and use the inferred centering type; otherwise, if cntr is neither explicitly set nor inferrable, a primitive setting is assumed.

source
Crystalline.pgirrepsFunction
pgirreps(iuclab::String, ::Val{D}=Val(3); mullikken::Bool=false) where D ∈ (1,2,3)
pgirreps(iuclab::String, D; mullikken::Bool=false)

Return the (crystallographic) point group irreps of the IUC label iuclab of dimension D as a Vector{PGIrrep{D}}.

See Crystalline.PG_IUC2NUM[D] for possible IUC labels in dimension D.

Notation

The irrep labelling follows the conventions of CDML [1] [which occasionally differ from those in e.g. Bradley and Cracknell, The Mathematical Theory of Symmetry in Solids (1972)].

To use Mulliken ("spectroscopist") irrep labels instead, set the keyword argument mulliken = true (default, false). See also mulliken.

Data sources

The data is sourced from the Bilbao Crystallographic Server [2]. If you are using this functionality in an explicit fashion, please cite the original reference [3].

References

source
Crystalline.pointgroupMethod
pointgroup(ops:AbstractVector{SymOperation{D}})
pointgroup(sg::AbstractGroup)

Computes the point group associated with a space group sg (characterized by a set of operators ops, which, jointly with lattice translations generate the space group), obtained by "taking away" any translational parts and then reducing to the resulting unique rotational operations. (technically, in the language of Bradley & Cracknell, this is the so-called isogonal point group of sg; see Sec. 1.5).

Returns a Vector of SymOperations.

source
Crystalline.pointgroupMethod
pointgroup(iuclab::String, ::Union{Val{D}, Integer}=Val(3))  -->  PointGroup{D}

Return the symmetry operations associated with the point group identified with label iuclab in dimension D as a PointGroup{D}.

source
Crystalline.pointgroupMethod
pointgroup(pgnum::Integer, ::Union{Val{D}, Integer}=Val(3), setting::Integer=1)
                                                                  -->  PointGroup{D}

Return the symmetry operations associated with the point group identfied with canonical number pgnum in dimension D as a PointGroup{D}. The connection between a point group's numbering and its IUC label is enumerated in Crystalline.PG_NUM2IUC[D] and Crystalline.IUC2NUM[D].

Certain point groups feature in multiple setting variants: e.g., IUC labels 321 and 312 both correspond to pgnum = 18 and correspond to the same group structure expressed in two different settings. The setting argument allows choosing between these setting variations.

source
Crystalline.realifyMethod
realify(lgirs::AbstractVector{<:LGIrrep}; verbose::Bool=false)
                                                    --> AbstractVector{<:LGIrrep}

From lgirs, a vector of LGIrreps, determine the associated (gray) co-representations, i.e. the "real", or "physical" irreps that are relevant in scenarios with time-reversal symmetry.

For LGIrrep that are REAL, or that characterize a k-point 𝐤 which is not equivalent to -𝐤 (i.e. its star does not include both 𝐤 and -𝐤; equivalently, the little group includes time-reversal symmetry), the associated co-representations are just the original irreps themselves. For PSEUDOREAL and COMPLEX LGIrreps where ±𝐤 are equivalent, the associated co-representations are built from pairs of irreps that "stick" together. This method computes this pairing and sets the LGIrrep field iscorep to true, to indicate that the resulting "paired irrep" (i.e. the co-representation) should be doubled with itself (PSEUDOREAL reality) or its complex conjugate (COMPLEX reality).

Background

For background, see p. 650-652 (and p. 622-626 for point groups) in Bradley & Cracknell's book. Their discussion is for magnetic groups (the "realified" irreps are, in fact, simply co-representations of the "gray" magnetic groups). Cornwell's book also explicates this at some length as does Inui et al. (p. 296-299).

Keyword arguments

  • verbose::Bool: if set to true, prints details about mapping from small irrep to small

corep for each LGIrrep (default: false).

source
Crystalline.realifyMethod
realify(pgirs::AbstractVector{T}) where T<:AbstractIrrep --> Vector{T}

Return physically real irreps (coreps) from a set of conventional irreps (as produced by e.g. pgirreps). Fallback method for point-group-like AbstractIrreps.

Example

julia> pgirs = pgirreps("4", Val(3));
julia> characters(pgirs)
CharacterTable{3}: ⋕9 (4)
───────┬────────────────────
       │ Γ₁  Γ₂    Γ₃    Γ₄ 
───────┼────────────────────
     1 │  1   1     1     1
  2₀₀₁ │  1   1    -1    -1
 4₀₀₁⁺ │  1  -1   1im  -1im
 4₀₀₁⁻ │  1  -1  -1im   1im
───────┴────────────────────

julia> characters(realify(pgirs))
CharacterTable{3}: ⋕9 (4)
───────┬──────────────
       │ Γ₁  Γ₂  Γ₃Γ₄ 
───────┼──────────────
     1 │  1   1     2
  2₀₀₁ │  1   1    -2
 4₀₀₁⁺ │  1  -1     0
 4₀₀₁⁻ │  1  -1     0
───────┴──────────────
source
Crystalline.reduce_opsMethod
reduce_ops(ops::AbstractVector{SymOperation{D}},
           cntr::Char,
           conv_or_prim::Bool=true,
           modw::Bool=true) --> Vector{SymOperation{D}}

Reduce the operations ops, removing operations that are identical in the primitive basis associated with the centering cntr.

If conv_or_prim = false, the reduced operations are returned in the primitive basis associated with cntr; otherwise, in the conventional. If modw = true, the comparison in the primitive basis is done modulo unit primitive lattice vectors; otherwise not. A final argument of type ::Val{P} can be specified to indicate a subperiodic group of periodicity dimension P, different from the spatial embedding dimension D.

source
Crystalline.rotationMethod
rotation(op::AbstractOperation{D}) --> SMatrix{D, D, Float64}

Return the D×Drotation part ofop`.

source
Crystalline.seitzMethod
seitz(op::SymOperation) --> String

Computes the correponding Seitz notation for a symmetry operation in triplet/xyzt form.

Implementation based on ITA5 Table 11.2.1.1, with 3D point group parts inferred from the trace and determinant of the matrix $\mathb{W}$ in the triplet $\{\mathbf{W}|\mathbf{w}\}$.

detW/trW-3-2-10123
123461
-1-1-6-4-3m

with the elements of the table giving the type of symmetry operation in in Hermann-Mauguin notation. The rotation axis and the rotation sense are computed following the rules in ITA6 Sec. 1.2.2.4(1)(b-c). See also .

Note that the orientation of the axis (i.e. its sign) does not necessarily match the orientation picked in Tables 1.4.2.1-5 of ITA6; it is a matter of (arbitrary) convention, and the conventions have not been explicated in ITA.

2D operations are treated by the same procedure, by elevation in a third dimension; 1D operations by a simple inspection of sign.

source
Crystalline.sitegroupMethod
sitegroup(
    sg::SpaceGroup{D},
    wp::WyckoffPosition{D}
) -> SiteGroup

Return the site symmetry group g::SiteGroup for a Wyckoff position wp in space group sg (or with space group number sgnum; in this case, the dimensionality is inferred from wp).

g is a group of operations that are isomorphic to the those listed in sg (in the sense that they might differ by lattice vectors) and that leave the Wyckoff position wp invariant, such that all(op -> wp == op*wp, g) == true.

The returned SiteGroup also contains the coset representatives of the Wyckoff position (that are again isomorphic to those featured in sg), accessible via cosets, which e.g. generate the orbit of the Wyckoff position (see orbit(::SiteGroup)) and define a left-coset decomposition of sg jointly with the elements in g.

Example

julia> sgnum = 16;

julia> D = 2;

julia> wp = wyckoffs(sgnum, D)[3] # pick a Wyckoff position
2b: [1/3, 2/3]

julia> sg = spacegroup(sgnum, D);

julia> g  = sitegroup(sg, wp)
SiteGroup{2} ⋕16 (p6) at 2b = [1/3, 2/3] with 3 operations:
 1
 {3⁺|1,1}
 {3⁻|0,1}

The group structure of a SiteGroup can be inspected with MultTable:

julia> MultTable(g)
3×3 MultTable{SymOperation{2}}:
──────────┬──────────────────────────────
          │        1  {3⁺|1,1}  {3⁻|0,1} 
──────────┼──────────────────────────────
        1 │        1  {3⁺|1,1}  {3⁻|0,1} 
 {3⁺|1,1} │ {3⁺|1,1}  {3⁻|0,1}         1
 {3⁻|0,1} │ {3⁻|0,1}         1  {3⁺|1,1}
──────────┴──────────────────────────────

The original space group can be reconstructed from a left-coset decomposition, using the operations and cosets contained in a SiteGroup:

julia> ops = [opʰ*opᵍ for opʰ in cosets(g) for opᵍ in g];

julia> Set(sg) == Set(ops)
true

Terminology

Mathematically, the site symmetry group is a stabilizer group for a Wyckoff position, in the same sense that the little group of k is a stabilizer group for a k-point.

source
Crystalline.spacegroupMethod
spacegroup(sgnum::Integer, ::Val{D}=Val(3))
spacegroup(sgnum::Integer, D::Integer)          --> SpaceGroup{D}

Return the space group symmetry operations for a given space group number sgnum and dimensionality D as a SpaceGroup{D}. The returned symmetry operations are specified relative to the conventional basis vectors, i.e. are not necessarily primitive (see centering). If desired, operations for the primitive unit cell can subsequently be generated using primitivize or Crystalline.reduce_ops.

The default choices for the conventional basis vectors follow the conventions of the Bilbao Crystallographic Server (or, equivalently, the International Tables of Crystallography), which are:

  • Unique axis b (cell choice 1) for monoclinic space groups.
  • Obverse triple hexagonal unit cell for rhombohedral space groups.
  • Origin choice 2: inversion centers are placed at (0,0,0). (relevant for certain centrosymmetric space groups with two possible choices; e.g., in the orthorhombic, tetragonal or cubic crystal systems).

See also directbasis.

Data sources

The symmetry operations returned by this function were originally retrieved from the Bilbao Crystallographic Server, SPACEGROUP GENPOS. The associated citation is: (Aroyo et al., Z. Kristallogr. Cryst. Mater. 221, 15 (2006).).

source
Crystalline.subduction_countMethod
subduction_count(Dᴳᵢ, Dᴴⱼ[, αβγᴴⱼ]) --> Int

For two groups $G$ and $H$, where $H$ is a subgroup of $G$, i.e. $H<G$, with associated irreducible representations Dᴳᵢ = $D^G_i(g)$ and Dᴴⱼ = $D^H_j(g)$ over operations $g∈G$ and $h∈H<G$, compute the compatibility relation between the two irreps from the subduction reduction formula (or "magic" formula/Schur orthogonality relation), returning how many times $n^{GH}_{ij}$ the subduced representation $D^G_i↓H$ contains the irrep $D^H_j$; in other words, this gives the compatibility between the two irreps.

Optionally, a vector αβγᴴⱼ may be provided, to evaluate the characters/irreps of Dᴳᵢ at a concrete value of $(α,β,γ)$. This is e.g. meaningful for LGIrreps at non-special k-vectors. Defaults to nothing.

The result is computed using the reduction formula [see e.g. Eq. (15) of arXiv:1706.09272v2]:

$n^{GH}_{ij} = |H|^{-1} \sum_h \chi^G_i(h)\chi^H_j(h)^*$

Example

Consider the two compatible k-vectors Γ (a point) and Σ (a line) in space group 207:

lgirsd  = lgirreps(207, Val(3));
Γ_lgirs = lgirsd["Γ"]; # at Γ ≡ [0.0, 0.0, 0.0]
Σ_lgirs = lgirsd["Σ"]; # at Σ ≡ [α, α, 0.0]

We can test their compatibility via:

[[subduction_count(Γi, Σj) for Γi in Γ_lgirs] for Σj in Σ_lgirs]
> # Γ₁ Γ₂ Γ₃ Γ₄ Γ₅
>  [ 1, 0, 1, 1, 2] # Σ₁
>  [ 0, 1, 1, 2, 1] # Σ₂

With following enterpretatation for compatibility relations between irreps at Γ and Σ:

Compatibility relationDegeneracies
Γ₁ → Σ₁1 → 1
Γ₂ → Σ₂1 → 1
Γ₃ → Σ₁ + Σ₂2 → 1 + 1
Γ₄ → Σ₁ + 2Σ₂3 → 1 + 2
Γ₅ → 2Σ₁ + Σ₂3 → 2 + 1

where, in this case, all the small irreps are one-dimensional.

source
Crystalline.subperiodicgroupMethod
subperiodicgroup(num::Integer, ::Val{D}=Val(3), ::Val{P}=Val(2))
subperiodicgroup(num::Integer, D::Integer, P::Integer)
                                                        --> ::SubperiodicGroup{D,P}

Return the operations of the subperiodic group num of embedding dimension D and periodicity dimension P as a SubperiodicGroup{D,P}.

The setting choices are those of the International Tables for Crystallography, Volume E.

Allowed combinations of D and P and their associated group names are:

  • D = 3, P = 2: Layer groups (num = 1, …, 80).
  • D = 3, P = 1: Rod groups (num = 1, …, 75).
  • D = 2, P = 1: Frieze groups (num = 1, …, 7).

Example

julia> subperiodicgroup(7, Val(2), Val(1))
SubperiodicGroup{2, 1} ⋕7 (𝓅2mg) with 4 operations:
 1
 2
 {m₁₀|½,0}
 {m₀₁|½,0}

Data sources

The symmetry operations returned by this function were originally retrieved from the Bilbao Crystallographic Database, SUBPERIODIC GENPOS.

source
Crystalline.translationMethod
translation(op::AbstractOperation{D}) --> SVector{D, Float64}

Return the D-dimensional translation part of op.

source
Crystalline.wyckoffsMethod
wyckoffs(sgnum::Integer) -> Vector{WyckoffPosition{3}}
wyckoffs(sgnum::Integer, ::Val{D}) -> Any

Return the Wyckoff positions of space group sgnum in dimension D as a Vector{WyckoffPosition{D}.

The positions are given in the conventional basis setting, following the conventions of the Bilbao Crystallographic Server (from which the underlying data is sourced [1]).

Example

julia> wps = wyckoffs(16, 2)
4-element Vector{WyckoffPosition{2}}:
 6d: [α, β]
 3c: [1/2, 0]
 2b: [1/3, 2/3]
 1a: [0, 0]

References

source

Exported macros

Crystalline.@S_strMacro
@S_str -> SymOperation

Construct a SymOperation from a triplet form given as a string.

Example

julia> S"-y,x"
4⁺ ──────────────────────────────── (-y,x)
 ┌ 0 -1 ╷ 0 ┐
 └ 1  0 ╵ 0 ┘
source

Exported constants

Crystalline.ENANTIOMORPHIC_PAIRSConstant
ENANTIOMORPHIC_PAIRS :: NTuple{11, Pair{Int,Int}}

Return the space group numbers of the 11 enantiomorphic space group pairs in 3D.

The space group types associated with each such pair (sgnum, sgnum') are related by a mirror transformation: i.e. there exists a transformation $\mathbb{P} = \{\mathbf{P}|\mathbf{p}\}$ between the two groups $G = \{g\}$ and $G' = \{g'\}$ such that $G' = \mathbb{P}^{-1}G\mathbb{P}$ where $\mathbf{P}$ is improper (i.e. $\mathrm{det}\mathbf{P} < 0$).

We define distinct space group types as those that cannot be related by a proper transformation (i.e. with $\mathrm{det}\mathbf{P} > 0$). With that view, there are 230 space group types. If the condition is relaxed to allow improper rotations, there are $230-11 = 219$ distinct affine space group types. See e.g. ITA5 Section 8.2.2.

The enantiomorphic space group types are also chiral space group types in 3D. There are no enantiomorphic pairs in lower dimensions; in 3D all enantiomorphic pairs involve screw symmetries, whose direction is inverted between pairs (i.e. have opposite handedness).

source
Crystalline.MAX_MSGNUMConstant
MAX_MSGNUM :: Tuple{Int,Int,Int}

Analogous to MAX_SUBGNUM, but for the number of magnetic space groups.

source
Crystalline.MAX_MSUBGNUMConstant
MAX_MSUBGNUM :: Tuple{Int,Int,Int}

Analogous to MAX_SUBGNUM, but for the number of magnetic subperiodic groups.

source
Crystalline.MAX_SGNUMConstant
MAX_SGNUM :: Tuple{Int,Int,Int}

Return the number of distinct space group types across dimensions 1, 2, and 3 (indexable by dimensionality).

source
Crystalline.MAX_SUBGNUMConstant
MAX_SUBGNUM :: ImmutableDict

An immutable dictionary with values v::Int and keys k::Tuple{Int,Int}, where v is the number of distinct subperiodic group types for a given key k = (D,P) describing a subperiodic group of dimensionality D and periodicity P:

  • layer groups: (D,P) = (3,2)
  • rod groups: (D,P) = (3,1)
  • frieze groups: (D,P) = (2,1)
source