Documentation

Miscelaneous

Python reference manual

class nrgljubljana_interface.Flat(half_bandwidth)[source]

The Hilbert transform of a flat density of states, with cut-off

\[g(z) = \int \frac{A(\omega)}{z-\omega} d\omega\]

where \(A(\omega) = \theta( D^2 - \omega^2)/(2D)\).

(Only works in combination with frequency Green’s functions.)

class nrgljubljana_interface.SemiCircular(half_bandwidth, chem_potential=0.0)[source]

Hilbert transform of a semicircular density of states, i.e.

\[g(z) = \int \frac{A(\omega)}{z-\omega} d\omega\]

where \(A(\omega) = \theta( D - |\omega|) 2 \sqrt{ D^2 - \omega^2}/(\pi D^2)\).

(Only works in combination with frequency Green’s functions.)

class nrgljubljana_interface.Solver(**params_kw)[source]
solve(**params_kw)[source]

Solve the impurity problem.

Parameters:params_kw : dict {‘param’:value} that is passed to the core solver.
class nrgljubljana_interface.SolverCore

The Solver class

A_w

The spectral function

B_w

The spectral function of the auxiliary correlator F_w

Delta_w

The hybridization function in real frequencies

F_w

The auxiliary Green function F_w = Sigma_w * G_w

G_w

The retarded Greens function

Sigma_w

The retarded Self energy

check_model_params()

Signature : (nrgljubljana_interface::solve_params_t sp) -> None

chi_NN_w

Charge susceptibility

chi_SS_w

Spin susceptibility

chi_struct

The susceptibility structure object

create_tempdir()

Signature : (str tempdir_) -> str

expv

Expectation values of local impurity operators

generate_param_file()

Signature : (float z) -> None

gf_struct

The Green function structure object

static hdf5_format()

Signature : () -> str

instantiate()

Signature : (float z, str taskdir) -> None

log_mesh

Logarithmic mesh

nrg_params

Low-level NRG parameters

readA()

Signature : (str name, std::optional<g_w_t> A_w, triqs::hilbert_space::gf_struct_t _gf_struct) -> None Read a block spectral function name-block-ij.dat; here we assume that the

spectral function is purely real.
readGF()

Signature : (str name, std::optional<g_w_t> G_w, triqs::hilbert_space::gf_struct_t _gf_struct) -> None Read a block Green’s function (im/re)name-block-ij.dat

read_structure()

Signature : (str filename, bool mandatory) -> triqs::hilbert_space::gf_struct_t

readexpv()

Signature : (int Nz) -> None Read expectation values

set_nrg_params()
Parameter Name Type Default Documentation
dmnrg bool false Perform DMNRG (density-matrix NRG) calculation
cfs bool false Perform CFS (complete Fock space) calculation
fdm bool true Perform FDM (full-density-matrix) calculation
fdmexpv bool true Calculate expectation values using FDM algorithm
dmnrgmats bool false DMNRG calculation on Matsubara axis
fdmmats bool false FDM calculation on Matsubara axis
mats size_t 100 Number of Matsubara points to collect
specgt std::string “” Conductance curves to compute
speci1t std::string “” I_1 curves to compute
speci2t std::string “” I_2 curves to compute
v3mm bool false Compute 3-leg vertex on matsubara/matsubara axis?
mMAX int -1 Number of sites in the star representation
Nmax int -1 Number of sites in the Wilson chain
xmax double -1.0 Largest x in the discretization ODE solver
discretization std::string “Z” Discretization scheme
z double 1.0 Parameter z in the logarithmic discretization
tri std::string “old” Tridiagonalisation approach
preccpp size_t 2000 Precision for tridiagonalisation
diag std::string “default” Eigensolver routine (dsyev|dsyevr|zheev|zheevr|default)
diagratio double 1.0 Ratio of eigenstates computed in partial diagonalisation
dsyevrlimit size_t 100 Minimal matrix size for dsyevr
zheevrlimit size_t 100 Minimal matrix size for zheevr
restart bool true Restart calculation to achieve truncation goal?
restartfactor double 2.0 Rescale factor for restart=true
safeguard double 1e-5 Additional states to keep in case of a near degeneracy
safeguardmax size_t 200 Maximal number of additional states
fixeps double 1e-15 Threshold value for eigenvalue splitting corrections
betabar double 1.0 Parameter bar{beta} for thermodynamics
gtp double 0.7 Parameter p for G(T) calculations
chitp double 1.0 Parameter p for chi(T) calculations
finite bool false Perform Costi-Hewson-Zlatic finite-T calculation
cfsgt bool false CFS greater correlation function
cfsls bool false CFS lesser correlation function
fdmgt bool false FDM greater correlation function?
fdmls bool false FDM lesser correlation function?
fdmexpvn size_t 0 Iteration where we evaluate the expectation values
finitemats bool false T>0 calculation on Matsubara axis
dm bool false Compute density matrixes?
broaden_min_ratio double 3.0 Auto-tune broaden_min parameter
omega0 double -1.0 Smallest energy scale in the problem
omega0_ratio double 1.0 omega0 = omega0_ratio x T
diagth int 1 Number of diagonalisation threads
substeps bool false Interleaved diagonalization scheme
strategy std::string “kept” Recalculation strategy
Ninit size_t 0 Initial Wilson chain ops
reim bool false Output imaginary parts of correlators?
dumpannotated size_t 0 Number of eigenvalues to dump
dumpabs bool false Dump in terms of absolute energies
dumpscaled bool true Dump using omega_N energy units
dumpprecision size_t 8 Dump with # digits of precision
dumpgroups bool true Dump by grouping degenerate states
grouptol double 1e-6 Energy tolerance for considering two states as degenerate
dumpdiagonal size_t 0 Dump diagonal matrix elements
savebins bool true Save binned (unbroadened) data
broaden bool false Enable broadening of spectra
emin double -1.0 Lower binning limit
emax double -1.0 Upper binning limit
bins size_t 1000 bins/decade for spectral data
accumulation double 0.0 Shift of the accumulation points for binning
linstep double 0 Bin width for linear mesh
discard_trim double 1e-16 Peak clipping at the end of the run
discard_immediately double 1e-16 Peak clipping on the fly
goodE double 2.0 Energy window parameter for patching
NN1 bool false Do N/N+1 patching?
NN2even bool true Use even iterations in N/N+2 patching
NN2avg bool false Average over even and odd N/N+2 spectra
NNtanh double 0.0 a in tanh[a(x-0.5)] window function
width_td size_t 16 Width of columns in ‘td’ output file
width_custom size_t 16 Width of columns in ‘custom’ output file
prec_td size_t 10 Precision of columns in ‘td’ output file
prec_custom size_t 10 Precision of columns in ‘custom’ output file
prec_xy size_t 10 Precision of spectral function output
resume bool false Attempt restart?
log std::string “” List of tokens to define what to log
logall bool false Log everything
done bool true Create DONE file?
calc0 bool true Perform calculations at 0-th iteration?
lastall bool false Keep all states in the last iteratio for DMNRG
lastalloverride bool false Override automatic lastall setting
dumpsubspaces bool false Save detailed subspace info
dump_f bool false Dump <f> matrix elements
dumpenergies bool false Dump (all) energies to file?
logenumber size_t 10 # of eigenvalues to show for log=e
stopafter std::string “” Stop calculation at some point?
forcestop int -1 Stop iteration?
removefiles bool true Remove temporary data files?
noimag bool true Do not output imaginary parts of expvs
checksumrules bool false Check operator sumrules
checkdiag bool false Test diag results
checkrho bool false Test tr(rho)=1
set_verbosity()

Signature : (bool v) -> None

solve()

Solve method that performs NRGLJUBLJANA_INTERFACE calculation

Parameter Name Type Default Documentation
Lambda double 2.0 Logarithmic discretization parameter
Nz int 1 Number of discretization meshes
Tmin double 1e-4 Lowest scale on the Wilson chain
keep size_t 100 Maximum number of states to keep at each step
keepenergy double -1.0 Cut-off energy for truncation
keepmin size_t 0 Minimum number of states to keep at each step
T double 0.001 Temperature, k_B T/D,
alpha double 0.3 Width of logarithmic gaussian
gamma double 0.2 Parameter for Gaussian convolution step
method std::string “fdm” Method for calculating the dynamical quantities
bandrescale double -1.0 Band rescaling factor (half-width of the support of the hybridisation function)
model_parameters std::map<std::string, double> Model parameters
solve_one()

Signature : (str taskdir) -> None

tdfdm

Thermodynamic variables (FDM algorithm)