green_mbtools.pesto¶
green_mbtools.pesto.ac_utils¶
- green_mbtools.pesto.ac_utils.dump_input_caratheodory_data(wsample, X_iw, ifile)[source]¶
Function to dump input data to files for caratheodory analytic continuation.
- Parameters:
wsample (numpy.ndarray) – 1D array of Matsubara frequencies
X_iw (numpy.ndarray) – 5D array of shape (nw, ns, nk, nao, nao)
ifile (string) –
filename for storing the Matsubara data. The data is stored in the format:
iw_1 iw_2 iw_3 ... iw_n Real X(iw_1) Matrix Imag X(iw_1) Matrix Real X(iw_2) Matrix Imag X(iw_2) Matrix ... and so on.
- green_mbtools.pesto.ac_utils.load_caratheodory_data(matrix_file, spectral_file, X_dims)[source]¶
Loads output data from Caratheodory analytic continuation.
- Parameters:
matrix_file (string) –
Path to output matrix file. The output file has the following format:
w1 Re.Xc[11](w1 + i eta) Im.Xc[11](w1 + i eta) Re.Xc[12](w1 + i eta) w2 Re.Xc[11](w2 + i eta) Im.Xc[11](w2 + i eta) Re.Xc[12](w2 + i eta) ... and so on
spectral_file (string) –
Path to output spectral file, which has the following format:
w1 XA(w1 + i eta) w2 XA(w2 + i eta) ... and so on
X_dims (tuple) – shape of tensor for reading and storing the output data into
- Returns:
numpy.ndarray – real frequencies on which analytically continued data is returned
numpy.ndarray – Complex valued output tensor with full matrix form of continued data
numpy.ndarray – Real valued spectral function corresponding to the complex continued data
green_mbtools.pesto.analyt_cont¶
green_mbtools.pesto.cheby¶
- green_mbtools.pesto.cheby.get_tau_mesh(beta, ncheb)[source]¶
Generate tau grid for a given beta in the chebishev basis.
green_mbtools.pesto.dyson¶
- green_mbtools.pesto.dyson.solve_dyson(fock, S=None, sigma=None, mu=0, ir=None)[source]¶
Solve Dyson equation and construct Green’s function. If overlap is None, it is set to identity.
- Parameters:
fock (numpy.ndarray) – Fock matrix of shape (ns, nk, nao, nao)
S (numpy.ndarray, optional) – Overlap matrix of shape (ns, nk, nao, nao), by default None
sigma (numpy.ndarray, optional) – Imaginary-time self-energy of shape (ntau, ns, nk, nao, nao), by default None (zero)
mu (float, optional) – chemical potential, by default 0
ir (green_mbtools.pesto.ir.IR_factory, optional) – Instance of IR factory for transforming between imaginary time and Matsubara frequencies, by default None
- Returns:
Green’s function on imaginary time axis of shape (ntau, ns, nk, nao, nao)
- Return type:
numpy.ndarray
- Raises:
ValueError – if ir is not provided
green_mbtools.pesto.ft¶
- green_mbtools.pesto.ft.compute_fourier_coefficients(kmesh, rmesh)[source]¶
Compute Fourier coefficients for direct and inverse transforms between k and real space
- Parameters:
kmesh (numpy.ndarray) – k-mesh of shape (nk, 3)
rmesh (numpy.ndarray) – real-space mesh
- Returns:
numpy.ndarray – coefficients to transform from r to k space
numpy.ndarray – coefficients to transform from k to r space
- green_mbtools.pesto.ft.construct_rmesh(nkx, nky, nkz)[source]¶
Generate real-space mesh, in the units of lattice translation vectors. e.g.,
nkx = 6 L = (nkx - 1) // 2 # = 2 rx = -2, -1, 0, 1, 2, 3
- Parameters:
nkx (int) – number of k-points in the x direction
nky (int) – number of k-points in the y direction
nkz (int) – number of k-points in the z direction
- Returns:
real-space mesh
- Return type:
numpy.ndarray
- green_mbtools.pesto.ft.construct_symmetric_rmesh(nkx, nky, nkz)[source]¶
Generate a real-space mesh that is symmetric along coordinate axes. e.g.,
nkx = 6 L = (nkx - 1) // 2 # = 2 rx = -3, -2, -1, 0, 1, 2, 3
- Parameters:
nkx (int) – number of k-points in the x direction
nky (int) – number of k-points in the y direction
nkz (int) – number of k-points in the z direction
- Returns:
unique / symmetric real-space mesh
- Return type:
numpy.ndarray
- green_mbtools.pesto.ft.k_to_real(frk, obj_k, weights)[source]¶
Perform Fourier transform from reciprocal to real space
- Parameters:
frk (numpy.ndarray) – coefficients to transform from k to r space
obj_k (numpy.ndarray) – Reciprocal space object
weights (numpy.ndarray) – Fourier coefficient’s degeneracy weights
- Returns:
Real space (inverse Fourier transformed) object
- Return type:
numpy.ndarray
- green_mbtools.pesto.ft.real_to_k(fkr, obj_i)[source]¶
Perform Fourier transform from real to reciprocal space
- Parameters:
fkr (numpy.ndarray) – coefficients to transform from r to k space
obj_i (numpy.ndarray) – Real space (inverse Fourier transformed) object
- Returns:
Reciprocal space object
- Return type:
numpy.ndarray
green_mbtools.pesto.ir¶
- class green_mbtools.pesto.ir.IR_factory(beta, ir_file=None, legacy_ir=False)[source]¶
Bases:
object
Class to handle Fourier transform between imaginary time and Matsubara frequency using the intermediate representation (IR) grids.
- tau_to_w(X_t)[source]¶
Transform from imaginary time to Matsubara axis
- Parameters:
X_t (numpy.ndarray) – Data on imaginary time axis
- Returns:
Data on Matsubara axis
- Return type:
numpy.ndarray
- tau_to_w_other(X_t, wsample)[source]¶
Transform from imaginary time axis to user-provided Matsubara frequency points. Typically, we deal with time and Matsubara grids that are pre-determined by intermediate representation grid files. This function allows transformation to other Matsubarar frequencies.
- Parameters:
X_t (numpy.ndarray) – Data on imaginary time axis
wsample (numpy.ndarray) – Frequency points on which Matsubara data is required
- Returns:
Data on Matsubara axis
- Return type:
numpy.ndarray
- tauf_to_wb(X_t)[source]¶
Transform from fermionic imaginary time to bosonic Matsubara axis. Imaginary time axis is the same for fermions and bosons, but Matsubara frequencies differ. E.g., Polarization from tau to i Omega axis.
- Parameters:
X_t (numpy.ndarray) – Data on imaginary time axis
- Returns:
Data on bosonic Matsubara axis
- Return type:
numpy.ndarray
- update(beta=None, ir_file=None)[source]¶
Update IR grid information for the IR_facory object
- Parameters:
beta (float, optional) – inverse temperature of the calculation, by default None
ir_file (string, required) – IR-grid file, by default None
- green_mbtools.pesto.ir.legacy_read_IR_matrices(ir_path, beta, ptype='fermi')[source]¶
Read IR file with legacy data format
- Parameters:
ir_path (string) – path to IR file
beta (float) – inverse temperature
ptype (str, optional) – parity type (‘fermi’ or ‘bose’), by default ‘fermi’
- Returns:
numpy.ndarray – Imaginary time grid
numpy.ndarray – Matsubara grid
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
- green_mbtools.pesto.ir.new_read_IR_matrices(ir_path, beta, ptype='fermi')[source]¶
Read IR file with new data format
- Parameters:
ir_path (string) – path to IR file
beta (float) – inverse temperature
ptype (str, optional) – parity type (‘fermi’ or ‘bose’), by default ‘fermi’
- Returns:
numpy.ndarray – Imaginary time grid
numpy.ndarray – Matsubara grid
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
numpy.ndarray – Transformation coefficients
green_mbtools.pesto.mb¶
green_mbtools.pesto.orth¶
- green_mbtools.pesto.orth.SAO_matrices(S)[source]¶
Löwdin’s symmetric orthogonalization transformation matrix
- Parameters:
S (numpy.ndarray) – overlap matrix
- Returns:
Transformation matrix for symmetric orthogonalization
- Return type:
numpy.ndarray
- green_mbtools.pesto.orth.canonical_matrices(S, thr=1e-07, type='f')[source]¶
Löwdin’s canonical orthogonalization transformation matrix
- Parameters:
S (numpy.ndarray) – overlap matrix
thr (float, optional) – threshold for orthogonalization and inverses, by default 1e-7
type (str, optional) – ‘f’ for Fock and ‘g’ for Green’s function or density matrix, by default ‘f’
- Returns:
Transformation matrix for canonical orthogonalization
- Return type:
numpy.ndarray
- Raises:
ValueError – if type of input is anything other than ‘g’ or ‘f’.
- green_mbtools.pesto.orth.canonical_orth(H, S, thr=1e-07, type='f')[source]¶
Löwdin’s canonical orthogonalization
- Parameters:
H (numpy.ndarray) – tensor/matrix to be orthogonalized (e.g., Fock, density or Green’s function). It has shape (…, ns, nk, nao, nao)
S (numpy.ndarray) – overlap matrix
thr (float, optional) – threshold for orthogonalization and inverses, by default 1e-7
type (str, optional) – ‘f’ for Fock and ‘g’ for Green’s function or density matrix, by default ‘f’
- Returns:
Orthogonalized matrix / tensor
- Return type:
numpy.ndarray
- green_mbtools.pesto.orth.sao_orth(X, S, type=None)[source]¶
Löwdin’s symmetric orthogonalization
- Parameters:
X (numpy.ndarray) – tensor/matrix to be orthogonalized (e.g., Fock, density or Green’s function). It has shape (…, ns, nk, nao, nao)
S (numpy.ndarray) – overlap matrix
type (str, optional) – ‘f’ for Fock and ‘g’ for Green’s function or density matrix, by default ‘f’
- Returns:
Orthogonalized matrix / tensor
- Return type:
numpy.ndarray
- Raises:
ValueError – if type of input is anything other than ‘g’ or ‘f’.
- green_mbtools.pesto.orth.transform(Z, X, X_inv)[source]¶
Transform Z into X basis using Z_X = X^* Z X
- Parameters:
Z (numpy.ndarray) – Object to be transformed
X (numpy.ndarray) – Transformation matrix
X_inv (numpy.ndarray) – Inverse transformation matrix (used to check the sanity of transformation)
- Returns:
Z in new basis
- Return type:
numpy.ndarray
- green_mbtools.pesto.orth.transform_per_k(Z, X, X_inv)[source]¶
Transform Z into X basis as Z_X = X^* Z X for each k-point
- Parameters:
Z (numpy.ndarray) – Object to be transformed
X (numpy.ndarray) – Transformation matrix
X_inv (numpy.ndarray) – Inverse transformation matrix (used to check the sanity of transformation)
- Returns:
Z in new basis for the given k-point
- Return type:
numpy.ndarray
green_mbtools.pesto.pes_utils¶
- green_mbtools.pesto.pes_utils.cvx_diag_projection(zM, GM_diag, w_cut=10, n_real=201, ofile='Giw', solver='SCS', **kwargs)[source]¶
Projection of noisy Green’s function diagonal data on to Nevanlinna manifold. Once this is performed, analytic continuation either using Nevanlinna or ES becomes easier. The main idea is to obtain a projected Matsubara Green’s function of the form
G(iw_m) = sum_n P_n / (iw_m - w_n)
- Parameters:
zM (numpy.ndarray) – n_imag number of imaginary frequency values
GM_matrix (numpy.ndarray) – Exact Matsubara Green’s function data in the shape (n_imag, n_orb, n_orb)
w_cut (float, optional) – cut-off to consider for real axis, by default 10.0
n_real (int, optional) – number of real frequencies to consider, by default 201
ofile (str, optional) – Output file in which the green’s function data will be dumped, by default ‘Giw’
solver (str, optional) – CVXPy solver to use, by default ‘SCS’
**kwargs (dict, optional) – dictionary of keyword arguments for CVXPy solver
- Returns:
Projected Green’s function
- Return type:
numpy.ndarray
- green_mbtools.pesto.pes_utils.cvx_gradient(poles, GM_matrix, zM)[source]¶
Computes the gradient for top-level optimization of pole positions in ES analytic continuation i.e., d Err / d poles where,
Err (poles) = min_{X} || Gimag (iw) - Gapprox [poles, X] (iw) ||
Specifically, G_approx is defined as
G_approx = sum_m X_vec(:,:,m)/(1j*zM - poles(m)) X_vec is shape (N_orb, N_orb, num_poles) and is just an array of optimized X_l matrices, for given vallue of pole positions poles
- Parameters:
poles (numpy.ndarray) – pole positions
GM_matrix (numpy.ndarray) – Exact Matsubara Green’s function data in the shape (num_imag, num_orb, num_orb)
zM (numpy.ndarray) – Matsubara frequencies
- Returns:
Gradient of E with respect to poles
- Return type:
numpy.ndarray
- green_mbtools.pesto.pes_utils.cvx_matrix_projection(zM, GM_matrix, w_cut=10, n_real=201, ofile='Giw', solver='SCS', **kwargs)[source]¶
Projection of noisy Green’s function matrix data on to Nevanlinna manifold. Once this is performed, analytic continuation either using Nevanlinna or ES becomes easier. The main idea is to obtain a projected Matsubara Green’s function of the form
G(iw_m) = sum_n P_n / (iw_m - w_n)
- Parameters:
zM (numpy.ndarray) – n_imag number of imaginary frequency values
GM_matrix (numpy.ndarray) – Exact Matsubara Green’s function data in the shape (n_imag, n_orb, n_orb)
w_cut (float, optional) – cut-off to consider for real axis, by default 10.0
n_real (int, optional) – number of real frequencies to consider, by default 201
ofile (str, optional) – Output file in which the green’s function data will be dumped, by default ‘Giw’
solver (str, optional) – CVXPy solver to use, by default ‘SCS’
**kwargs (dict, optional) – dictionary of keyword arguments for CVXPy solver
- Returns:
Projected Green’s function
- Return type:
numpy.ndarray
- green_mbtools.pesto.pes_utils.cvx_optimize(poles, GM_matrix, zM, solver='SCS', **kwargs)[source]¶
Top level target error function in the ES approach for matrix analytic continuation. It returns
Err (poles) = min_{X} || Gimag (iw) - Gapprox [poles, X] (iw) ||
- Parameters:
poles (numpy.ndarray) – 1D array of pole positions
GM_matrix (numpy.ndarray) – Exact Matsubara Green’s function data in the shape (num_imag, num_orb, num_orb)
zM (numpy.ndarray) – 1D Matsubara frequencies
solver (str, optional) – CVXPy solver for semi-definite relaxation in the ES continuation, by default ‘SCS’
**kwargs (dict, optional) – dictionary of keyword arguments for CVXPy solver
- Returns:
float – Optimal value of error function Err (poles) for fixed poles
numpy.ndarray – optimal X_vec in the shape (num_poles, num_orb, num_orb)
numpy.ndarray – approximated optimal Green’s function for the specified Matsubara frequencies. The shape of the array is (num_imag, num_orb, num_orb)
- green_mbtools.pesto.pes_utils.cvx_optimize_spectral(poles, GM_diags, zM, solver='SCS', **kwargs)[source]¶
Top level target error function in the ES approach for spectral analytic continuation. It returns
Err (poles) = min_{X} || Diag(Gimag (iw) - Gapprox [poles, X] (iw)) ||
- Parameters:
poles (numpy.ndarray) – 1D array of pole positions
GM_diags (numpy.ndarray) – Exact Matsubara Green’s function diagonal data in the shape (num_imag, num_orb)
zM (numpy.ndarray) – 1D Matsubara frequencies
solver (str, optional) – CVXPy solver for semi-definite relaxation in the ES continuation, by default ‘SCS’
**kwargs (dict, optional) – dictionary of keyword arguments for CVXPy solver
- Returns:
float – Optimal value of error function Err (poles) for fixed poles
numpy.ndarray – optimal X_vec in the shape (num_poles, num_orb)
numpy.ndarray – approximated optimal Green’s function diagonals for the specified Matsubara frequencies. The shape of the array is (num_imag, num_orb)
- green_mbtools.pesto.pes_utils.run_es(iw_vals, G_iw, re_w_vals, diag=True, eta=0.01, eps_pol=1, ofile='Xw_real.txt', solver='SCS', **kwargs)[source]¶
Pole estimation and semi-definite (PES) relaation algorithm based on Huang, Gull and Lin, 10.1103/PhysRevB.107.075151. The output is saved in ofile.
- Parameters:
iw_vals (numpy.ndarray) – 1D array of imaginary frequency values (value only, without i)
G_iw (numpy.ndarray) – matsubara quantity to be analytically continued
re_w_vals (numpy.ndarray) – real frequency grid to perform continuation on
diag (bool, optional) – perform continuation for diagonal entries only if set to True, by default True
eta (float, optional) – broadening parameter, by default 0.01
eps_pol (float, optional) – maximum imag part of the poles to be considered, by default 1
ofile (str, optional) – Path or name of output file for storing the continued data
solver (str, optional) – CVXPy solver to use (options: SCS, MOSEK, CLARABEL, etc.), by default ‘SCS’
**kwargs (dict, optional) – dictionary of keyword arguments for CVXPy solver
green_mbtools.pesto.quasiparticle¶
green_mbtools.pesto.spectral¶
- green_mbtools.pesto.spectral.compute_mo(F, S=None)[source]¶
Solve the generalized eigen problem: FC = SCE
- Parameters:
F (numpy.ndarray) – Fock matrix, dim = (ns, nk, nao, nao)
S (numpy.ndarray, optional) – Overlap matrix, dim = (ns, nk, nao, nao), by default None
- Returns:
numpy.ndarray – molecular orbital energies
numpy.ndarray – molecular orbital coefficient in AO basis
- green_mbtools.pesto.spectral.compute_no(dm, S=None)[source]¶
Compute natural orbitals by diagonalizing density matrix
- Parameters:
dm (numpy.ndarray) – density matrix of shape (ns, nk, nao, nao)
- Returns:
numpy.ndarray – natural orbital occupations
numpy.ndarray – natural orbital coefficients / vectors
green_mbtools.pesto.winter¶
- green_mbtools.pesto.winter.interpolate(obj_k, kmesh, kpts_inter, dim=3, hermi=False, debug=False)[source]¶
Perform Wannier interpolation for static quantities
- Parameters:
obj_k (numpy.ndarray) – input object on full Brillouin zone k-mesh, of shape (ns, nk, nao, nao)
kmesh (numpy.ndarray) – original k-mesh (in scaled units)
kpts_inter (numpy.ndarray) – k-points to interpolate on (in scaled units)
dim (int, optional) – dimensionality of crystal (2D or 3D)
hermi (bool, optional) – set to True if the interpolated matrix Hermitian, by default False
debug (bool, optional) – set to True if debug information needs to be printed, by default False
- Returns:
Interpolated tensor
- Return type:
numpy.ndarray
- Raises:
ValueError – if dim other than 2 or 3 is specified
- green_mbtools.pesto.winter.interpolate_G(Fk, Sigma_tk, mu, Sk, kmesh, kpts_inter, ir, dim=3, hermi=False, debug=False)[source]¶
Interpolate Green’s function from full BZ k-mesh to specified k-points
- Parameters:
Fk (numpy.ndarray) – Fock matrix
Sigma_tk (numpy.ndarray) – Self-energy on imaginar time axis
mu (float) – Chemical potential
Sk (numpy.ndarray) – overlap matrix
kmesh (numpy.ndarray) – input k-mesh in the Brillouin zone
kpts_inter (numpy.ndarray) – interpolation k-points
ir (IR_factory) – handler for Fourier transforms between imaginary time and Matsubara frequencies
dim (int, optional) – dimensionality of latticej, by default 3
hermi (bool, optional) – Is Green’s function expected to be Hermitian, by default False
debug (bool, optional) – print extra messages for debugging, by default False
- Returns:
Green’s function interpolated on k-points
- Return type:
numpy.ndarray
- Raises:
ValueError – if dim other than 2 or 3 is specified
- green_mbtools.pesto.winter.interpolate_tk_object(obj_tk, kmesh, kpts_inter, dim=3, hermi=False, debug=False)[source]¶
Perform Wannier interpolation for dynamic quantities
- Parameters:
obj_k (numpy.ndarray) – input object on full Brillouin zone k-mesh, of shape (ntau, ns, nk, nao, nao)
kmesh (numpy.ndarray) – original k-mesh (in scaled units)
kpts_inter (numpy.ndarray) – k-points to interpolate on (in scaled units)
dim (int, optional) – dimensionality of crystal (2D or 3D)
hermi (bool, optional) – set to True if the interpolated matrix Hermitian, by default False
debug (bool, optional) – set to True if debug information needs to be printed, by default False
- Returns:
Interpolated tensor
- Return type:
numpy.ndarray
- Raises:
ValueError – if dim other than 2 or 3 is specified