gw_embedding.iaft

class gw_embedding.iaft.IAFT(beta: float, lmbda: float, prec: float = 1e-15, verbal: bool = True)[source]

Driver for FT on the imaginary axis. Given inverse temperature, lambda and precision, the IAFT class evaluate the corresponding IR basis and sparse sampling points on-the-fly.

Dependency:

sparse-ir with xprec supports (https://sparse-ir.readthedocs.io/en/latest/) To install sparse-ir with xprec supports: “pip install sparse-ir[xprec]”.

Attributes: beta: float

Inverse temperature (a.u.)

lmbda: float

Dimensionless lambda parameter for constructing the IR basis

prec: float

Precision for IR basis

bases: sparse-ir.FiniteTempBasisSet

IR basis instance

tau_mesh_f: numpy.ndarray(dim=1)

Fermionic tau sampling points

tau_mesh_b: numpy.ndarray(dim=1)

Bosonic tau sampling points

wn_mesh_f: numpy.ndarray(dim=1)

Fermionic Matsubara “indices” sampling points. NOT PHYSICAL FREQUENCIES. Physical Matsubara frequencies are wn_mesh_f * numpy.pi / beta

wn_mesh_b: numpy.ndarray(dim=1)

Bosonic Matsubara “indices” sampling points. NOT PHYSICAL FREQUENCIES. Physical Matsubara frequencies are wn_mesh_f * numpy.pi / beta

nt_f: int

Number of fermionic tau sampling points

nt_b: int

Number of bosonic tau sampling points

nw_f: int

Number of fermionic frequency sampling points

nw_b: int

Number of bosonic frequency sampling points

Methods

check_leakage(Ot, stats[, name, w_input])

Check decay of the IR coefficients to assess the quality of IR basis for the beta and lambda. The coefficients should decay exponentially, and the leakage is defined as: leakage = the smallest coefficients / the largest coefficients :param Ot: :param stats: :param name: :param w_input: :return:.

tau_interpolate(Ot, tau_mesh_interp, stats)

Interpolate a dynamic object to arbitrary points on the imaginary-time axis.

tau_interpolate_phsym(Ot, tau_mesh_interp, stats)

Interpolate a dynamic object to arbitrary points on the imaginary-time axis.

tau_to_w(Ot, stats)

Fourier transform from imaginary-time axis to Matsubara-frequency axis :param Ot: numpy.ndarray imaginary-time object with dimensions (nts, ...) :param stats: str statistics: 'f' for fermions and 'b' for bosons

tau_to_w_phsym(Ot, stats)

Fourier transform from imaginary-time axis to Matsubara-frequency axis w/ particle-hole symmetry :param Ot: numpy.ndarray imaginary-time object with dimensions (nts, ...) :param stats: str statistics: 'f' for fermions and 'b' for bosons

w_interpolate(Ow, wn_mesh_interp, stats[, ...])

Interpolate a dynamic object to arbitrary points on the Matsubara axis.

w_interpolate_phsym(Ow, wn_mesh_interp, stats)

Interpolate a dynamic object to arbitrary points on the Matsubara axis.

w_to_tau(Ow, stats)

Fourier transform from Matsubara-frequency axis to imaginary-time axis.

w_to_tau_phsym(Ow, stats)

Fourier transform from Matsubara-frequency axis to imaginary-time axis w/ particle-hole symmetry.

wn_mesh(stats[, ir_notation])

Return Matsubara frequency indices. :param stats: str statistics: 'f' for fermions and 'b' for bosons :param ir_notation: bool Whether wn_mesh_interp is in sparse_ir notation where iwn = n*pi/beta for both fermions and bosons. Otherwise, iwn = (2n+1)*pi/beta for fermions and 2n*pi/beta for bosons.

check_leakage(Ot, stats: str, name: str = '', w_input: bool = False)[source]

Check decay of the IR coefficients to assess the quality of IR basis for the beta and lambda. The coefficients should decay exponentially, and the leakage is defined as:

leakage = the smallest coefficients / the largest coefficients

Parameters:
  • Ot

  • stats

  • name

  • w_input

Returns:

tau_interpolate(Ot, tau_mesh_interp, stats: str)[source]

Interpolate a dynamic object to arbitrary points on the imaginary-time axis.

Parameters:
  • Ot – numpy.ndarray Dynamic object on the imaginary-time sampling points, self.tau_mesh.

  • tau_mesh_interp – numpy.ndarray(dim=1, dtype=float) Target tau points.

  • stats – str Statistics, ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Imaginary-time object with dimensions (nt_interp, …)

tau_interpolate_phsym(Ot, tau_mesh_interp, stats: str)[source]

Interpolate a dynamic object to arbitrary points on the imaginary-time axis.

Parameters:
  • Ot – numpy.ndarray Dynamic object on the imaginary-time sampling points, self.tau_mesh.

  • tau_mesh_interp – numpy.ndarray(dim=1, dtype=float) Target tau points.

  • stats – str Statistics, ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Imaginary-time object with dimensions (nt_interp, …)

tau_to_w(Ot, stats: str)[source]

Fourier transform from imaginary-time axis to Matsubara-frequency axis :param Ot: numpy.ndarray

imaginary-time object with dimensions (nts, …)

Parameters:

stats – str statistics: ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Matsubara-frequency object with dimensions (nw, …)

tau_to_w_phsym(Ot, stats: str)[source]

Fourier transform from imaginary-time axis to Matsubara-frequency axis w/ particle-hole symmetry :param Ot: numpy.ndarray

imaginary-time object with dimensions (nts, …)

Parameters:

stats – str statistics: ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Matsubara-frequency object with dimensions (nw, …)

w_interpolate(Ow, wn_mesh_interp, stats: str, ir_notation: bool = True)[source]

Interpolate a dynamic object to arbitrary points on the Matsubara axis.

Parameters:
  • Ow – numpy.ndarray Dynamic object on the Matsubara sampling points, self.wn_mesh.

  • wn_mesh_interp – numpy.ndarray(dim=1, dtype=int) Target frequencies “INDICES”. The physical Matsubara frequencies are wn_mesh_interp * pi/beta.

  • stats – str Statistics, ‘f’ for fermions and ‘b’ for bosons.

  • ir_notation – bool Whether wn_mesh_interp is in sparse_ir notation where iwn = n*pi/beta for both fermions and bosons. Otherwise, iwn = (2n+1)*pi/beta for fermions and 2n*pi/beta for bosons.

Returns:

numpy.ndarray Matsubara-frequency object with dimensions (nw_interp, …)

w_interpolate_phsym(Ow, wn_mesh_interp, stats: str, ir_notation: bool = True)[source]

Interpolate a dynamic object to arbitrary points on the Matsubara axis.

Parameters:
  • Ow – numpy.ndarray Dynamic object on the Matsubara sampling points, self.wn_mesh.

  • wn_mesh_interp – numpy.ndarray(dim=1, dtype=int) Target frequencies “INDICES”. The physical Matsubara frequencies are wn_mesh_interp * pi/beta.

  • stats – str Statistics, ‘f’ for fermions and ‘b’ for bosons.

  • ir_notation – bool Whether wn_mesh_interp is in sparse_ir notation where iwn = n*pi/beta for both fermions and bosons. Otherwise, iwn = (2n+1)*pi/beta for fermions and 2n*pi/beta for bosons.

Returns:

numpy.ndarray Matsubara-frequency object with dimensions (nw_interp, …)

w_to_tau(Ow, stats)[source]

Fourier transform from Matsubara-frequency axis to imaginary-time axis.

Parameters:
  • Ow – numpy.ndarray Matsubara-frequency object with dimensions (nw, …)

  • stats – str statistics, ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Imaginary-time object with dimensions (nt, …)

w_to_tau_phsym(Ow, stats)[source]

Fourier transform from Matsubara-frequency axis to imaginary-time axis w/ particle-hole symmetry.

Parameters:
  • Ow – numpy.ndarray Matsubara-frequency object with dimensions (nw, …)

  • stats – str statistics, ‘f’ for fermions and ‘b’ for bosons

Returns:

numpy.ndarray Imaginary-time object with dimensions (nt, …)

wn_mesh(stats: str, ir_notation: bool = True)[source]

Return Matsubara frequency indices. :param stats: str

statistics: ‘f’ for fermions and ‘b’ for bosons

Parameters:

ir_notation – bool Whether wn_mesh_interp is in sparse_ir notation where iwn = n*pi/beta for both fermions and bosons. Otherwise, iwn = (2n+1)*pi/beta for fermions and 2n*pi/beta for bosons.

Returns:

numpy.ndarray(dim=1) Matsubara frequency indices

Classes

IAFT(beta, lmbda[, prec, verbal])

Driver for FT on the imaginary axis.