Documentation

Index

Lattices

Plotting

Symmetry

Utilites

Functions

Lattices

SymmetryReduceBZ.Lattices.check_reducedMethod
check_reduced(basis)

Verify a lattice basis is Minkowski reduced

Arguments

  • basis::AbstractMatrix{<:Real}: the lattice basis given by the columns of a 2x2 or 3x3 matrix.

Returns

  • Bool: a boolean that indicates if the lattice basis is reduced.

Examples

using SymmetryReduceBZ
basis = [1 0; 0 1]
SymmetryReduceBZ.Lattices.check_reduced(basis)
# output
true
source
SymmetryReduceBZ.Lattices.genlat_BCCFunction
genlat_BCC(a)

Generate a body-centered cubic lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
SymmetryReduceBZ.Lattices.genlat_BCC(a)
# output
3×3 Matrix{Float64}:
 -0.5   0.5   0.5
  0.5  -0.5   0.5
  0.5   0.5  -0.5
source
SymmetryReduceBZ.Lattices.genlat_BCTFunction
genlat_BCT(a,c)

Generate a body-centered tetragonal lattice.

Arguments

  • a::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
c=1.2;
SymmetryReduceBZ.Lattices.genlat_BCT(a,c)
# output
3×3 Matrix{Float64}:
 -0.5   0.5   0.5
  0.5  -0.5   0.5
  0.6   0.6  -0.6
source
SymmetryReduceBZ.Lattices.genlat_CUBFunction
genlat_CUB(a)

Generate a simple cubic lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
SymmetryReduceBZ.Lattices.genlat_CUB(a)
# output
3×3 Matrix{Int64}:
 1  0  0
 0  1  0
 0  0  1
source
SymmetryReduceBZ.Lattices.genlat_FCCFunction
genlat_FCC(a)

Generate a face-centered cubic lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
SymmetryReduceBZ.Lattices.genlat_FCC(a)
# output
3×3 Matrix{Float64}:
 0.0  0.5  0.5
 0.5  0.0  0.5
 0.5  0.5  0.0
source
SymmetryReduceBZ.Lattices.genlat_HEXFunction
genlat_HEX(a,c)

Generate a hexagonal lattice.

Arguments

  • a::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
c=1.2;
SymmetryReduceBZ.Lattices.genlat_HEX(a,c)
# output
3×3 Matrix{Float64}:
  0.5       0.5       0.0
 -0.866025  0.866025  0.0
  0.0       0.0       1.2
source
SymmetryReduceBZ.Lattices.genlat_HXGFunction
genlat_HXG(a)

Generate a 2D hexagonal lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
SymmetryReduceBZ.Lattices.genlat_HXG(a)
# output
2×2 Matrix{Float64}:
 1.0  -0.5
 0.0   0.866025
source
SymmetryReduceBZ.Lattices.genlat_MCLFunction
genlat_MCL(a,b,c,α)

Generate a monoclinic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • α::Real: a lattice angle in radians less than π/2
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
c=1.4
α=π/6;
SymmetryReduceBZ.Lattices.genlat_MCL(a,b,c,α)
# output
3×3 Matrix{Float64}:
 1.0  0.0  0.0
 0.0  1.2  1.21244
 0.0  0.0  0.7
source
SymmetryReduceBZ.Lattices.genlat_MCLCFunction
genlat_MCLC(a,b,c,α)

Generate a base-centered monoclinic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • α::Real: a lattice angle in radians less than π/2
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
c=1.4
α=π/6;
SymmetryReduceBZ.Lattices.genlat_MCLC(a,b,c,α)
# output
3×3 Matrix{Float64}:
 0.5  -0.5  0.0
 0.6   0.6  1.21244
 0.0   0.0  0.7
source
SymmetryReduceBZ.Lattices.genlat_OBLFunction
genlat_OBL(a,b,θ)

Generate an oblique lattice.

Arguments

  • a::Real: the lattice constant
  • b::Real: the lattice constant
  • θ::Real: the lattice angle
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
θ=π/3
SymmetryReduceBZ.Lattices.genlat_OBL(a,b,θ)
# output
2×2 Matrix{Float64}:
 1.0  0.6
 0.0  1.03923
source
SymmetryReduceBZ.Lattices.genlat_ORCFunction
genlat_ORC(a,b,c)

Generate an orthorhombic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.4;
c=1.2;
SymmetryReduceBZ.Lattices.genlat_ORC(a,b,c)
# output
3×3 Matrix{Float64}:
 1.0  0.0  0.0
 0.0  1.4  0.0
 0.0  0.0  1.2
source
SymmetryReduceBZ.Lattices.genlat_ORCCFunction
genlat_ORCC(a,b,c)

Generate a base-centered orthorhombic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2;
c=1.4;
SymmetryReduceBZ.Lattices.genlat_ORCC(a,b,c)
# output
3×3 Matrix{Float64}:
  0.5  0.5  0.0
 -0.6  0.6  0.0
  0.0  0.0  1.4
source
SymmetryReduceBZ.Lattices.genlat_ORCFFunction
genlat_ORCF(a,b,c)

Generate a face-centered orthorhombic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.4;
c=1.2;
SymmetryReduceBZ.Lattices.genlat_ORCF(a,b,c)
# output
3×3 Matrix{Float64}:
 0.0  0.5  0.5
 0.7  0.0  0.7
 0.6  0.6  0.0
source
SymmetryReduceBZ.Lattices.genlat_ORCIFunction
genlat_ORCI(a,b,c)

Generate a body-centered orthorhombic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.4;
c=1.2;
SymmetryReduceBZ.Lattices.genlat_ORCI(a,b,c)
# output
3×3 Matrix{Float64}:
 -0.5   0.5   0.5
  0.7  -0.7   0.7
  0.6   0.6  -0.6
source
SymmetryReduceBZ.Lattices.genlat_RECFunction
genlat_REC(a,b)

Generate a rectangular lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
SymmetryReduceBZ.Lattices.genlat_REC(a,b)
# output
2×2 Matrix{Float64}:
 1.0  0.0
 0.0  1.2
source
SymmetryReduceBZ.Lattices.genlat_RECIFunction
genlat_RECI(a,b)

Generate a body-centered rectangular lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
SymmetryReduceBZ.Lattices.genlat_RECI(a,b)
# output
2×2 Matrix{Float64}:
  0.5  0.5
 -0.6  0.6
source
SymmetryReduceBZ.Lattices.genlat_RHLFunction
genlat_RHL(a,α)

Generate a rhombohedral lattice.

Arguments

  • a::Real: a lattice constant
  • α::Real: a lattice angle in radians
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
α=π/6;
SymmetryReduceBZ.Lattices.genlat_RHL(a,α)
# output
3×3 Matrix{Float64}:
  0.965926  0.965926  0.896575
 -0.258819  0.258819  0.0
  0.0       0.0       0.442891
source
SymmetryReduceBZ.Lattices.genlat_SQRFunction
genlat_SQR(a)

Generate a square lattice.

Arguments

  • a::Real: the lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
SymmetryReduceBZ.Lattices.genlat_SQR(a)
# output
2×2 Matrix{Int64}:
 1  0
 0  1
source
SymmetryReduceBZ.Lattices.genlat_TETFunction
genlat_TET(a,c)

Generate a simple tetragonal lattice.

Arguments

  • a::Real: a lattice constant
  • c::Real: a lattice constant
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
c=1.2;
SymmetryReduceBZ.Lattices.genlat_TET(a,c)
# output
3×3 Matrix{Float64}:
 1.0  0.0  0.0
 0.0  1.0  0.0
 0.0  0.0  1.2
source
SymmetryReduceBZ.Lattices.genlat_TRIFunction
genlat_TRI(a,b,c,α,β,γ)

Generate a triclinic lattice.

Arguments

  • a::Real: a lattice constant
  • b::Real: a lattice constant
  • c::Real: a lattice constant
  • α::Real: a lattice angle in radians
  • β::Real: a lattice angle in radians
  • γ::Real: a lattice angle in radians
  • type::String="primitive": the lattice type: "primitive" or "conventional".

Returns

  • AbstractMatrix{<:Real}: the basis of the lattice as columns of an array.

Examples

using SymmetryReduceBZ
a=1
b=1.2
c=1.4
α=π/6;
β=π/3;
γ=π/4;
SymmetryReduceBZ.Lattices.genlat_TRI(a,b,c,α,β,γ)
# output
3×3 Matrix{Float64}:
 1.0  0.848528  0.7
 0.0  0.848528  1.01464
 0.0  0.0       0.663702
source
SymmetryReduceBZ.Lattices.get_latparamsMethod
get_latparams(latvecs)

Calculate the lattice constants and angles of a lattice basis.

Arguments

  • latvecs::AbstractMatrix{<:Real}: the lattice basis as columns of an array.

Returns

  • A list where the first element is a list lattice constants (a,b,c) and second lattice angles in radians (α,β,γ).

Examples

using SymmetryReduceBZ
latvecs = [1 0; 0 1]
SymmetryReduceBZ.Lattices.get_latparams(latvecs)
# output
2-element Vector{Vector{Float64}}:
 [1.0, 1.0]
 [1.5707963267948966, 1.5707963267948966]
source
SymmetryReduceBZ.Lattices.get_recip_latvecsFunction
get_recip_latvecs(real_latvecs, convention)

Calculate the reciprocal lattice vectors.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the real-space lattice vectors or primitive translation vectors as columns of a 2x2 or 3x3 array.
  • convention::String="ordinary": the convention used to go between real and reciprocal space. The two conventions are ordinary (temporal) frequency and angular frequency. The transformation from real to reciprocal space is unitary if the convention is ordinary.

Returns

  • recip_latvecs::Array{<:Real,2} the reciprocal lattice vectors (reciprocal primitive translation vectors) as columns of a 2x2 or 3x3 array.

Examples

using SymmetryReduceBZ
real_latvecs=[1 0 0; 0 1 0; 0 0 1]
convention="angular"
SymmetryReduceBZ.Lattices.get_recip_latvecs(real_latvecs,convention)
# output
3×3 Matrix{Float64}:
 6.28319  0.0      0.0
 0.0      6.28319  0.0
 0.0      0.0      6.28319
source
SymmetryReduceBZ.Lattices.minkowski_reduceMethod
minkowski_reduce(basis;rtol,atol)

Minkowski reduce a lattice basis. Follows the logic of Fig. 4 in "Low-Dimensional Lattice Basis Reduction Revisited" by Nguyen, 2009.

Arguments

  • basis::AbstractMatrix{<:Real}: the lattice basis given by the columns of a 2x2 or 3x3 array.
  • rtol::Real=sqrt(eps(float(maximum(basis)))): a relative tolerance.
  • atol::Real=1e-9: an absolute tolerance.

Returns

  • rbasis:: the Minkowski reduced lattice basis as columns of an array.

Examples

using SymmetryReduceBZ
basis = [1 2 0; 0 1 0; 0 0 1]
SymmetryReduceBZ.Lattices.minkowski_reduce(basis)
# output
3×3 Matrix{Int64}:
 0  1  0
 0  0  1
 1  0  0
source
SymmetryReduceBZ.Lattices.reduce_basis!Method
reduce_basis!(basis,k;rtol,atol)

Reduces the kth lattice vector. This is accomplished by locating the lattice point closest to the projection of the kth lattice vector onto the line or plane given by the other lattice vector(s), subtracting the closest lattice point from the kth lattice vector, and reordering the lattice vectors by increasing Euclidean norms.

Arguments

  • basis::AbstractMatrix{<:Real}: the lattice basis as columns of an array.
  • k::Int: Keeps track of which lattice vector needs to be reduced.
  • rtol::Real=sqrt(eps(float(maximum(basis)))): a relative tolerance.
  • atol::Real=1e-9: an absolute tolerance.

Returns

  • basis::AbstractMatrix{<:Real}: the partially reduced lattice basis as columns of an array.
  • k::Int: The index of the lattice vector that needs to be reduced next.

Examples

using SymmetryReduceBZ
basis = Array([1 2 0; 0 1 0; 3 2 1]')
k=2
SymmetryReduceBZ.Lattices.reduce_basis!(basis,k)
basis
# output
3×3 Matrix{Int64}:
 0  1  3
 1  2  2
 0  0  1
source
SymmetryReduceBZ.Lattices.lattice_typesConstant

A list of lattice types. Follows the naming convention in the article High-throughput electronic band structure calculations: Challenges and tools by Wahyu Setyawan and Stefano Curtarolo except triclinic lattices have "β" instead of "b" as subscripts.

source

Plotting

SymmetryReduceBZ.Plotting.plot_2DconvexhullFunction
plot_2Dconvexhull(convexhull, ax, color)

Plot a 2D convex hull

Julia 1.9 and above

This function is implemented in a package extension and requires using PyPlot.

Arguments

  • convexhull: a polyhedron or convex hull object.
  • ax::PyObject: an axes object from matplotlib.
  • facecolor::String="blue": the color of the area within the convex hull.
  • alpha::Real=0.3: the transparency of the convex hull.
  • linewidth::Real=3: the width of the edges.
  • edgecolor::String="black": the color of the edges.

Returns

  • ax::PyObject: updated ax that includes a plot of the convex hull.

Examples

ENV["MPLBACKEND"]="qt5agg"
using PyPlot
using SymmetryReduceBZ.Symmetry: calc_bz, calc_ibz
using SymmetryReduceBZ.Plotting: plot_2Dconvexhull
real_latvecs = [1 0; 0 1]
convention="ordinary"
atom_types=[0]
atom_pos = Array([0 0]')
coords = "Cartesian"
makeprim=false
bz = calc_bz(real_latvecs,atom_types,atom_pos,coords,makeprim,convention)
ibz = calc_ibz(real_latvecs,atom_types,atom_pos,coords,makeprim,convention)
ax = plot_2Dconvexhull(bz,facecolor="deepskyblue",linewidth=3,edgecolor="cyan",alpha=0.2)
ax = plot_2Dconvexhull(ibz,ax;facecolor="coral",linewidth=3,edgecolor="magenta",alpha=0.4)
# output
PyObject <AxesSubplot: >
source
SymmetryReduceBZ.Plotting.plot_3DconvexhullFunction
plot_3Dconvexhull(convexhull,ax;color)

Plot a 3D convex hull

Julia 1.9 and above

This function is implemented in a package extension and requires using PyPlot.

Arguments

  • convexhull: a polyhedron or convex hull object.
  • ax::PyObject: an axes object from matplotlib.
  • facecolors::String="blue": the color of the faces of the convex hull.
  • alpha::Real=0.3: the transparency of the faces of the convex hull.
  • linewidths::Real=1: the width of the edges of the convex hull.
  • edgecolors::String="black": the color of the edges of the convex hull.

Returns

  • ax::PyObject: updated ax that includes a plot of the convex hull.

Examples

ENV["MPLBACKEND"]="qt5agg"
using PyPlot
using SymmetryReduceBZ.Symmetry: calc_bz, calc_ibz
using SymmetryReduceBZ.Plotting: plot_3Dconvexhull
real_latvecs = [1 0 0; 0 1 0; 0 0 1]
convention="ordinary"
atom_types=[0]
atom_pos = Array([0 0 0]')
coords = "Cartesian"
makeprim=false
bz = calc_bz(real_latvecs,atom_types,atom_pos,coords,makeprim,convention)
ibz = calc_ibz(real_latvecs,atom_types,atom_pos,coords,makeprim,convention)
using3D()
fig = figure()
ax = fig.add_subplot(111, projection="3d")
ax = plot_3Dconvexhull(ibz,ax,facecolors="coral",alpha=1,edgecolors="black",linewidths = 1)
ax = plot_3Dconvexhull(bz,ax,facecolors="deepskyblue",edgecolors="white",linewidths=1,alpha=0.2)
# output
PyObject <Axes3DSubplot: >
source
SymmetryReduceBZ.Plotting.plot_convexhullsFunction
plot_convexhulls(real_latvecs,atom_types,atom_pos,coords,makeprim,convention;rtol,atol)

Plot the Brillouin and Irreducible Brillouin zone in 2D or 3D.

Julia 1.9 and above

This function is implemented in a package extension and requires using PyPlot.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the basis of a real-space lattice as columns of a matrix.
  • atom_types:AbstractVector{<:Int}: a list of atom types as integers.
  • atom_pos::AbstractMatrix{<:Real}: the positions of atoms in the crystal structure as columns of a matrix.
  • coords::String: indicates the positions of the atoms are in "lattice" or "Cartesian" coordinates.
  • makeprim::Bool=false: make the unit cell primitive before calculating the the IBZ if equal to true.
  • convention::String="ordinary": the convention used to go between real and reciprocal space. The two conventions are ordinary (temporal) frequency and angular frequency. The transformation from real to reciprocal space is unitary if the convention is ordinary.
  • library::Polyhedra.Library=CDDLib.Library(): a polyhedron manipulation library
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))) a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • ax::PyObject: an updated ax with plots of the BZ and IBZ.

Examples

ENV["MPLBACKEND"]="qt5agg"
using PyPlot
using SymmetryReduceBZ
real_latvecs = [1 0; .5 1]
atom_types=[0]
atom_pos = Array([0 0]')
coords = "Cartesian"
makeprim = true
convention = "ordinary"
ax=plot_convexhulls(real_latvecs,atom_types,atom_pos,coords,makeprim,convention)
# output
PyObject <AxesSubplot: >
source

Symmetry

SymmetryReduceBZ.Symmetry.calc_bzFunction
calc_bz(real_latvecs,atom_types,atom_pos,coordinates,bzformat,makeprim,
    convention,library;rtol,atol)

Calculate the Brillouin zone for the given real-space lattice basis.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the real-space lattice vectors or primitive translation vectors as columns of a 2x2 or 3x3 matrix.
  • atom_typesAbstractVector{<:Int}: a list of atom types as integers.
  • atom_pos::AbstractMatrix{<:Real}: the positions of atoms in the crystal structure as columns of a matrix.
  • coordinates::String: indicates the positions of the atoms are in "lattice" or "Cartesian" coordinates.
  • bzformat::String: the format of the Brillouin zone. Options include "convex hull" and "half-space".
  • makeprim::Bool=false: make the unit cell primitive before calculating the the BZ if equal to true.
  • convention::String="ordinary": the convention used to go between real and reciprocal space. The two conventions are ordinary (temporal) frequency and angular frequency. The transformation from real to reciprocal space is unitary if the convention is ordinary.
  • library::Polyhedra.Library=CDDLib.Library(): a polyhedron manipulation library
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))) a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • bz: a polyhedron conforming to the Polyhedra.jl interface that represents the convex hull of the Brillouin zone.

Examples

using SymmetryReduceBZ
real_latvecs = [1 0; 0 1]
convention="ordinary"
atom_types=[0]
atom_pos = Array([0 0]')
coordinates = "Cartesian"
makeprim=false
SymmetryReduceBZ.Symmetry.calc_bz(real_latvecs,atom_types,atom_pos,coordinates,
    makeprim,convention)
# output
Polyhedron CDDLib.Polyhedron{Float64}:
25-element iterator of Polyhedra.HalfSpace{Float64, Vector{Float64}}:
 HalfSpace([-2.0, -2.0], 4.0)
 HalfSpace([-2.0, -1.0], 2.5)
 HalfSpace([-2.0, 0.0], 2.0)
 HalfSpace([-2.0, 1.0], 2.5)
 HalfSpace([-2.0, 2.0], 4.0)
 HalfSpace([-1.0, -2.0], 2.5)
 HalfSpace([-1.0, -1.0], 1.0)
 HalfSpace([-1.0, 0.0], 0.5)
 HalfSpace([-1.0, 1.0], 1.0)
 HalfSpace([-1.0, 2.0], 2.5)
 HalfSpace([0.0, -2.0], 2.0)
 HalfSpace([0.0, -1.0], 0.5)
 HalfSpace([0.0, 0.0], 0.0)
 HalfSpace([0.0, 1.0], 0.5)
 HalfSpace([0.0, 2.0], 2.0)
 HalfSpace([1.0, -2.0], 2.5)
 HalfSpace([1.0, -1.0], 1.0)
 HalfSpace([1.0, 0.0], 0.5)
 HalfSpace([1.0, 1.0], 1.0)
 HalfSpace([1.0, 2.0], 2.5)
  ⋮:
4-element iterator of Vector{Float64}:
 [0.5, 0.5]
 [0.5, -0.5]
 [-0.5, -0.5]
 [-0.5, 0.5]
source
SymmetryReduceBZ.Symmetry.calc_ibzFunction
calc_ibz(real_latvecs,atom_types,atom_pos,coordinates,ibzformat,makeprim,
    convention,library;rtol,atol)

Calculate the irreducible Brillouin zone of a crystal structure in 2D or 3D.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the basis of a real-space lattice as columns of a matrix.
  • atom_types:AbstractVector{<:Int}: a list of atom types as integers.
  • atom_pos::AbstractMatrix{<:Real}: the positions of atoms in the crystal structure as columns of a matrix.
  • coordinates::String: indicates the positions of the atoms are in "lattice" or "Cartesian" coordinates.
  • ibzformat::String: the format of the irreducible Brillouin zone. Options include "convex hull" and "half-space".
  • makeprim::Bool=false: make the unit cell primitive before calculating the IBZ if true.
  • convention::String="ordinary": the convention used to go between real and reciprocal space. The two conventions are ordinary (temporal) frequency and angular frequency. The transformation from real to reciprocal space is unitary if the convention is ordinary.
  • library::Polyhedra.Library=CDDLib.Library(): a polyhedron manipulation library
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))) a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • ibz: the irreducible Brillouin zone as a polyhedron conforming to the Polyhedra.jl interface.

Examples

using SymmetryReduceBZ
real_latvecs = [1 0; 0 1]
convention="ordinary"
atom_types=[0]
atom_pos = Array([0 0]')
coordinates = "Cartesian"
makeprim=false
SymmetryReduceBZ.Symmetry.calc_ibz(real_latvecs,atom_types,atom_pos,coordinates,
    makeprim,convention)
# output
Polyhedron CDDLib.Polyhedron{Float64}:
32-element iterator of Polyhedra.HalfSpace{Float64, Vector{Float64}}:
 HalfSpace([-2.0, -2.0], 4.0)
 HalfSpace([-2.0, -1.0], 2.5)
 HalfSpace([-2.0, 0.0], 2.0)
 HalfSpace([-2.0, 1.0], 2.5)
 HalfSpace([-2.0, 2.0], 4.0)
 HalfSpace([-1.0, -2.0], 2.5)
 HalfSpace([-1.0, -1.0], 1.0)
 HalfSpace([-1.0, 0.0], 0.5)
 HalfSpace([-1.0, 1.0], 1.0)
 HalfSpace([-1.0, 2.0], 2.5)
 HalfSpace([0.0, -2.0], 2.0)
 HalfSpace([0.0, -1.0], 0.5)
 HalfSpace([0.0, 0.0], 0.0)
 HalfSpace([0.0, 1.0], 0.5)
 HalfSpace([0.0, 2.0], 2.0)
 HalfSpace([1.0, -2.0], 2.5)
 HalfSpace([1.0, -1.0], 1.0)
 HalfSpace([1.0, 0.0], 0.5)
 HalfSpace([1.0, 1.0], 1.0)
 HalfSpace([1.0, 2.0], 2.5)
  ⋮:
3-element iterator of Vector{Float64}:
 [0.0, 0.0]
 [0.5, 0.0]
 [0.5, 0.5]
source
SymmetryReduceBZ.Symmetry.calc_pointgroupMethod
calc_pointgroup(latvecs;rtol,atol)

Calculate the point group of a lattice in 2D or 3D.

Arguments

  • latvecs::AbstractMatrix{<:Real}: the basis of the lattice as columns of a matrix.
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))): a relative tolerance for floating point comparisons. It is used to compare lengths of vectors and the volumes of primitive cells.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons. It is used to compare lengths of vectors and the volumes of primitive cells.

Returns

  • pointgroup::AbstractVector{<:AbstractMatrix{<:Real}}: the point group of the lattice. The operators operate on points in Cartesian coordinates from the right.

Examples

using SymmetryReduceBZ
basis = [1 0; 0 1]
SymmetryReduceBZ.Symmetry.calc_pointgroup(basis)
# output
8-element Vector{Matrix{Float64}}:
 [0.0 -1.0; -1.0 0.0]
 [0.0 -1.0; 1.0 0.0]
 [-1.0 0.0; 0.0 -1.0]
 [1.0 0.0; 0.0 -1.0]
 [-1.0 0.0; 0.0 1.0]
 [1.0 0.0; 0.0 1.0]
 [0.0 1.0; -1.0 0.0]
 [0.0 1.0; 1.0 0.0]
source
SymmetryReduceBZ.Symmetry.calc_spacegroupMethod
calc_spacegroup(real_latvecs,atom_types,atom_pos,coordinates;rtol,atol)

Calculate the space group of a crystal structure.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the basis of the lattice as columns of a matrix.
  • atom_types::AbstractVector{<:Int}: a list of atom types as integers.
  • atom_pos::AbstractMatrix{<:Real}: the positions of atoms in the crystal structure as columns of a matrix.
  • coordinates::String: indicates the positions of the atoms are in "lattice" or "Cartesian" coordinates.
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))) a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • spacegroup::Tuple: the space group of the crystal structure. The first element of spacegroup is a list of fractional translations, and the second element is a list of point operators. The translations are in Cartesian coordinates, and the operators operate on points in Cartesian coordinates.

Examples

using SymmetryReduceBZ
real_latvecs = Array([1 0; 2 1]')
atom_types = [0, 1]
atom_pos = Array([0 0; 0.5 0.5]')
coordinates = "Cartesian"
SymmetryReduceBZ.Symmetry.calc_spacegroup(real_latvecs,atom_types,atom_pos,
    coordinates)
# output
([[0.0, 0.0], [0.0, 0.0], [0.0, 0.0], [0.0, 0.0], [0.0, 0.0], [0.0, 0.0], [0.0, 0.0], [0.0, 0.0]], [[0.0 -1.0; -1.0 0.0], [0.0 -1.0; 1.0 0.0], [-1.0 0.0; 0.0 -1.0], [1.0 0.0; 0.0 -1.0], [-1.0 0.0; 0.0 1.0], [1.0 0.0; 0.0 1.0], [0.0 1.0; -1.0 0.0], [0.0 1.0; 1.0 0.0]])
source
SymmetryReduceBZ.Symmetry.complete_orbitMethod
complete_orbit(pt,pointgroup,rtol=sqrt(eps(float(maximum(pt)))),atol=1e-9)

Complete the orbits of multiple points.

Arguments

  • pt::AbstractMatrix{<:Real}: the Cartesian coordinates of a points as columns of a matrix.
  • pointgroup::AbstractVector{<:AbstractMatrix{<:Real}}: the point group operators in a nested list. The operators operate on points in Cartesian coordinates from the right.
  • rtol::Real=sqrt(eps(float(maximum(pt)))): a relative tolerance.
  • atol::Real=1e-9: an absolute tolerance.

Returns

  • ::AbstractMatrix{<:Real}: the unique points of the orbits in Cartesian coordinates as columns of a matrix.

Examples

import SymmetryReduceBZ.Symmetry: complete_orbit
pts = [0.0 0.05 0.1; 0.0 0.0 0.0]
pointgroup = [[0.0 -1.0; -1.0 0.0], [0.0 -1.0; 1.0 0.0], [-1.0 0.0; 0.0 -1.0], [1.0 0.0; 0.0 -1.0], [-1.0 0.0; 0.0 1.0], [1.0 0.0; 0.0 1.0], [0.0 1.0; -1.0 0.0], [0.0 1.0; 1.0 0.0]]
complete_orbit(pts,pointgroup)
# output
2×9 Matrix{Float64}:
 0.0   0.0   0.0   -0.05  0.05   0.0  0.0  -0.1  0.1
 0.0  -0.05  0.05   0.0   0.0   -0.1  0.1   0.0  0.0
source
SymmetryReduceBZ.Symmetry.complete_orbitMethod
complete_orbit(pt,pointgroup,rtol=sqrt(eps(float(maximum(pt)))),atol=1e-9)

Complete the orbit of a point.

Arguments

  • pt::AbstractVector{<:Real}: the Cartesian coordinates of a point.
  • pointgroup::AbstractVector{<:AbstractMatrix{<:Real}}: the point group operators in a nested list. The operators operate on points in Cartesian coordinates from the right.
  • rtol::Real=sqrt(eps(float(maximum(pt)))): a relative tolerance.
  • atol::Real=1e-9: an absolute tolerance.

Returns

  • ::AbstractMatrix{<:Real}: the points of the orbit in Cartesian coordinates as columns of a matrix.

Examples

import SymmetryReduceBZ.Symmetry: complete_orbit
pt = [0.05, 0.0]
pointgroup = [[0.0 -1.0; -1.0 0.0], [0.0 -1.0; 1.0 0.0], [-1.0 0.0; 0.0 -1.0], [1.0 0.0; 0.0 -1.0], [-1.0 0.0; 0.0 1.0], [1.0 0.0; 0.0 1.0], [0.0 1.0; -1.0 0.0], [0.0 1.0; 1.0 0.0]]
complete_orbit(pt,pointgroup)
# output
2×4 Matrix{Float64}:
  0.0   0.0   -0.05  0.05
 -0.05  0.05   0.0   0.0
source
SymmetryReduceBZ.Symmetry.inhullMethod
inhull(point, chull; rtol, atol)

Check if a point lies within a convex hull (including the boundaries).

Arguments

  • point::AbstractVector{<:Real}: a point in Cartesian coordinates.
  • chull::Chull: a convex hull in 2D or 3D.
  • rtol::Real=sqrt(eps(float(maximum(flatten(chull.points))))): a relative tolerance for floating point comparisons. Needed when a point is on the boundary of the convex hull.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • inside::Bool: if true, the point lies within the convex hull.

Examples

import QHull: chull
import SymmetryReduceBZ.Symmetry: inhull
pts = [0.5 0.25; 0.5 -0.25; -0.5 -0.25; -0.5 0.25]
pt = [0,0]
convexhull = chull(pts)
inhull(pt,convexhull)
# output
true
source
SymmetryReduceBZ.Symmetry.make_primitiveMethod
make_primitive(real_latvecs,atom_types,atom_pos,coordinates;rtol,atol)

Make a given unit cell primitive.

This is a Julia translation of the function by the same in https://github.com/msg-byu/symlib.

Arguments

  • real_latvecs::AbstractMatrix{<:Real}: the basis of the lattice as columns of a matrix.
  • atom_types::AbstractVector{<:Int}: a list of atom types as integers.
  • atom_pos::AbstractMatrix{<:Real}: the positions of atoms in the crystal structure as columns of a matrix.
  • coordinates::String: indicates the positions of the atoms are in "lattice" or "Cartesian" coordinates.
  • rtol::Real=sqrt(eps(float(maximum(real_latvecs)))) a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • prim_latvecs::AbstractMatrix{<:Real}: the primitive lattice vectors as columns of a matrix.
  • prim_types::AbstractVector{<:Int}: a list of atom types as integers in the primitive unit cell.
  • prim_pos::AbstractMatrix{<:Real}: the positions of the atoms in in the crystal structure as columns of a matrix in Cartesian coordinates.

Examples

import SymmetryReduceBZ.Lattices: genlat_CUB
import SymmetryReduceBZ.Symmetry: make_primitive
a = 1.0
real_latvecs = genlat_CUB(a)
atom_types = [0,0]
atom_pos = Array([0 0 0; 0.5 0.5 0.5]')
ibzformat = "convex hull"
coordinates = "Cartesian"
convention = "ordinary"
make_primitive(real_latvecs, atom_types, atom_pos, coordinates)
# output
([1.0 0.0 0.5; 0.0 1.0 0.5; 0.0 0.0 0.5], [0], [0.0; 0.0; 0.0])
source
SymmetryReduceBZ.Symmetry.mapto_bzMethod
mapto_bz(kpoint,recip_latvecs,inv_latvecs,coordinates;rtol,atol)

Map a k-point to a translationally equivalent point within the Brillouin zone.

Arguments

  • kpoint::AbstractVector{<:Real}: a single k-point in lattice or Cartesian coordinates.
  • recip_latvecs::AbstractMatrix{<:Real}: the reciprocal lattice vectors as columns of a matrix.
  • inv_latvecs::AbstractMatrix{<:Real}: the inverse matrix of the reciprocal lattice vectors.
  • coordinates::String: the coordinates of the given point, either "lattice" or "Cartesian". The point returned will be in the same coordinates.
  • rtol::Real=sqrt(eps(float(maximum(recip_latvecs)))): a relative tolerance for floating point comparisons. Finite precision errors creep in when pt is transformed to lattice coordinates because the transformation requires calculating a matrix inverse. The components of the k-point in lattice coordinates are checked to ensure that values close to 1 are equal to 1.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • bz_point::AbstractVector{<:Real}: the symmetrically equivalent k-point within the Brillouin zone in either lattice or Cartesian coordinates, depending on the coordinates specified.

Examples

import SymmetryReduceBZ.Symmetry: mapto_bz
import LinearAlgebra: inv
recip_latvecs = [1 0 0; 0 1 0; 0 0 1]
inv_latvecs = inv(recip_latvecs)
kpoint = [2, 3, 2]
coordinates = "Cartesian"
mapto_bz(kpoint, recip_latvecs, inv_latvecs, coordinates)
# output
3-element Vector{Float64}:
 0.0
 0.0
 0.0
source
SymmetryReduceBZ.Symmetry.mapto_ibzFunction
mapto_ibz(kpoints,recip_latvecs,inv_rlatvecs,ibz,pointgroup,coordinates;
    rtol,atol)

Map points as columns of a matrix to the IBZ and then remove duplicate points.

source
SymmetryReduceBZ.Symmetry.mapto_ibzMethod
mapto_ibz(kpoint,recip_latvecs,inv_rlatvecs,ibz,pointgroup,coordinates;rtol,
    atol)

Map a point to a symmetrically equivalent point within the IBZ.

Arguments

  • kpoint::AbstractVector{<:Real}: a k-point in 2D or 3D in Cartesian coordinates.
  • recip_latvecs::AbstractMatrix{<:Real}: the reciprocal lattice vectors as columns of a a matrix.
  • inv_rlatvecs::AbstractMatrix{<:Real}: the inverse of the square matrix recip_latvecs.
  • ibz::Chull: the irreducible Brillouin zone as as a convex hull objects from QHull.
  • pointgroup::AbstractVector{<:AbstractMatrix{<:Real}}: a list of point symmetry operators in matrix form that operate on points from the left.
  • coordinates::String: the coordinates the k-point is in. Options are "lattice" and "Cartesian". The k-point within the IBZ is returned in the same coordinates.
  • rtol::Real=sqrt(eps(float(maximum(recip_latvecs)))): a relative tolerance for floating point comparisons. The k-point is first mapped the unit cell and rtol is used when comparing components of the k-point to 1. It is also used for comparing floats to zero when checking if the point lies within ibz.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons. This is used everywhere rtol is used.

Returns

  • rot_point::AbstractVector{<:Real}: a symmetrically equivalent k-point to kpoint within the irreducible Brillouin zone in the same coordinates as coordinates.

Examples

import SymmetryReduceBZ.Lattices: get_recip_latvecs
import SymmetryReduceBZ.Symmetry: calc_spacegroup, mapto_ibz
import LinearAlgebra: inv
import QHull: chull
real_latvecs = [1 0; 0 2]
atom_types=[0]
atom_pos=Array([0 0]')
coordinates="Cartesian"
convention="ordinary"
recip_latvecs = get_recip_latvecs(real_latvecs,convention)
inv_rlatvecs = inv(recip_latvecs)
(ftrans,pg) = calc_spacegroup(real_latvecs,atom_types,atom_pos,coordinates)
ibz = chull([0.0 0.25; 0.0 0.0; 0.5 0.0; 0.5 0.25])
kpoint = [2,3]
ibz_point = mapto_ibz(kpoint,recip_latvecs,inv_rlatvecs,ibz,pg,coordinates)
# output
2-element Vector{Float64}:
 0.0
 0.0
source
SymmetryReduceBZ.Symmetry.mapto_unitcellMethod
mapto_unitcell(pt,latvecs,inv_latvecs,coordinates;rtol,atol)

Map a point to the first unit cell.

Arguments

  • pt::AbstractVector{<:Real}: a point in lattice or Cartesian coordinates.
  • latvecs::AbstractMatrix{<:Real}: the basis vectors of the lattice as columns of a matrix.
  • inv_latvecs::AbstractMatrix{<:Real}: the inverse of the matrix of that contains the lattice vectors.
  • coordinates::String: indicates whether pt is in "Cartesian" or "lattice" coordinates.
  • rtol::Real=sqrt(eps(float(maximum(inv_latvecs)))): a relative tolerance for floating point comparisons. Finite precision errors creep up when pt is transformed to lattice coordinates because the transformation requires calculating a matrix inverse. The components of the point in lattice coordinates are checked to ensure that values close to 1 are equal to 1.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • AbstractVector{<:Real}: a translationally equivalent point to pt in the first unit cell in the same coordinates.

Examples

using SymmetryReduceBZ
real_latvecs = [0 1 2; 0 -1 1; 1 0 0]
inv_latvecs=inv(real_latvecs)
pt=[1,2,3.2]
coordinates = "Cartesian"
SymmetryReduceBZ.Symmetry.mapto_unitcell(pt,real_latvecs,inv_latvecs,
    coordinates)
# output
3-element Vector{Float64}:
 0.0
 0.0
 0.20000000000000018
source

Utilities

SymmetryReduceBZ.Utilities.affine_transMethod
affine_trans(pts)

Calculate the affine transformation that maps the points to the xy-plane.

Arguments

  • pts::AbstractMatrix{<:Real}: Cartesian points as the columns of a matrix. The points must all lie on a plane in 3D.

Returns

  • M::AbstractMatrix{<:Real}: the affine transformation matrix that operates on points in homogeneous coordinates from the left.

Examples

using SymmetryReduceBZ
pts = [0.5 0.5 0.5; 0.5 -0.5 0.5; -0.5 0.5 0.5; -0.5 -0.5 0.5]'
SymmetryReduceBZ.Utilities.affine_trans(pts)
# output
4×4 Matrix{Float64}:
  0.0  -1.0   0.0  0.5
 -1.0   0.0   0.0  0.5
  0.0   0.0  -1.0  0.5
  0.0   0.0   0.0  1.0
source
SymmetryReduceBZ.Utilities.containsMethod
contains(array,arrays;rtol,atol)

Check if an array of arrays contains an array.

Arguments

  • array::AbstractArray: an array of reals of arbitrary dimension.
  • arrays::AbstractArray: a nested array of arrays of arbitrary dimension.
  • rtol::Real=sqrt(eps(float(maximum(pts)))): a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • Bool: a boolean that indicates the presence of absence of array in arrays.

Examples

using SymmetryReduceBZ.Utilities: contains
arrays = [[1 2; 2 3], [2 3; 4 5]]
array = [1 2; 2 3]
contains(array, arrays)
# output
true
source
SymmetryReduceBZ.Utilities.containsMethod
contains(pt,pts;rtol,atol)

Check if a point is contained in a matrix of points as columns.

Arguments

  • pt::AbstractVector{<:Real}: a point whose coordinates are the components of a vector.
  • pts::AbstractMatrix{<:Real}: coordinates of points as columns of a matrix.
  • rtol::Real=sqrt(eps(float(maximum(pts)))): a relative tolerance for floating point comparisons
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • Bool: a boolean that indicates the presence or absence of pt in pts.

Examples

using SymmetryReduceBZ.Utilities: contains
pts = Array([1 2; 2 3; 3 4; 4 5]')
pt = [1,2]
contains(pt,pts)
# output
true
source
SymmetryReduceBZ.Utilities.edgelengthsMethod
edgelengths(basis,radius;rtol,atol)

Calculate the edge lengths of a parallelepiped circumscribed by a sphere.

Arguments

  • basis::AbstractMatrix{<:Real}: a 2x2 or 3x3 matrix whose columns give the parallelogram or parallelepiped directions, respectively.
  • radius::Real: the radius of the sphere.
  • rtol::Real=sqrt(eps(float(radius))): a relative tolerace for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • [la,lb,lc]::AbstractVector{<:Real}: a list of parallelepiped lengths.

Examples

using SymmetryReduceBZ
basis=Array([1. 0. 0.; 0. 1. 0.; 0. 0. 1.])
radius=3.0
SymmetryReduceBZ.Utilities.edgelengths(basis,radius)
# output
3-element Vector{Float64}:
 3.0
 3.0
 3.0
source
SymmetryReduceBZ.Utilities.get_uniquefacetsMethod
get_uniquefacets(ch)

Calculate the unique facets of a convex hull.

QHull.jl package extension

This function is compatible with QHull.jl convex hulls through a package extension. After installing Python, SciPy, and QHull.jl, do using QHull to load it.

Arguments

  • ch::Polyhedron: a polyhedron from Polyhedra.jl or convex hull in 3D from QHull.jl.

Returns

  • unique_facets::Vector{<:Vector}: a nested list of the vertices of points that lie on each face.

Examples

using SymmetryReduceBZ.Utilities: get_uniquefacets
using SymmetryReduceBZ.Symmetry: calc_bz
real_latvecs = [1 0 0; 0 1 0; 0 0 1]
atom_types = [0]
atom_pos = Array([0 0 0]')
coordinates = "Cartesian"
makeprim = false
convention = "ordinary"
bz = calc_bz(real_latvecs,atom_types,atom_pos,coordinates,makeprim,convention)
get_uniquefacets(bz)
# output
6-element Vector{Vector{Vector{Float64}}}:
 [[-0.5, 0.5, -0.5], [-0.5, 0.5, 0.5], [-0.5, -0.5, 0.5], [-0.5, -0.5, -0.5]]
 [[0.5, -0.5, 0.5], [0.5, -0.5, -0.5], [-0.5, -0.5, -0.5], [-0.5, -0.5, 0.5]]
 [[0.5, -0.5, -0.5], [0.5, 0.5, -0.5], [-0.5, 0.5, -0.5], [-0.5, -0.5, -0.5]]
 [[0.5, -0.5, 0.5], [0.5, 0.5, 0.5], [-0.5, 0.5, 0.5], [-0.5, -0.5, 0.5]]
 [[0.5, 0.5, -0.5], [0.5, 0.5, 0.5], [-0.5, 0.5, 0.5], [-0.5, 0.5, -0.5]]
 [[0.5, -0.5, 0.5], [0.5, -0.5, -0.5], [0.5, 0.5, -0.5], [0.5, 0.5, 0.5]]
source
SymmetryReduceBZ.Utilities.get_uniquefacetsindicesMethod
get_uniquefacetsindices(ch)

Calculate the indices of unique facets of a convex hull.

QHull.jl package extension

This function is compatible with QHull.jl convex hulls through a package extension. After installing Python, SciPy, and QHull.jl, do using QHull to load it.

Arguments

  • ch::Polyhedron: a polyhedron from Polyhedra.jl or convex hull in 3D from QHull.jl.

Returns

  • unique_facets::Vector{<:Vector}: a nested list of the indices of points that lie on each face. See get_uniquefacets if you want the facets' points.

Examples

using QHull: chull
using SymmetryReduceBZ.Utilities: get_uniquefacetsindices, vertices
using SymmetryReduceBZ.Symmetry: calc_bz
real_latvecs = [1 0 0; 0 1 0; 0 0 1]
atom_types = [0]
atom_pos = Array([0 0 0]')
coordinates = "Cartesian"
makeprim = false
convention = "ordinary"
bz = calc_bz(real_latvecs,atom_types,atom_pos,coordinates,makeprim,convention)
ch = chull(permutedims(reduce(hcat, vertices(bz))))
get_uniquefacetsindices(ch)
# output
6-element Vector{Vector{Int32}}:
 [1, 2, 3, 4]
 [7, 2, 3, 5]
 [6, 4, 3, 5]
 [7, 2, 1, 8]
 [6, 4, 1, 8]
 [8, 7, 5, 6]
source
SymmetryReduceBZ.Utilities.mapto_xyplaneMethod
function mapto_xyplane(pts)

Map Cartesian points embedded in 3D on a plane to the xy-plane embedded in 2D.

Arguments

  • pts::AbstractMatrix{<:Real}: Cartesian points embedded in 3D as columns of a matrix.

Returns

  • AbstractMatrix{<:Real}: Cartesian points in 2D as columns of a matrix.

Examples

using SymmetryReduceBZ
pts = [0.5 -0.5 0.5; 0.5 -0.5 -0.5; 0.5 0.5 -0.5; 0.5 0.5 0.5]'
SymmetryReduceBZ.Utilities.mapto_xyplane(pts)
# output
2×4 Matrix{Float64}:
 0.0  1.0  1.0  0.0
 0.0  0.0  1.0  1.0
source
SymmetryReduceBZ.Utilities.points_in_ballMethod
points_in_ball(points,radius,offset,rtol=sqrt(eps(float(radius))),atol=1e-9)

Calculate the points within a ball (circle, sphere, ...).

Arguments

  • points::AbstractMatrix{<:Real}: points in Cartesian coordinates as columns of a matrix.
  • radius::Real: the radius of the ball.
  • offset::AbstractVector{<:Real}: the location of the center of the ball in Cartesian coordinates.
  • rtol::Real=sqrt(eps(float(radius))): a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance.

Returns

  • ball_points::AbstractVector{<:Int}: the indices of points in points within the ball.

Examples

using SymmetryReduceBZ.Utilities: points_in_ball
points = [0.0 0.05 0.1 0.15 0.2 0.25 0.3 0.35 0.4 0.45 0.5 0.05 0.1 0.15 0.2 0.25 0.3 0.35 0.4 0.45 0.5 0.1 0.15 0.2 0.25 0.3 0.35 0.4 0.45 0.5 0.15 0.2 0.25 0.3 0.35 0.4 0.45 0.5 0.2 0.25 0.3 0.35 0.4 0.45 0.5 0.25 0.3 0.35 0.4 0.45 0.5 0.3 0.35 0.4 0.45 0.5 0.35 0.4 0.45 0.5 0.4 0.45 0.5 0.45 0.5 0.5; 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.05 0.05 0.05 0.05 0.05 0.05 0.05 0.05 0.05 0.05 0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.1 0.15 0.15 0.15 0.15 0.15 0.15 0.15 0.15 0.2 0.2 0.2 0.2 0.2 0.2 0.2 0.25 0.25 0.25 0.25 0.25 0.25 0.3 0.3 0.3 0.3 0.3 0.35 0.35 0.35 0.35 0.4 0.4 0.4 0.45 0.45 0.5]
radius = 0.1
offset = [0,0]
points_in_ball(points,radius,offset)
# output
4-element Vector{Int64}:
  1
  2
  3
 12
source
SymmetryReduceBZ.Utilities.remove_duplicatesMethod
remove_duplicates(points;rtol,atol)

Remove duplicates from an array.

Arguments

  • points::AbstractVector: items in a vector, which can be floats or arrays.
  • rtol::Real=sqrt(eps(float(maximum(points)))): relative tolerance.
  • atol::Real=1e-9: absolute tolerance.

Returns

  • uniquepts::AbstractVector: an vector with only unique elements.

Examples

using SymmetryReduceBZ.Utilities: remove_duplicates
test = [1.,1.,2,2,]
remove_duplicates(test)
# output
2-element Vector{Float64}:
 1.0
 2.0
source
SymmetryReduceBZ.Utilities.sample_circleFunction
sample_circle(basis,radius,offset;rtol,atol)

Sample uniformly within a circle centered about a point.

Arguments

  • basis::AbstractMatrix{<:Real}: a 2x2 matrix whose columns are the grid generating vectors.
  • radius::Real: the radius of the circle.
  • offset::AbstractVector{<:Real}=[0.,0.]: the xy-coordinates of the center of the circle.
  • rtol::Real=sqrt(eps(float(radius))): a relative tolerace for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • pts::AbstractMatrix{<:Real} a matrix whose columns are sample points in Cartesian coordinates.

Examples

using SymmetryReduceBZ
basis=Array([1. 0.; 0. 1.]')
radius=1.0
offset=[0.,0.]
SymmetryReduceBZ.Utilities.sample_circle(basis,radius,offset)
# output
2×5 Matrix{Float64}:
  0.0  -1.0  0.0  1.0  0.0
 -1.0   0.0  0.0  0.0  1.0
source
SymmetryReduceBZ.Utilities.sample_sphereFunction
sample_sphere(basis,radius,offset;rtol,atol)

Sample uniformly within a sphere centered about a point.

Arguments

  • basis::AbstractMatrix{<:Real}: a 3x3 matrix whose columns are the grid generating vectors.
  • radius::Real: the radius of the sphere.
  • offset::AbstractVector{<:Real}=[0.,0.]: the xy-coordinates of the center of the circle.
  • rtol::Real=sqrt(eps(float(radius))): a relative tolerace for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • pts::AbstractMatrix{<:Real} a matrix whose columns are sample points in Cartesian coordinates.

Examples

using SymmetryReduceBZ
basis=Array([1. 0. 0.; 0. 1. 0.; 0. 0. 1.])
radius=1.0
offset=[0.,0.,0.]
SymmetryReduceBZ.Utilities.sample_sphere(basis,radius,offset)
# output
3×7 Matrix{Float64}:
  0.0   0.0  -1.0  0.0  1.0  0.0  0.0
  0.0  -1.0   0.0  0.0  0.0  1.0  0.0
 -1.0   0.0   0.0  0.0  0.0  0.0  1.0
source
SymmetryReduceBZ.Utilities.shoelaceMethod
shoelace(vertices)

Calculate the area of a polygon with the shoelace algorithm.

Arguments

  • vertices::AbstractMatrix{<:Real}: the xy-coordinates of the vertices of the polygon as the columns of a matrix.

Returns

  • area::Real: the area of the polygon.

Examples

using SymmetryReduceBZ.Utilities: shoelace
pts = [0 0 1; -1 1 0]
shoelace(pts)
# output
1.0
source
SymmetryReduceBZ.Utilities.sortpts2DMethod
function sortpts2D(pts)

Calculate the permutation vector that sorts 2D Cartesian points counterclockwise with respect to the average of the points.

Arguments

  • pts::AbstractMatrix{<:Real}: Cartesian points in 2D.

Returns

  • perm::AbstractVector{<:Real}: the permutation vector that orders the points clockwise or counterclockwise.

```

source
SymmetryReduceBZ.Utilities.sortpts_permMethod
function sortpts_perm(pts)

Calculate the permutation vector that sorts Cartesian points embedded in 3D that lie on a plane (counter)clockwise with respect to the average of all points.

Arguments

  • pts::AbstractMatrix{<:Real}: Cartesian points embedded in 3D that all lie on a plane. The points are columns of a matrix.

Returns

  • ::AbstractVector{<:Real}: the permutation vector that orders the points clockwise or counterclockwise.

Examples

using SymmetryReduceBZ.Utilities: sortpts_perm
pts = [0.5 -0.5 0.5; 0.5 -0.5 -0.5; 0.5 0.5 -0.5; 0.5 0.5 0.5]'
perm=sortpts_perm(pts)
pts[:,perm]
# output
3×4 Matrix{Float64}:
  0.5   0.5   0.5  0.5
 -0.5  -0.5   0.5  0.5
  0.5  -0.5  -0.5  0.5
source
SymmetryReduceBZ.Utilities.unique_pointsMethod
unique_points(points;rtol,atol)

Remove duplicate points.

Arguments

  • points::AbstractMatrix{<:Real}: the points are columns of a matrix.
  • rtol::Real=sqrt(eps(float(maximum(flatten(points))))): a relative tolerance for floating point comparisons.
  • atol::Real=1e-9: an absolute tolerance for floating point comparisons.

Returns

  • uniquepts::AbstractMatrix{<:Real}: the unique points as columns of a matrix.

Examples

using SymmetryReduceBZ
points=Array([1 2; 2 3; 3 4; 1 2]')
SymmetryReduceBZ.Utilities.unique_points(points)
# output
2×3 Matrix{Int64}:
 1  2  3
 2  3  4
source
SymmetryReduceBZ.Utilities.verticesMethod
vertices(hull)

Arguments

  • hull::Polyhedron: a convex hull of a polytope

Returns

  • pts: a iterator of points (vectors) at the vertices of the convex hull

Examples

using SymmetryReduceBZ.Utilities: vertices
using SymmetryReduceBZ.Symmetry: calc_bz
real_latvecs = [1 0 0; 0 1 0; 0 0 1]
atom_types = [0]
atom_pos = Array([0 0 0]')
coordinates = "Cartesian"
makeprim = false
convention = "ordinary"
bz = calc_bz(real_latvecs,atom_types,atom_pos,coordinates,makeprim,convention)
vertices(bz)
# output
8-element iterator of Vector{Float64}:
 [0.5, -0.5, 0.5]
 [0.5, -0.5, -0.5]
 [0.5, 0.5, -0.5]
 [0.5, 0.5, 0.5]
 [-0.5, 0.5, -0.5]
 [-0.5, 0.5, 0.5]
 [-0.5, -0.5, -0.5]
 [-0.5, -0.5, 0.5]
source
SymmetryReduceBZ.Utilities.volumeMethod
volume(hull)

Arguments

  • hull::Polyhedron

Returns

  • vol::Real

Examples

using SymmetryReduceBZ.Utilities: volume
using SymmetryReduceBZ.Symmetry: calc_bz
real_latvecs = [1 0 0; 0 1 0; 0 0 1]
atom_types = [0]
atom_pos = Array([0 0 0]')
coordinates = "Cartesian"
makeprim = false
convention = "ordinary"
bz = calc_bz(real_latvecs,atom_types,atom_pos,coordinates,makeprim,convention)
volume(bz)
# output
1.0
source