Cluster Dynamical Mean Field Theory

This submodule provides functions that perform the CDMFT algorithm. The class CDMFT contains the CDMFT algorithm per se and implements it. Some auxiliary classes are also defined:

  • convergence_manager defines a test for assessing convergence of the CDMFT procedure

  • frequency_grid defines a grid of imaginary frequencies for computing the CDMFT distance function.

  • hybridization defines the hybridization function, mostly for import/export.

  • general_bath defines general cluster models with bath, with automatic bath parameter definitions, etc.

frequency grids

All frequency grids used in CDMFT are defined along the imaginary axis. The grid is built by passing a frequency_grid instance to CDMFT() via its grid argument; if none is supplied a default grid is constructed. The frequency_grid constructor takes a grid_type and a specs tuple. The available grid types are:

  • grid_type='legendre' (default): a non-uniform Gauss-Legendre grid, defined by the tuple (w1, w2, n1, n2, n3). A Gauss-Legendre grid is defined in the interval \((0,\omega_1)\) with \(n_1\) points, another one in the interval \((\omega_1,\omega_2)\) with \(n_2\) points, and a final one as a function of \(1/\omega\) in the interval \((0,1/\omega_2)\) with \(n_3\) points.

  • grid_type='matsubara' : the set of fermionic Matsubara frequencies \(\omega_n\) at inverse temperature beta between 0 and a maximum frequency wc, defined by the tuple (wc, beta). All frequencies carry the same integration weight.

  • grid_type='regular' : a uniform grid defined by the tuple (wc, n1, n2), with \(n_1\) points in \((0, \omega_c)\) and \(n_2\) points in the high-frequency complement.

In addition, the constructor takes an opt argument that modifies the CDMFT distance weights:

  • opt='self' : the CDMFT weights are scaled like the norm of the self-energy at each frequency (the Hartree-Fock part is subtracted), giving more weight to those frequencies at which correlations are more important.

  • opt='ifreq' : the CDMFT weights are scaled like \(1/\omega_n\).

Subbaths

It is possible to associate more than one bath system to a given cluster. In defining the cluster, one must simply provide a liste of cluster_model objects instead of one. This feature is useless unless a bath system is defined in the cluster_model object. If \(N_s\) cluster models are associated to a cluster (we call them systems), then we are faced with \(N_s\) sets of bath parameters and as many different self-energies \(\Sigma_i\) and hybridization functions \(\Gamma_i\) . The default behavior is to average these self-energies and hybridization functions when computing the host, from the projected Green function. This allows in principle a tighter correspondence between the host and the set of impurity models. This is newer feature (version > 2.19).

Description of functions

class CDMFT(model, varia, grid=None, maxiter=32, miniter=0, convergence='parameters', depth=2, accur_bath=1e-06, accur=0.0001, accur_dist=1e-08, hartree=None, converge_with_stdev=False, iteration='broyden', alpha=0.0, anderson_depth=5, anderson_beta=1.0, method='Nelder-Mead', lm_max_nfev=2000, file='cdmft.tsv', iter_file='cdmft_iter.tsv', eps_algo=0, initial_step=0.1, SEF=False, check_ground_state=False, max_function_eval=1000000, compute_potential_energy=False, host_function=None, pre_host=None, max_value=100, bias=None, post_min=None)

class containing the elements of a CDMFT computation. The constructor executes the computation.

Parameters:
  • varia ([str]) – list of variational parameters.

  • grid (str) – frequency_grid object. If None, constructs a default frequency_grid object.

  • maxiter (int) – maximum number of CDMFT iterations

  • miniter (int) – minimum number of CDMFT iterations

  • accur_bath (float) – the x-tolerance for distance function optimization

  • convergence ([str]) – the convergence tests (sequence of strings or single string)

  • accur ([float]) – the tolerance for the various convergence tests (sequence of floats or single float)

  • accur_dist ([float]) – relative tolerance of the distance function when minimizing it

  • hartree ([hartree]) – mean-field hartree couplings to incorportate in the convergence procedure

  • converge_with_stdev (bool) – If True, checks convergence using the standard deviation of the convergence tests, not the difference

  • iteration (str) – method of iteration of parameters (‘fixed_point’, ‘broyden’, or ‘anderson’)

  • alpha (float) – if iteration=’fixed_point’, damping parameter (fraction of the previous iteration in the new one). If iteration=’broyden’, 1+alpha is the inverse initial Jacobian (or alpha can literally be a matrix, the inverse Jacobian from a previous run). Unused for ‘anderson’.

  • anderson_depth (int) – history depth m for Anderson mixing (number of previous steps to mix); only used when iteration=’anderson’

  • anderson_beta (float) – mixing parameter for Anderson mixing (1.0 = full step, < 1.0 = damped); only used when iteration=’anderson’

  • method (str) – minimization method. Derivative-free choices: ‘Nelder-Mead’ (default), ‘Powell’, ‘CG’, ‘ANNEAL’, NLopt methods ‘NELDERMEAD’, ‘COBYLA’, ‘BOBYQA’, ‘PRAXIS’, ‘SUBPLEX’. Analytical-Jacobian choices (the Jacobian qcm.CDMFT_gradient is activated automatically): ‘trf’ (Trust Region Reflective via scipy.least_squares), ‘BFGS’, ‘L-BFGS-B’. The finite-difference step for the Jacobian is cdmft_jacobian_delta (default 1e-5, tunable via pyqcm.set_global_parameter).

  • lm_max_nfev (int) – maximum number of function/gradient evaluations for the jac-capable methods (default 2000; ignored for derivative-free methods)

  • file (str) – name of the file where the solution is written

  • iter_file (str) – name of the file where the CDMFT iterations are recorded

  • eps_algo (int) – number of elements in the epsilon algorithm convergence accelerator = 2*eps_algo + 1 (0 = no acceleration)

  • initial_step (float) – initial step in the minimization routine

  • SEF (bool) – if True, computes the Potthoff functional at the end

  • check_ground_state (bool) – if True, checks the ground state consistency and raises exception if inconsistent

  • max_function_eval (int) – maximum number of distance function evaluations when minimizing distance

  • compute_potential_energy (bool) – If True, computes Tr(Sigma*G) along with the averages

  • host_function (ndarray) – if not None, function that computes the host array and passes it to qcm

  • pre_host (function) – function to be executed before computing the host. Takes a model instance as argument

  • max_value (float) – maximum absolute value of variational parameters

  • bias (bias_field) – bias field (for spontaneous symmetry breaking) that decreases with iterations

  • post_min (function) – function to be executed before writing to the progress file

Variables:
  • model (lattice_model) – (unique) model on which the computation is based

  • Hyb (ndarray) – host function

  • Hyb_down (ndarray) – host function for the spin down component in the case of mixing=4

  • I – current model instance (changes in the course of the computation)

CDMFT_step()

Performs a CDMFT step that brings the system from the current values (self.x) of the bath and Hartree parameters and updates it to the next set of values

check_convergence()

Performs the CDMFT convergence tests returns True if all are converged, otherwise False

set_Hyb()

Computes the hybridization function, i.e. an array of arrays of matrices. Hyb[i], for cluster #i, is a (nw,d,d) Numpy array. with nw frequencies, and d sites Hyb_down[i] is the same, for the spin-down part, if mixing=4

:returns None

diff_matrix(X, Y)

Computes a difference in hybridization functions between the current iteration (Hyb) and the previous one (Hyb0) :param object X: the current test object :param object Y: any past test object (same structure as X) :returns float: the difference in hybridization arrays

set_sigma()

Computes the self-energy on the frequency grid an array of arrays of matrices. sigma[i], for cluster #i, is a (nw,d,d) Numpy array. with nw frequencies, and d sites

:returns None

plot_iterations()

Plots the variational parameters as a function of iteration, for diagnostics purposes

class convergence_manager(name, diff_func, tol, depth=2, stdev=False)

class managing a convergence test. One or more instances of this class must be used by the CDMFT procedure. One instance corresponds to a quantity that is evaluated at each iteration and then compared with the corresponding quantities at a number (depth) of previous iterations. A difference function (diff_func) is applied to the current value of the quantity and each of its predecessors, and convergence is declared when the result of this difference (a float) is smaller than the tolerance (tol) for each of the depth previous iterations. The possible values of the ‘name’ argument are: * ‘parameters’ : the set of bath parameters is the quantity used. * ‘GS energy’ : the ground state energy of the impurity problem is used. In the case of more than one cluster, the sum is used. * ‘hybridization’ : the hybridization function (or set thereof, in the case of many clusters). The test is based on the norm of the matrix differences, summed over grid frequencies. * ‘self-energy’ : same as above, but using the impurity self-energy instead. * ‘distance’ : the distance function (difference between successive iterations) * any one-body or anomalous lattice operator name. The lattice average is used as test value.

test(x)

Tests for convergence and stores the test object x in an array for comparisons with ‘depth’ future iterations

Parameters:

x (object) – the object containing the current test quantity

Returns:

True if converged, False otherwise

stdev_test(x)

Tests for convergence using the standard deviation of the accumulated sequence of values. Appends x to the internal sequence and computes the standard deviation of the mean; convergence is declared when this standard deviation falls below self.tol.

Parameters:

x (float) – the new value to add to the convergence sequence

Returns:

True if the standard deviation of the mean is below the tolerance, False otherwise

Return type:

bool

print()

Prints the status of the convergence (differences are printed, with a number “depth” of previous iterations)

class frequency_grid(grid_type='legendre', specs=(1, 10, 5, 10, 5), opt='')

This class contains the imaginary frequency grid data, including weights

Parameters:
  • grid_type (str) – type of frequency grid along the imaginary axis : ‘legendre’, ‘matsubara’, ‘regular’

  • specs (tuple) – specific parameters for each grid type. For ‘legendre’, specs=(w1, w2, n1, n2, n3) where w1 and w2 define 3 regions along the imaginary frequency axis. Below w1, n1 grid points are used; between w1 and w2, n2 grid points are used, and n3 grid points are used between w2 and infinity. For ‘matsubara’, specs=(wc, beta), where wc is the frequency cutoff and beta the inverse effective temperature. For ‘regular’, specs=(wc, n1, n2), where wc is the boundary between the low- and high-frequency regions, and n1 and n2 the number of points in each region, respectively.

  • opt (str) – if opt=’self’, the cdmft weights (not the integration weights) are not the same as the integration weights, but scale like the norm of the self-energy at the corresponding frequency (the Hartree-Fock part of the self-energy is subtracted). If opt = ‘ifreq’, the cdmft weights are rather scaled like 1/frequency.

Variables:
  • wr – (real array) the frequencies along the imaginary axis

  • weight – the weight associated to each frequency in an integral

  • cdmft_weight – the weight associated to each frequency in the CDMFT distance function

  • name – the name of the frequency array scheme chosen (for reference)

udpate(I)

Updates the cdmft_weight array depending on the model_instance

Parameters:

I (model_instance) – current model instance

class hybridization(file)

Defines hybridization data from a data file. the frequency is purely imaginary; it is its imaginary part that appears in the file

Parameters:

file (str) – name of the file or string. Format : each line starts with a frequency and then has N*N columns for Delta_{ij}(w)

distance(I)

Evaluates the distance function bewteen the data and the hybridization function of an impurity model instance of label I

Parameters:

I ([model_instance]) – model instance

Returns:

the squared Frobenius norm of the difference matrix, averaged over frequencies

Return type:

float

optimize_bath(model, varia, x, method='Nelder-Mead', accur=0.0001, accur_dist=1e-08)

Optimizes the bath parameters to fit the hybridization function delta :param lattice_model model: the lattice model whose bath parameters are to be optimized :param [str] varia: list of variational parameters (bath parameters) from PyQCM :param [float] x: starting values of the parameters :param str method: minimization method passed to optimize() (default ‘Nelder-Mead’) :param float accur: requested accuracy in the parameters (x-tolerance) :param float accur_dist: requested accuracy in the distance function (f-tolerance)

class general_bath(ns, nb, name='clus', spin_dependent=False, spin_flip=False, singlet=False, triplet=False, complex=False, sites=None)

Class that construct a cluster model with nb bath orbitals, in the most general way possible, with possible restrictions.

Parameters:
  • name (str) – name of the cluster-bath model to be defined

  • ns (int) – number of sites in the cluster

  • nb (int) – number of bath orbitals in the cluster

  • spin_dependent (bool) – if True, the parameters are spin dependent

  • spin_flip (bool) – if True, spin-flip hybridizations are present

  • singlet (bool) – if True, defines anomalous singlet hybridizations

  • triplet (bool) – if True, defines anomalous triplet hybridizations

  • complex (bool) – if True, defines imaginary parts as well, when appropriate

  • sites (list[list[int]]) – 2-level list of sites to couple to the bath orbitals (labels from 1 to ns). Format resembles [[site labels to bind to orbital 1], …] .

Variables:
  • ns (int) – number of physical sites in the cluster

  • nb (int) – number of bath orbitals in the cluster

  • name (str) – name of the cluster model

  • var_E ([str]) – list of bath parameters of the type “energy level”

  • var_H ([str]) – list of bath parameters of the type “hybridization”

  • spin_dependent (bool) – True if the bath parameter conserve spin, but are different for spins up and down

  • spin_flip (bool) – True if we include spin-flip hybridizations terms

  • singlet (bool) – True if we include anomalous hybridizations of the singlet type

  • triplet (bool) – True if we include anomalous hybridizations of the triplet type

  • complex (bool) – True if the model is complex-valued and thus hybridization operators have both real and imaginary components

varia_E(c=1)

returns a list of parameter names from the bath energies with the suffix appropriate for cluster c

Parameters:

c (int) – label of the cluster (starts at 1)

Return [str]:

list of parameter names from the bath energies with the suffix appropriate for cluster c

varia_H(c=1)

returns a list of parameter names from the bath hybridization with the suffix appropriate for cluster c

Parameters:

c (int) – label of the cluster (starts at 1)

Return [str]:

list of parameter names from the bath hybridization with the suffix appropriate for cluster c

varia(H=None, E=None, c=1, spin_down=False)

creates a dict of variational parameters to values taken from the hybridization matrix H and the energies E, for cluster c

Parameters:
  • H (ndarray) – matrix of hybridization values, shape (nmixed*nb, nmixed*ns)

  • E (ndarray) – array of energy values for the bath orbitals

  • c (int) – label of the cluster (starts at 1), used as suffix in parameter names

  • spin_down (bool) – True for the spin-down values

Return {str,float}:

dict of variational parameters to values

starting_values(c=1, e=(0.5, 1.5), hyb=(0.5, 0.2), shyb=(0.1, 0.05), pr=False)

returns an initialization string for the bath parameters

Parameters:
  • c (int) – cluster label (starts at 1)

  • e ((float,float)) – bounds of the values for the bath energies (absolute value)

  • hyb ((float,float)) – average and deviation of the normal hybridization parameters

  • shyb ((float,float)) – average and deviation of the anomalous hybridization parameters

  • pr (bool) – prints the starting values if True

Return str:

initialization string

starting_values_PH(c=1, e=(1, 0.5), hyb=(0.5, 0.2), phi=None, pr=False)

returns an initialization string for the bath parameters, in the particle-hole symmetric case.

Parameters:
  • c (int) – cluster label

  • e ((float)) – range of bath energies

  • hyb ((float,float)) – average and deviation of the normal hybridization parameters

  • phi ([int]) – PH phases of the cluster sites proper

  • pr (bool) – if True, prints info

Return str:

initialization string