Input Variables
Description
Contains all global input variables which can be set by the user through the input file. A specific preocedure ed_read_input() should be called to read the input file using parse_input_variable() procedure from SciFortran. All variables are automatically set to a default, looked for and updated by reading into the file and, sequentially looked for and updated from command line (std.input) using the notation variable_name=variable_value(s) (case independent).
Quick access
- Variables:
a_ph,bath_type,beta,bfile,cg_ftol,cg_grad,cg_method,cg_minimize_hh,cg_minimize_ver,cg_niter,cg_norm,cg_pow,cg_scheme,cg_stop,cg_weight,chidens_flag,chiexct_flag,chipair_flag,chispin_flag,cutoff,deltasc,dmft_error,ed_all_g,ed_finite_temp,ed_hw_bath,ed_input_file,ed_mode,ed_obs_all,ed_offset_bath,ed_print_chidens,ed_print_chiexct,ed_print_chipair,ed_print_chispin,ed_print_g,ed_print_g0,ed_print_sigma,ed_read_umatrix,ed_sectors,ed_sectors_shift,ed_solve_offdiag_gf,ed_sparse_h,ed_total_ud,ed_twin,ed_use_kanamori,ed_verbose,eps,exc_field,g_ph,g_ph_diag,gphfile,gs_threshold,hfile,hfmode,hlocfile,jh,jp,jx,jz_basis,jz_max,jz_max_value,lanc_dim_threshold,lanc_method,lanc_ncv_add,lanc_ncv_factor,lanc_ngfiter,lanc_niter,lanc_nstates_sector,lanc_nstates_step,lanc_nstates_total,lanc_tolerance,lfit,lmats,logfile,lpos,lreal,ltau,nbath,ncoeff,ndelta,nerr,nloop,norb,nph,nread,nspin,nsuccess,pair_field,ph_type,print_input_vars,print_sector_eigenvalues,rdm_flag,sb_field,sectorfile,spin_field_x,spin_field_y,spin_field_z,uloc,umatrix_file,ust,w0_ph,wfin,wini,xmax,xmin,xmu- Routines:
Used modules
ed_version: Contains the EDIpack version (git commit hash)
External modules
-
iso_c_binding
Variables
- ed_input_vars/a_ph
Phonon forcing field coupled to displacement operator (constant)
- Type:
real
- Default:
0d0
- ed_input_vars/bath_type
Flag to select the bath geometry
normal: each impurity orbital has a set of bath levels in a star geometryhybrid: all impurity orbitals communicate with the same set of bath levels in a star geometryreplica: the impurity communicates with clusters of the same form via an hybridization term \(V\mathbb{I}\)general: extendsreplicaso that each orbital has a different hybridization strength
- Type:
character(len=8)
- Default:
normal
- ed_input_vars/beta
Inverse temperature, at zero temperature it is used as a IR cut-off.
- Type:
real
- Bindings:
Language = c, Name = beta
- Default:
1000d0
- ed_input_vars/bfile
File where to retrieve/store the bath parameters
- Type:
character(len=100)
- Default:
hbasis[.used/restart]
- ed_input_vars/cg_ftol
Tolerance in the conjugate-gradient fitting procedure
- Type:
real
- Default:
1d-5
- ed_input_vars/cg_grad
Flag to set the type of gradient evaluation (if
cg_method=0)0: analytic1: numeric
- Type:
integer
- Default:
0
- ed_input_vars/cg_method
Conjugate-gradient fitting routine to be used:
0: Numerical Recipes1:minimize (FORTRAN77 code)
- Type:
integer
- Default:
0
- ed_input_vars/cg_minimize_hh
Unknown parameter used in the CG minimize procedure ( for
cg_grad=1)- Type:
real
- Default:
1d-4
- ed_input_vars/cg_minimize_ver
If
cg_grad=1, select which version ofminimize.fto useT: Lichtenstein (newer)F: Krauth (older)
- Type:
logical
- Default:
F
- ed_input_vars/cg_niter
Maximum number of iterations in the bath fitting procedure
- Type:
integer
- Default:
500
- ed_input_vars/cg_norm
Which norm to use in the evaluation of the \(\chi^{2}\) for matrix quantities.
ELEMENTAL: \(\chi^{2}\) is the sum of each component's \(\chi^{2}\)FROBENIUS: \(\chi^{2}\) is calculated on the Frobenius norm (Matrix distance)
Warning
The Frobenius norm is currently implemented only for
ed_bath=replica, general. Also, forcg_pow\(\neq 2\) the Frobenius norm is ill-defined, at least with respect to its usual mathematical meaning. The behavior ofcg_norm=frobeniusmight be changed or removed in future versions of the code, breaking back-compatibility.- Type:
character(len=12)
- Default:
ELEMENTAL
- ed_input_vars/cg_pow
Power exponent in the \(\chi^{2}\) , according to \(\vert \mathcal{G}_{0}(i\omega_{n}) - G_{0}^{\mathrm{And}}(i\omega_{n}) \vert ^{\mathrm{cg\_pow}}\) or \(\vert \Delta(i\omega_{n})- \Delta^{\mathrm{And}}(i\omega_{n}) \vert ^{\mathrm{cg\_pow}}\)
- Type:
integer
- Default:
2
- ed_input_vars/cg_scheme
Which quantity to use in the bath fitting routine:
WEISS: the lattice Weiss field Green's function \(\mathcal{G}_{0}(i\omega_{n})\) is fittedDELTA: the lattice Hybridization function \(\Delta(i\omega_{n})\) is fitted
- Type:
character(len=6)
- Default:
Weiss
- ed_input_vars/cg_stop
Conjugate-gradient stopping condition
0:1 .and. 21: \(\vert F_{n-1} -F_{n} \vert < \mathrm{tol} \cdot (1+F_{n})\)2: \(\vert\vert x_{n-1} -x_{n} \vert\vert <\mathrm{tol} \cdot (1+ \vert\vert x_{n} \vert\vert\))
- Type:
integer
- Default:
0
- ed_input_vars/cg_weight
Weight assigned to the imaginary frequency data points in the calculation of the \(\chi^{2}\)
0: \(1\)1: \(\frac{1}{n}\)2: \(\frac{1}{\omega_{n}}\)
- Type:
integer
- Default:
0
- ed_input_vars/chidens_flag
Flag to activate charge susceptibility evaluation
- Type:
logical
- Bindings:
Language = c, Name = chidens_flag
- Default:
F
- ed_input_vars/chiexct_flag
Flag to activate excitonic susceptibility evaluation
- Type:
logical
- Bindings:
Language = c, Name = chiexct_flag
- Default:
F
- ed_input_vars/chipair_flag
Flag to activate pairing susceptibility evaluation
- Type:
logical
- Bindings:
Language = c, Name = chipair_flag
- Default:
F
- ed_input_vars/chispin_flag
Flag to activate spin susceptibility evaluation
- Type:
logical
- Bindings:
Language = c, Name = chispin_flag
- Default:
F
- ed_input_vars/cutoff
Spectrum cut-off, used to determine the number states to be retained
- Type:
real
- Default:
1d-9
- ed_input_vars/deltasc
Value of the SC symmetry breaking term (only used if
ed_mode=superc)- Type:
real
- Default:
2d-2
- ed_input_vars/dmft_error
Error threshold for DMFT convergence
- Type:
real
- Bindings:
Language = c, Name = dmft_error
- Default:
1d-5
- ed_input_vars/ed_all_g
Flag to evaluate all the components of the impurity Green`s functions irrespective of the symmetries
- Type:
logical
- Default:
F
- ed_input_vars/ed_finite_temp
Flag to set whether the problem is at finite temperature.If
Tthenlanc_nstates_totalmust be greater than 1.- Type:
logical
- Default:
F
- ed_input_vars/ed_hw_bath
Half-bandwidth for bath level initialization if
bath_type=normal, hybrid. The levels will be equispaced in the range[-hw,hw]- Type:
real
- Default:
2d0
- ed_input_vars/ed_input_file
Name of input file
- Type:
character(len=200)
- ed_input_vars/ed_mode
Flag to set the ED mode
normal: normalsuperc: s-wave superconductive (singlet pairing)nonsu2: broken SU(2) symmetry
- Type:
character(len=7)
- Default:
normal
- ed_input_vars/ed_obs_all
Flag to print observables for every loop
- Type:
logical
- Default:
T
- ed_input_vars/ed_offset_bath
Offset for the initialization of diagonal terms if
bath_type=replica, general. The replicas will be equally spaced in the range[-offset,offset]- Type:
real
- Default:
1d-1
- ed_input_vars/ed_print_chidens
Flag to print impurity dens susceptibility, if calculated
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_chiexct
Flag to print impurity exct susceptibility, if calculated
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_chipair
Flag to print impurity pair susceptibility, if calculated
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_chispin
Flag to print impurity spin susceptibility, if calculated
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_g
Flag to print impurity Green`s functions
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_g0
Flag to print impurity non-interacting Green`s functions
- Type:
logical
- Default:
T
- ed_input_vars/ed_print_sigma
Flag to print impurity Self-energies
- Type:
logical
- Default:
T
- ed_input_vars/ed_read_umatrix
Flag to enable (
T) or not (F) reading the two-body terms from an external filedefined byumatrix_file- Type:
logical
- Bindings:
Language = c, Name = ed_read_umatrix
- Default:
F
- ed_input_vars/ed_sectors
Flag to reduce sector scan for the spectrum to specific sectors
- Type:
logical
- Default:
F
- ed_input_vars/ed_sectors_shift
Additional sectors to be evaluated if
ed_sectorsis set. These are sectors with all the quantum numbers varying of at most byed_sectors_shiftaround the sectors listed ined_sectors.- Type:
integer
- Default:
1
- ed_input_vars/ed_solve_offdiag_gf
Flag to select the calculation of the off-diagonal impurity GF. Set to
Tby default ifbath_typeis notnormal- Type:
logical
- Default:
F
- ed_input_vars/ed_sparse_h
Flag to select storage of the Fock space Hamiltonian as a sparse matrix
T: \(H\) is stored. CPU intensiveF: on-the-fly \(H \cdot v\) product is stored. Memory intensive
- Type:
logical
- Default:
T
- ed_input_vars/ed_total_ud
Flag to select which type of quantum numbers have to be considered (if
ed_mode=normal)T: blocks have different total \(N_{\uparrow}\) and \(N_{\downarrow}\)F: blocks have different total \(N^{\alpha}_{\uparrow}\) and \(N^{\alpha}_{\downarrow}\)
where \(\alpha\) is the orbital index. Speeds up calculation in the case where orbitals are not hybridized
- Type:
logical
- Bindings:
Language = c, Name = ed_total_ud
- Default:
T
- ed_input_vars/ed_twin
Flag to reduce (
T) or not (F) the number of visited sector using twin symmetry- Type:
logical
- Bindings:
Language = c, Name = ed_twin
- Default:
F
- ed_input_vars/ed_use_kanamori
Flag to enable (
T) or not (F) the use of the input variablesULOC,UST,JH,JX,JPfor Hubbard-Kanamori coefficients.- Type:
logical
- Bindings:
Language = c, Name = ed_use_kanamori
- Default:
T
- ed_input_vars/ed_verbose
Verbosity level:
0: almost nothing…
3: most of the verbose output…
5: everything. Really, everything
- Type:
integer
- Default:
3
- ed_input_vars/eps
Broadening on the real frequency axis for Green's function and Susceptibility calculations.
- Type:
real
- Bindings:
Language = c, Name = eps
- Default:
1d-2
- ed_input_vars/exc_field
External field coupling to exciton order parameter
- Type:
real(4)
- Default:
zero
- ed_input_vars/g_ph
Electron-phonon coupling constant all
- Type:
complex(•, •)
- Attributes:
allocatable
- Default:
zero
- ed_input_vars/g_ph_diag
- Type:
real(•)
- Attributes:
allocatable
- ed_input_vars/gphfile
File of Phonon couplings. Set to NONE to use only density couplings.
- Type:
character(len=100)
- Default:
NONE
- ed_input_vars/gs_threshold
Energy threshold for ground state degeneracy loop up
- Type:
real
- Default:
1d-9
- ed_input_vars/hfile
File where to retrieve/store the bath parameters
- Type:
character(len=100)
- Default:
hamiltonian[.used/restart]
- ed_input_vars/hfmode
Flag to set the form of the Hubbard-Kanamori interaction
T: \(U(n_{\uparrow}-\frac{1}{2})(n_{\downarrow}-\frac{1}{2})\)F: \(Un_{\uparrow}n_{\downarrow}\)
- Type:
logical
- Default:
T
- ed_input_vars/hlocfile
File to read the input local H from
- Type:
character(len=100)
- Default:
inputHLOC.in
- ed_input_vars/jh
Hund's coupling constant
- Type:
real
- Bindings:
Language = c, Name = jh
- Default:
0.d0
- ed_input_vars/jp
Coupling constant for the Pair-hopping interaction term
- Type:
real
- Bindings:
Language = c, Name = jp
- Default:
0.d0
- ed_input_vars/jx
Coupling constant for the spin-eXchange interaction term
- Type:
real
- Bindings:
Language = c, Name = jx
- Default:
0.d0
- ed_input_vars/jz_basis
Flag to enable the \(J_{z}\) basis in SOC calculations
- Type:
logical
- Default:
F
- ed_input_vars/jz_max
Flag to enable a maximum value for \(J_{z}\)
- Type:
logical
- Default:
F
- ed_input_vars/jz_max_value
Maximum value for Jz
- Type:
real
- Default:
1000d0
- ed_input_vars/lanc_dim_threshold
Minimal sector dimension for Lanczos diagonalization. Smaller sectors will be solved with Exact Diagonalization provided by Lapack
- Type:
integer
- Default:
1024
- ed_input_vars/lanc_method
Flag to select the Lanczos method to be used in the determination of the spectrum.
ARPACK: uses the Arnoldi algorithmLANCZOS: uses an in-house Lanczos algorithm (limited to zero temperature)
- Type:
character(len=12)
- Default:
ARPACK
- ed_input_vars/lanc_ncv_add
Offset to add to the size of the block to prevent it to become too small according to \(N_{cv}=\mathrm{lanc\_ncv\_factor} \cdot \mathrm{Neigen} + \mathrm{lanc\_ncv\_add}\)
- Type:
integer
- Default:
0
- ed_input_vars/lanc_ncv_factor
Size of the block used in Lanczos-Arpack by multiplying the required
Neigenaccording to \(N_{cv}=\mathrm{lanc\_ncv\_factor} \cdot \mathrm{Neigen} + \mathrm{lanc\_ncv\_add}\)- Type:
integer
- Default:
10
- ed_input_vars/lanc_ngfiter
Number of Lanczos iteration in GF determination. Number of moments.
- Type:
integer
- Default:
200
- ed_input_vars/lanc_niter
Max number of Lanczos iterations
- Type:
integer
- Default:
512
- ed_input_vars/lanc_nstates_sector
Max number of required eigenvalues per sector
- Type:
integer
- Default:
2
- ed_input_vars/lanc_nstates_step
Number of states added at each step for finite-temperature calculations: if the latest state included in thepartition function has a Boltzmann weight higher than
cutoff,lanc_nstates_totalwill be increased by this amount at the next DMFT iteration.- Type:
integer
- Default:
2
- ed_input_vars/lanc_nstates_total
Max number of states contributing to the partition function for finite-temperature calculations.It must be set to
1for zero-temperature calculations, and to a greater value for finite-temperature calculations.- Type:
integer
- Default:
1
- ed_input_vars/lanc_tolerance
Tolerance for the Lanczos iterations as used in Arpack and plain Lanczos
- Type:
real
- Default:
1d-18
- ed_input_vars/lfit
Number of frequencies for bath fitting
- Type:
integer
- Bindings:
Language = c, Name = lfit
- Default:
1000
- ed_input_vars/lmats
Number of Matsubara frequencies
- Type:
integer
- Bindings:
Language = c, Name = lmats
- Default:
4096
- ed_input_vars/logfile
Logfile unit
- Type:
integer
- Attributes:
save
- Bindings:
Language = c, Name = logfile
- Default:
6
- ed_input_vars/lpos
Number of points in Probability Distribution Function lattice for phonons
- Type:
integer
- Bindings:
Language = c, Name = lpos
- Default:
100
- ed_input_vars/lreal
Number of real-axis frequencies
- Type:
integer
- Bindings:
Language = c, Name = lreal
- Default:
5000
- ed_input_vars/ltau
Number of imaginary time points
- Type:
integer
- Bindings:
Language = c, Name = ltau
- Default:
1024
- ed_input_vars/nbath
Number of bath sites:
- Type:
integer
- Bindings:
Language = c, Name = nbath
- Default:
6
- ed_input_vars/ncoeff
Multiplier for the
ndeltavalue ifxmuand its error are read from a file ( \(\mathrm{ndelta} \rightarrow \mathrm{ndelta} \cdot \mathrm{ncoeff}\) )- Type:
real
- Default:
1d0
- ed_input_vars/ndelta
Initial chemical potential variation for fixed-density calculations
- Type:
real
- Default:
1d-1
- ed_input_vars/nerr
Error threshold for fixed-density calculations
- Type:
real
- Default:
1d-4
- ed_input_vars/nloop
Maximum number of DMFT loops
- Type:
integer
- Bindings:
Language = c, Name = nloop
- Default:
100
- ed_input_vars/norb
Number of impurity orbitals (max
5)- Type:
integer
- Bindings:
Language = c, Name = norb
- Default:
1
- ed_input_vars/nph
Max number of phonons allowed (cut off)
- Type:
integer
- Bindings:
Language = c, Name = nph
- Default:
0
- ed_input_vars/nread
Target occupation value for fixed-density calculations. If set to
0.0the calculation is assumed to be at fixedxmu- Type:
real
- Bindings:
Language = c, Name = nread
- Default:
0d0
- ed_input_vars/nspin
If
1, assume \(H_{\downarrow}\) = \(H_{\uparrow}\) . If2, the Hamiltonian needs to explicitly include spin-up and spin-down blocks- Type:
integer
- Bindings:
Language = c, Name = nspin
- Default:
1
- ed_input_vars/nsuccess
Number of repeated success to fall below convergence threshold
- Type:
integer
- Bindings:
Language = c, Name = nsuccess
- Default:
1
- ed_input_vars/pair_field
Pair field per orbital coupling to s-wave order parameter component which explicitly appears in the impurity Hamiltonian
- Type:
real(•)
- Attributes:
allocatable
- Default:
zero
- ed_input_vars/ph_type
Shape of the e part of the e-ph interaction:
1= orbital occupation2= orbital hybridization
- Type:
integer
- Default:
1
- ed_input_vars/print_input_vars
Flag to toggle the printing on the terminal output of a list of input variables and their values
- Type:
logical
- Default:
T
- ed_input_vars/print_sector_eigenvalues
Flag to toggle the printing of the eigenvalues of each sectors to a file.
- Type:
logical
- Default:
T
- ed_input_vars/rdm_flag
Flag to activate Reduced Density Matrix evaluation
- Type:
logical
- Default:
F
- ed_input_vars/sb_field
Value of a symmetry breaking field for magnetic solutions
- Type:
real
- Bindings:
Language = c, Name = sb_field
- Default:
1d-1
- ed_input_vars/sectorfile
File where to retrieve/store the sectors contributing to the spectrum
- Type:
character(len=100)
- Default:
sectors[.used/restart]
- ed_input_vars/spin_field_x
Magnetic field per orbital coupling to X-spin component
- Type:
real(•)
- Attributes:
allocatable
- Default:
zero
- ed_input_vars/spin_field_y
Magnetic field per orbital coupling to Y-spin component
- Type:
real(•)
- Attributes:
allocatable
- Default:
zero
- ed_input_vars/spin_field_z
Magnetic field per orbital coupling to Z-spin component
- Type:
real(•)
- Attributes:
allocatable
- Default:
zero
- ed_input_vars/uloc
Values of the local interaction per orbital (max
5)- Type:
real(5)
- Bindings:
Language = c, Name = uloc
- Default:
(/( 2d0,i=1,Norb )/)
- ed_input_vars/umatrix_file
File containing the list of two-body operators
- Type:
character(len=100)
- Default:
umatrix[.used/restart]
- ed_input_vars/ust
Value of the inter-orbital interaction term
- Type:
real
- Bindings:
Language = c, Name = ust
- Default:
0d0
- ed_input_vars/w0_ph
Phonon frequency
- Type:
real
- Default:
0d0
- ed_input_vars/wfin
Maximum value of the real frequency range
- Type:
real
- Bindings:
Language = c, Name = wfin
- Default:
5d0
- ed_input_vars/wini
Minimum value of the real frequency range
- Type:
real
- Bindings:
Language = c, Name = wini
- Default:
-5d0
- ed_input_vars/xmax
Maximum of the x-range for the local lattice probability distribution function (phonons)
- Type:
real
- Bindings:
Language = c, Name = xmax
- Default:
3d0
- ed_input_vars/xmin
Minimum of the x-range for the local lattice probability distribution function (phonons)
- Type:
real
- Bindings:
Language = c, Name = xmin
- Default:
-3d0
- ed_input_vars/xmu
Chemical potential. If
HFMODE=T,xmu=0.d0indicates half-filling condition.- Type:
real
- Bindings:
Language = c, Name = xmu
- Default:
0d0
Subroutines and functions
- subroutine ed_input_vars/ed_read_input(inputunit)
This functions reads the input file provided by
INPUTunitand sets the global variables accordingly- Parameters:
inputunit [character(len=*)]
- Use :
mpi,sf_mpi
- subroutine ed_input_vars/ed_update_input(name, vals)
This functions updates some variables in the input file, namely
exc_field,pair_field,exc_field,spin_field_x,spin_field_y, andspin_field_z.- Parameters:
name [character(len=*)] – the name of the variable to update
vals (•) [real] – the new value of the variable