Bravais.jl

Bravais types, basis systems, and transformations between conventional and primitive settings.

API

Types

Bravais.AbstractBasisType
AbstractBasis <: StaticVector{D, SVector{D,Float64}}

Abstract supertype of a D-dimensional basis in D-dimensional space.

source
Bravais.DirectBasisType
DirectBasis{D} <: AbstractBasis{D}

A wrapper type over D distinct D-dimensional vectors (given as a SVector{D, SVector{D,Float64}}), defining a lattice basis in direct space.

source
Bravais.ReciprocalBasisType
ReciprocalBasis{D} <: AbstractBasis{D}

A wrapper type over D distinct D-dimensional vectors (given as a SVector{D, SVector{D,Float64}}), defining a lattice basis in reciprocal space.

source
Bravais.AbstractPointType
AbstractPoint{D, T} <: StaticVector{D, T}

Abstract supertype of a D-dimensional point with elements of type T.

source
Bravais.DirectPointType
DirectPoint{D} <: AbstractPoint{D}

A wrapper type over an SVector{D, Float64}, defining a single point in D-dimensional direct space.

The coordinates of a DirectPoint are generally assumed specified relative to an associated DirectBasis. To convert to Cartesian coordinates, see cartesianize.

source
Bravais.ReciprocalPointType
ReciprocalPoint{D} <: AbstractPoint{D}

A wrapper type over an SVector{D, Float64}, defining a single point in D-dimensional reciprocal space.

The coordinates of a ReciprocalPoint are generally assumed specified relative to an associated ReciprocalBasis. To convert to Cartesian coordinates, see cartesianize.

source

Crystal systems & Bravais types

Bravais.crystalsystemFunction
crystalsystem(Rs::DirectBasis{D})  -->  String

Determine the crystal system of a point lattice Rs, assuming the conventional setting choice defined in the International Tables of Crystallography [ITA6].

There are 4 crystal systems in 2D and 7 in 3D (see Section 2.1.2(iii) of [ITA5]):

DSystemConditionsFree parameters
1Dlinearnonea
2Dsquarea=b & γ=90°a
rectangularγ=90°a,b
hexagonala=b & γ=120°a
obliquenonea,b,γ
3Dcubica=b=c & α=β=γ=90°a
hexagonala=b & α=β=90° & γ=120°a,c
trigonala=b & α=β=90° & γ=120°a,c (a,α for hR)
tetragonala=b & α=β=γ=90°a,c
orthorhombicα=β=γ=90°a,b,c
monoclinicα=γ=90°a,b,c,β≥90°
triclinicnonea,b,c,α,β,γ

The Rs must specify a set of conventional basis vectors, i.e., not generally primitive. For primitive basis vectors, the crystal system can be further reduced into 5 Bravais types in 2D and 14 in 3D (see bravaistype).

source
Bravais.bravaistypeFunction
bravaistype(sgnum::Integer, D::Integer=3; normalize::Bool=false)  -->  String

Return the Bravais type of sgnum in dimension D as a string (as the concatenation of the single-character crystal abbreviation and the centering type).

Keyword arguments

  • normalize: If the centering type associated with sgnum is 'A', we can choose (depending on the keyword argument normalize, defaulting to false) to "normalize" to the centering type 'C', since the difference between 'A' and 'C' centering only amounts to a basis change. With normalize=true we then have only the canonical 14 Bravais type, i.e. unique(bravaistype.(1:230, 3), normalize=true) returns only 14 distinct types, rather than 15.

    This only affects space groups 38-41 (normalizing their conventional Bravais types from "oA" to "oC").

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

Return the conventional centering type of a group.

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

source
centering(sgnum::Integer, D::Integer=3)  -->  Char

Return the conventional centering type cntr of the space group with number sgnum and dimension D.

The centering type is equal to the first letter of the Hermann-Mauguin notation's label, i.e., centering(sgnum, D) == first(Crystalline.iuc(sgnum, D)). Equivalently, the centering type is the second and last letter of the Bravais type (bravaistype), i.e., centering(sgnum, D) == bravaistype(sgnum, D).

Possible values of cntr, depending on dimensionality D, are (see ITA Sec. 9.1.4):

  • D = 1:
    • cntr = 'p': no centering (primitive)
  • D = 2:
    • cntr = 'p': no centring (primitive)
    • cntr = 'c': face centered
  • D = 3:
    • cntr = 'P': no centring (primitive)
    • cntr = 'I': body centred (innenzentriert)
    • cntr = 'F': all-face centred
    • cntr = 'A' or 'C': one-face centred, (b,c) or (a,b)
    • cntr = 'R': hexagonal cell rhombohedrally centred
source

Basis construction

Bravais.crystalFunction
crystal(a, b, c, α, β, γ)  -->  DirectBasis{3}

Calculate basis vectors $\mathbf{R}_1$, $\mathbf{R}_2$, $\mathbf{R}_3$ in a 3D Cartesian basis for a right-handed coordinate system with specified basis vector lengths a, b, c (associated with $\mathbf{R}_1$, $\mathbf{R}_2$, & $\mathbf{R}_3$, respectively) and specified interaxial angles α $= ∠(\mathbf{R}_2,\mathbf{R}_3)$, β $= ∠(\mathbf{R}_3,\mathbf{R}_1)$, γ $= ∠(\mathbf{R}_1,\mathbf{R}_2)$, with $∠$ denoting the angle between two vectors.

For definiteness, the $\mathbf{R}_1$ basis vector is oriented along the $x$-axis of the Cartesian coordinate system, and the $\mathbf{R}_2$ axis is placed in the $xy$-plane.

source
crystal(a, b, γ)  -->  DirectBasis{2}

Calculate basis vectors $\mathbf{R}_1$, $\mathbf{R}_2$ in a 2D Cartesian basis for a right-handed coordinate system with specified basis vector lengths a, b (associated with $\mathbf{R}_1$ & $\mathbf{R}_2$, respectively) and specified interaxial angle γ $= ∠(\mathbf{R}_1,\mathbf{R}_2)$.

For definiteness, the $\mathbf{R}_1$ basis vector is oriented along the $x$-axis of the Cartesian coordinate system.

source
crystal(a)  -->  DirectBasis{1}

Return a one-dimensional crystal with lattice period a.

source
Bravais.directbasisFunction
directbasis(sgnum, D=3;    abclims, αβγlims)
directbasis(sgnum, Val(D); abclims, αβγlims) --> DirectBasis{D}

Return a random (conventional) DirectBasis for a crystal compatible with the space group number sgnum and dimensionality D. Free parameters in the lattice vectors are chosen randomly, with limits optionally supplied in abclims (lengths) and αβγlims (angles). By convention, the length of the first lattice vector (a) is set to unity, such that the second and third (b and c) lattice vectors' lengths are relative to the first.

Limits on the relative uniform distribution of lengths b and c can be specified as 2-tuple kwarg abclims; similarly, limits on the angles α, β, γ can be set via αβγlims (only affects oblique, monoclinic, & triclinic lattices).

source
Bravais.reciprocalbasisFunction
reciprocalbasis(Rs)  -->  ::ReciprocalBasis{D}

Return the reciprocal basis of a direct basis Rs in D dimensions, provided as a StaticVector of AbstractVectors (e.g., a DirectBasis{D}) or a D-dimensional NTuple of AbstractVectors, or a (or, type-unstably, as any iterable of AbstractVectors).

source

Transformations

Bravais.primitivebasismatrixFunction
primitivebasismatrix(cntr::Char, ::Val{D}=Val(3)) --> SMatrix{D,D,Float64}
primitivebasismatrix(cntr::Char, ::Val{D}, ::Val{P}) --> SMatrix{D,D,Float64}

Return the transformation matrix $\mathbf{P}$ that transforms a conventional unit cell with centering cntr to the corresponding primitive unit cell (in dimension D and periodicity P) in CDML setting.

If P is not provided, it default to D (as e.g., applicable to Crystalline.jl's spacegroup). If D and P differ, a subperiodic group setting is assumed (as e.g., applicable to Crystalline.jl's subperiodicgroup).

Transformations in direct and reciprocal space

Bases

The returned transformation matrix $\mathbf{P}$ transforms a direct-space conventional basis $(\mathbf{a}\ \mathbf{b}\ \mathbf{c})$ to the direct-space primitive basis

\[ (\mathbf{a}'\ \mathbf{b}'\ \mathbf{c}') = (\mathbf{a}\ \mathbf{b}\ \mathbf{c})\mathbf{P}.\]

Analogously, $\mathbf{P}$ transforms a reciprocal-space conventional basis $(\mathbf{a}^*\ \mathbf{b}^*\ \mathbf{c}^*)$ to

\[(\mathbf{a}^{*\prime}\ \mathbf{b}^{*\prime}\ \mathbf{c}^{*\prime}) = (\mathbf{a}^*\ \mathbf{b}^*\ \mathbf{c}^*)(\mathbf{P}^{-1})^{\text{T}}.\]

see also transform(::DirectBasis, ::AbstractMatrix{<:Real}) and transform(::ReciprocalBasis, ::AbstractMatrix{<:Real})).

Coordinates

The coordinates of a point in either direct or reciprocal space, each referred to a basis, also transform under $\mathbf{P}$. Concretely, direct- and reciprocal-space conventional points $\mathbf{r} = (r_1, r_2, r_3)^{\text{T}}$ and $\mathbf{k} = (k_1, k_2, k_3)^{\text{T}}$, respectively, transform to a primitive setting under $\mathbf{P}$ according to:

\[\mathbf{r}' = \mathbf{P}^{-1}\mathbf{r},\\ \mathbf{k}' = \mathbf{P}^{\text{T}}\mathbf{k}.\]

See also transform(::DirectPoint, ::AbstractMatrix{<:Real}) and transform(::ReciprocalPoint, ::AbstractMatrix{<:Real})).

Setting conventions

The setting choice for the primitive cell implied by the returned $\mathbf{P}$ follows the widely adopted Cracknell-Davies-Miller-Love (CDML) convention.[CDML] This convention is explicated e.g. in Table 2 of [Aroyo] (or, alternatively, can be inferred from Tables 1.5.4.1 and 1.5.4.2 of [ITB2]) and is followed e.g. on the Bilbao Crystallographic Server[BCS], in the CDML reference work on space group irreps[CDML], and in the C library spglib.[spglib]

Note that this setting choice is not what is frequently referred to as the "ITA primitive setting", from which it differs for hP, hR, and oA Bravais types.

The setting choice is usually referred to as the CDML primitive setting, or, less frequently and more ambiguously, as the crystallographic primitive setting.

source
Bravais.transformFunction
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
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
transform(Rs::DirectBasis, P::AbstractMatrix{<:Real})

Transform a direct basis Rs $= (\mathbf{a}\ \mathbf{b}\ \mathbf{c})$ under the transformation matrix P $= \mathbf{P}$, returning Rs′ $= (\mathbf{a}'\ \mathbf{b}'\ \mathbf{c}') = (\mathbf{a}\ \mathbf{b}\ \mathbf{c}) \mathbf{P}$.

source
transform(Gs::ReciprocalBasis, P::AbstractMatrix{<:Real})

Transform a reciprocal basis Gs $= (\mathbf{a}^* \mathbf{b}^* \mathbf{c}^*)$ under the transformation matrix P $= \mathbf{P}$, returning Gs′ $= (\mathbf{a}^{*\prime}\ \mathbf{b}^{*\prime}\ \mathbf{c}^{*\prime}) = (\mathbf{a}^*\ \mathbf{b}^*\ \mathbf{c}^*)(\mathbf{P}^{-1})^{\text{T}}$.

source
transform(r::DirectPoint, P::AbstractMatrix{<:Real})  -->  r′::typeof(r)

Transform a point in direct space r $= (r_1, r_2, r_3)^{\text{T}}$ under the transformation matrix P $= \mathbf{P}$, returning r′ $= (r_1', r_2', r_3')^{\text{T}} = \mathbf{P}^{-1}(r_1', r_2', r_3')^{\text{T}}$.

source
transform(k::ReciprocalPoint, P::AbstractMatrix{<:Real})  -->  k′::typeof(k)

Transform a point in reciprocal space k $= (k_1, k_2, k_3)^{\text{T}}$ under the transformation matrix P $= \mathbf{P}$, returning k′ $= (k_1', k_2', k_3')^{\text{T}} = \mathbf{P}^{\text{T}}(k_1', k_2', k_3')^{\text{T}}$.

source
Bravais.primitivizeFunction
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
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
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
primitivize(V::Union{AbstractBasis, AbstractPoint}, 
            cntr_or_sgnum::Union{Char, <:Integer})   -->  V′::typeof(V)

Return the primitive basis or point V′ associated with the input conventional AbstractBasis or AbstractPoint V.

The assumed centering type is specified by cntr_or_sgnum, given either as a centering character (::Char) or inferred from a space group number (::Integer) and the dimensionality of V (see also centering(::Integer, ::Integer)).

source
primitivize(Rs::DirectBasis, cntr_or_sgnum::Union{Char, <:Integer}) --> Rs′::typeof(Rs)

Return the primitive direct basis Rs′ corresponding to the input conventional direct basis Rs.

source
primitivize(Gs::ReciprocalBasis, cntr_or_sgnum::Union{Char, <:Integer}) --> Gs′::typeof(Gs)

Return the primitive reciprocal basis Gs′ corresponding to the input conventional reciprocal basis Gs.

source
primitivize(r::DirectPoint, cntr_or_sgnum::Union{Char, <:Integer}) --> r′::typeof(r)

Return the direct point r′ with coordinates in a primitive basis, corresponding to the input point r with coordinates in a conventional basis.

source
primitivize(k::ReciprocalPoint, cntr_or_sgnum::Union{Char, <:Integer}) --> k′::typeof(k)

Return the reciprocal point k′ with coordinates in a primitive basis, corresponding to the input point k with coordinates in a conventional basis.

source
Bravais.conventionalizeFunction
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
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
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
conventionalize(V′::Union{AbstractBasis, AbstractPoint}, 
                cntr_or_sgnum::Union{Char, <:Integer})    -->  V::typeof(V′)

Return the conventional basis or point V associated with the input primitive AbstractBasis or AbstractPoint V′.

The assumed centering type is specified by cntr_or_sgnum, given either as a centering character (::Char) or inferred from a space group number (::Integer) and the dimensionality of V′ (see also centering(::Integer, ::Integer)).

source
conventionalize(Rs′::DirectBasis, cntr_or_sgnum::Union{Char, <:Integer}) --> Rs::typeof(Rs′)

Return the conventional direct basis Rs corresponding to the input primitive direct basis Rs′.

source
conventionalize(Gs′::ReciprocalBasis, cntr_or_sgnum::Union{Char, <:Integer}) --> Gs::typeof(Gs′)

Return the conventional reciprocal basis Gs corresponding to the input primitive reciprocal basis Gs′.

source
conventionalize(r′::DirectPoint, cntr_or_sgnum::Union{Char, <:Integer}) --> r::typeof(r′)

Return the direct point r with coordinates in a conventional basis, corresponding to the input point r′ with coordinates in a primitive basis.

source
conventionalize(k′::ReciprocalPoint, cntr_or_sgnum::Union{Char, <:Integer}) --> k::typeof(k′)

Return the reciprocal point k with coordinates in a conventional basis, corresponding to the input point k′ with coordinates in a primitive basis.

source
Bravais.cartesianizeFunction
cartesianize

Transform an object with coordinates in an lattice basis to an object with coordinates in a Cartesian basis.

Depending on the object, the basis may be inferrable directly from the object; if not, it must be supplied explicitly. @doc

source
Bravais.cartesianize!Function
cartesianize!

In-place transform an object with coordinates in an lattice basis to an object with coordinates in a Cartesian basis.

Depending on the object, the basis may be inferrable directly from the object; if not, it must be supplied explicitly.

source
Bravais.latticizeFunction
latticize

Transform an object with coordinates in a Cartesian basis to an object with coordinates in a lattice basis.

Depending on the object, the basis may be inferrable directly from the object; if not, it must be supplied explicitly.

source
Bravais.latticize!Function
latticize!

In-place transform object with coordinates in a Cartesian basis to an object with coordinates in a lattice basis.

Depending on the object, the basis may be inferrable directly from the object; if not, it must be supplied explicitly.

source

Miscellaneous

Bravais.volumeFunction
volume(Vs::AbstractBasis)

Return the volume $V$ of the unit cell associated with the basis Vs::AbstractBasis{D}.

The volume is computed as $V = \sqrt{\mathrm{det}\mathbf{G}}$ with with $\mathbf{G}$ denoting the metric matrix of Vs (cf. the International Tables of Crystallography, Volume A, Section 5.2.2.3).

See also metricmatrix.

source
Bravais.metricmatrixFunction
metricmatrix(Vs::AbstractBasis)

Return the (real, symmetric) metric matrix of a basis Vs, i.e., the matrix with elements $G_{ij} =$ dot(Vs[i], Vs[j]), as defined in the International Tables of Crystallography, Volume A, Section 5.2.2.3.

Equivalently, this is the Gram matrix of Vs, and so can also be expressed as Vm' * Vm with Vm denoting the columnwise concatenation of the basis vectors in Vs.

See also volume.

source
Bravais.stackFunction
stack(Vs::AbstractBasis)

Return a matrix [Vs[1] Vs[2] .. Vs[D]] from Vs::AbstractBasis{D}, i.e., the matrix whose columns are the basis vectors of Vs.

source

Crystalline.jl extensions of Bravais.jl functions

SymOperation

Bravais.transformFunction
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
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.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

AbstractVec

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.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.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

AbstractFourierLattice

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.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
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th ed. (2016): Tables 3.1.2.1 and 3.1.2.2 (or Tables 2.1.2.1, 9.1.7.1, and 9.1.7.2 of [ITA5]).
  • ITA5T. Hahn, International Tables of Crystallography, Vol. A, 5th ed. (2005).
  • CDMLCracknell, Davies, Miller, & Love, Kroenecker Product Tables, Vol. 1 (1979).
  • BCSBilbao Crystallographic Server, KVEC.
  • AroyoAroyo et al., Acta Cryst. A70, 126 (2014): Table 2 gives $(\mathbf{P}^{-1})^{\text{T}}$.
  • ITB2Hahn, International Tables of Crystallography, Vol. B, 2nd edition (2001).
  • spglibSpglib documentation: Transformation to the primitive setting. Thus, Bravais.jl and Spglib.jl transform to identical primitive settings and are hence mutually compatible.
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th edition (2016): Secs. 1.5.1.2 and 1.5.2.1.
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th ed. (2016).
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th edition (2016): Secs. 1.5.1.2 and 1.5.2.1.
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th ed. (2016).
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th ed. (2016).
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th ed. (2016).
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th edition (2016): Secs. 1.5.1.2 and 1.5.2.1.
  • ITA6M.I. Aroyo, International Tables of Crystallography, Vol. A, 6th edition (2016): Secs. 1.5.1.2 and 1.5.2.1.