Package reference#

Interactive problem definition#

To define the interactive problem, we use the following class:

class meanfi.model.Model(h_0, h_int, filling)[source]#

Data class which defines the interacting tight-binding problem.

Parameters:
  • h_0 (dict[tuple[int, ...], ndarray]) – Non-interacting hermitian Hamiltonian tight-binding dictionary.

  • h_int (dict[tuple[int, ...], ndarray]) – Interaction hermitian Hamiltonian tight-binding dictionary.

  • filling (float) – Number of particles in a unit cell. Used to determine the Fermi level.

Notes

The interaction h_int must be of density-density type. For example, h_int[(1,)][i, j] = V means a repulsive interaction of strength V between two particles with internal degrees of freedom i and j separated by 1 lattice vector.

density_matrix(rho, nk=20)[source]#

Computes the density matrix from a given initial density matrix.

Parameters:
  • rho (dict[tuple[int, ...], ndarray]) – Initial density matrix tight-binding dictionary.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

Returns:

Density matrix tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

mfield(mf, nk=20)[source]#

Computes a new mean-field correction from a given one.

Parameters:
  • mf (dict[tuple[int, ...], ndarray]) – Initial mean-field correction tight-binding dictionary.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

Returns:

new mean-field correction tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

Mean-field and density matrix#

meanfi.mf.density_matrix(h, filling, nk)[source]#

Compute the real-space density matrix tight-binding dictionary.

Parameters:
  • h (dict[tuple[int, ...], ndarray]) – Hamiltonian tight-binding dictionary from which to construct the density matrix.

  • filling (float) – Number of particles in a unit cell. Used to determine the Fermi level.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

Returns:

Density matrix tight-binding dictionary and Fermi energy.

Return type:

Tuple[dict[tuple[int, …], ndarray], float]

meanfi.mf.density_matrix_kgrid(kham, filling)[source]#

Calculate density matrix on a k-space grid.

Parameters:
  • kham (ndarray) – Hamiltonian from which to construct the density matrix. The hamiltonian is sampled on a grid of k-points and has shape (nk, nk, …, ndof, ndof), where ndof is number of internal degrees of freedom.

  • filling (float) – Number of particles in a unit cell. Used to determine the Fermi level.

Returns:

Density matrix on a k-space grid with shape (nk, nk, …, ndof, ndof) and Fermi energy.

Return type:

Tuple[ndarray, float]

meanfi.mf.fermi_on_kgrid(vals, filling)[source]#

Compute the Fermi energy on a grid of k-points.

Parameters:
  • vals (ndarray) – Eigenvalues of a hamiltonian sampled on a k-point grid with shape (nk, nk, …, ndof, ndof), where ndof is number of internal degrees of freedom.

  • filling (float) – Number of particles in a unit cell. Used to determine the Fermi level.

Returns:

Fermi energy

Return type:

float

meanfi.mf.meanfield(density_matrix, h_int)[source]#

Compute the mean-field correction from the density matrix.

Parameters:
  • density_matrix (dict[tuple[int, ...], ndarray]) – Density matrix tight-binding dictionary.

  • h_int (dict[tuple[int, ...], ndarray]) – Interaction hermitian Hamiltonian tight-binding dictionary.

Returns:

Mean-field correction tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

Notes

The interaction h_int must be of density-density type. For example, h_int[(1,)][i, j] = V means a repulsive interaction of strength V between two particles with internal degrees of freedom i and j separated by 1 lattice vector.

Observables#

meanfi.observables.expectation_value(density_matrix, observable)[source]#

Compute the expectation value of an observable with respect to a density matrix.

Parameters:
Returns:

Expectation value.

Return type:

complex

Solvers#

meanfi.solvers.cost_density(rho_params, model, nk=20)[source]#

Defines the cost function for root solver.

The cost function is the difference between the computed and inputted density matrix reduced to the hoppings only present in the h_int.

Parameters:
  • rho_params (ndarray) – 1D real array that parametrises the density matrix reduced to the hoppings (keys) present in h_int.

  • Model – Interacting tight-binding problem definition.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

  • model (Model)

Returns:

1D real array that is the difference between the computed and inputted density matrix parametrisations reduced to the hoppings present in h_int.

Return type:

ndarray

meanfi.solvers.cost_mf(mf_param, model, nk=20)[source]#

Defines the cost function for root solver.

The cost function is the difference between the computed and inputted mean-field.

Parameters:
  • mf_param (ndarray) – 1D real array that parametrises the mean-field correction.

  • Model – Interacting tight-binding problem definition.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

  • model (Model)

Returns:

1D real array that is the difference between the computed and inputted mean-field parametrisations

Return type:

ndarray

meanfi.solvers.solver(model, mf_guess, nk=20, optimizer=<function anderson>, optimizer_kwargs={'M': 0, 'line_search': 'wolfe'})#

Solve for the mean-field correction through self-consistent root finding by finding the density matrix fixed point.

Parameters:
  • model (Model) – Interacting tight-binding problem definition.

  • mf_guess (dict[tuple[int, ...], ndarray]) – The initial guess for the mean-field correction in the tight-binding dictionary format.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

  • optimizer (Callable | None) – The solver used to solve the fixed point iteration. Default uses scipy.optimize.anderson.

  • optimizer_kwargs (dict[str, str] | None) – The keyword arguments to pass to the optimizer.

Returns:

Mean-field correction solution in the tight-binding dictionary format.

Return type:

dict[tuple[int, …], ndarray]

meanfi.solvers.solver_mf(model, mf_guess, nk=20, optimizer=<function anderson>, optimizer_kwargs={'M': 0})[source]#

Solve for the mean-field correction through self-consistent root finding by finding the mean-field correction fixed point.

Parameters:
  • model (Model) – Interacting tight-binding problem definition.

  • mf_guess (ndarray) – The initial guess for the mean-field correction in the tight-binding dictionary format.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

  • optimizer (Callable | None) – The solver used to solve the fixed point iteration. Default uses scipy.optimize.anderson.

  • optimizer_kwargs (dict[str, str] | None) – The keyword arguments to pass to the optimizer.

Returns:

Mean-field correction solution in the tight-binding dictionary format.

Return type:

dict[tuple[int, …], ndarray]

Tight-binding dictionary#

Manipulation#

meanfi.tb.tb.add_tb(tb1, tb2)[source]#

Add up two tight-binding dictionaries together.

Parameters:
Returns:

Sum of the two tight-binding dictionaries.

Return type:

dict[tuple[int, …], ndarray]

meanfi.tb.tb.scale_tb(tb, scale)[source]#

Scale a tight-binding dictionary by a constant.

Parameters:
  • tb (dict[tuple[int, ...], ndarray]) – Tight-binding dictionary.

  • scale (float) – Constant to scale the tight-binding dictionary by.

Returns:

Scaled tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

Brillouin zone transformations#

meanfi.tb.transforms.ifftn_to_tb(ifft_array)[source]#

Convert the result of scipy.fft.ifftn to a tight-binding dictionary.

Parameters:

ifft_array (ndarray) – Result of scipy.fft.ifftn to convert to a tight-binding dictionary. The input to scipy.fft.ifftn should be from tb_to_khamvector.

Returns:

Tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

meanfi.tb.transforms.kgrid_to_tb(kgrid_array)[source]#

Convert a k-space grid array to a tight-binding dictionary.

Parameters:

kgrid_array (ndarray) – K-space grid array to convert to a tight-binding dictionary. The array should be of shape (nk, nk, …, ndof, ndof), where ndof is number of internal degrees of freedom.

Returns:

Tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

meanfi.tb.transforms.tb_to_kfunc(tb)[source]#

Fourier transforms a real-space tight-binding model to a k-space function.

Parameters:

tb (dict[tuple[int, ...], ndarray]) – Tight-binding dictionary.

Returns:

A function that takes a k-space vector and returns a complex np.array.

Return type:

Callable

Notes

Function doesn’t work for zero dimensions.

meanfi.tb.transforms.tb_to_kgrid(tb, nk)[source]#

Evaluate a tight-binding dictionary on a k-space grid.

Parameters:
  • tb (dict[tuple[int, ...], ndarray]) – Tight-binding dictionary to evaluate on a k-space grid.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

Returns:

Tight-binding dictionary evaluated on a k-space grid. Has shape (nk, nk, …, ndof, ndof), where ndof is number of internal degrees of freedom.

Return type:

ndarray

Parametrisation#

meanfi.params.rparams.complex_to_real(z)[source]#

Split and concatenate real and imaginary parts of a complex array.

Parameters:

z (ndarray) – Complex array.

Returns:

Real array that concatenates the real and imaginary parts of the input array.

Return type:

ndarray

meanfi.params.rparams.real_to_complex(z)[source]#

Undo complex_to_real.

Parameters:

z (ndarray)

Return type:

ndarray

meanfi.params.rparams.rparams_to_tb(tb_params, tb_keys, ndof)[source]#

Extract a hermitian tight-binding dictionary from a real vector parametrisation.

Parameters:
  • tb_params (ndarray) – 1D real array that parametrises the tight-binding dictionary.

  • tb_keys (list[tuple[None] | tuple[int, ...]]) – List of keys of the tight-binding dictionary.

  • ndof (int) – Number internal degrees of freedom within the unit cell.

Returns:

Tight-biding dictionary.

Return type:

dict[tuple[int, …], ndarray]

meanfi.params.rparams.tb_to_rparams(tb)[source]#

Parametrise a hermitian tight-binding dictionary by a real vector.

Parameters:

tb (dict[tuple[int, ...], ndarray]) – tight-binding dictionary.

Returns:

1D real vector that parametrises the tight-binding dictionary.

Return type:

ndarray

Utility functions#

meanfi.tb.utils.fermi_energy(tb, filling, nk=100)[source]#

Calculate the Fermi energy of a given tight-binding dictionary.

Parameters:
  • tb (dict[tuple[int, ...], ndarray]) – Tight-binding dictionary.

  • filling (float) – Number of particles in a unit cell. Used to determine the Fermi level.

  • nk (int) – Number of k-points in a grid to sample the Brillouin zone along each dimension. If the system is 0-dimensional (finite), this parameter is ignored.

Returns:

Fermi energy.

meanfi.tb.utils.generate_tb_keys(cutoff, dim)[source]#

Generate tight-binding dictionary keys up to a cutoff.

Parameters:
  • cutoff (int) – Maximum distance along each dimension to generate tight-bindign dictionary keys for.

  • dim (int) – Dimension of the tight-binding dictionary.

Returns:

List of generated tight-binding dictionary keys up to a cutoff.

Return type:

list[tuple[None] | tuple[int, …]]

meanfi.tb.utils.guess_tb(tb_keys, ndof, scale=1)[source]#

Generate hermitian guess tight-binding dictionary.

Parameters:
  • tb_keys (list[tuple[None] | tuple[int, ...]]) – List of hopping vectors (tight-binding dictionary keys) the guess contains.

  • ndof (int) – Number internal degrees of freedom within the unit cell.

  • scale (float) – Scale of the random guess.

Returns:

Hermitian guess tight-binding dictionary.

Return type:

dict[tuple[int, …], ndarray]

kwant interface#

meanfi.kwant_helper.utils.build_interacting_syst(builder, lattice, func_onsite, func_hop=None, max_neighbor=1)[source]#

Construct an auxiliary kwant system that encodes the interactions.

Parameters:
  • builder (Builder) – Non-interacting kwant.builder.Builder system.

  • lattice (Polyatomic) – Lattice of the system.

  • func_onsite (Callable) – Onsite interactions function.

  • func_hop (Callable | None) – Hopping/inter unit cell interactions function.

  • max_neighbor (int) – The maximal number of neighbouring unit cells (along a lattice vector) connected by interaction. Interaction goes to zero after this distance.

Returns:

Auxiliary kwant.builder.Builder that encodes the interactions of the system.

Return type:

Builder

meanfi.kwant_helper.utils.builder_to_tb(builder, params={}, return_data=False)[source]#

Construct a tight-binding dictionary from a kwant.builder.Builder system.

Parameters:
  • builder (Builder) – system to convert to tight-binding dictionary.

  • params (dict) – Dictionary of parameters to evaluate the builder on.

  • return_data (bool) – Returns dictionary with sites and number of orbitals per site.

Returns:

  • Tight-binding dictionary that corresponds to the builder.

  • Data with sites and number of orbitals. Only if return_data=True.

Return type:

dict[tuple[int, …], ndarray]

meanfi.kwant_helper.utils.tb_to_builder(h_0, sites_list, periods)[source]#

Construct a kwant.builder.Builder from a tight-binding dictionary.

Parameters:
  • h_0 (dict[tuple[int, ...], ndarray]) – Tight-binding dictionary.

  • sites_list (list[Site, ...]) – List of sites in the builder’s unit cell.

  • periods (ndarray) – 2d array with periods of the translational symmetry.

Returns:

kwant.builder.Builder that corresponds to the tight-binding dictionary.

Return type:

Builder