Public
TopologicalNumbers.BP — TypeFukui-Hatsugai-Suzuki method [1]
Definition
The Berry phase of the $n$th band $\nu_{n}$ is defined by
\[\nu_{n}=\frac{1}{\pi}\sum_{k\in\mathrm{BZ}}U_{n}(k)\]
The range $\mathrm{BZ}$(Brillouin Zone) is $k\in[0,2\pi]$. $U_{n,i}(k)$ is the link variable at wavenumber $k$. $e_{1}$ is the unit vector.
\[U_{n}(k)=\braket{\Psi_{n}(k)|\Psi_{n}(k+e_{1})}\]
\[\ket{\Psi_{n}(k)}\]
is the wave function of the $n$th band.
TopologicalNumbers.BPProblem — Typestruct BPProblem{T1<:Function,T2<:Union{Tuple,AbstractVector,Int},T3<:Real,T4<:Bool} <: TopologicalNumbersProblemsThe BPProblem struct represents a problem for calculating Berry phase.
Fields
H::T1: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 1.N::T2: The number of points for one direction in the Brillouin zone. Default is 51.gapless::T3: The threshold for considering a band as gapless. Default is 0.0.rounds::T4: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.BPProblem — MethodBPProblem(H, N)Constructs a Berry phase problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 1.N: The number of points for one direction in the Brillouin zone.
Returns
A BPProblem object.
Example
julia> TopologicalNumbers.BPProblem — MethodBPProblem(H)Constructs a Berry phase problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 1.
Returns
A BPProblem object.
Example
julia> TopologicalNumbers.BPSolution — TypeBPSolution{T1,T2} <: TopologicalNumbersSolutionsThe BPSolution struct represents a solution for calculating Berry phase.
Fields
TopologicalNumber::T1: The Berry phase for each energy bands.Total::T2: The total Berry phase.
TopologicalNumbers.Evar — TypeTopologicalNumbers.FCProblem — TypeFCProblem{T1<:Function,T2<:Union{Tuple,AbstractVector,Int},T3<:Real,T4<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the first Chern number.
Fields
H::T1: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.N::T2: The number of points for one direction in the Brillouin zone. Default is 51.gapless::T3: The threshold for considering a band as gapless. Default is 0.0.rounds::T4: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.FCProblem — MethodFCProblem(H, N)Constructs a first Chern number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.N: The number of points for one direction in the Brillouin zone.
Returns
A FCProblem object.
Example
julia> TopologicalNumbers.FCProblem — MethodFCProblem(H)Constructs a first Chern number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.
Returns
A FCProblem object.
Example
julia> TopologicalNumbers.FCSolution — TypeFCSolution{T1,T2} <: TopologicalNumbersSolutionsThe FCSolution struct represents a solution for calculating the first Chern number.
Fields
TopologicalNumber::T1: The first Chern number for each energy bands.Total::T2: The total of the first Chern number.
TopologicalNumbers.FHS — TypeFukui-Hatsugai-Suzuki method [1]
Definition
The first Chern number of the $n$th band $\nu_{n}$ is defined by
\[\nu_{n}=\frac{1}{2\pi}\sum_{\bm{k}\in\mathrm{BZ}}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,2}^{*}(\bm{k})\right]\right]\]
The range $\mathrm{BZ}$(Brillouin Zone) is $\bm{k}\in[0,2\pi]^{2}$. $U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $n$th band.
TopologicalNumbers.FHS2 — TypeFukui-Hatsugai-Suzuki method [2]
TopologicalNumbers.FHSlocal2 — TypeFukui-Hatsugai-Suzuki method [1]
Definition
The Berry flux at the wavenumber $\bm{k}$ of the $n$th band $F_{n}(\bm{k})$ is defined by
\[F_{n}(\bm{k})=f_{n}(\bm{k})-df_{n}(\bm{k})\]
\[f_{n}(\bm{k})=\frac{1}{2\pi}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,1}^{*}(\bm{k})\right]\right]\]
\[df_{n}(\bm{k})=\frac{1}{2\pi}\mathrm{Im}\left[\mathrm{Log}[U_{n,1}(\bm{k})]+\mathrm{Log}[U_{n,2}(\bm{k}+\bm{e}_{1})]-\mathrm{Log}[U_{n,1}(\bm{k}+\bm{e}_{2})]-\mathrm{Log}[U_{n,1}(\bm{k})]\right]\]
$U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $n$th band.
TopologicalNumbers.FHSlocal3 — TypeTopologicalNumbers.FHSsurface — TypeTopologicalNumbers.LBFProblem — TypeLBFProblem{T1<:Function,T2<:AbstractVector,T3<:Union{Tuple,AbstractVector,Int},T4<:Real,T5<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the $k$-local value of Berry flux.
Fields
H::T1: The Hamiltonian functionH=H(k)that defines the system.kis anAbstractVector(or aTuple) of the wavenumber vector. Dimension ofkmust be 2.n::T2: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$) when calculating Berry flux. Dimension ofnmust be 2.N::T3: The number of points for one direction in the Brillouin zone. Default is 51.gapless::T4: The threshold for considering a band as gapless. Default is 0.0.rounds::T5: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.LBFProblem — MethodLBFProblem(H, n, N)Constructs a local Berry flux problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 2.n: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$) when calculating Berry flux. Dimension ofnmust be 2.N: The number of points for one direction in the Brillouin zone.
Returns
A LBFProblem object.
Example
julia> TopologicalNumbers.LBFProblem — MethodLBFProblem(H, n)Constructs a local Berry flux problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 2.n: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$) when calculating Berry flux. Dimension ofnmust be 2.
Returns
A LBFProblem object.
Example
julia> TopologicalNumbers.LBFSolution — TypeLBFSolution{T1,T2} <: TopologicalNumbersSolutionsThe LBFSolution struct represents a solution for calculating the $k$-local value of Berry flux.
Fields
TopologicalNumber::T1: The local Berry flux for each energy bands.n::T2: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$) when calculating Berry flux.
TopologicalNumbers.SCProblem — TypeSCProblem{T1<:Function,T2<:Union{Int,Nothing},T3<:Union{Tuple,AbstractVector,Int},T4<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the second Chern number.
Fields
H::T1: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 4.Nfill::T2: The number of filled bands.nothingmeans the half-filling. Default isnothing.N::T3: The number of points in the Brillouin zone. Type ofIntmeans the uniform mesh. You can also specify the mesh by giving a tuple or a vector. Default is 30.RV::T4: A boolean indicating whether to return arealvalue. Default is true.
Example
julia> TopologicalNumbers.SCProblem — MethodSCProblem(H, N, Nf)Constructs a second Chern number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 4.N: The number of points in the Brillouin zone. Type ofIntmeans the uniform mesh. You can also specify the mesh by giving a tuple or a vector.Nf: The number of filled bands.nothingmeans the half-filling.
Returns
A SCProblem object.
Example
julia> TopologicalNumbers.SCProblem — MethodSCProblem(H, N)Constructs a second Chern number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 4.N: The number of points in the Brillouin zone. Type ofIntmeans the uniform mesh. You can also specify the mesh by giving a tuple or a vector.
Returns
A SCProblem object.
Example
julia> TopologicalNumbers.SCProblem — MethodSCProblem(H)Constructs a second Chern number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 4.
Returns
A SCProblem object.
Example
julia> TopologicalNumbers.SCSolution — TypeSCSolution{T} <: TopologicalNumbersSolutionsThe SCSolution struct represents a solution for calculating the second Chern number.
Fields
TopologicalNumber::T: The second Chern number.
TopologicalNumbers.Shio — TypeShiozaki method [4]
Definition
The $\mathbb{Z}_{2}$ number of the $2n$th (and $2n-1$th) band $\nu_{n}$ is defined by
\[\nu_{n}=F_{n}-\left(P_{n}(0)-P_{n}(\pi)\right)\]
$F_{n}$ is the Berry flux of the $n$th band in the $\mathrm{BZ}'$. The range $\mathrm{BZ}'$ is $\bm{k}\in[0,2\pi]\times[0,\pi]$ half of BZ(Brillouin Zone).
\[F_{n}=\frac{1}{2\pi}\sum_{\bm{k}\in\mathrm{BZ}'}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,1}^{*}(\bm{k})\right]\right]\]
$P_{n}(k_{2})$ is the time-reversal polarization at wavenumber $k_{2}$.
\[P_{n}(k_{2})=\frac{1}{2\pi}\frac{\mathrm{Pf}[\omega(0,k_{2})]}{\mathrm{Pf}[\omega(\pi,k_{2})]}\sum_{k_{1}=0}^{\pi-e_{1}}U_{n,1}(\bm{k})\]
$U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $2n$th (and $2n-1$th) band. $\omega(\bm{k})$ is the unitary matrix given by
\[\omega(\bm{k})=\bra{\Psi(-\bm{k})}T\ket{\Psi(\bm{k})}\]
$T$ is the time-reversal operator.
TopologicalNumbers.UseMPI — TypeUseMPI{T}(MPI::T) <: TopologicalNumbersMultiProcessA struct representing the use of MPI for parallel computing in the TopologicalNumbers package.
Arguments
MPI: An object representing the MPI library.
TopologicalNumbers.UseSingleThread — TypeUseSingleThread <: TopologicalNumbersSingleProcessA struct representing the use of a single thread for parallel processing in the TopologicalNumbers module.
TopologicalNumbers.WCSProblem — TypeWCSProblem{T1<:Function,T2<:String,T3<:Int,T4<:Union{Tuple,AbstractVector,Int},T5<:Real,T6<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for finding and calculating the Weyl points.
Fields
H::T1: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.kn::T2: Compute the Chern number of the plane perpendicular to the"kn"direction in Brillouin zone ("k1","k2","k3").kn_mesh::T3: Number of mesh in"kn"direction. Default is 51.N::T4: The number of points for one direction in the Brillouin zone. Default is 51.gapless::T5: The threshold for considering a band as gapless. Default is 0.0.rounds::T6: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.WCSProblem — MethodWCSProblem(H, kn, N1, N2)Constructs a problem for finding and calculating the Weyl points with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.kn: Compute the Chern number of the plane perpendicular to the"kn"direction in Brillouin zone ("k1","k2","k3").N1: Number of mesh in"kn"direction.N2: The number of points for one direction in the Brillouin zone.
Returns
A WCSProblem object.
Example
julia> TopologicalNumbers.WCSProblem — MethodWCSProblem(H, kn)Constructs a problem for finding and calculating the Weyl points with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.kn: Compute the Chern number of the plane perpendicular to the"kn"direction in Brillouin zone ("k1","k2","k3").
Returns
A WCSProblem object.
Example
julia> TopologicalNumbers.WCSProblem — MethodWCSProblem(H, kn, N::T) where {T<:Int}Constructs a problem for finding and calculating the Weyl points with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.kn: Compute the Chern number of the plane perpendicular to the"kn"direction in Brillouin zone ("k1","k2","k3").N: The number of points for one direction in the Brillouin zone.
Returns
A WCSProblem object.
Example
julia> TopologicalNumbers.WCSSolution — TypeWCSSolution{T1,T2,T3} <: TopologicalNumbersSolutionsThe WCSSolution struct represents a solution for finding and calculating the Weyl points.
Fields
kn::T1: Compute the Chern number of the plane perpendicular to the"kn"direction in Brillouin zone ("k1","k2","k3").param::T2: Wavenumber parameters in"kn"direction.nums::T3: The Chern number for each energy bands and each wavenumber parameters.
TopologicalNumbers.WNProblem — TypeWNProblem{T1<:Function,T2<:AbstractVector,T3<:Int,T4<:Real,T5<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the Weyl nodes.
Fields
H::T1: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.n::T2: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$). Dimension ofnmust be 3.N::T3: The number of points for one direction in the Brillouin zone. Default is 51.gapless::T4: The threshold for considering a band as gapless. Default is 0.0.rounds::T5: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.WNProblem — MethodWNProblem(H, n, N)Constructs a problem for calculating the Weyl nodes with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.n: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$). Dimension ofnmust be 3.N: The number of points for one direction in the Brillouin zone.
Returns
A WNProblem object.
Example
julia> TopologicalNumbers.WNProblem — MethodWNProblem(H, n)Constructs a problem for calculating the Weyl nodes with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.n: AnAbstractVector(or aTuple) including two elements ofInt, which represents wavenumber ($2\pi n/N$). Dimension ofnmust be 3.
Returns
A WNProblem object.
Example
julia> TopologicalNumbers.WNSolution — TypeWNSolution{T1,T2,T3} <: TopologicalNumbersSolutionsThe WNSolution struct represents a solution for calculating the Weyl nodes.
Fields
TopologicalNumber::T1:n::T2:N::T3:
TopologicalNumbers.WPProblem — TypeWPProblem{T1<:Function,T2<:Int,T3<:AbstractVector,T4<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the Weyl points.
Fields
H::T1: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.N::T2: The number of meshes when discretizing the Brillouin Zone. The $n$th iteration divides the Brillouin zone into $1/N^n$. Default is 10.gapless::T3The threshold that determines the state to be degenerate. The $n$th iteration adopts the threshold value of the $n$th value of the vector. The number of iterations can be varied by the length of the vector. Default is[1e-1, 1e-2, 1e-3, 1e-4].rounds::T4: A boolean indicating whether to round a returned variable. Default is true.
Example
julia> TopologicalNumbers.WPProblem — MethodWPProblem(H, N)Constructs a problem for calculating the Weyl points with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.N: The number of meshes when discretizing the Brillouin Zone. The $n$th iteration divides the Brillouin zone into $1/N^n$.
Returns
A WPProblem object.
Example
julia> TopologicalNumbers.WPProblem — MethodWPProblem(H)Constructs a problem for calculating the Weyl points with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector. Dimension ofkmust be 3.
Returns
A WPProblem object.
Example
julia> TopologicalNumbers.WPSolution — TypeWPSolution{T1,T2,T3} <: TopologicalNumbersSolutionsThe WPSolution struct represents a solution for calculating the Weyl points.
Fields
WeylPoint::T1:N::T2:Nodes::T3:
TopologicalNumbers.Z2Problem — TypeZ2Problem{T1<:Function,T2<:Union{Int,Nothing},T3<:Union{Tuple,AbstractVector,Int},T4<:Bool} <: TopologicalNumbersProblemsA struct representing a problem for calculating the Z2 number.
Fields
H::T1: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.Nfill::T2: The number of filled bands.nothingmeans the half-filling. Default isnothing.N::T3: The number of points for one direction in the Brillouin zone. Default is 50.rounds::T4: A boolean indicating whether to round a returned variable. Default is true.TR::T4: A boolean indicating whether to calculate the remaining part of the Brillouin zone. Iftrue,solvereturns an additional resultTRTopologicalNumber. If the calculation is done nomally,TRTopologicalNumberis equal toTopologicalNumber. Default is false.
Example
julia> TopologicalNumbers.Z2Problem — MethodZ2Problem(H, Nf, N)Constructs a Z2 number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.Nf: The number of filled bands.nothingmeans the half-filling.N: The number of points for one direction in the Brillouin zone.
Returns
A Z2Problem object.
Example
julia> TopologicalNumbers.Z2Problem — MethodZ2Problem(H, Nf)Constructs a Z2 number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.Nf: The number of filled bands.nothingmeans the half-filling.
Returns
A Z2Problem object.
Example
julia> TopologicalNumbers.Z2Problem — MethodZ2Problem(H)
Constructs a Z2 number problem with the default parameters.
Arguments
H: The Hamiltonian functionH=H(k, p)that defines the system.kis a abstract vector (or a tuple) of the wavenumber vector andpcontains parameter. Dimension ofkmust be 2.
Returns
A Z2Problem object.
Example
julia> TopologicalNumbers.Z2Solution — TypeZ2Solution{T1,T2,T3} <: TopologicalNumbersSolutionsThe Z2Solution struct represents a solution for calculating Z2 number.
Fields
TopologicalNumber::T1: The Z2 number for each pair of energy bands.TRTopologicalNumber::T2: The Z2 number for the remaining part of the Brillouin zone. This field is only returned whenTRistrueinZ2Problem.Total::T3: The total Z2 number.
TopologicalNumbers.BHZ — MethodHamiltonian of the Bernevig–Hughes–Zhang (BHZ) model.
BHZ(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Union{AbstractVector,Tuple}}Arguments
k::T1: two-dimensional wavenumberk.p::T2: parameters defined as below.
Definition
Hamiltonian of the BHZ model is defined as
\[H(k)=\begin{pmatrix} h_{0}(\bm{k})+h_{3}(\bm{k}) & 0 & h_{5}(\bm{k})-ih_{4}(\bm{k}) & 0 \\ 0 & h_{0}(\bm{k})-h_{3}(\bm{k}) & 0 & h_{5}(\bm{k})-ih_{4}(\bm{k}) \\ h_{5}(\bm{k})+ih_{4}(\bm{k}) & 0 & h_{0}(\bm{k})-h_{3}(\bm{k}) & 0 \\ 0 & h_{5}(\bm{k})+ih_{4}(\bm{k}) & 0 & h_{0}(\bm{k})+h_{3}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{0}(\bm{k})=&-(t_{ss}-t_{pp})(\cos(k_{1})+\cos(k_{2}))+\frac{(\epsilon_{s}+\epsilon_{p})}{2} \\ h_{3}(\bm{k})=&2t_{sp}\sin(k_{2}) \\ h_{4}(\bm{k})=&2t_{sp}\sin(k_{1}) \\ h_{5}(\bm{k})=&-(t_{ss}+t_{pp})(\cos(k_{1})+\cos(k_{2}))+\frac{(\epsilon_{s}-\epsilon_{p})}{2} \end{align*}\]
where $t_{sp}$ is the amplitude of the mixed hopping between s- and p-orbitals in the tight binding model. $t_{ss}$ is the hopping between s-orbitals, $t_{pp}$ is the hopping between p-orbitals, and $t_{1}=t_{s}-t_{p}$, $t_{2}=t_{s}+t_{p}$. $\epsilon_{s}$ and $\epsilon_{p}$ are the site potential terms for the s- and p-orbitals, $\epsilon_{1}\epsilon_{s}s+\epsilon_{p}$, $\epsilon_{2}=\epsilon_{s}-\epsilon_{p}. Nondimensionalize with $t_{sp}=t_{1}=\epsilon_{1}=1$ and $\epsilon_{2},t_{2}=$p.
Examples
julia> BHZ([0, π/3], [0.5, 0.5])
4×4 Matrix{ComplexF64}:
0.732051+0.0im 0.0+0.0im -0.5+0.0im 0.0+0.0im
0.0+0.0im -2.73205+0.0im 0.0+0.0im -0.5+0.0im
-0.5+0.0im 0.0+0.0im -2.73205+0.0im 0.0+0.0im
0.0+0.0im -0.5+0.0im 0.0+0.0im 0.732051+0.0imTopologicalNumbers.Flux2d — MethodHamiltonian of the two-dimensional square lattice with flux model.
Flux2d(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Union{AbstractVector,Tuple}}Arguments
k::T1: one-dimensional wavenumberk.p::T2: parameters defined as below.
Definition
Hamiltonian of the two-dimensional square lattice with flux model is defined as
\[H(k)=\begin{pmatrix} -2t\cos(k_{2}-2\pi m\frac{1}{n}) & -t & 0 & 0 & \cdots & e^{-ik_{1}} \\ -t & -2t\cos(k_{2}-2\pi m\frac{2}{n}) & -t & 0 & \cdots & 0 \\ 0 & -t & \ddots & -t & \cdots & \vdots \\ \vdots & \vdots & & -t & -2t\cos(k_{2}-2\pi m\frac{(n-1)}{n}) & -t \\ e^{ik_{1}} & 0 & \cdots & 0 & -t & -2t\cos(k_{2}-2\pi m) \end{pmatrix}\]
where $t$ is the amplitude of the nearest neighbor hopping in the tight binding model. $n$ and $m$ are the size of the unit cell of the phase produced by the application of the static magnetic field. Nondimensionalize with $t=1$ and $n,m=$p.
Examples
julia> Flux2d([0, π/3], [4, 1])
4×4 Matrix{ComplexF64}:
-1.73205+0.0im -1.0+0.0im 0.0+0.0im -1.0+0.0im
-1.0+0.0im 1.0+0.0im -1.0+0.0im 0.0+0.0im
0.0+0.0im -1.0+0.0im 1.73205+0.0im -1.0+0.0im
-1.0-0.0im 0.0+0.0im -1.0+0.0im -1.0+0.0imTopologicalNumbers.Haldane — MethodHamiltonian of the Haldane model.
Haldane(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Union{AbstractVector,Tuple}}Arguments
k::T1: two-dimensional wavenumberk.p::T2: parameters defined as below.
Definition
Hamiltonian of the Haldane model is defined as
\[H(\bm{k})=\begin{pmatrix} h_{0}(\bm{k})+h_{3}(\bm{k}) & h_{1}(\bm{k})-ih_{2}(\bm{k}) \\ h_{1}(\bm{k})+ih_{2}(\bm{k}) & h_{0}(\bm{k})-h_{3}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{0}(\bm{k})=&2t_{2}\cos(\phi)(\sin(k_{1})+\sin(k_{2})+\sin(k_{1}+k_{2})) \\ h_{1}(\bm{k})=&-t_{1}(1+\cos(k_{1})+\cos(k_{2})) \\ h_{2}(\bm{k})=&-t_{1}(-\sin(k_{1})+\sin(k_{2})) \\ h_{3}(\bm{k})=&m+2t_{2}\sin(\phi)(\sin(k_{1})+\sin(k_{2})-\sin(k_{1}+k_{2})) \end{align*}\]
where $t_{1}$ and $t_{2}$ are the amplitudes of the nearest and the next nearest neighbor hopping in the tight binding model. $\phi$ is the phase produced by the application of the static magnetic field. Nondimensionalize with $t_{1}=1$ and $t_{2},\phi,m=$p.
Examples
julia> Haldane([0, π/3], [0.5, π/2, 0.3])
2×2 Matrix{ComplexF64}:
0.3-0.0im -2.5+0.866025im
-2.5-0.866025im -0.3-0.0imTopologicalNumbers.KaneMele — MethodHamiltonian of the Kane–Mele model.
KaneMele(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Real}Arguments
k::T1: two-dimensional wavenumberk.p::T2: parameter defined as below.
Definition
Hamiltonian of the Kane–Mele model is defined as
\[H(k)=\begin{pmatrix} h_{3}(\bm{k}) & 0 & h_{5}(\bm{k})-ih_{4}(\bm{k}) & 0 \\ 0 & -h_{3}(\bm{k}) & 0 & h_{5}(\bm{k})-ih_{4}(\bm{k}) \\ h_{5}(\bm{k})+ih_{4}(\bm{k}) & 0 & -h_{3}(\bm{k}) & 0 \\ 0 & h_{5}(\bm{k})+ih_{4}(\bm{k}) & 0 & h_{3}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{3}(\bm{k})=&2\lambda_{\mathrm{SO}}(\sin(k_{1})-\sin(k_{2})-\sin(k_{1}-k_{2})) \\ h_{4}(\bm{k})=&-t(\sin(k_{1})+\sin(k_{2})) \\ h_{5}(\bm{k})=&-t(\cos(k_{1})+\cos(k_{2})+1) \end{align*}\]
where $t$ is the amplitude of the nearest neighbor hopping in the tight binding model. $\lambda_{\mathrm{SO}}$ is the magnitude of spin-orbit interaction. Nondimensionalize with $t=1$ and $\lambda_{\mathrm{SO}}=$p.
Examples
julia> KaneMele([0, π/3], 0.5)
4×4 Matrix{ComplexF64}:
0.0+0.0im 0.0+0.0im -2.5+0.866025im 0.0+0.0im
0.0+0.0im 0.0+0.0im 0.0+0.0im -2.5+0.866025im
-2.5-0.866025im 0.0+0.0im 0.0+0.0im 0.0+0.0im
0.0+0.0im -2.5-0.866025im 0.0+0.0im 0.0+0.0imTopologicalNumbers.KitaevChain — MethodHamiltonian of the Kitaev chain model.
KitaevChain(k::T1, p::T2) where {T1<:Real, T2<:Union{AbstractVector,Tuple}}Arguments
k::T1: one-dimensional wavenumberk.p::T2: parameters defined as below.
Definition
Hamiltonian of the Kitaev chain model is defined as
\[H(k)=\begin{pmatrix} -\mu-2t\cos(k) & 2i\Delta\sin(k) \\ -2i\Delta\sin(k) & \mu+2t\cos(k) \end{pmatrix}\]
where $t$ is the amplitude of the nearest neighbor hopping in the tight binding model. $\Delta$ is the pairing with the nearest neighbor site, $\mu$ is the chemical potential. Nondimensionalize with $t=1$ and $\mu,\Delta=$p.
Examples
julia> SSH(π/3, 0.5)
2×2 Matrix{ComplexF64}:
-0.5+0.0im 0.0+0.866025im
0.0-0.866025im 0.5+0.0imTopologicalNumbers.KitaevHoneycomb — MethodHamiltonian of the Kitaev honeycomb model.
KitaevHoneycomb(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Real}Arguments
k::T1: two-dimensional wavenumberk.p::T2: parameter defined as below.
Definition
Hamiltonian of the Kitaev honeycomb model is defined as
\[H(\bm{k})=\begin{pmatrix} h_{3}(\bm{k}) & h_{1}(\bm{k})-ih_{2}(\bm{k}) \\ h_{1}(\bm{k})+ih_{2}(\bm{k}) & -h_{3}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{1}(\bm{k})=&-K_{2}(\sin(k_{2})-\sin(k_{1})+\sin(k_{1}-k_{2})) \\ h_{2}(\bm{k})=&-K_{1}(\sin(k_{1})+\sin(k_{2})) \\ h_{3}(\bm{k})=&K_{1}(\cos(k_{1})+\cos(k_{2})+1) \end{align*}\]
where $K_{1}$ is the magnitude of Kitaev interaction. $K_{2}$ is the magnitude of spin triples term due to magnetic field. Nondimensionalize with $K_{1}=1$ and $K_{2}=$p.
Examples
julia> KitaevHoneycomb([0, π/3], 0.5)
2×2 Matrix{ComplexF64}:
2.5-0.0im 0.0+0.866025im
0.0-0.866025im -2.5-0.0imTopologicalNumbers.LatticeDirac — MethodHamiltonian of the lattice Dirac model.
LatticeDirac(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple},T2<:Real}Arguments
k::T1: four-dimensional wavenumberk.p::T2: parameter defined as below.
Definition
Hamiltonian of the lattice Dirac model is defined as
\[H(k)=\begin{pmatrix} h_{4}(\bm{k}) & h_{2}(\bm{k})-ih_{3}(\bm{k}) & h_{0}(\bm{k})-ih_{1}(\bm{k}) & 0 \\ h_{2}(\bm{k})+ih_{3}(\bm{k}) & -h_{4}(\bm{k}) & 0 & h_{0}(\bm{k})-ih_{1}(\bm{k}) \\ h_{0}(\bm{k})+ih_{1}(\bm{k}) & 0 & -h_{4}(\bm{k}) & -h_{2}(\bm{k})+ih_{3}(\bm{k}) \\ 0 & h_{0}(\bm{k})+ih_{1}(\bm{k}) & -h_{2}(\bm{k})-ih_{3}(\bm{k}) & h_{4}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{0}(\bm{k})=&m+c(cos(k_{1})+cos(k_{2})+cos(k_{3})+cos(k_{4})) \\ h_{1}(\bm{k})=&sin(k_{1}) \\ h_{2}(\bm{k})=&sin(k_{2}) \\ h_{3}(\bm{k})=&sin(k_{3}) \\ h_{4}(\bm{k})=&sin(k_{4}) \end{align*}\]
where,,,
Examples
```julia julia> LatticeDirac([0, 0, 0, π/2], 0.5) 4×4 Matrix{ComplexF64}: 1.0+0.0im 0.0+0.0im 3.5+0.0im 0.0+0.0im 0.0+0.0im -1.0+0.0im 0.0+0.0im 3.5+0.0im 3.5+0.0im 0.0+0.0im -1.0+0.0im 0.0+0.0im 0.0+0.0im 3.5+0.0im 0.0+0.0im 1.0+0.0im
TopologicalNumbers.SSH — MethodHamiltonian of the Su-Schrieffer-Heeger model.
SSH(k::T1, p::T2) where {T1<:Real,T2<:Real}Arguments
k::T1: one-dimensional wavenumberk.p::T2: parameter defined as below.
Definition
Hamiltonian of the Su-Schrieffer-Heeger model is defined as
\[H(k)=\begin{pmatrix} 0 & t_{1}+t_{2}e^{-ik} \\ t_{1}+t_{2}e^{ik} & 0 \end{pmatrix}\]
where $t_{1}$ and $t_{2}$ are the amplitudes of the nearest neighbor hopping in the tight binding model. Nondimensionalize with $t_{1}=1$ and $t_{2}=$p.
Examples
julia> SSH(π/3, 0.5)
2×2 Matrix{ComplexF64}:
0.0+0.0im 1.25-0.433013im
1.25+0.433013im 0.0+0.0imTopologicalNumbers.ThoulessPump — MethodHamiltonian of the Thouless pumping model.
ThoulessPump(k::T1, p::T2) where {T1<:Union{AbstractVector,Tuple}, T2<:Union{AbstractVector,Tuple}}Arguments
k::T1: when one-dimensional wavenumberk₁and timet,k=[k₁, t].p::T2: parameters defined as below.
Definition
Hamiltonian of the Kane–Mele model is defined as
\[H(k)=\begin{pmatrix} h_{3}(\bm{k}) & h_{5}(\bm{k})-ih_{4}(\bm{k}) & 0 & -h_{5}(\bm{k})-ih_{1}(\bm{k}) \\ h_{5}(\bm{k})+ih_{4}(\bm{k}) & -h_{3}(\bm{k}) & h_{5}(\bm{k})-ih_{1}(\bm{k}) & 0 \\ 0 & h_{5}(\bm{k})+ih_{1}(\bm{k}) & -h_{3}(\bm{k}) & h_{5}(\bm{k})-ih_{4}(\bm{k}) \\ -h_{5}(\bm{k})+ih_{1}(\bm{k}) & 0 & h_{5}(\bm{k})+ih_{4}(\bm{k}) & h_{3}(\bm{k}) \end{pmatrix}\]
\[\begin{align*} h_{1}(\bm{k})=&-\Delta\sin(k_{1}) \\ h_{2}(\bm{k})=&-\Delta(1-\cos(k_{1})) \\ h_{3}(\bm{k})=&m_{0}\sin(t) \\ h_{4}(\bm{k})=&-t_{0}(1+\cos(t))\sin(k_{1}) \\ h_{5}(\bm{k})=&-t_{0}(1-\cos(t))-t_{0}(1+\cos(t)))\cos(k_{1}) \end{align*}\]
where $t_{0}$ is the amplitude of the nearest neighbor hopping in the tight binding model. $m_{0}$ is the magnitude of the alternating magnetic field. $\Delta$ is the magnitude of spin flip interaction. Nondimensionalize with $t_{0}=1$ and $m_{0},\Delta=$p.
Examples
julia> ThoulessPump([0, π/3], (-1.0, 0.5))
4×4 Matrix{ComplexF64}:
-0.866025-0.0im -2.0+0.0im -0.0-0.0im 0.0+0.0im
-2.0-0.0im 0.866025-0.0im -0.0+0.0im -0.0-0.0im
-0.0-0.0im -0.0-0.0im 0.866025-0.0im -2.0+0.0im
0.0-0.0im -0.0-0.0im -2.0-0.0im -0.866025-0.0imTopologicalNumbers.calcBerryFlux — MethodCalculate the Berry flux in the two-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
calcBerryFlux(Hamiltonian::Function, n::Vector{Int64}; N::Int=51, gapless::Real=0.0, rounds::Bool=true)Arguments
- Hamiltionian::Function: the Hamiltonian matrix with one-dimensional wavenumber
kas an argument. - n::Vector{Int64}: The wavenumber($2\pi n/N$) when calculating Berry flux.
- N::Int=51: The number of meshes when discretizing the Brillouin Zone. It is preferable for
Nto be an odd number to increase the accuracy of the calculation. - gapless::Real: The threshold that determines the state to be degenerate. Coarsening the mesh(
N) but increasinggaplesswill increase the accuracy of the calculation. - rounds::Bool=true: An option to round the value of the topological number to an integer value. The topological number returns a value of type
Intwhentrue, and a value of typeFloatwhenfalse.
Definition
The Berry flux at the wavenumber $\bm{k}$ of the $n$th band $F_{n}(\bm{k})$ is defined by
\[F_{n}(\bm{k})=f_{n}(\bm{k})-df_{n}(\bm{k})\]
\[f_{n}(\bm{k})=\frac{1}{2\pi}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,1}^{*}(\bm{k})\right]\right]\]
\[df_{n}(\bm{k})=\frac{1}{2\pi}\mathrm{Im}\left[\mathrm{Log}[U_{n,1}(\bm{k})]+\mathrm{Log}[U_{n,2}(\bm{k}+\bm{e}_{1})]-\mathrm{Log}[U_{n,1}(\bm{k}+\bm{e}_{2})]-\mathrm{Log}[U_{n,1}(\bm{k})]\right]\]
$U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $n$th band.
TopologicalNumbers.calcBerryPhase — MethodCalculate the winding numbers in the one-dimensional case.
calcBerryPhase(Hamiltonian::Function; N::Int=51, gapless::Real=0.0, rounds::Bool=true)Arguments
Hamiltonian::Function: the Hamiltonian matrix function with one-dimensional wavenumberkas an argument.N::Int: the number of meshes when discretizing the Brillouin Zone. It is preferable forNto be an odd number to increase the accuracy of the calculation.gapless::Real: the threshold that determines the state to be degenerate. Coarsening the mesh(N) but increasinggaplesswill increase the accuracy of the calculation.rounds::Bool: an option to round the value of the topological number to an integer value. The topological number returns a value of typeIntwhentrue, and a value of typeFloatwhenfalse.
Definition
The Berry phase of the $n$th band $\nu_{n}$ is defined by
\[\nu_{n}=\frac{1}{\pi}\sum_{k\in\mathrm{BZ}}U_{n}(k)\]
The range $\mathrm{BZ}$(Brillouin Zone) is $k\in[0,2\pi]$. $U_{n,i}(k)$ is the link variable at wavenumber $k$. $e_{1}$ is the unit vector.
\[U_{n}(k)=\braket{\Psi_{n}(k)|\Psi_{n}(k+e_{1})}\]
$\ket{\Psi_{n}(k)}$ is the wave function of the $n$th band.
TopologicalNumbers.calcChern — MethodCalculate the first Chern numbers in the two-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
calcChern(Hamiltonian::Function; N::Int=51, gapless::Real=0.0, rounds::Bool=true)Arguments
- Hamiltionian::Function: The Hamiltonian matrix with two-dimensional wavenumber
kas an argument. - N::Int=51: The number of meshes when discretizing the Brillouin Zone. It is preferable for
Nto be an odd number to increase the accuracy of the calculation. - gapless::Real: The threshold that determines the state to be degenerate. Coarsening the mesh(
N) but increasinggaplesswill increase the accuracy of the calculation. - rounds::Bool=true: An option to round the value of the topological number to an integer value. The topological number returns a value of type
Intwhentrue, and a value of typeFloatwhenfalse.
Definition
The first Chern number of the $n$th band $\nu_{n}$ is defined by
\[\nu_{n}=\frac{1}{2\pi}\sum_{\bm{k}\in\mathrm{BZ}}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,2}^{*}(\bm{k})\right]\right]\]
The range $\mathrm{BZ}$(Brillouin Zone) is $\bm{k}\in[0,2\pi]^{2}$. $U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $n$th band.
TopologicalNumbers.calcChernSurface — MethodcalcChernSurface(H::Function, kn::String; kn_mesh::Int=51, N::Int=51, gapless::Real=0.0, rounds::Bool=true, plot::Bool=false)Arguments
- Hamiltionian::Function: The Hamiltonian matrix with three-dimensional wavenumber
kas an argument. - kn::String: Compute the Chern number of the plane perpendicular to the
"kn"direction in Brillouin zone ("k1","k2","k3"). - kn_mesh::T: Number of mesh in
"kn"direction. - N::Int=51: The number of meshes when discretizing the Brillouin Zone. It is preferable for
Nto be an odd number to increase the accuracy of the calculation. - gapless::Real: The threshold that determines the state to be degenerate. Coarsening the mesh(
N) but increasinggaplesswill increase the accuracy of the calculation. - rounds::Bool=true: An option to round the value of the topological number to an integer value. The topological number returns a value of type
Intwhentrue, and a value of typeFloatwhenfalse.
TopologicalNumbers.calcPhaseDiagram — MethodCalculate the one-dimensional phase diagram for the Berry phase over a specified parameter range.
calcPhaseDiagram(prob::BPProblem, param_range::T1, alg::T2=BP(); parallel::T3=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:BerryPhaseAlgorithms,T3<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the Berry phase for a one-dimensional system described by BPProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::BPProblem: A structure that contains the problem settings for Berry phase calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the Berry phase.alg::BerryPhaseAlgorithms=BP(): The algorithm used to compute the Berry phase. The default isBP(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the Berry phase calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param::AbstractVector: A copy of the inputparam_range.nums::AbstractMatrix: A matrix containing the computed Berry phases for each parameter. The size of this matrix is (number of bands × number of parameter values).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = SSH(k, p)
julia> # Set up the Berry phase problem
prob = BPProblem(H₀)
julia> # Define a parameter range from -2.0 to 2.0 with 9 divisions
param_range = range(-2.0, 2.0, length=9)
julia> # Calculate the phase diagram (using the default BP algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param_range; plot=true, progress=true)
(param = -2.0:0.5:2.0, nums = [1 1; 1 1; 0 0; 0 0; 0 0; 0 0; 0 0; 1 1; 1 1])
julia> result.param
-2.0:0.5:2.0
julia> result.nums
9×2 Matrix{Int64}:
1 1
1 1
0 0
0 0
0 0
0 0
0 0
1 1
1 1TopologicalNumbers.calcPhaseDiagram — MethodCalculate the one-dimensional phase diagram for the first Chern number over a specified parameter range.
calcPhaseDiagram(prob::FCProblem, param_range::T1, alg::T2=FHS(); parallel::T3=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:FirstChernAlgorithms,T3<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the first Chern number for a one-dimensional system described by FCProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::FCProblem: A structure that contains the problem settings for first Chern number calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the first Chern number.alg::FirstChernAlgorithms=FHS(): The algorithm used to compute the first Chern number. The default isFHS(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the first Chern number calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param::AbstractVector: A copy of the inputparam_range.nums::AbstractMatrix: A matrix containing the computed first Chern numbers for each parameter. The size of this matrix is (number of bands × number of parameter values).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = Haldane(k, p)
H(k, p) = H₀(k, (1, p, 2.5))
julia> # Set up the first Chern number problem
prob = FCProblem(H)
julia> # Define a parameter range
param = range(-π, π, length=10);
julia> # Calculate the phase diagram (using the default FHS algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param; plot=true, progress=true)
(param = -3.141592653589793:0.6981317007977318:3.141592653589793, nums = [0 0; -1 1; -1 1; -1 1; 0 0; 0 0; 1 -1; 1 -1; 1 -1; 0 0])
julia> result.param
-3.141592653589793:0.6981317007977318:3.141592653589793
julia> result.nums
10×2 Matrix{Int64}:
0 0
-1 1
-1 1
-1 1
0 0
0 0
1 -1
1 -1
1 -1
0 0TopologicalNumbers.calcPhaseDiagram — MethodCalculate the one-dimensional phase diagram for the second Chern number over a specified parameter range.
calcPhaseDiagram(prob::SCProblem, param_range::T1, alg::T2=FHS2(); parallel::T3=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:SecondChernAlgorithms,T3<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the second Chern number for a one-dimensional system described by SCProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::SCProblem: A structure that contains the problem settings for second Chern number calculations, including the Hamiltonian, mesh size, number of filled bands, gap information, and other parameters. The default is a half-filled band case.param_range::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the second Chern number.alg::SecondChernAlgorithms=FHS2(): The algorithm used to compute the second Chern number. The default isFHS2(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly.plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param::AbstractVector: A copy of the inputparam_range.nums::AbstractVector: A vector containing the computed second Chern numbers for each parameter.
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = LatticeDirac(k, p)
julia> N = 10; # number of k-points in each direction
julia> # Set up the second Chern number problem
prob = SCProblem(H₀, N)
julia> # Define a parameter range
param = range(-4.9, 4.9, length=4)
julia> # Calculate the phase diagram (using the default FHS2 algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param; plot=true, progress=true)
(param = -4.9:3.2666666666666666:4.9, nums = [0.0010237313095167136, -2.0667333080974726, 2.1572606447321445, -0.0009805850180973081])
julia> result.param
-4.9:3.2666666666666666:4.9
julia> result.nums
4-element Vector{Float64}:
0.0010237313095167136
-2.0667333080974726
2.1572606447321445
-0.0009805850180973081TopologicalNumbers.calcPhaseDiagram — MethodCalculate the one-dimensional phase diagram for the $\mathbb{Z}_2$ number over a specified parameter range.
calcPhaseDiagram(prob::Z2Problem, param_range::T1, alg::T2=Shio(); parallel::T3=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:Z2Algorithms,T3<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the $\mathbb{Z}_2$ number for a one-dimensional system described by Z2Problem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::Z2Problem: A structure that contains the problem settings for $\mathbb{Z}_2$ number calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the $\mathbb{Z}_2$ number.alg::Z2Algorithms=Shio(): The algorithm used to compute the $\mathbb{Z}_2$ number. The default isShio(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the $\mathbb{Z}_2$ number calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param::AbstractVector: A copy of the inputparam_range.nums::AbstractMatrix: A matrix containing the computed $\mathbb{Z}_2$ numbers for each parameter. The size of this matrix is (number of bands × number of parameter values).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = BHZ(k, p)
H(k, p) = H₀(k, (p, 0.25))
julia> # Set up the $\mathbb{Z}_2$ number problem
prob = Z2Problem(H)
julia> # Define a parameter range
param = range(-2, 2, length=10)
julia> # Calculate the phase diagram (using the default Shio algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param; plot=true, progress=true)
(param = -2.0:0.4444444444444444:2.0, nums = [0 0; 0 0; 0 0; 1 1; 1 1; 1 1; 1 1; 0 0; 0 0; 0 0])
julia> result.param
-2.0:0.4444444444444444:2.0
julia> result.nums
10×2 Matrix{Int64}:
0 0
0 0
0 0
1 1
1 1
1 1
1 1
0 0
0 0
0 0TopologicalNumbers.calcPhaseDiagram — MethodCalculate the two-dimensional phase diagram for the Berry phase over a specified parameter range.
calcPhaseDiagram(prob::BPProblem, param_range1::T1, param_range2::T2, alg::T3=BP(); parallel::T4=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:AbstractVector,T3<:BerryPhaseAlgorithms,T4<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the Berry phase for a two-dimensional system described by BPProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::BPProblem: A structure that contains the problem settings for Berry phase calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range1::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the Berry phase.param_range2::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the Berry phase.alg::BerryPhaseAlgorithms=BP(): The algorithm used to compute the Berry phase. The default isBP(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the Berry phase calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param1::AbstractVector: A copy of the inputparam_range1.param2::AbstractVector: A copy of the inputparam_range2.nums::AbstractArray: An array containing the computed Berry phases for each parameter. The size of this array is (number of bands × number of param1 × number of param2).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = KitaevChain(k, p)
julia> # Set up the Berry phase problem
prob = BPProblem(H₀)
julia> # Define a parameter range
param1 = range(-3.0, 3.0, length=4)
param2 = range(-1.0, 1.0, length=4)
julia> # Calculate the phase diagram (using the default BP algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param1, param2; plot=true, progress=true)
(param1 = -3.0:2.0:3.0, param2 = -1.0:0.6666666666666666:1.0, nums = [0 1 1 0; 0 1 1 0;;; 0 1 1 0; 0 1 1 0;;; 0 1 1 0; 0 1 1 0;;; 0 1 1 0; 0 1 1 0])
julia> result.param1
-3.0:2.0:3.0
julia> result.param2
-1.0:0.6666666666666666:1.0
julia> result.nums
2×4×4 Array{Int64, 3}:
[:, :, 1] =
0 1 1 0
0 1 1 0
[:, :, 2] =
0 1 1 0
0 1 1 0
[:, :, 3] =
0 1 1 0
0 1 1 0
[:, :, 4] =
0 1 1 0
0 1 1 0TopologicalNumbers.calcPhaseDiagram — MethodCalculate the two-dimensional phase diagram for the first Chern number over a specified parameter range.
calcPhaseDiagram(prob::FCProblem, param_range1::T1, param_range2::T2, alg::T3=FHS(); parallel::T4=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:AbstractVector,T3<:FirstChernAlgorithms,T4<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the first Chern number for a two-dimensional system described by FCProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::FCProblem: A structure that contains the problem settings for first Chern number calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range1::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the first Chern number.param_range2::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the first Chern number.alg::FirstChernAlgorithms=FHS(): The algorithm used to compute the first Chern number. The default isFHS(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the first Chern number calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param1::AbstractVector: A copy of the inputparam_range1.param2::AbstractVector: A copy of the inputparam_range2.nums::AbstractArray: An array containing the computed first Chern numbers for each parameter. The size of this array is (number of bands × number of param1 × number of param2).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = Haldane(k, p)
H(k, p) = H₀(k, (1, p[1], p[2]))
julia> # Set up the first Chern number problem
prob = FCProblem(H₀)
julia> # Define a parameter range
param1 = range(-π, π, length=4)
param2 = range(-6.0, 6.0, length=4)
julia> # Calculate the phase diagram (using the default FHS algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param1, param2; plot=true, progress=true)
(param1 = -3.141592653589793:2.0943951023931953:3.141592653589793, param2 = -6.0:4.0:6.0, nums = [0 0 0 0; 0 0 0 0;;; 0 -1 1 0; 0 1 -1 0;;; 0 -1 1 0; 0 1 -1 0;;; 0 0 0 0; 0 0 0 0])
julia> result.param1
-3.141592653589793:2.0943951023931953:3.141592653589793
julia> result.param2
-6.0:4.0:6.0
julia> result.nums
2×4×4 Array{Int64, 3}:
[:, :, 1] =
0 0 0 0
0 0 0 0
[:, :, 2] =
0 -1 1 0
0 1 -1 0
[:, :, 3] =
0 -1 1 0
0 1 -1 0
[:, :, 4] =
0 0 0 0
0 0 0 0TopologicalNumbers.calcPhaseDiagram — MethodcalcPhaseDiagram(H::Function, param_range::T1, alg::T2; N::T3=51, parallel::T4=UseSingleThread(), gapless::Real=0.0, rounds::Bool=true, plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:Union{String,TopologicalNumbersAlgorithms},T3<:Union{Int,Tuple,AbstractVector},T4<:TopologicalNumbersParallel}
TopologicalNumbers.calcPhaseDiagram — MethodCalculate the two-dimensional phase diagram for the second Chern number over a specified parameter range.
calcPhaseDiagram(prob::SCProblem, param_range1::T1, param_range2::T2, alg::T3=FHS2(); parallel::T4=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:AbstractVector,T3<:SecondChernAlgorithms,T4<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the second Chern number for a two-dimensional system described by SCProblem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::SCProblem: A structure that contains the problem settings for second Chern number calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range1::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the second Chern number.param_range2::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the second Chern number.alg::SecondChernAlgorithms=FHS2(): The algorithm used to compute the second Chern number. The default isFHS2(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the second Chern number calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param1::AbstractVector: A copy of the inputparam_range1.param2::AbstractVector: A copy of the inputparam_range2.nums::AbstractArray: An array containing the computed second Chern numbers for each parameter. The size of this array is (number of bands × number of param1 × number of param2).
TopologicalNumbers.calcPhaseDiagram — MethodCalculate the two-dimensional phase diagram for the $\mathbb{Z}_2$ number over a specified parameter range.
calcPhaseDiagram(prob::Z2Problem, param_range1::T1, param_range2::T2, alg::T3=Shio(); parallel::T4=UseSingleThread(), plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:AbstractVector,T3<:Z2Algorithms,T4<:TopologicalNumbersParallel}Description
calcPhaseDiagram calculates the $\mathbb{Z}_2$ number for a two-dimensional system described by Z2Problem over a range of parameters (param_range), constructing a phase diagram. By changing the alg or parallel options, you can flexibly switch between different calculation methods and parallelization strategies. Moreover, setting plot=true generates a simple plot of the results upon completion of the calculation.
Arguments
prob::Z2Problem: A structure that contains the problem settings for $\mathbb{Z}_2$ number calculations, including the Hamiltonian, mesh size, gap information, and other parameters.param_range1::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the $\mathbb{Z}_2$ number.param_range2::AbstractVector: A list or range of parameters (e.g.,range(-2.0, 2.0, length=50)) over which to scan. The Hamiltonian is evaluated at each parameter value in order to compute the $\mathbb{Z}_2$ number.alg::Z2Algorithms=Shio(): The algorithm used to compute the $\mathbb{Z}_2$ number. The default isShio(). You can implement and specify other algorithms as needed.parallel::TopologicalNumbersParallel=UseSingleThread(): The parallelization strategy. The default is single-threaded (UseSingleThread()). If you wish to run in a distributed environment, you can specify it accordingly. (Note that, at present, the $\mathbb{Z}_2$ number calculation does not support multithreading.)plot::Bool=false: If set totrue, a simple plot of the calculated results is shown after the computation finishes. This is helpful for visually confirming topological phase transitions.progress::Bool=false: If set totrue, progress is displayed during loop calculations, which can be useful for large-scale computations.
Returns
NamedTuple{(:param, :nums)}:param1::AbstractVector: A copy of the inputparam_range1.param2::AbstractVector: A copy of the inputparam_range2.nums::AbstractArray: An array containing the computed $\mathbb{Z}_2$ numbers for each parameter. The size of this array is (number of bands × number of param1 × number of param2).
Examples
julia> using TopologicalNumbers
julia> # Define the Hamiltonian
H₀(k, p) = BHZ(k, p)
julia> # Set up the $\mathbb{Z}_2$ number problem
prob = Z2Problem(H₀)
julia> # Define a parameter range
param1 = range(-2, 2, length=4);
param2 = range(-0.4, 0.4, length=4)
julia> # Calculate the phase diagram (using the default Shio algorithm, single-threaded,
# with plotting enabled and progress display)
result = calcPhaseDiagram(prob, param1, param2; plot=true, progress=true)
(param1 = -2.0:1.3333333333333333:2.0, param2 = -0.4:0.26666666666666666:0.4, nums = [0 1 1 0; 0 1 1 0;;; 0 0 0 0; 0 0 0 0;;; 0 0 0 0; 0 0 0 0;;; 0 1 1 0; 0 1 1 0])
julia> result.param1
-2.0:1.3333333333333333:2.0
julia> result.param2
-0.4:0.26666666666666666:0.4
julia> result.nums
2×4×4 Array{Int64, 3}:
[:, :, 1] =
0 1 1 0
0 1 1 0
[:, :, 2] =
0 0 0 0
0 0 0 0
[:, :, 3] =
0 0 0 0
0 0 0 0
[:, :, 4] =
0 1 1 0
0 1 1 0TopologicalNumbers.calcPhaseDiagram — MethodcalcPhaseDiagram(H::Function, param_range1::T1, param_range2::T2, alg::T3; N::T4=51, parallel::T5=UseSingleThread(), gapless::Real=0.0, rounds::Bool=true, plot::Bool=false, progress::Bool=false) where {T1<:AbstractVector,T2<:AbstractVector,T3<:Union{String,TopologicalNumbersAlgorithms},T4<:Union{Int,Tuple,AbstractVector},T5<:TopologicalNumbersParallel}TopologicalNumbers.calcSecondChern — MethodCalculate the second Chern numbers in the four-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1, 2].
calcSecondChern(Hamiltonian::Function; Nfill::T1=nothing, N::T2=(30, 30, 30, 30), returnRealValue::Bool=true, parallel::T3=UseSingleThread()) where {T1<:Union{Int,Nothing},T2<:Union{AbstractVector,Tuple},T3<:TopologicalNumbersParallel}Arguments
Hamiltionian: The Hamiltonian matrix with two-dimensional wavenumberkas an argument.Nfill::T1: The filling number. The default value isHs ÷ 2, whereHsis the size of the Hamiltonian matrix.N::T2: The numbers of meshes when discretizing the Brillouin Zone. Each element ofNis the number of meshes in the x, y, z, and w directions, respectively.returnRealValue::Bool: An option to return the value of the topological number by an real value. The topological number returns a value of typeFloat64whentrue, and a value of typeComplexF64whenfalse.
Examples
TopologicalNumbers.calcWeylNode — MethodCalculate the Weyl node in the three-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
calcWeylNode(Hamiltonian::Function, n::Vector{Int64}; N::Int=51, gapless::Real=0.0, rounds::Bool=true)Arguments
- Hamiltionian::Function: the Hamiltonian matrix with three-dimensional wavenumber
kas an argument. - n::Vector{Int64}: The wavenumber($2\pi\bm{n}/N$) when calculating Weyl node.
- N::Int=51: The number of meshes when discretizing the Brillouin Zone. It is preferable for
Nto be an odd number to increase the accuracy of the calculation. - gapless::Real: The threshold that determines the state to be degenerate. Coarsening the mesh(
N) but increasinggaplesswill increase the accuracy of the calculation. - rounds::Bool=true: An option to round the value of the topological number to an integer value. The topological number returns a value of type
Intwhentrue, and a value of typeFloatwhenfalse.
TopologicalNumbers.calcZ2 — MethodCalculate the $\mathbb{Z}_2$ numbers in the two-dimensional case with reference to Shiozaki method [3, 4].
calcZ2(Hamiltonian::Function; Nfill::T1=nothing, N::Int=50, rounds::Bool=true, TR::Bool=false) where {T1<:Union{Int,Nothing}}Arguments
Hamiltonian::Functionis a matrix with one-dimensional wavenumberkas an argument.Nfill::T1: The filling number. The default value isHs ÷ 2, whereHsis the size of the Hamiltonian matrix.N::Intis the number of meshes when discretizing the Brillouin Zone. It is preferable forNto be an odd number to increase the accuracy of the calculation.rounds::Boolis an option to round the value of the topological number to an integer value. The topological number returns a value of typeIntwhentrue, and a value of typeFloatwhenfalse.
Definition
The $\mathbb{Z}_{2}$ number of the $2n$th (and $2n-1$th) band $\nu_{n}$ is defined by
\[\nu_{n}=F_{n}-\left(P_{n}(0)-P_{n}(\pi)\right)\]
$F_{n}$ is the Berry flux of the $n$th band in the $\mathrm{BZ}'$. The range $\mathrm{BZ}'$ is $\bm{k}\in[0,2\pi]\times[0,\pi]$ half of BZ(Brillouin Zone).
\[F_{n}=\frac{1}{2\pi}\sum_{\bm{k}\in\mathrm{BZ}'}\mathrm{Im}\left[\mathrm{Log}\left[U_{n,1}(\bm{k})U_{n,2}(\bm{k}+\bm{e}_{1})U_{n,1}^{*}(\bm{k}+\bm{e}_{2})U_{n,1}^{*}(\bm{k})\right]\right]\]
$P_{n}(k_{2})$ is the time-reversal polarization at wavenumber $k_{2}$.
\[P_{n}(k_{2})=\frac{1}{2\pi}\frac{\mathrm{Pf}[\omega(0,k_{2})]}{\mathrm{Pf}[\omega(\pi,k_{2})]}\sum_{k_{1}=0}^{\pi-e_{1}}U_{n,1}(\bm{k})\]
$U_{n,i}(\bm{k})$ is the link variable at wavenumber $\bm{k}$. $\bm{e}_{i}$ is the unit vector.
\[U_{n,i}(\bm{k})=\braket{\Psi_{n}(\bm{k})|\Psi_{n}(\bm{k}+\bm{e}_{i})}\]
$\ket{\Psi_{n}(\bm{k})}$ is the wave function of the $2n$th (and $2n-1$th) band. $\omega(\bm{k})$ is the unitary matrix given by
\[\omega(\bm{k})=\bra{\Psi(-\bm{k})}T\ket{\Psi(\bm{k})}\]
$T$ is the time-reversal operator.
TopologicalNumbers.findWeylPoint — MethodfindWeylPoint(Hamiltonian::Function; N::Int=10, gapless::T=[1e-1, 1e-2, 1e-3, 1e-4], rounds::Bool=true) where {T<:AbstractVector{Float64}}Arguments
- Hamiltionian::Function: The Hamiltonian matrix with three-dimensional wavenumber
kas an argument. - N::Int=10: The number of meshes when discretizing the Brillouin Zone. The $n$th iteration divides the Brillouin zone into $1/N^n$.
- gapless<:AbstractVector{Float64}: The threshold that determines the state to be degenerate. The $n$th iteration adopts the threshold value of the $n$th value of the vector. The number of iterations can be varied by the length of the vector.
- rounds::Bool=true: An option to round the value of the topological number to an integer value. The topological number returns a value of type
Intwhentrue, and a value of typeFloatwhenfalse.
TopologicalNumbers.pfaffian — Methodpfaffian(A::AbstractMatrix{T}, overwrite_a=false, method='P') where {T<:Number}Compute the Pfaffian of a real or complex skew-symmetric matrix A (A=-A^T). If overwrite_a=true, the matrix A is overwritten in the process. This function uses either the Parlett-Reid algorithm (method='P', default), or the Householder tridiagonalization (method='H')
TopologicalNumbers.pfaffian_LTL — Methodpfaffian_LTL(A::AbstractMatrix{T}; overwrite_a=false) where {T<:Number}Compute the Pfaffian of a real or complex skew-symmetric matrix A (A=-A^T). If overwrite_a=true, the matrix A is overwritten in the process. This function uses the Parlett-Reid algorithm.
TopologicalNumbers.pfaffian_householder — Methodpfaffian_householder(A::AbstractMatrix{T}; overwrite_a=false) where {T<:Complex}Compute the Pfaffian of a real or complex skew-symmetric matrix A (A=-A^T). If overwrite_a=true, the matrix A is overwritten in the process. This function uses the Householder tridiagonalization.
Note that the function pfaffianschur() can also be used in the real case. That function does not make use of the skew-symmetry and is only slightly slower than pfaffianhouseholder().
TopologicalNumbers.pfaffian_householder — Methodpfaffian_householder(A::AbstractMatrix{T}; overwrite_a=false) where {T<:Real}Compute the Pfaffian of a real or complex skew-symmetric matrix A (A=-A^T). If overwrite_a=true, the matrix A is overwritten in the process. This function uses the Householder tridiagonalization.
Note that the function pfaffianschur() can also be used in the real case. That function does not make use of the skew-symmetry and is only slightly slower than pfaffianhouseholder().
TopologicalNumbers.pfaffian_schur — Methodpfaffian_schur(A::AbstractMatrix{T}; overwrite_a=false) where {T<:Real}Calculate Pfaffian of a real antisymmetric matrix using the Schur decomposition. (Hessenberg would in principle be faster, but scipy-0.8 messed up the performance for scipy.linalg.hessenberg()).
This function does not make use of the skew-symmetry of the matrix A, but uses a LAPACK routine that is coded in FORTRAN and hence faster than python. As a consequence, pfaffian_schur is only slightly slower than pfaffian().
TopologicalNumbers.plot1D — Methodplot1D(result::NamedTuple; labels::Bool=true, disp::Bool=true, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="phaseDiagram")Plot a 1D phase diagram.
Arguments
result::NamedTuple: A named tuple containing the result data.labels::Bool=true: Whether to display axis labels.disp::Bool=true: Whether to display the plot.png::Bool=false: Whether to save the plot as a PNG file.pdf::Bool=false: Whether to save the plot as a PDF file.svg::Bool=false: Whether to save the plot as an SVG file.filename::String="phaseDiagram": The filename for the saved plot.
Returns
fig: The matplotlib figure object.
Example
julia>TopologicalNumbers.plot1D — Methodplot1D(nums::T1, param_range::T2; labels::Bool=true, disp::Bool=true, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="phaseDiagram") where {T1<:AbstractMatrix,T2<:AbstractVector}Plot a 1D phase diagram.
Arguments
nums::T1: A matrix or vector containing the data to be plotted.param_range::T2: A vector representing the parameter range.labels::Bool=true: Whether to display axis labels. Default istrue.disp::Bool=true: Whether to display the plot. Default istrue.png::Bool=false: Whether to save the plot as a PNG file. Default isfalse.pdf::Bool=false: Whether to save the plot as a PDF file. Default isfalse.svg::Bool=false: Whether to save the plot as an SVG file. Default isfalse.filename::String="phaseDiagram": The filename for the saved plot. Default is "phaseDiagram".
Examples
nums = [1 2 3; 4 5 6]
param_range = [0.1, 0.2, 0.3]
plot1D(nums, param_range)Returns
fig: The matplotlib figure object representing the plot.
TopologicalNumbers.plot1D — Methodplot1D(::T1, nums::T2, param_range::T3; labels::Bool=true, disp::Bool=true, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="phaseDiagram") where {T1<:SecondChernAlgorithms,T2<:AbstractVector,T3<:AbstractVector}Plot a 1D phase diagram.
Arguments
::T1: An object of typeT1that implements theSecondChernAlgorithmsinterface.nums::T2: A vector of numbers representing the values to be plotted.param_range::T3: A vector representing the parameter range.labels::Bool=true: Whether to display labels on the plot. Default istrue.disp::Bool=true: Whether to display the plot. Default istrue.png::Bool=false: Whether to save the plot as a PNG file. Default isfalse.pdf::Bool=false: Whether to save the plot as a PDF file. Default isfalse.svg::Bool=false: Whether to save the plot as an SVG file. Default isfalse.filename::String="phaseDiagram": The filename for the saved plot. Default is"phaseDiagram".
Examples
plot1D(obj, nums, param_range)Returns
fig: The matplotlib figure object representing the plot.
TopologicalNumbers.plot2D — Methodplot2D(result::NamedTuple; labels::Bool=true, disp::Bool=true, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="phaseDiagram")Plot a 2D phase diagram.
Arguments
result::NamedTuple: A named tuple containing the result data.labels::Bool=true: Whether to display axis labels.disp::Bool=true: Whether to display the plot.png::Bool=false: Whether to save the plot as a PNG file.pdf::Bool=false: Whether to save the plot as a PDF file.svg::Bool=false: Whether to save the plot as an SVG file.filename::String="phaseDiagram": The filename for the saved plot.
Returns
fig: The matplotlib figure object.
Example
julia>TopologicalNumbers.plot2D — Methodplot2D(nums::T1, param_range1::T2, param_range2::T3; labels::Bool=true, disp::Bool=true, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="phaseDiagram") where {T1<:AbstractArray,T2<:AbstractVector,T3<:AbstractVector}Plot a 2D phase diagram.
Arguments
nums::T1: An array of numbers representing the phase diagram.param_range1::T2: A vector representing the range of the first parameter.param_range2::T3: A vector representing the range of the second parameter.labels::Bool=true: Whether to display axis labels. Default istrue.disp::Bool=true: Whether to display the plot. Default istrue.png::Bool=false: Whether to save the plot as a PNG file. Default isfalse.pdf::Bool=false: Whether to save the plot as a PDF file. Default isfalse.svg::Bool=false: Whether to save the plot as an SVG file. Default isfalse.filename::String="phaseDiagram": The filename for the saved plot. Default is "phaseDiagram".
Returns
fig: The matplotlib figure object.
Example
julia>TopologicalNumbers.showBand — MethodshowBand(Hamiltonian::Function; N::Int=51, labels::Bool=true, value::Bool=true, disp::Bool=false, png::Bool=false, pdf::Bool=false, svg::Bool=false, filename::String="Band")This function generates a band structure plot for a given Hamiltonian.
Arguments
Hamiltonian::Function: The Hamiltonian function that takes a wave number parameterkand returns the corresponding Hamiltonian matrix.N::Int: The number of points in the Brillouin zone. Default is 51.labels::Bool: Whether to display the labels of the figure. Default istrue.value::Bool: Whether to output the values of the wave number and the energy in the matrix form. Default istrue.disp::Bool: Whether to display the plot. Default isfalse.png::Bool: Whether to save the plot as a PNG file. Default isfalse.pdf::Bool: Whether to save the plot as a PDF file. Default isfalse.svg::Bool: Whether to save the plot as an SVG file. Default isfalse.filename::String: The filename for the saved plot. Default is "Band".
Examples
julia> H(k) = SSH(k, (0.9, 1.0)) # $N \times N$ Hamiltonian matrix with a wavenumber parameter k
julia> showBand(H)
(k = -3.141592653589793:0.12566370614359174:3.141592653589793, Ene = [-0.09999999999999998 0.09999999999999998; -0.15554271964299698 0.15554271964299698; … ; -0.15554271964299698 0.15554271964299698; -0.09999999999999998 0.09999999999999998])TopologicalNumbers.skew_LTL — MethodT, L, P = skew_LTL(A::AbstractMatrix{T}; overwrite_a=false, calc_L=true, calc_P=true) where {T<:Number}Bring a real or complex skew-symmetric matrix (A=-A^T) into tridiagonal form T (with zero diagonal) with a lower unit triangular matrix L such that P A P^T= L T L^T
A is overwritten if overwritea=true (default: false), L and P only calculated if calcL=true or calc_P=true, respectively (default: true).
TopologicalNumbers.skew_tridiagonalize — MethodT, Q = skew_tridiagonalize(A::AbstractMatrix{T}; overwrite_a=false, calc_q=true) where {T<:Number}or
T = skew_tridiagonalize(A::AbstractMatrix{T}; overwrite_a=false, calc_q=false) where {T<:Number}Bring a real or complex skew-symmetric matrix (A=-A^T) into tridiagonal form T (with zero diagonal) with a orthogonal (real case) or unitary (complex case) matrix U such that A = Q T Q^T (Note that Q^T and not Q^dagger also in the complex case)
A is overwritten if overwritea=true (default: false), and Q only calculated if calcq=true (default: true)
TopologicalNumbers.solve — MethodCalculate the winding numbers in the one-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
solve(prob::BPProblem, alg::T1=BP(); parallel::T2=UseSingleThread()) where {T1<:BerryPhaseAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::BPProblem: The BPProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=BP(): The algorithm to use for calculating the Berry phases. Default isBPalgorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
BPSolution: A struct that contains the calculated Berry phases.
Examples
julia> H(k) = SSH(k, 1.1)
julia> prob = BPProblem(H)
julia> result = solve(prob)
BPSolution{Vector{Int64}, Int64}([1, 1], 0)
julia> result.TopologicalNumber
2-element Vector{Int64}:
1
1TopologicalNumbers.solve — MethodCalculate the first Chern numbers in the two-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
solve(prob::FCProblem, alg::T1=FHS(); parallel::T2=UseSingleThread()) where {T1<:FirstChernAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::FCProblem: The FCProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=FHS(): The algorithm to use for calculating the first Chern numbers. Default isFHSalgorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
FCSolution: A struct that contains the calculated first Chern numbers.
Examples
julia> H(k) = Haldane(k, (1, 0.5, 1.0))
julia> prob = FCProblem(H)
julia> result = solve(prob)
FCSolution{Vector{Int64}, Int64}([1, -1], 0)
julia> result.TopologicalNumber
2-element Vector{Int64}:
1
-1TopologicalNumbers.solve — MethodCalculate the Berry flux in the two-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
solve(prob::LBFProblem, alg::T1=FHSlocal2(); parallel::T2=UseSingleThread()) where {T1<:BerryFluxAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::LBFProblem: The LBFProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=FHSlocal2(): The algorithm to use for calculating the second Chern numbers. Default isFHSlocal2algorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
LBFSolution: A struct that contains the calculated Berry flux.
Examples
julia> H(k) = Flux2d(k, (6, 1))
julia> n = [0, 0]
julia> prob = LBFProblem(H, n)
julia> result = solve(prob)
LBFSolution{Vector{Int64}, Vector{Int64}}([0, 0, 0, 0, -1, 0], [0, 0])
julia> result.TopologicalNumber
6-element Vector{Int64}:
0
0
0
0
-1
0TopologicalNumbers.solve — MethodCalculate the second Chern numbers in the four-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1, 2].
solve(prob::SCProblem, alg::T1=FHS2(); parallel::T2=UseSingleThread()) where {T1<:SecondChernAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::SCProblem: The SCProblem struct that contains the Hamiltonian matrix function in the wave number space, filling number, mesh numbers, and other parameters.alg::T1=FHS2(): The algorithm to use for calculating the second Chern numbers. Default isFHS2algorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
SCSolution: A struct that contains the calculated second Chern numbers.
Examples
julia> H(k) = LatticeDirac(k, -3.0)
julia> prob = SCProblem(H)
julia> result = solve(prob)
SCSolution{Float64}(0.9793607631927381)
julia> result.TopologicalNumber
0.9793607631927381TopologicalNumbers.solve — MethodCalculate the sliced first Chern numbers in the three-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
solve(prob::WCSProblem, alg::T1=FHSsurface(); parallel::T2=UseSingleThread()) where {T1<:WeylPointsAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::WCSProblem: The WCSProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=FHSsurface(): The algorithm to use for calculating the sliced first Chern numbers. Default isFHSsurfacealgorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
WCSSolution: A struct that contains the calculated sliced first Chern numbers.
Examples
julia> function H₀(k, p) # Weyl
k1, k2, k3 = k
t1, t2, t3, m, k0 = p
h0 = 0
hx = 2t1*(cos(k1) - cos(k0)) + m*(2 - cos(k2) - cos(k3))
hy = 2t2*sin(k2)
hz = 2t3*sin(k3)
s0 = [1 0; 0 1]
sx = [0 1; 1 0]
sy = [0 -im; im 0]
sz = [1 0; 0 -1]
h0 .* s0 .+ hx .* sx .+ hy .* sy .+ hz .* sz
end
julia> p0 = (1, 1, 1, 2, 2pi*2/5);
julia> H(k) = H₀(k, p0);
julia> prob = WCSProblem(H, "k1");
julia> sol = solve(prob)
WCSSolution{String, StepRangeLen{Float64, Base.TwicePrecision{Float64}, Base.TwicePrecision{Float64}, Int64}, Matrix{Int64}}("k1", 6.283185307179587e-5:0.12319971190548208:6.160048427127176, [0 0; 0 0; … ; 0 0; 0 0])
julia> sol.nums
51×2 Matrix{Int64}:
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
0 0
0 0
0 0
0 0
1 -1
1 -1
1 -1
1 -1
1 -1
1 -1
1 -1
1 -1
1 -1
1 -1
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
0 0
0 0
0 0TopologicalNumbers.solve — MethodCalculate the Weyl node in the three-dimensional case with reference to Fukui-Hatsugai-Suzuki method [1].
solve(prob::WNProblem, alg::T1=FHSlocal3(); parallel::T2=UseSingleThread()) where {T1<:WeylPointsAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::WNProblem: The WNProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=FHSlocal3(): The algorithm to use for calculating the Weyl nodes. Default isFHSlocal3algorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
WNSolution: A struct that contains the calculated Weyl nodes.
Examples
julia> function H₀(k, p) # Weyl
k1, k2, k3 = k
t1, t2, t3, m, k0 = p
h0 = 0
hx = 2t1*(cos(k1) - cos(k0)) + m*(2 - cos(k2) - cos(k3))
hy = 2t2*sin(k2)
hz = 2t3*sin(k3)
s0 = [1 0; 0 1]
sx = [0 1; 1 0]
sy = [0 -im; im 0]
sz = [1 0; 0 -1]
h0 .* s0 .+ hx .* sx .+ hy .* sy .+ hz .* sz
end
julia> p0 = (1, 1, 1, 2, 2pi*2/5);
julia> H(k) = H₀(k .- 2pi*1e-5, p0);
julia> weylPoint = [4000, 0, 0];
julia> N = 10000;
julia> prob = WNProblem(H, weylPoint, N);
julia> sol = solve(prob)
WNSolution{Vector{Int64}, Vector{Int64}, Int64}([1, -1], [4000, 0, 0], 10000)
julia> sol.TopologicalNumber
2-element Vector{Int64}:
1
-1TopologicalNumbers.solve — MethodCalculate the Weyl points in the three-dimensional case using energy variational method.
solve(prob::WPProblem, alg::T1=Evar(); parallel::T2=UseSingleThread()) where {T1<:WeylPointsAlgorithms,T2<:TopologicalNumbersParallel}Arguments
prob::WPProblem: The WPProblem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=Evar(): The algorithm to use for calculating the Weyl points. Default isEvaralgorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
WPSolution: A struct that contains the calculated Weyl points.
Examples
julia> function H₀(k, p) # Weyl
k1, k2, k3 = k
t1, t2, t3, m, k0 = p
h0 = 0
hx = 2t1*(cos(k1) - cos(k0)) + m*(2 - cos(k2) - cos(k3))
hy = 2t2*sin(k2)
hz = 2t3*sin(k3)
s0 = [1 0; 0 1]
sx = [0 1; 1 0]
sy = [0 -im; im 0]
sz = [1 0; 0 -1]
h0 .* s0 .+ hx .* sx .+ hy .* sy .+ hz .* sz
end
julia> p0 = (1, 1, 1, 2, 2pi*2/5);
julia> H(k) = H₀(k, p0);
julia> prob = WPProblem(H)
julia> result = solve(prob)
WPSolution{Vector{Vector{Vector{Int64}}}, Int64, Vector{Vector{Int64}}}([[[4000, 0, 0], [6000, 0, 0]], [[4000, 0, 0], [6000, 0, 0]]], 10000, [[1, -1], [-1, 1]])
julia> 2pi*result.WeylPoint[1] / result.N .- pi*[ones(3), ones(3)]
2-element Vector{Vector{Float64}}:
[-0.6283185307179586, -3.141592653589793, -3.141592653589793]
[0.6283185307179586, -3.141592653589793, -3.141592653589793]TopologicalNumbers.solve — MethodCalculate the $\mathbb{Z}_2$ numbers in the two-dimensional case with reference to Shiozaki method [3, 4].
solve(prob::Z2Problem, alg::T1=Shio(); parallel::T2=UseSingleThread()) where {T1<:Z2Algorithms,T2<:TopologicalNumbersParallel}Arguments
prob::Z2Problem: The Z2Problem struct that contains the Hamiltonian matrix function in the wave number space and other parameters.alg::T1=Shio(): The algorithm to use for calculating the Z2 numbers. Default isShioalgorithm.parallel::T2=UseSingleThread(): The parallelization strategy to use. Default is to use a single thread.
Returns
Z2Solution: A struct that contains the calculated Z2 numbers.
Examples
julia> H(k) = KaneMele(k, 1.0)
julia> prob = Z2Problem(H)
julia> result = solve(prob)
Z2Solution{Vector{Int64}, Nothing, Int64}([1, 1], nothing, 0)
julia> result.TopologicalNumber
2-element Vector{Int64}:
1
1