Wave functions

class gpaw.wavefunctions.base.WaveFunctions(gd, nvalence, setups, bd, dtype, collinear, world, kd, kptband_comm, timer)[source]

setups:

List of setup objects.

symmetry:

Symmetry object.

kpt_u:

List of k-point objects.

nbands: int

Number of bands.

nspins: int

Number of spins.

dtype: dtype

Data type of wave functions (float or complex).

bzk_kc: ndarray

Scaled k-points used for sampling the whole Brillouin zone - values scaled to [-0.5, 0.5).

ibzk_kc: ndarray

Scaled k-points in the irreducible part of the Brillouin zone.

weight_k: ndarray

Weights of the k-points in the irreducible part of the Brillouin zone (summing up to 1).

kpt_comm:

MPI-communicator for parallelization over k-points.

calculate_atomic_density_matrices(D_asp)[source]

Calculate atomic density matrices from projections.

calculate_atomic_density_matrices_with_occupation(D_asp, f_un)[source]

Calculate atomic density matrices from projections with custom occupation f_un.

calculate_density_contribution(nt_sG)[source]

Calculate contribution to pseudo density from wave functions.

Array entries are written to (not added to).

collect_array(name, k, s, subset=None)[source]

Helper method for collect_eigenvalues and collect_occupations.

For the parallel case find the rank in kpt_comm that contains the (k,s) pair, for this rank, collect on the corresponding domain a full array on the domain master and send this to the global master.

collect_auxiliary(value, k, s, shape=1, dtype=<class 'float'>)[source]

Helper method for collecting band-independent scalars/arrays.

For the parallel case find the rank in kpt_comm that contains the (k,s) pair, for this rank, collect on the corresponding domain a full array on the domain master and send this to the global master.

collect_projections(k, s)[source]

Helper method for collecting projector overlaps across domains.

For the parallel case find the rank in kpt_comm that contains the (k,s) pair, for this rank, send to the global master.

get_homo_lumo(spin=None)[source]

Return HOMO and LUMO eigenvalues.

get_orbital_density_matrix(a, kpt, n)[source]

Add the nth band density from kpt to density matrix D_sp

get_wave_function_array(n, k, s, realspace=True, periodic=False)[source]

Return pseudo-wave-function array on master.

n: int

Global band index.

k: int

Global IBZ k-point index.

s: int

Spin index (0 or 1).

realspace: bool

Transform plane wave or LCAO expansion coefficients to real-space.

For the parallel case find the ranks in kd.comm and bd.comm that contains to (n, k, s), and collect on the corresponding domain a full array on the domain master and send this to the global master.

class gpaw.wavefunctions.fd.FDWaveFunctions(stencil, parallel, initksl, gd, nvalence, setups, bd, dtype, world, kd, kptband_comm, timer, reuse_wfs_method=None, collinear=True)[source]
ibz2bz(atoms)[source]

Transform wave functions in IBZ to the full BZ.

random_wave_functions(nao)[source]

Generate random wave functions.

class gpaw.wavefunctions.pw.PWWaveFunctions(ecut, fftwflags, dedepsilon, parallel, initksl, reuse_wfs_method, collinear, gd, nvalence, setups, bd, dtype, world, kd, kptband_comm, timer)[source]
apply_pseudo_hamiltonian(kpt, ham, psit_xG, Htpsit_xG)[source]

Apply the pseudo Hamiltonian i.e. without PAW corrections.

get_wave_function_array(n, k, s, realspace=True, cut=True, periodic=False)[source]

Return pseudo-wave-function array on master.

n: int

Global band index.

k: int

Global IBZ k-point index.

s: int

Spin index (0 or 1).

realspace: bool

Transform plane wave or LCAO expansion coefficients to real-space.

For the parallel case find the ranks in kd.comm and bd.comm that contains to (n, k, s), and collect on the corresponding domain a full array on the domain master and send this to the global master.

initialize_from_lcao_coefficients(basis_functions, block_size=10)[source]

Convert from LCAO to PW coefficients.

class gpaw.wavefunctions.pw.PW(ecut=340, *, fftwflags=0, cell=None, gammacentered=False, pulay_stress=None, dedecut=None, force_complex_dtype=False, interpolation='fft')[source]

Plane-wave basis mode.

Only one of dedecut and pulay_stress can be used.

Parameters:
  • ecut (float) – Plane-wave cutoff in eV.

  • gammacentered (bool) – Center the grid of chosen plane waves around the gamma point or q/k-vector

  • dedecut (float | str | None) – Estimate of derivative of total energy with respect to plane-wave cutoff. Used to calculate the Pulay-stress.

  • pulay_stress (float or None) – Pulay-stress correction.

  • fftwflags (int) –

    Flags for making an FFTW plan. There are 4 possibilities (default is MEASURE):

    from gpaw.fftw import ESTIMATE, MEASURE, PATIENT, EXHAUSTIVE
    

  • cell (np.ndarray) – Use this unit cell to chose the plane-waves.

  • interpolation (str or int) – Interpolation scheme to construct the density on the fine grid. Default is 'fft' and alternatively a stencil (integer) can be given to perform an explicit real-space interpolation.

class gpaw.wavefunctions.lcao.LCAOWaveFunctions(ksl, gd, nvalence, setups, bd, dtype, world, kd, kptband_comm, timer, atomic_correction=None, collinear=True)[source]
add_to_density_from_k_point_with_occupation(nt_sG, kpt, f_n)[source]

Add contribution to pseudo electron-density. Do not use the standard occupation numbers, but ones given with argument f_n.

calculate_atomic_density_matrices_with_occupation(D_asp, f_un)[source]

Calculate atomic density matrices from projections with custom occupation f_un.

initialize_wave_functions_from_lcao()[source]

Fill the calc.wfs.kpt_[u].psit_nG arrays with useful data.

Normally psit_nG is NOT used in lcao mode, but some extensions (like ase.dft.wannier) want to have it. This code is adapted from fd.py / initialize_from_lcao_coefficients() and fills psit_nG with data constructed from the current lcao coefficients (kpt.C_nM).

(This may or may not work in band-parallel case!)

initialize_wave_functions_from_restart_file()[source]

Dummy function to ensure compatibility to fd mode

class gpaw.kpoint.KPoint(weightk, weight, s, k, q, phase_cd=None)[source]

Class for a single k-point.

The KPoint class takes care of all wave functions for a certain k-point and a certain spin.

Attributes:

eps_n: float ndarray

Eigenvalues.

f_n: float ndarray

Occupation numbers. The occupation numbers already include the k-point weight in supercell calculations.

Construct k-point object.

Parameters:

weightk:

Weight of this k-point.

weight:

Old confusing weight.

s: int

Spin index: up or down (0 or 1).

k: int

k-point IBZ-index.

q: int

local k-point index.

phase_cd: complex ndarray

Bloch phase-factors for translations - axis c=0,1,2 and direction d=0,1.

class gpaw.eigensolvers.eigensolver.Eigensolver(keep_htpsit=True, blocksize=1)[source]
calculate_residuals(kpt, wfs, ham, psit, P, eps_n, R, C, n_x=None, calculate_change=False)[source]

Calculate residual.

From R=Ht*psit calculate R=H*psit-eps*S*psit.

iterate(ham, wfs)[source]

Solves eigenvalue problem iteratively

This method is inherited by the actual eigensolver which should implement iterate_one_k_point method for a single iteration of a single kpoint.

iterate_one_k_point(ham, kpt)[source]

Implemented in subclasses.

subspace_diagonalize(ham, wfs, kpt, rotate_psi=True)[source]

Diagonalize the Hamiltonian in the subspace of kpt.psit_nG

Htpsit_nG is a work array of same size as psit_nG which contains the local part of the Hamiltonian times psit on exit

First, the Hamiltonian (defined by kin, vt_sG, and dH_asp) is applied to the wave functions, then the H_nn matrix is calculated and diagonalized, and finally, the wave functions (and also Htpsit_nG are rotated. Also the projections P_ani are rotated.

It is assumed that the wave functions psit_nG are orthonormal and that the integrals of projector functions and wave functions P_ani are already calculated.

Return rotated wave functions and H applied to the rotated wave functions if self.keep_htpsit is True.

weights(wfs)[source]

Calculate convergence weights for all eigenstates.