Fuzzifino

Fuzzifino is a module for exact diagonalisation (ED) calculation on the fuzzy sphere for systems with both bosons and fermions. The usage is similar to FuzzifiED, with new types SQNDiag, SQNOffd, SConf, SBasis, STerm and SOperator defined. To use the module, include also at the start of your Julia script

using FuzzifiED.Fuzzifino

Environment parameter

FuzzifiED.Fuzzifino.LibpathinoConstant
Libpathino :: String = FuzzifiED_jll.LibpathFuzzifino

define path of the Fortran library libfuzzifino.so. You do not need to modify that by yourself. However, if you compile the Fortran codes by yourself, you need to point this to your compiled library.

source

Quantum numbers

The diagonal and off-diagonal quantum numbers are implemented as

FuzzifiED.Fuzzifino.SQNDiagType
SQNDiag

The mutable type SQNDiag records the information of a diagonal $\mathrm{U}(1)$ or $ℤ_p$ quantum number in the form of a symmetry charge

\[Q=∑_{o=1}^{N_{of}}q_{f,o}n_{f,o}+∑_{o=1}^{N_{ob}}q_{b,o}n_{b,o}\]

or

\[Q=∑_{o=1}^{N_{of}}q_{f,o}n_{f,o}+∑_{o=1}^{N_{ob}}q_{b,o}n_{b,o}\ \mathrm{mod}\ p\]

where $i=1,…,N_U$ is the index of quantum number, $o$ is the index of site, $N_{of}$ and $N_{ob}$ are the number of fermionic and bosonic sites, $n_{f,o}=f^†_of_o$, $n_{b,o}=b^†_ob_o$, and $q_{f,o},q_{b,o}$ are a set of symmetry charges that must be integer valued.

Fields

  • name :: String is the name of the diagonal quantum number
  • chargef :: Vector{Int64} is the symmetry charge $q_{f,o}$ of each site
  • chargeb :: Vector{Int64} is the symmetry charge $q_{b,o}$ of each site
  • modul :: Vector{Int64} is the modulus $p$, set to 1 for $\mathrm{U}(1)$ SQNDiags.

Initialisation

It can be initialised by the following method

SQNDiag([name :: String, ]chargef :: Vector{Int64}, chargeb :: Vector{Int64}[, modul :: Int64]) :: SQNDiag

The arguments name and modul are facultative. By default name is set to "QN" and modul is set to 1.

source
FuzzifiED.Fuzzifino.SQNOffdType
SQNOffd

The mutable type SQNOffd records the information of an off-diagonal $ℤ_p$ quantum number in the form of a discrete transformation

\[𝒵:\ f_o↦ α_{f,o}^* f^{(p_{f,o})}_{π_{f,o}},  f_o^†↦α_{f,o} c^{(1-p_{f,o})}_{π_{f,o}},  b_o^†↦α_{b,o} b^†_{π_{b,o}}\]

where we use a notation $c^{(1)}=c^†$ and $c^{0}=c$ for convenience, $π_{f,o},π_{b,o}$ are permutations of $1,…,N_{of}$ or $N_{ob}$, $α_{f,o},α_{b,o}$ are coefficients, and $p_{f,o}$ specified whether or not particle-hole transformation is performed for the fermionic site. Note that one must guarentee that all these transformations commute with each other and also commute with the diagonal QNs.

Arguments

  • permf :: Vector{Int64} is a length-$N_{of}$ vector that records the fermion permutation $π_{f,o}$.
  • permb :: Vector{Int64} is a length-$N_{ob}$ vector that records the boson permutation $π_{b,o}$.
  • phf :: Vector{Int64} is a length-$N_{of}$ vector that records $p_{f,o}$ to determine whether or not to perform a particle-hole transformation
  • facf :: Vector{ComplexF64} is a length-$N_{of}$ vector that records the factor $α_{f,o}$ in the transformation.
  • facb :: Vector{ComplexF64} is a length-$N_{ob}$ vector that records the factor $α_{b,o}$ in the transformation.
  • cyc :: Int64 is the cycle $p$.

Initialisation

It can be initialised by the following method

SQNOffd(permf :: Vector{Int64}, permb :: Vector{Int64}[, phf :: Vector{Int64}][, facf :: Vector{ComplexF64}, facb :: Vector{ComplexF64}][, cyc :: Int64]) :: SQNOffd
SQNOffd(permf :: Vector{Int64}, permb :: Vector{Int64}, phf_q :: Bool[, fac :: Vector{ComplexF64}, facb :: Vector{ComplexF64}]) :: SQNOffd

The arguments phf, facf, facb and cyc are facultative. By default ph is set all 0, facf, facb is set to all 1 and cyc is set to 2. If phf_q is a bool and true, then ph is set to all 1.

source

Configurations

FuzzifiED.Fuzzifino.SConfsType
SConfs

The mutable type SConfs stores all the configurations that respects the diagonal quantum numbers (SQNDiag) and also a table to inversely look up the index from the configuration.

Fields

  • nof :: Int64 is the number of fermionic sites.
  • nof :: Int64 is the number of bosonic sites.
  • nebm :: Int64 is the maximal number of boson occupation.
  • ncf :: Int64 is the number of configurations.
  • conff :: Vector{Int64} is an array of length ncf containing all the fermion configurations. Each configuration is expressed in a binary number. If the o-1-th bit of conf[i] is 1, then the o-th site in the i-th configuration is occupied ; if the bit is 0, then the site is empty.
  • confb :: Vector{Int64} is an array of length ncf containing all the boson configurations. Each configuration is expressed in a binary number that has $N_{b,o}$ 1's and $N_{eb}$ 0's and the number of 0's following each 1 records the number of bosons in that site.
  • norf :: Int64, nobf :: Int64, lid :: Vector{Int64} and rid :: Vector{Int64} contain the information of Lin table that is used to inversely look up the index i from the configuration.
source

It can be generated from the QNDiags.

FuzzifiED.Fuzzifino.SConfsMethod
SConfs(nof :: Int64, nob :: Int64, nebm :: Int64, secd :: Vector{Int64}, qnd :: Vector{SQNDiag} ; , num_th :: Int64, disp_std :: Bool) :: Confs

generates the configurations from the list of QNDiags.

Arguments

  • nof :: Int64 is the number of fermionic sites $N_{of}$.
  • nob :: Int64 is the number of bosonic sites $N_{ob}$.
  • nebm :: Int64 is the maximal number of total bosons.
  • secd :: Vector{Int64} is the set of $Q_i$ for the selected configurations in the sector.
  • qnd :: Vector{SQNDiag} is the set of SQNDiags.
  • norf :: Int64 and norb :: Int64 are the number of less significant bits used to generate the Lin table. Facultative, $N_{of}/2$ and $N_{ob}/2$ by default.
  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.

Output

  • cfs :: SConfs is a SConfs object.

Note

If your qnd has negative entries, QNDiags must contain the total number of particles (i.e., bosons plus fermions).

source

Basis

FuzzifiED.Fuzzifino.SBasisType
SBasis

The mutable type SBasis stores the information of the SBasis that respects both diagonal and off-diagonal quantum numbers. The states in the SBasis is in the form

\[|I⟩=λ_{i_{I1}}|i_{I1}⟩+λ_{i_{I2}}|i_{I2}⟩+⋯+λ_{i_{Im_I}}|i_{Im_I}⟩\]

where $|i⟩$ is a direct product state, i.e., the configurations $|i_{Ik}⟩$ are grouped into a state $|I⟩$.

Fields

  • cfs :: SConfs stores the configurations that respect the QNDiags ;
  • dim :: Int64 is the dimension of the SBasis ;
  • szz :: Int64 records the maximum size $\max m_g$ of groups;
  • cfgr :: Vector{Int64} is a vector of length cfs.ncf and records which group $|I⟩$ each configuration $|i⟩$ belong to ;
  • cffac :: Vector{ComplexF64} is a vector of length cfs.ncf and records the coefficients $λ_i$ of each configuration ;
  • grel :: Matrix{Int64} is a szz*dim matrix that records the configurations in each group $|i_{Ik}⟩ (k = 1,…,m_I)$
  • grsz :: Vector{Int64} is a vector of length dim that records the size $m_I$ of each group.
source

It can be generated by the following methods.

FuzzifiED.Fuzzifino.SBasisMethod
SBasis(cfs :: SConfs, secf :: Vector{ComplexF64}, qnf :: Vector{SQNOffd} ; num_th :: Int64, disp_std :: Bool) :: SBasis

generates the SBasis that respects the off-diagonal $ℤ_p$ quantum numbers (secfQNOffd)

Arguments

  • cfs :: SConfs is the diagonal QN–preserving configurations ;
  • secf :: Vector{ComplexF64} is a vector of length the same as the number of discrete symmetries that records the eigenvalue of each transformation in the sector ;
  • qnf :: Vector{SQNOffd} is a vector of off-diagonal quantum numbers.
  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.

Output

  • bs :: SBasis is the resulting SBasis object
source
FuzzifiED.Fuzzifino.SBasisMethod
SBasis(cfs :: SConfs) :: SBasis

Generate a SBasis from the configurations without off-diagonal $ℤ_n$ symmetries.

Arguments

  • cfs :: SConfs is the diagonal QN–preserving configurations ;

Output

  • bs :: SBasis is the resulting SBasis object
source

Term

FuzzifiED.Fuzzifino.STermType
STerm

The mutable type STerm records a STerm that looks like $Ua^{(p_1)}_{o_1}a^{(p_2)}_{o_2}… a^{(p_l)}_{o_l}$ in an operator, where positive $o$ denotes fermions and negative $o$ denotes bosons

\[ a^{(0)}_o=f_o, a^{(1)}_o=f_o^†, a^{(0)}_{-o}=b_o, a^{(1)}_{-o}=b_o^†\]

Fields

  • coeff :: ComplexF64 records the coefficient $U$
  • cstr :: Vector{Int64} is a length-$2l$ vector $(p_1,o_1,p_2,o_2,… p_l,o_l)$ recording the operator string

Method

It can be generated by the function

STerm(coeff :: ComplexF64, cstr :: Vector{Int64}) :: ComplexF64
source
FuzzifiED.Fuzzifino.STermsType

STerms is an alias for Vector{STerm} for convenience

Initialisation

STerms(coeff :: Number, cstr :: Vector{Int64})

Gives a STerms with a single STerm.

source

The product of terms with a number, the sum and product of terms, adjoint and particle-hole transformation are defined

Base.:*Method
*(fac :: Number, tms :: STerms) :: STerms
-(tms :: STerms) :: STerms
*(tms :: STerms, fac :: Number) :: STerms
/(tms :: STerms, fac :: Number) :: STerms

Return the product of a collection of STerms with a number.

source
Base.:+Method
+(tms1 :: STerms, tms2 :: STerms) :: STerms
-(tms1 :: STerms, tms2 :: STerms) :: STerms

Return the naive sum of two series of STerms by taking their union.

source
Base.:*Method
*(tms1 :: STerms, tms2 :: STerms) :: STerms
^(tms :: STerms, pow :: Int64) :: STerms

Return the naive product of two series of STerms or the power of one STerms. The number of STerms equals the product of the number of STerms in tms1 and tms2. For each STerm in tms1 $Ua^{(p_1)}_{o_1}…$ and tms2 $U'a^{(p'_1)}_{o'_1}…$, a new STerm is formed by taking $UU'a^{(p_1)}_{o_1}… a^{(p'_1)}_{o'_1}…$

source
Base.adjointMethod
adjoint(tm :: STerm) :: STerm
adjoint(tms :: STerms) :: STerms

Return the Hermitian conjugate of a series of STerms. For each STerm $Ua^{(p_1)}_{o_1}a^{(p_2)}_{o_2}… a^{(p_l)}_{o_l}$, the adjoint is $\bar{U}a^{(1-p_l)}_{o_l}… a^{(1-p_2)}_{o_2}a^{(1-p_1)}_{o_1}$

source

The terms can be simplified by

FuzzifiED.NormalOrderMethod
NormalOrder(tm :: STerm) :: STerms

rearrange a STerm such that

  • the creation operators must be commuted in front of the annihilation operator
  • the site index of the creation operators are in ascending order and the annihilation operators in descending order.

return a list of STerms whose result is equal to the original STerm.

source
FuzzifiED.SimplifyTermsMethod
SimplifyTerms(tms :: STerms ; cutoff :: Float64 = eps(Float64)) :: STerms

simplifies the sum of STerms such that

  • each STerm is normal ordered,
  • like STerms are combined, and STerms with zero coefficients are removed.

Argument

  • cutoff :: Float64 is the cutoff such that STerms with smaller absolute value of coefficients will be neglected. Facultative, eps(Float64) by default.
source

Operator

FuzzifiED.Fuzzifino.SOperatorType
SOperator

The mutable type SOperator records the sum of terms together with information about its symmetry and the basis of the state it acts on and the basis of the resulting state.

Fields

  • bsd :: SBasis is the basis of the initial state ;
  • bsf :: SBasis is the basis of the final state ;
  • red_q :: Int64 is a flag that records whether or not the conversion to a sparse martrix can be simplified : if bsd and bsf have exactly the same set of quantum numbers, and the operator fully respects the symmetries, and all the elements in bsd.cffac and bsf.cffac has the same absolute value, then red_q = 1 ; otherwise red_q = 0 ;
  • sym_q :: Int64 records the symmetry of the operator : if the matrix is Hermitian, then sym_q = 1 ; if it is symmetric, then sym_q = 2 ; otherwise sym_q = 0 ;
  • ntm :: Int64 is the number of terms ;
  • nc :: Int64 is the maximum number of operators in an operator string
  • cstrs :: Matrix{Int64} is a matrix recording the operator string of each term. Each column corresponds to a term and is padded to the maximum length with -1's.
  • coeffs :: Vector{ComplexF64} corresponds to the coefficients in each term.
source

It can be generated by the following methods.

FuzzifiED.Fuzzifino.SOperatorMethod
SOperator(bsd :: SBasis[, bsf :: SBasis], terms :: STerms ; red_q :: Int64, sym_q :: Int64, num_th :: Int64, disp_std :: Bool) :: SOperator

generates an operator object from a series of terms.

Arguments

  • bsd :: SBasis is the basis of the initial state ;
  • bsf :: SBasis is the basis of the final state. Facultative, the same as bsd by default.
  • terms :: STerms records the terms ;
  • red_q :: Int64 is a flag that records whether or not the conversion to a sparse martrix can be simplified : if bsd and bsf have exactly the same set of quantum numbers, and the operator fully respects the symmetries, and all the elements in bsd.cffac and bsf.cffac has the same absolute value, then red_q = 1 ; otherwise red_q = 0 ; Facultative, if bsf is not given, 1 by default, otherwise 0 by default.
  • sym_q :: Int64 records the symmetry of the operator : if the matrix is Hermitian, then sym_q = 1 ; if it is symmetric, then sym_q = 2 ; otherwise sym_q = 0. Facultative, if bsf is not given, 1 by default, otherwise 0 by default.
  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.
source

The product of an operator on a state and the inner product of a final state, an operator and an initial state can be calculated by

Base.:*Method
*(op :: SOperator, st_d :: Vector{ComplexF64} ; num_th :: Int64, disp_std :: Bool) :: Vector{ComplexF64}
*(op :: SOperator, st_d :: Vector{Float64} ; num_th :: Int64, disp_std :: Bool) :: Vector{Float64}

Measure the action of an operator on a state. st_d must be of length op.bsd.dim. Returns a vector of length op.bsf.dim that represents the final state.

Facultative arguments

  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.
source
Base.:*Method
*(st_fp :: LinearAlgebra.Adjoint{ComplexF64, Vector{ComplexF64}}, op :: SOperator, st_d :: Vector{ComplexF64} ; num_th :: Int64, disp_std :: Bool) :: ComplexF64
*(st_fp :: LinearAlgebra.Adjoint{Float64, Vector{Float64}}, op :: SOperator, st_d :: Vector{Float64} ; num_th :: Int64, disp_std :: Bool) :: Float64

Measuring the inner product between two states and an operator. st_d must be of length op.bsd.dim and st_fp must be of length op.bsf.dim, and st_fp must be an adjoint.

Facultative arguments

  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.
source

Sparse matrix

The OpMat can be generated from SOperator by the following methods.

FuzzifiED.OpMatMethod
OpMat[{type}](op :: SOperator ; num_th :: Int64, disp_std :: Bool) :: OpMat{type}

Generates the sparse matrix from the operator. The parameter type is either Float64 or ComplexF64 ; it is facultative, given by ElementType by default.

Arguments

  • op :: SOperator is the operator ;
  • type :: DataType specifies the type of the matrix. It can either be ComplexF64 or Float64. Facultative, the same as ElementType by default
  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.
source

After the generation of sparse matrix, the diagonalisation can be condicted with FuzzifiED.

  • test_boson.jl tests the nearest-neighbour tight-binding model $H=\sum_i(b^\dagger_ib_{i+1}+f^\dagger_if_{i+1}+\mathrm{h.c.})$. The example diagonalises the sector with the number of bosons and fermions both $N_o/2$, and even under the reflection with respect to a bond center $i\mapsto N_o+1-i$, and measures the total particle number squared $\left[\sum_i(b_i^\dagger b_i+f^\dagger_if_i)\right]^2$.
  • ising_frac_boson.jl calculates the spectrum of 3d Ising model on fuzzy sphere for bosons at fractional filling $ν = 1/2$. This example reproduces Figure 12a,b in Voinea 2024.