Core Functions
Environment Parameters
The following environment parameters defines how FuzzifiED works, viz. whether it outputs logs, how many threads it uses and where it reads its libraries. In general, you can keep it at default.
FuzzifiED.SilentStd — ConstantFuzzifiED.SilentStd :: Bool = falsea flag to determine whether logs of the FuzzifiED functions should be turned off. False by default. If you want to evaluate without log, put FuzzifiED.SilentStd = true. This parameter can be defined for each process separately.
FuzzifiED.NumThreads — ConstantFuzzifiED.NumThreads :: Int = Threads.nthreads()an integer to define how many threads OpenMP uses. By default, it is the same as the number of threads in Julia. If you use Jupyter notebooks, which by default uses one core only, you may need to define this by hand, e. g., FuzzifiED.NumThreads = 8. This parameter can be defined for each process separately.
FuzzifiED.Libpath — ConstantFuzzifiED.Libpath :: String = FuzzifiED_jll.LibpathFuzzifiEDdefine path of the Fortran library libfuzzified.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.
FuzzifiED.ElementType — TypeFuzzifiED.ElementType :: DataType = ComplexF64set the default type of the operator elements, either ComplexF64 or Float64. ComplexF64 by default.
FuzzifiED.OpenHelp! — MethodFuzzifiED.OpenHelp!()A shortcut to open the link for documentation docs.fuzzified.world in the system browser.
Quantum Numbers
FuzzifiED implements diagonal and off-diagonal quantum numbers. They are defined as
FuzzifiED.QNDiag — TypeQNDiagThe mutable type QNDiag records the information of a diagonal $\mathrm{U}(1)$ or $ℤ_p$ quantum number in the form of a symmetry charge
\[Q=∑_{o=1}^{N_o}q_on_o\]
or
\[Q=∑_{o=1}^{N_o}q_on_o\ \mathrm{mod}\ p\]
where $i=1,…,N_U$ is the index of quantum number, $o$ is the index of site, $n_o=c^†_oc_o$, and $q_o$ is a set of coefficients that must be integer valued.
Fields
name :: Stringis the name of the diagonal quantum number.charge :: Vector{Int64}is the symmetry charge $q_o$ of each site.modul :: Int64is the modulus $p$, set to 1 for $\mathrm{U}(1)$ QNDiags.
Initialisation
It can be initialised by the following method
QNDiag([name :: String, ]charge :: Vector{Int64}[, modul :: Int64])The arguments name and modul are facultative. By default name is set to "QN" and modul is set to 1.
FuzzifiED.QNOffd — TypeQNOffdThe mutable type QNOffd records the information of an off-diagonal $ℤ_p$ quantum number in the form of a discrete transformation
\[𝒵:\ c_o↦ α_o^* c^{(p_o)}_{π_o}, c_o^†↦ α_o c^{(1-p_o)}_{π_o}\]
where we use a notation $c^{(1)}=c^†$ and $c^{0}=c$ for convenience, $π_o$ is a permutation of $1,…,N_o$, $α_o$ is a coefficient, and $p_o$ specified whether or not particle-hole transformation is performed for the site. Note that one must guarentee that all these transformations commute with each other and also commute with the diagonal QNs.
Arguments
perm :: Vector{Int64}is a length-$N_o$ vector that records the permutation $π_o$.ph :: Vector{Int64}is a length-$N_o$ vector that records $p_o$ to determine whether or not to perform a particle-hole transformation.fac :: Vector{ComplexF64}is a length-$N_o$ vector that records the factor $α_o$ in the transformation.cyc :: Int64is the cycle $p$.
Initialisation
It can be initialised by the following method
QNOffd(perm :: Vector{Int64}[, ph :: Vector{Int64}][, fac :: Vector{ComplexF64}][, cyc :: Int64])
QNOffd(perm :: Vector{Int64}, ph_q :: Bool[, fac :: Vector{ComplexF64}])The arguments ph, fac and cyc are facultative. By default ph is set all 0, fac is set to all 1 and cyc is set to 2. If ph_q is a bool and true, then ph is set to all 1.
The QNDiag can be added or multiplied by a number
Base.:+ — Method+(qnd1 :: QNDiag, qnd2 :: QNDiag) :: QNDiag
-(qnd1 :: QNDiag, qnd2 :: QNDiag) :: QNDiagreturns the sum or substraction of two QNDiags, 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 QNDiag will be returned.
Base.:* — Method*(fac :: Int64, qnd :: QNDiag) :: QNDiag
*(qnd :: QNDiag, fac :: Int64) :: QNDiag
÷(qnd :: QNDiag, fac :: Int64) :: QNDiag
-(qnd :: QNDiag) :: QNDiagreturns the QNDiag 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 QNDiag will be returned.
The QNOffd can be composed
Base.:* — Method*(qnf1 :: QNOffd, qnf2 :: QNOffd) :: QNOffdreturns the composition of two QNOffd transformations. The cycle is set to be the LCM of two QNOffds.
Empty orbitals can be added to the left or right for QNDiag and QNOffd
FuzzifiED.PadQNDiag — FunctionPadQNDiag(qnd :: QNDiag, nol :: Int64, nor :: Int64)adds nol empty orbitals to the left and nor empty orbitals to the right, implemented as
QNDiag(qnd.name, [ fill(0, nol) ; qnd.charge ; fill(0, nor) ], qnd.modul)FuzzifiED.PadQNOffd — FunctionPadQNOffd(qnf :: QNOffd, nol :: Int64, nor :: Int64)adds nol empty orbitals to the left and nor empty orbitals to the right, implemented as
QNOffd(
[ collect(1 : nol) ; qnf.perm .+ nol ; collect(1 : nor) .+ (nol + length(qnf.perm)) ],
[ fill(0, nol) ; qnf.ph ; fill(0, nor) ],
[ fill(ComplexF64(1), nol) ; qnf.fac ; fill(ComplexF64(1), nor) ], qnf.cyc)Configurations
FuzzifiED.Confs — TypeConfsThe mutable type Confs stores all the configurations that respects the diagonal quantum numbers (QNDiag) and also a table to inversely look up the index from the configuration.
Fields
no :: Int64is the number of sites.ncf :: Int64is the number of configurations.conf :: Vector{Int64}is an array of lengthncfcontaining all the 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.nor :: 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.Confs — MethodConfs(no :: Int64, secd :: Vector{Int64}, qnd :: Vector{QNDiag} ; nor :: Int64 = div(no, 2), num_th :: Int64, disp_std :: Bool)generates the configurations from the list of QNDiags.
Arguments
no :: Int64is the number of sites $N_o$.secd :: Vector{Int64}is the set of $Q_i$ for the selected configurations in the sector.qnd :: Vector{QNDiag}is the set of QNDiags.nor :: Int64is the number of less significant bits used to generate the Lin table. Facultative, $N_o/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 :: Confsis a Confs object.
Note
If your qnd has negative entries, QNDiags must contain the number of electrons.
The inverse look-back for a configuration can be done by
FuzzifiED.GetConfId — MethodGetConfId(cfs :: Confs, cf :: Int64) :: Int64inversely look up the index from the configuration
Arguments
cfs :: Confsstores the configurations.cf :: Int64stores the configuration to be looked-up 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.
Output
id :: Int64is the id of the configuration such thatcfs.conf[id] == cf.
Basis
FuzzifiED.Basis — TypeBasisThe mutable type Basis stores the information of the basis that respects both diagonal and off-diagonal quantum numbers. The states in the basis 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 :: Confsstores the configurations that respect the QNDiags.dim :: Int64is the dimension of the basis.szz :: Int64records the maximum size $\max m_I$ 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.Basis — MethodBasis(cfs :: Confs, secf :: Vector{ComplexF64}, qnf :: Vector{QNOffd} ; num_th :: Int64, disp_std :: Bool)generates the basis that respects the off-diagonal $ℤ_p$ quantum numbers (QNOffd)
Arguments
cfs :: Confsis 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{QNOffd}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 :: Basisis the resulting Basis object.
FuzzifiED.Basis — MethodBasis(cfs :: Confs)Generate a basis from the configurations without off-diagonal $ℤ_n$ symmetries.
Arguments
cfs :: Confsis the diagonal QN–preserving configurations.
Output
bs :: Basisis the resultingBasisobject.
The look-back of the weight of a configuration in a state can be done by
FuzzifiED.GetConfWeight — MethodGetConfWeight(bs :: Basis, st :: Vector{<:Number}, cf :: Int64) :: ComplexF64looks up a the weight of a configuration in a state.
Arguments
bs :: Basisis the basis of the state.st :: Vector{ComplexF64}orst :: Vector{Float64}is a vector of lengthbs.dimthat stores the state.cf :: Int64stores the configuration to be looked-up 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.
Output
- The weight of the configuration in the state.
Term
FuzzifiED.Term — TypeTermThe mutable type Term records a term that looks like $Uc^{(p_1)}_{o_1}c^{(p_2)}_{o_2}… c^{(p_l)}_{o_l}$ in an operator
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.
Initialisation
It can be generated by the function
Term(coeff :: Number, cstr :: Vector{Int64})FuzzifiED.Terms — TypeTermsTerms is an alias for Vector{Term} for convenience
Initialisation
Terms(coeff :: Number, cstr :: Vector{Int64})Gives a Terms with a single Term.
Special elements
The zero and identity terms are defined.
zero(Terms) = Term[]
one(Terms) = Terms(1, [-1, -1])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 :: Terms) :: Terms
-(tms :: Terms) :: Terms
*(tms :: Terms, fac :: Number) :: Terms
/(tms :: Terms, fac :: Number) :: TermsReturn the product of a collection of terms with a number.
Base.:+ — Method+(tms1 :: Terms, tms2 :: Terms) :: Terms
-(tms1 :: Terms, tms2 :: Terms) :: TermsReturn the naive sum of two series of terms by taking their union.
Base.:* — Method*(tms1 :: Terms, tms2 :: Terms) :: Terms
^(tms :: Terms, pow :: Int64) :: TermsReturn the naive product of two series of terms or the power of one terms. The number of terms equals the product of the number of terms in tms1 and tms2. For each term in tms1 $Uc^{(p_1)}_{o_1}…$ and tms2 $U'c^{(p'_1)}_{o'_1}…$, a new term is formed by taking $UU'c^{(p_1)}_{o_1}… c^{(p'_1)}_{o'_1}…$.
Base.adjoint — Methodadjoint(tm :: Term) :: Term
adjoint(tms :: Terms) :: TermsReturn the Hermitian conjugate of a series of terms. For each term $Uc^{(p_1)}_{o_1}c^{(p_2)}_{o_2}… c^{(p_l)}_{o_l}$, the adjoint is $\bar{U}c^{(1-p_l)}_{o_l}… c^{(1-p_2)}_{o_2}c^{(1-p_1)}_{o_1}$.
FuzzifiED.ParticleHole — MethodParticleHole(tm :: Term) :: Term
ParticleHole(tms :: Terms) :: TermsReturn the particle-hole transformation of a series of terms. For each term $Uc^{(p_1)}_{o_1}c^{(p_2)}_{o_2}… c^{(p_l)}_{o_l}$, the transformation results in $Uc^{(1-p_1)}_{o_1}c^{(1-p_2)}_{o_2}…c^{(1-p_l)}_{o_l}$.
The terms can be simplified by
FuzzifiED.NormalOrder — MethodNormalOrder(tm :: Term) :: Termsrearrange a term 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 terms whose result is equal to the original term.
FuzzifiED.SimplifyTerms — MethodSimplifyTerms(tms :: Terms ; cutoff :: Float64 = eps(Float64)) :: Termssimplifies the sum of terms such that
- each term is normal ordered,
- like terms are combined, and terms with zero coefficients are removed.
Argument
cutoff :: Float64is the cutoff such that terms with smaller absolute value of coefficients will be neglected. Facultative,eps(Float64)by default.
We allow the removal or relabelling of orbitals by
FuzzifiED.RemoveOrbs — MethodRemoveOrbs(tms :: Terms, o_rm :: Vector{Int64}) :: Termsremove all the terms that has creation or annihilation operators of sites in the set o_rm.
FuzzifiED.RelabelOrbs — MethodRelabelOrbs(tms :: Terms, 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.PadTerms — FunctionPadTerms(tms :: Terms, nol :: Int64)adds nol empty orbitals to the left by shifting each orbital index.
Operator
FuzzifiED.Operator — TypeOperatorThe mutable type Operator 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 :: Basisis the basis of the initial state.bsf :: Basisis 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 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.
It can be generated by the following methods.
FuzzifiED.Operator — MethodOperator(bsd :: Basis[, bsf :: Basis], terms :: Terms ; red_q :: Int64, sym_q :: Int64, num_th :: Int64, disp_std :: Bool)generates an operator object from a series of terms.
Arguments
bsd :: Basisis the basis of the initial state.bsf :: Basisis the basis of the final state. Facultative, the same asbsdby default.terms :: Termsrecords 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 :: Operator, st_d :: Vector{ComplexF64} ; num_th :: Int64, disp_std :: Bool) :: Vector{ComplexF64}
*(op :: Operator, 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 :: Operator, st_d :: Vector{ComplexF64} ; num_th :: Int64, disp_std :: Bool) :: ComplexF64
*(st_fp :: LinearAlgebra.Adjoint{Float64, Vector{Float64}}, op :: Operator, 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
FuzzifiED.OpMat — TypeOpMat{ComplexF64}
OpMat{Float64}The mutable type OpMat{T} where the element type T can be Float64 and ComplexF64 stores a sparse matrix in the same form as SparseMatrixCSC in SparseArrays. If the matrix is Hermitian or symmetric, only the lower triangle is stored.
Fields
dimd :: Int64anddimf :: Int64are the number of columns and rows of the matrix.symq :: Int64records whether or not the matrix is Hermitian or symmetric.nel :: Int64records the number of elements.colptr :: Vector{Int64},rowid :: Vector{Int64}andelval :: Vector{ComplexF64}records the elements of the sparse matrix as in theSparseMatrixCSCelements of Julia.
It can be generated by the following methods.
FuzzifiED.OpMat — MethodOpMat[{T}](op :: Operator ; num_th :: Int64, disp_std :: Bool)Generates the sparse matrix from the operator. The parametric type T is either Float64 or ComplexF64 ; it is facultative, given by ElementType by default.
Arguments
op :: Operatoris the operator.T :: DataTypespecifies the type of the matrix. It can either beComplexF64orFloat64. Facultative, the same asElementTypeby 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 lowest eigenstates of the sparse matrix can be calculated by
FuzzifiED.GetEigensystem — MethodGetEigensystem(mat :: OpMat{ComplexF64}, nst :: Int64 ; tol :: Float64, ncv :: Int64, initvec :: Vector{ComplexF64}, num_th :: Int64, disp_std :: Bool) :: Tuple{Vector{ComplexF64}, Matrix{ComplexF64}}
GetEigensystem(mat :: OpMat{Float64}, nst :: Int64 ; tol :: Float64, ncv :: Int64, initvec :: Vector{Float64}, num_th :: Int64, disp_std :: Bool) :: Tuple{Vector{Float64}, Matrix{Float64}}calls the Arpack package to calculate the lowest eigenstates of sparse matrix.
Arguments
mat :: OpMat{ComplexF64}ormat :: OpMat{Float64}is the matrix.nst :: Int64is the number of eigenstates to be calculated.tol :: Float64is the tolerence for the Arpack process. The default value is1E-8.ncv :: Int64is an auxiliary parameter needed in the Arpack process. The default value ismax(2 * nst, nst + 10).initvec :: Vector{ComplexF64}orinitvec :: Vector{Float64}is the initial vector. If empty, a random initialisation shall be used. Facultative, empty 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
- A length-
nstarray that has the same type asmatrecording the eigenvalues, and - A
dimd×nstmatrix that has the same type asmatwhere every column records an eigenstate.
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*(mat :: OpMat{ComplexF64}, st_d :: Vector{ComplexF64} ; num_th :: Int64) :: Vector{ComplexF64}
*(mat :: OpMat{Float64}, st_d :: Vector{Float64} ; num_th :: Int64) :: Vector{Float64}Measure the action of a sparse matrix on a state. st_d must be of length mat.dimd. Returns a vector of length mat.dimf that represents the final state.
Facultative argument
num_th :: Int64, the number of threads. Facultative,NumThreadsby default.
Base.:* — Method*(st_fp :: LinearAlgebra.Adjoint{ComplexF64, Vector{ComplexF64}}, mat :: OpMat{ComplexF64}, st_d :: Vector{ComplexF64} ; num_th :: Int64) :: ComplexF64
*(st_fp :: LinearAlgebra.Adjoint{Float64, Vector{Float64}}, mat :: OpMat{Float64}, st_d :: Vector{Float64} ; num_th :: Int64) :: Float64Measuring the inner product between two states and a sparse matrix. st_d must be of length mat.dimd and st_fp must be of length mat.dimf, and st_fp must be an adjoint.
Facultative argument
num_th :: Int64, the number of threads. Facultative,NumThreadsby default.
Note that sometimes it is needed to transform a state from one basis to another. This can be done by constructing an identity operator.
stf = Operator(bsd, bsf, one(Terms)) * stdTransformation
FuzzifiED.Transf — TypeTransfThe mutable type Transf records a transformation in the same form as a QNOffd
\[𝒵:\ 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 :: Basisis the basis of the initial state.bsf :: Basisis the basis of the final state.perm :: Vector{Int64},ph :: Vector{Int64}andfac :: Vector{ComplexF64}records the transformation in the same form as a QNOffd.
FuzzifiED.Transf — MethodTransf(bsd :: Basis[, bsf :: Basis], qnf :: QNOffd)generates a transformation object from a QNOffd.
Arguments
bsd :: Basisis the basis of the initial state.bsf :: Basisis the basis of the final state. Facultative, the same asbsdby default.qnf :: QNOffdrecords the transformation.
Base.:* — Method*(trs :: Transf, st_d :: Vector{ComplexF64} ; num_th = NumThreads) :: Vector{ComplexF64}
*(trs :: Transf, 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.
Entanglement
FuzzifiED.StateDecompMat — MethodStateDecompMat(st :: Vector{<:Number}, bs0 :: Basis, bsa :: Basis, bsb :: Basis, amp_oa :: Vector{ComplexF64}, amp_ob :: Vector{ComplexF64}) :: 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 :: Basisis the total basis.bsa :: Basisis the basis for the subsystem A.bsb :: Basisis the basis for the subsystem B.amp_oa :: Vector{ComplexF64}is a complex list of lengthnothat specifies the 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_ob :: Vector{ComplexF64}is a complex list of lengthnothat specifies the 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 :: Basis, secd_lst :: Vector{Vector{Vector{Int64}}}, secf_lst :: Vector{Vector{Vector{<:Number}}} ; qnd_a :: Vector{QNDiag}, qnd_b :: Vector{QNDiag} = qnd_a, qnf_a :: Vector{QNOffd}, qnf_b :: Vector{QNOffd} = 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 :: Basisis 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{QNDiag}, qnd_b :: Vector{QNDiag} = 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.