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.FuzzifinoCore Functions
Environment Parameter
FuzzifiED.Fuzzifino.Libpathino — ConstantFuzzifino.Libpathino :: String = FuzzifiED_jll.LibpathFuzzifinodefine 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.
Quantum Numbers
The diagonal and off-diagonal quantum numbers are implemented as
FuzzifiED.Fuzzifino.SQNDiag — TypeSQNDiagThe 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 :: Stringis the name of the diagonal quantum numberchargef :: Vector{Int64}is the symmetry charge $q_{f,o}$ of each sitechargeb :: Vector{Int64}is the symmetry charge $q_{b,o}$ of each sitemodul :: 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]) :: SQNDiagThe arguments name and modul are facultative. By default name is set to "QN" and modul is set to 1.
FuzzifiED.Fuzzifino.SQNOffd — TypeSQNOffdThe 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 transformationfacf :: 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 :: Int64is 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}]) :: SQNOffdThe 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.
The SQNs can be converted from the fermionic QNs by setting the bosonic part to identity.
FuzzifiED.Fuzzifino.SQNDiag — MethodSQNDiag(qnd :: QNDiag, nob :: Int64) :: SQNDiagconverts a pure fermionic QNDiag to a SQNDiag with the boson charges set to empty, implemented as
SQNDiag(qnd.name, qnd.charge, fill(0, nob), qnd.modul)FuzzifiED.Fuzzifino.SQNOffd — MethodSQNOffd(qnf :: QNOffd, nob :: Int64) :: SQNOffdconverts a pure fermionic QNOffd to a SQNOffd with the boson transformations set to identity, implemented as
SQNOffd(qnf.perm, collect(1 : nob), qnf.ph, qnf.fac, fill(ComplexF64(1), nob), qnf.cyc)Several operations of the SQNs are supported.
Base.:* — Method*(fac :: Int64, qnd :: SQNDiag) :: SQNDiag
*(qnd :: SQNDiag, fac :: Int64) :: SQNDiag
÷(qnd :: SQNDiag, fac :: Int64) :: SQNDiag
-(qnd :: SQNDiag) :: SQNDiagreturns the SQNDiag multiplied or divided by an integer factor, where the charge is multiplied or integer-divided by the factor. For $ℤ_p$ quantum numbers, their modulus will be multiplied or integer-divided by the absolute value. If qnd.modul ÷ abs(fac) ≤ 1, a trivial SQNDiag will be returned.
Base.:+ — Method+(qnd1 :: SQNDiag, qnd2 :: SQNDiag) :: SQNDiag
-(qnd1 :: SQNDiag, qnd2 :: SQNDiag) :: SQNDiagreturns the sum or substraction of two SQNDiags, whose name is the samea as qnd1, charge is the same as qnd1 ± qnd2, and modulus is the GCD of qnd1 and qnd2. If qnd1 and qnd2 are both $ℤ_p$ quantum numbers and their modulus are coprime, a trivial SQNDiag will be returned.
Base.:* — Method*(qnf1 :: SQNOffd, qnf2 :: SQNOffd) :: SQNOffdreturns the composition of two SQNOffd transformations. The cycle is set to be the LCM of two QNOffds.
FuzzifiED.Fuzzifino.PadSQNDiag — FunctionPadSQNDiag(qnd :: SQNDiag, nofl :: Int64, nobl :: Int64, nofr :: Int64, nobr :: Int64)adds nol empty orbitals to the left and nor empty orbitals to the right, implemented as
SQNDiag(qnd.name,
[ fill(0, nofl) ; qnd.chargef ; fill(0, nofr) ],
[ fill(0, nobl) ; qnd.chargeb ; fill(0, nobr) ],
qnd.modul)FuzzifiED.Fuzzifino.PadSQNOffd — FunctionPadSQNOffd(qnf :: SQNOffd, nofl :: Int64, nobl :: Int64, nofr :: Int64, nobr :: Int64)adds nofl fermionic and nobl bosonic empty orbitals to the left and nofr fermionic and nobr bosonic empty orbitals to the right, implemented as
SQNOffd(
[ collect(1 : nofl) ; qnf.permf .+ nofl ; collect(1 : nofr) .+ (nofl + length(qnf.permf)) ],
[ collect(1 : nobl) ; qnf.permb .+ nobl ; collect(1 : nobr) .+ (nobl + length(qnf.permb)) ],
[ fill(0, nofl) ; qnf.phf ; fill(0, nofr) ],
[ fill(ComplexF64(1), nofl) ; qnf.facf ; fill(ComplexF64(1), nofr) ],
[ fill(ComplexF64(1), nobl) ; qnf.facb ; fill(ComplexF64(1), nobr) ],
qnf.cyc)Configurations
FuzzifiED.Fuzzifino.SConfs — TypeSConfsThe 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 :: Int64is the number of fermionic sites.nof :: Int64is the number of bosonic sites.nebm :: Int64is the maximal number of boson occupation.ncf :: Int64is the number of configurations.conff :: Vector{Int64}is an array of lengthncfcontaining all the fermion configurations. Each configuration is expressed in a binary number. If theo-1-th bit ofconf[i]is 1, then theo-th site in thei-th configuration is occupied ; if the bit is 0, then the site is empty.confb :: Vector{Int64}is an array of lengthncfcontaining 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}andrid :: Vector{Int64}contain the information of Lin table that is used to inversely look up the indexifrom the configuration.
It can be generated from the QNDiags.
FuzzifiED.Fuzzifino.SConfs — MethodSConfs(nof :: Int64, nob :: Int64, nebm :: Int64, secd :: Vector{Int64}, qnd :: Vector{SQNDiag} ; , num_th :: Int64, disp_std :: Bool) :: Confsgenerates the configurations from the list of QNDiags.
Arguments
nof :: Int64is the number of fermionic sites $N_{of}$.nob :: Int64is the number of bosonic sites $N_{ob}$.nebm :: Int64is 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 :: Int64andnorb :: Int64are 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,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
Output
cfs :: SConfsis a SConfs object.
Note
If your qnd has negative entries, QNDiags must contain the total number of particles (i. e., bosons plus fermions).
Basis
FuzzifiED.Fuzzifino.SBasis — TypeSBasisThe 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 :: SConfsstores the configurations that respect the QNDiags.dim :: Int64is the dimension of the SBasis.szz :: Int64records the maximum size $\max m_g$ of groups.cfgr :: Vector{Int64}is a vector of lengthcfs.ncfand records which group $|I⟩$ each configuration $|i⟩$ belong to.cffac :: Vector{ComplexF64}is a vector of lengthcfs.ncfand records the coefficients $λ_i$ of each configuration.grel :: Matrix{Int64}is aszz×dimmatrix that records the configurations in each group $|i_{Ik}⟩ (k = 1,…,m_I)$grsz :: Vector{Int64}is a vector of lengthdimthat records the size $m_I$ of each group.
It can be generated by the following methods.
FuzzifiED.Fuzzifino.SBasis — MethodSBasis(cfs :: SConfs, secf :: Vector{ComplexF64}, qnf :: Vector{SQNOffd} ; num_th :: Int64, disp_std :: Bool)generates the SBasis that respects the off-diagonal $ℤ_p$ quantum numbers (secfQNOffd)
Arguments
cfs :: SConfsis 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,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
Output
bs :: SBasisis the resulting SBasis object.
FuzzifiED.Fuzzifino.SBasis — MethodSBasis(cfs :: SConfs)Generate a SBasis from the configurations without off-diagonal $ℤ_n$ symmetries.
Arguments
cfs :: SConfsis the diagonal QN–preserving configurations.
Output
bs :: SBasisis the resultingSBasisobject.
Term
FuzzifiED.Fuzzifino.STerm — TypeSTermThe 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 :: ComplexF64records 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 or directly from a term
STerm(coeff :: ComplexF64, cstr :: Vector{Int64})
STerm(tm :: Term)FuzzifiED.Fuzzifino.STerms — TypeSTerms is an alias for Vector{STerm} for convenience
Initialisation
STerms(coeff :: Number, cstr :: Vector{Int64})Gives a STerms with a single STerm.
STerms(tms :: Terms)Generates STerms from Terms
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) :: STermsReturn the product of a collection of STerms with a number.
Base.:+ — Method+(tms1 :: STerms, tms2 :: STerms) :: STerms
-(tms1 :: STerms, tms2 :: STerms) :: STermsReturn the naive sum of two series of STerms by taking their union.
Base.:* — Method*(tms1 :: STerms, tms2 :: STerms) :: STerms
^(tms :: STerms, pow :: Int64) :: STermsReturn 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}…$
Base.adjoint — Methodadjoint(tm :: STerm) :: STerm
adjoint(tms :: STerms) :: STermsReturn 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}$
The terms can be simplified and manipulated by
FuzzifiED.NormalOrder — MethodNormalOrder(tm :: STerm) :: STermsrearrange 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.
FuzzifiED.SimplifyTerms — MethodSimplifyTerms(tms :: STerms ; cutoff :: Float64 = eps(Float64)) :: STermssimplifies the sum of STerms such that
- each STerm is normal ordered,
- like STerms are combined, and STerms with zero coefficients are removed.
Argument
cutoff :: Float64is the cutoff such that STerms with smaller absolute value of coefficients will be neglected. Facultative,eps(Float64)by default.
FuzzifiED.RemoveOrbs — MethodRemoveOrbs(tms :: STerms, o_rm :: Vector{Int64}) :: STermsremove all the terms that has creation or annihilation operators of sites in the set o_rm.
FuzzifiED.RelabelOrbs — MethodRelabelOrbs(tms :: STerms, dict_o :: Dict{Int64, Int64})Relabel all the orbitals. Terms that has operators on the sites that are not in the keys of dict_o are removed.
FuzzifiED.Fuzzifino.PadSTerms — FunctionPadSTerms(tms :: STerms, nol :: Int64)adds nofl fermionic and nobl bosonic empty orbitals to the left by shifting each orbital index.
Operator
FuzzifiED.Fuzzifino.SOperator — TypeSOperatorThe 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 :: SBasisis the basis of the initial state.bsf :: SBasisis the basis of the final state.red_q :: Int64is a flag that records whether or not the conversion to a sparse martrix can be simplified : ifbsdandbsfhave exactly the same set of quantum numbers, and the operator fully respects the symmetries, and all the elements inbsd.cffacandbsf.cffachas the same absolute value, thenred_q = 1; otherwisered_q = 0.sym_q :: Int64records the symmetry of the operator : if the matrix is Hermitian, thensym_q = 1; if it is symmetric, thensym_q = 2; otherwisesym_q = 0.ntm :: Int64is the number of terms.nc :: Int64is the maximum number of operators in an operator stringcstrs :: 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.
It can be generated by the following methods.
FuzzifiED.Fuzzifino.SOperator — MethodSOperator(bsd :: SBasis[, bsf :: SBasis], terms :: STerms ; red_q :: Int64, sym_q :: Int64, num_th :: Int64, disp_std :: Bool) :: SOperatorgenerates an operator object from a series of terms.
Arguments
bsd :: SBasisis the basis of the initial state.bsf :: SBasisis the basis of the final state. Facultative, the same asbsdby default.terms :: STermsrecords the terms.red_q :: Int64is a flag that records whether or not the conversion to a sparse martrix can be simplified : ifbsdandbsfhave exactly the same set of quantum numbers, and the operator fully respects the symmetries, and all the elements inbsd.cffacandbsf.cffachas the same absolute value, thenred_q = 1; otherwisered_q = 0; Facultative, ifbsfis not given, 1 by default, otherwise 0 by default.sym_q :: Int64records the symmetry of the operator : if the matrix is Hermitian, thensym_q = 1; if it is symmetric, thensym_q = 2; otherwisesym_q = 0. Facultative, ifbsfis not given, 1 by default, otherwise 0 by default.num_th :: Int64, the number of threads. Facultative,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
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,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
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) :: Float64Measuring 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,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
Sparse Matrix
The OpMat can be generated from SOperator by the following methods.
FuzzifiED.OpMat — MethodOpMat[{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 :: SOperatoris the operator.type :: DataTypespecifies the type of the matrix. It can either beComplexF64orFloat64. Facultative, the same asElementTypeby defaultnum_th :: Int64, the number of threads. Facultative,NumThreadsby default.disp_std :: Bool, whether or not the log shall be displayed. Facultative,!SilentStdby default.
After the generation of sparse matrix, the diagonalisation can be condicted with FuzzifiED.
Entanglement
FuzzifiED.StateDecompMat — MethodStateDecompMat(st :: Vector{<:Number}, bs0 :: SBasis, bsa :: SBasis, bsb :: SBasis, amp_ofa :: Vector{<:Number}, amp_oba :: Vector{<:Number}, amp_ofb :: Vector{<:Number}, amp_obb :: Vector{<:Number}) :: Matrix{ComplexF64}Decompose a state $|ψ⟩=v_I|I⟩$ into a direct-product basis of two subsystems $|ψ⟩=M_{JI}|I_A⟩|J_B⟩$
Arguments
st :: Vector{<:Number}is the state to be decomposed into direct-product basis of two subsystems.bs0 :: SBasisis the total basis.bsa :: SBasisis the basis for the subsystem A.bsb :: SBasisis the basis for the subsystem B.amp_ofa :: Vector{ComplexF64}is a complex list of lengthnothat specifies the fermionic amplitute of each orbital in the subsystem A. For a non-local basis, we decompose each electron into creation operators in two subsystems $c^†_o=a_{o,A}c^†_{o,A}+a_{o,B}c^†_{o,B}$ and this list specifies $a_{o,A}$. This is equivalent to $√{ℱ_{m,A}}$ in PRB 85, 125308 (2012) with an extra phase factor.amp_oba :: Vector{ComplexF64}is a complex list of lengthnothat specifies the bosonic amplitute of each orbital in the subsystem A.amp_ofb :: Vector{ComplexF64}is a complex list of lengthnothat specifies the fermionic amplitute of each orbital in the subsystem B.amp_obb :: Vector{ComplexF64}is a complex list of lengthnothat specifies the bosonic amplitute of each orbital in the subsystem B.
Output
A complex matrix of dimension bsb.dim * bsa.dim that corresponds to the state in the decomposed basis $|ψ⟩=M_{JI}|I_A⟩|J_B⟩$. This is equivalent to $R_{μν}^A/√p$ in PRB 85, 125308 (2012). After calculating all the sectors, the reduced density matrix will be $ρ_B=⊕𝐌𝐌^†$.
FuzzifiED.GetEntSpec — MethodGetEntSpec(st :: Vector{<:Number}, bs0 :: SBasis, secd_lst :: Vector{Vector{Vector{Int64}}}, secf_lst :: Vector{Vector{Vector{<:Number}}} ; qnd_a :: Vector{SQNDiag}, qnd_b :: Vector{SQNDiag} = qnd_a, qnf_a :: Vector{SQNOffd}, qnf_b :: Vector{SQNOffd} = qnf_a, amp_oa :: Vector{<:Number}, amp_ob :: Vector{<:Number} = sqrt.(1 .- abs.(amp_oa .^ 2))) :: Dict{@NamedTuple{secd_a, secf_a, secd_b, secf_b}, Vector{Float64}}Arguments
st :: Vector{<:Number}is the state to be decomposed into direct-product basis of two subsystems.bs0 :: SBasisis the total basis.secd_lst :: Vector{Vector{Vector{Int64}}}gives the list of QNDiag sectors of subsystems to be calculated. Each of its elements is a two element vector ; the first specifies the sector for subsystem A, and the second specifies the sector for subsystem B.secf_lst :: Vector{Vector{Vector{ComplexF64}}}gives the list of QNOffd sectors of subsystems to be calculated. Each of its elements is a two element vector ; the first specifies the sector for subsystem A, and the second specifies the sector for subsystem B.qnd_a :: Vector{SQNDiag}, qnd_b :: Vector{SQNDiag} = qnd_a, qnf_a :: Vector{QNOffd}, qnf_b :: Vector{QNOffd}specifies the diagonal and off-diagonal quantum numbers of the subsystems A and B.qnd_bandqnf_bare facultative and the same asqnd_aandqnf_aby default.amp_oa :: Vector{ComplexF64}andamp_ob :: Vector{ComplexF64}are complex lists of lengthnothat specify the amplitute of each orbital in the subsystems A and B. For a non-local basis, we decompose each electron into creation operators in two subsystems $c^†_o=a_{o,A}c^†_{o,A}+a_{o,B}c^†_{o,B}$ and this list specifies $a_{o,A}$. This is equivalent to $√{ℱ_{m,A}}$ in PRB 85, 125308 (2012) with an extra phase factor.
Output
A dictionary whose keys are named tuples that specify the sector containing entries secd_a, secf_a, secd_b, secf_b and values are lists of eigenvalues of the density matrix in those sectors.
Transformation
FuzzifiED.Fuzzifino.STransf — TypeSTransfThe mutable type STransf records a transformation in the same form as a SQNOffd
\[𝒵:\ c_o↦ α_o^* c^{(p_o)}_{π_o}, c_o^†↦ α_o c^{(1-p_o)}_{π_o}\]
together with information about its symmetry and the basis of the state it acts on and the basis of the resulting state.
Fields
bsd :: SBasisis the basis of the initial state.bsf :: SBasisis the basis of the final state.permf :: Vector{Int64}, permb :: Vector{Int64},phf :: Vector{Int64}andfacf :: Vector{ComplexF64},facb :: Vector{ComplexF64}` records the transformation in the same form as a SQNOffd.
FuzzifiED.Fuzzifino.STransf — MethodSTransf(bsd :: SBasis, bsf :: SBasis, qnf :: SQNOffd)generates a transformation object from a SQNOffd.
Arguments
bsd :: SBasisis the basis of the initial state.bsf :: SBasisis the basis of the final state. Facultative, the same asbsdby default.qnf :: SQNOffdrecords the transformation ;
Base.:* — Method*(trs :: STransf, st_d :: Vector{ComplexF64} ; num_th = NumThreads) :: Vector{ComplexF64}
*(trs :: STransf, st_d :: Vector{Float64} ; num_th = NumThreads) :: Vector{Float64}Act a transformation on a state. st_d must be of length trs.bsd.dim. Returns a vector of length trs.bsf.dim that represents the final state.
Facultative arguments
num_th :: Int64, the number of threads. Facultative,NumThreadsby default.
Models
Several quantum numbers and operator terms for the pure bosonic models are built-in.
FuzzifiED.Fuzzifino.GetNeSQNDiag — FunctionGetNeSQNDiag(nof :: Int64, nob :: Int64) :: SQNDiagReturn the SQNDiag of the number of particles, implemented as
SQNDiag("Ne", fill(1, nof), fill(1, nob))FuzzifiED.Fuzzifino.GetBosonLz2SQNDiag — FunctionGetBosonLz2SQNDiag(nof :: Int64, nm :: Int64, nf :: Int64) :: SQNDiagReturn the SQNDiag of twice the angular momentum $2L_z$ for pure bosons, implemented as
SQNDiag("Lz", fill(0, nof), collect(0 : nm * nf - 1) .÷ nf .* 2 .- (nm - 1))FuzzifiED.Fuzzifino.GetBosonFlavSQNDiag — FunctionGetBosonFlavSQNDiag(nof :: Int64, nm :: Int64, nf :: Int64, qf :: Dict{Int64, Int64}[, id :: Int64 = 1, modul :: Int64 = 1]) :: QNDiag
GetBosonFlavSQNDiag(nof :: Int64, nm :: Int64, nf :: Int64, qf :: Vector{Int64}[, id :: Int64 = 1, modul :: Int64 = 1]) :: QNDiagReturn the QNDiag of linear combination of number of particles in each flavour for pure bosons,
\[ Q = ∑_{f}q_fn_f\]
the factor $q_f$ can either be given by a length-$N_f$ vector or a dictionary containing non-zero terms. E.g., for $Q=n_{f=1}-n_{f=3}$ in a 4-flavour system, qf = [1, 0, -1, 0] or qf = Dict(1 => 1, 3 => -3). id is an index to be put in the name to distinguish. For qf given as vector, the function is implemented as
SQNDiag("Sz$id", fill(0, nof), qf[collect(0 : nm * nf - 1) .% nf .+ 1], modul)FuzzifiED.Fuzzifino.GetBosonFlavPermSQNOffd — FunctionGetBosonFlavPermSQNOffd(nof :: Int64, nm :: Int64, nf :: Int64, permf[, fac][, cyc :: Int64])Return the flavour permutaiton transformation for pure bosons
\[ 𝒵: b^†_{mf}↦α_fb^†_{mπ_f}\]
Arguments
nm :: Int64andnf :: Int64are the number of orbitals and the flavours.permf :: Dict{Int64, Int64},permf :: Vector{Vector{Int64}}orVector{Int64}gives the flavour permutation $π_f$. It is either a vector of the cycles, a vector of the target flavours, or a dictionary of the changed elements. E.g., a permutation $1↦4,2↦5,3↦3,4↦1,5↦6,6↦2$ can be expressed as[4,5,3,1,6,2],[[1,4],[2,5,6]]orDict(1=>4,2=>5,4=>1,5=>6,6=>2). Facultative, identity by default.fac :: Dict{Int64, <: Number}orVector{<: Number}gives the factor $α_f$. It is either a vector of all vectors, or a dictionary of all non-unity elements. Facultative, all unity by default.cyc :: Int64is the period of the permutation.
FuzzifiED.Fuzzifino.GetBosonRotySQNOffd — FunctionGetBosonRotySQNOffd(nof :: Int64, nm :: Int64, nf :: Int64)Return the $π$-rotation with respect to the $y$-axis for pure bosons.
\[ ℛ_y: b^†_{mf}↦(-)^{m+s}b^†_{-mf}\]
FuzzifiED.Fuzzifino.GetBosonDenIntSTerms — FunctionGetBosonDenIntSTerms(nm :: Int64, nf :: Int64[, ps_pot :: Vector{<:Number}][, mat_a :: Matrix{<:Number}[, mat_b :: Matrix{<:Number}]]) :: TermsReturn the normal-ordered density-density term in the Hamiltonian for bosons
\[∑_{\{m_i,f_i\}}U_{m_1m_2m_3m_4}M^A_{f_1f_4}M^B_{f_2f_3}b^{†}_{m_1f_1}b^{†}_{m_2f_2}b_{m_3f_3}b_{m_4f_4}.\]
Arguments
nm :: Int64is the number of orbitals.nf :: Int64is the number of flavours.ps_pot :: Vector{<:Number}is a list of numbers specifying the pseudopotentials for the interacting matrix $U_{m_1m_2m_3m_4}$. Facultative,[1.0]by default.mat_a :: Matrix{<:Number}is anf×nfmatrix specifying $M^A_{ff'}$. Facultative, $I_{N_f}$ by default.mat_b :: Matrix{<:Number}is anf×nfmatrix specifying $M^B_{ff'}$. Facultative, the Hermitian conjugate ofmat_aby default.
GetBosonDenIntSTerms(nm :: Int64, nf :: Int64, ps_pot :: Vector{<:Number}, mat_a :: Vector{<:AbstractMatrix{<:Number}}[, mat_b :: Vector{<:AbstractMatrix{<:Number}}]) :: TermsReturn the sum of a series of normal-ordered density-density term in the Hamiltonian for bosons
\[∑_{\{m_i,f_i,α\}}U_{m_1m_2m_3m_4}(M^A_{α})_{f_1f_4}(M^B_{α})_{f_2f_3}b^{†}_{m_1f_1}b^{†}_{m_2f_2}b_{m_3f_3}b_{m_4f_4}.\]
Arguments
nm :: Int64is the number of orbitals.nf :: Int64is the number of flavours.ps_pot :: Vector{<:Number}is a list of numbers specifying the pseudopotentials for the interacting matrix $U_{m_1m_2m_3m_4}$. Facultative,[1.0]by default.mat_a :: Vector{<:AbstractMatrix{<:Number}}is a vector ofnf×nfmatrix specifying $(M^A_{α})_{ff'}$. Facultative, $I_{N_f}$ by default.mat_b :: Vector{<:AbstractMatrix{<:Number}}is a vector ofnf×nfmatrix specifying $(M^B_{α})_{ff'}$. Facultative, the Hermitian conjugate ofmat_aby default.
FuzzifiED.Fuzzifino.GetBosonPairIntSTerms — FunctionGetBosonPairIntSTerms(nm :: Int64, nf :: Int64, ps_pot :: Vector{<:Number}, mat_a :: Matrix{<:Number}[, mat_b :: Matrix{<:Number}]) :: TermsReturn the normal-ordered pair-pair interaction term in the Hamiltonian for bosons
\[∑_{\{m_i,f_i\}}U_{m_1m_2m_3m_4}M^A_{f_1f_2}M^B_{f_3f_4}b^{†}_{m_1f_1}b^{†}_{m_2f_2}b_{m_3f_3}b_{m_4f_4}.\]
Arguments
nm :: Int64is the number of orbitals.nf :: Int64is the number of flavours.ps_pot :: Vector{<:Number}is a list of numbers specifying the pseudopotentials for the interacting matrix $U_{m_1m_2m_3m_4}$.mat_a :: Matrix{<:Number}is anf×nfmatrix specifying $M^A_{ff'}$. Facultative, $I_{N_f}$ by default.mat_b :: Matrix{<:Number}is anf×nfmatrix specifying $M^B_{ff'}$. Facultative, the Hermitian conjugate ofmat_aby default.
FuzzifiED.Fuzzifino.GetBosonPolSTerms — FunctionGetBosonPolSTerms(nm :: Int64, nf :: Int64[, mat :: Matrix{<:Number}][ ; fld_m :: Vector{<:Number}]) :: STermsReturn the polarisation term in the Hamiltonian for bosons
\[∑_{mff'}c^{†}_{mf}M_{ff'}c_{mf'}.\]
Arguments
nm :: Int64is the number of orbitals.nf :: Int64is the number of flavours.mat :: Matrix{<:Number}is anf×nfmatrix specifying $M_{ff'}$. Facultative, $I_{N_f}$ by default.fld_m :: Vector{<:Number}gives an orbital dependent polarisation
\[∑_{mff'}h_mc^{†}_{mf}M_{ff'}c_{mf'}\]
Facultative.
FuzzifiED.Fuzzifino.GetBosonLpLzSTerms — FunctionGetBosonLpLzSTerms(nm :: Int64, nf :: Int64) :: Tuple{Terms, Terms}Return the a pair of terms containing $L^z$ and $L^+$ for bosons.
Arguments
nm :: Int64is the number of bosonic orbitals.nf :: Int64is the number of bosonic flavours.
FuzzifiED.Fuzzifino.GetL2STerms — MethodGetL2Terms(tms_lzlp :: Tuple{STerms, STerms}) :: STermsReturn the terms for the total angular momentum from a tuple of terms containing $L^z$ and $L^+$, implemented as
tms_lz, tms_lp = tms_lzlp
return SimplifyTerms(tms_lz * tms_lz - tms_lz + tms_lp * tms_lp')FuzzifiED.Fuzzifino.GetBosonC2STerms — FunctionGetBosonC2STerms(nm :: Int64, nf :: Int64, mat_gen :: Vector{Matrix{<:Number}}[, mat_tr :: Vector{Matrix{<:Number}}]) :: STermsReturn the terms for the quadratic Casimir of the flavour symmetry for bosons.
\[ C_2=∑_{imm'}\frac{(b^†_{mf_1}G_{i,f_1f_2}b_{mf_2})(b^†_{m'f_3}G^†_{i,f_3f_4}b_{m'f_4})}{2\operatorname{tr}G_i^†G_i}-∑_{imm'}\frac{(b^†_{mf_1}T_{i,f_1f_2}b_{mf_2})(b^†_{m'f_3}T^†_{i,f_3f_4}b_{m'f_4})}{2\operatorname{tr}T_i^†T_i}\]
where $G_i$ are the generator matrices, and $T_i$ are the trace matrices.
Arguments
nm :: Int64is the number of orbitals.nf :: Int64is the number of flavours.mat_gen :: Vector{Matrix{Number}}is a list of the matrices that gives the generators. It will automatically be normalised such that its square traces to $1/2$.mat_tr :: Vector{Matrix{Number}}is a list of trace matrices that will be normalised automatically and substracted. Facultative.
The types SSphereObs and SAngModes are similarly defined as their fermionic counterparts
FuzzifiED.Fuzzifino.SSphereObs — TypeSSphereObsThe mutable type SSphereObs stores the information of a local observable (or local operator) $𝒪$ that can be decomposed into angular components.
\[ 𝒪(\Omega)=∑_{lm}𝒪_{lm}Y^{(s)}_{lm}\]
Fields
s2 :: Int64is twice the spin $2s$ of the observable.l2m :: Int64is twice the maximal angular momentum $2l_{\max}$ of the components of the observable.get_comp :: Functionis a functionget_comp(l2 :: Int64, m2 :: Int64) :: STermsthat sends the component specified by a tuple of integers $(2l,2m)$ where $|s|\leq l\leq l_{\max}, -l\leq m\leq l$ to a list of terms that specifies the expression of the component.stored_q :: Boolis a boolean that specifies whether or not each component of the observable is stored.comps :: Dict{Tuple{Int64, Int64}, STerms}stores each component of the observable in the format of a dictionary whose keys are the tuples of integers $(2l,2m)$ and values are the lists of terms that specifies the expression of the component.
Methods
The methods for this type is similarly defined as in SphereObs.
SSphereObs(s2 :: Int64, l2m :: Int64, get_comp :: Function) :: SSphereObs
SSphereObs(s2 :: Int64, l2m :: Int64, comps :: Dict{Tuple{Int64, Int64}, STerms}) :: SSphereObs
StoreComps(obs :: SSphereObs) :: SSphereObs
*(fac :: Number, obs :: SSphereObs) :: SSphereObs
+(obs1 :: SSphereObs, obs2 :: SSphereObs) :: SSphereObs
adjoint(obs :: SSphereObs) :: SSphereObs
*(obs1 :: SSphereObs, obs2 :: SSphereObs) :: SSphereObs
PadSSphereObs(obs :: SSphereObs, nofl :: Int64, nobl :: Int64) :: SSphereObs
Laplacian(obs :: SSphereObs[ ; norm_r2 :: Float64]) :: SSphereObs
GetIntegral(obs :: SSphereObs[ ; norm_r2 :: Float64]) :: STerms
GetComponent(obs :: SSphereObs, l :: Number, m :: Number) :: STerms
FilterComponent(obs :: SSphereObs, flt) :: AngModes
GetPointValue(obs :: SSphereObs, θ :: Float64, ϕ :: Float64) :: STerms
GetFermionSObs(nm :: Int64, nf :: Int64, f :: Int64[ ; norm_r2 :: Float64]) :: SSphereObs
GetBosonSObs(nm :: Int64, nf :: Int64, f :: Int64[ ; norm_r2 :: Float64]) :: SSphereObs
GetFerDensitySObs(nm :: Int64, nf :: Int64[, mat :: Matrix{<:Number}][ ; norm_r2 :: Float64]) :: SSphereObs
GetBosDensitySObs(nm :: Int64, nf :: Int64[, mat :: Matrix{<:Number}][ ; norm_r2 :: Float64]) :: SSphereObsFuzzifiED.Fuzzifino.SAngModes — TypeSAngModesThe mutable type SAngModes stores angular momentum components of an operator on the sphere $Φ_{lm}$ and superposes in the rule of Clebsch-Gordan coefficients. The usage is similar to the spherical observables, except that SphereObs superposes in the rule of spherical harmonics and has the notion of locality
Fields
l2m :: Int64is twice the maximal angular momentum $2l_{\max}$ of the components of the modes object.get_comp :: Functionis a functionget_comp(l2 :: Int64, m2 :: Int64) :: STermsthat sends the component specified by a tuple of integers $(2l,2m)$ where $|s|\leq l\leq l_{\max}, -l\leq m\leq l$ to a list of terms that specifies the expression of the component.stored_q :: Boolis a boolean that specifies whether or not each component of the modes object is stored.comps :: Dict{Tuple{Int64, Int64}, STerms}stores each component of the modes object in the format of a dictionary whose keys are the tuples of integers $(2l,2m)$ and values are the lists of terms that specifies the expression of the component.
Methods
The methods for this type is similarly defined as in AngModes.
SAngModes(l2m :: Int64, get_comp :: Function) :: SAngModes
SAngModes(l2m :: Int64, cmps :: Dict{Tuple{Int64, Int64}, STerms}) :: SAngModes
StoreComps!(amd :: SAngModes)
StoreComps(amd :: SAngModes) :: SAngModes
*(fac :: Number, amd :: SAngModes) :: SAngModes
+(amd1 :: SAngModes, amd2 :: SAngModes) :: SAngModes
adjoint(amd :: SAngModes) :: SAngModes
*(amd1 :: SAngModes, amd2 :: SAngModes) :: SAngModes
PadSAngModes(amd :: SAngModes, nofl :: Int64, nobl :: Int64) :: SAngModes
GetComponent(amd :: SAngModes, l :: Number, m :: Number) :: STerms
ContractMod(amd1 :: SAngModes, amd2 :: amd2, comps :: Dict) :: STerms
ContractMod(amd1 :: SAngModes, amd2 :: SAngModes, l0 :: Number) :: STerms
FilterComponent(amd :: SAngModes, flt) :: SAngModes
FilterL2(amd :: SAngModes, l :: Number) :: SAngModes
GetFermionSMod(nm :: Int64, nf :: Int64, f :: Int64) :: SAngModes
GetBosonSMod(nm :: Int64, nf :: Int64, f :: Int64) :: SAngModes
GetFerPairingSMod(nm :: Int64, nf :: Int64, mat :: Matrix{<:Number}) :: SAngModes
GetBosPairingSMod(nm :: Int64, nf :: Int64, mat :: Matrix{<:Number}) :: SAngModes
GetFerDensitySMod(nm :: Int64, nf :: Int64, mat :: Matrix{<:Number}) :: SAngModes
GetBosDensitySMod(nm :: Int64, nf :: Int64, mat :: Matrix{<:Number}) :: SAngModesRelated Examples
ising_frac_boson.jlcalculates the spectrum of 3d Ising model on the fuzzy sphere for bosons at fractional filling $ν = 1/2$. This example reproduces Figure 12a,b in Voinea 2024.su2_1_scal_spectrum.jlcalculates the spectrum of the $\mathrm{SU}(2)_1$ coupled to a complex scalar Chern-Simons matter CFT on the fuzzy sphere. This example reproduces Figs. 2d and 3 in Zhou 2025Jul.free_majorana_spectrum.jlcalculates the spectrum of free Majorana fermion. This example reproduces Figure 5 in Zhou 2025Sep.free_majorana_correlator.jlcalculates the 2-pt correlator of free Majorana fermion. This example reproduces Figure 7 in Zhou 2025Sep.bpf_220_spectrum.jlcalculates the spectrum of the transition between bosonic Pfaffian and Halperin 220 described by gauged Majorana fermion CFT # This example reproduces Figure 4a in Voinea 2025