Miscellaneous objects and functions

class gpaw.lfc.LocalizedFunctionsCollection(gd, spline_aj, kd=None, cut=False, dtype=<class 'float'>, integral=None, forces=None, xp=<module 'numpy' from '/home/docs/checkouts/readthedocs.org/user_builds/gpaw/envs/latest/lib/python3.12/site-packages/numpy/__init__.py'>)[source]

Utilizes that localized functions can be stored on a spherical subset of the uniform grid, as opposed to LocalizedFunctionsCollection which is just a wrapper around the old localized_functions which use rectangular grids.

add(a_xG, c_axi=1.0, q=-1)[source]

Add localized functions to extended arrays.

         --  a     a
a (G) += >  c   Phi (G)
 x       --  xi    i
         a,i
add_derivative(a, v, a_xG, c_axi=1.0, q=-1)[source]

Add derivative of localized functions on atom to extended arrays.

Parameters:

a: int

Atomic index of the derivative

v: int

Cartesian coordinate of the derivative (0, 1 or 2)

This function adds the following sum to the extended arrays:

         --  a      a
a (G) += >  c   dPhi  (G)
 x       --  xi     iv
          i

where:

    a        d     a
dPhi  (G) =  -- Phi (g)
    iv       dv    i

is the derivative of the Phi^a and v is either x, y, or z.

derivative(a_xG, c_axiv, q=-1)[source]

Calculate x-, y-, and z-derivatives of localized function integrals.

          /              a*
c_axiv =  | dG a (G) dPhi  (G)
          /     x        iv

where:

    a        d     a
dPhi  (G) =  -- Phi (g)
    iv       dv    i

and v is either x, y, or z, and R^a_v is the center of Phi^a.

Notice that d Phi^a_i / dR^a_v == - d Phi^a_i / d v.

griditer()[source]

Iterate over grid points.

integrate(a_xG, c_axi, q=-1, add_to=False)[source]

Calculate integrals of arrays times localized functions.

         /             a*
c_axi =  | dG a (G) Phi  (G)
         /     x       i
second_derivative(a_xG, c_axivv, q=-1)[source]

Calculate second derivatives.

Works only for this type of input for now:

second_derivative(self, a_G, c_avv, q=-1)
                    2 a _ _a
         /  _   _  d f (r-R )
c_avv =  | dr a(r) ----------
         /             a  a
                     dR dR
                       i  j
set_positions(spos_ac, atom_partition=None)[source]

Set positions and return True if any atoms have migrated to another rank.

class gpaw.lfc.BasisFunctions(gd, spline_aj, kd=None, cut=False, dtype=<class 'float'>, integral=None, forces=None, xp=<module 'numpy' from '/home/docs/checkouts/readthedocs.org/user_builds/gpaw/envs/latest/lib/python3.12/site-packages/numpy/__init__.py'>)[source]
add_to_density(nt_sG, f_asi)[source]

Add linear combination of squared localized functions to density.

\[\tilde{n}_{s} (r) += \sum^{}_{a,i} f^{a}_{si} [Φ^{a}_{i}(r)]^{2}\]
calculate_force_contribution(vt_G, rhoT_MM, q)[source]

Calculate derivatives of potential matrix elements.

            /     *  _
           |   Phi  (r)
 ~c        |      mu    ~ _        _   _
DV      += |   -------- v(r) Phi  (r) dr
  mu nu    |     dr             nu
          /        c

Results are added to DVt_vMM.

calculate_potential_matrices(vt_G)[source]

Calculate lower part of potential matrix.

          /
~         |     *  _  ~ _        _   _
V      =  |  Phi  (r) v(r) Phi  (r) dr    for  mu >= nu
 mu nu    |     mu            nu
          /

Overwrites the elements of the target matrix Vt_MM.

calculate_potential_matrix(vt_G, Vt_MM, q)[source]

Calculate lower part of potential matrix.

          /
~         |     *  _  ~ _        _   _
V      =  |  Phi  (r) v(r) Phi  (r) dr    for  mu >= nu
 mu nu    |     mu            nu
          /

Overwrites the elements of the target matrix Vt_MM.

calculate_potential_matrix_derivative(vt_G, DVt_vMM, q)[source]

Calculate derivatives of potential matrix elements.

            /     *  _
           |   Phi  (r)
 ~c        |      mu    ~ _        _   _
DV      += |   -------- v(r) Phi  (r) dr
  mu nu    |     dr             nu
          /        c

Results are added to DVt_vMM.

construct_density(rho_MM, nt_G, q)[source]

Calculate electron density from density matrix.

rho_MM: ndarray

Density matrix.

nt_G: ndarray

Pseudo electron density.

~        --      *
n(r) +=  >    Phi (r) rho     Phi (r)
         --     M1       M1M2   M2
        M1,M2
integrate2(a_xG, c_xM, q=-1)[source]

Calculate integrals of arrays times localized functions.

                     /       *
c_xM += <Phi | a > = | dG Phi (G) a (G)
            M   x    /       M     x
lcao_to_grid(C_xM, psit_xG, q, block_size=10)[source]

Deploy basis functions onto grids according to coefficients.

           ----
 ~   _     \                 _
psi (r) +=  )    C     Phi  (r)
   n       /      n mu    mu
           ----
            mu
set_positions(spos_ac)[source]

Set positions and return True if any atoms have migrated to another rank.

class gpaw.spline.Spline(spline)[source]

Spline object

classmethod from_data(l, rmax, f_g)[source]

The integer l gives the angular momentum quantum number and the list contains the spline values from r=0 to r=rcut.

The array f_g gives the radial part of the function on the grid. The radial function is multiplied by a real solid spherical harmonics (r^l * Y_lm).

get_angular_momentum_number()[source]

Return the angular momentum quantum number.

get_cutoff()[source]

Return the radial cutoff.

get_value_and_derivative(r)[source]

Return the value and derivative.

map(r_x)[source]

Map f(r) onto a given radial grid.

class gpaw.poisson.FDPoissonSolver(nn=3, relax='J', eps=2e-10, maxiter=1000, remove_moment=None, use_charge_center=False, metallic_electrodes=False, use_charged_periodic_corrections=False, **kwargs)[source]
create_laplace(gd, scale=1.0, n=1, dtype=<class 'float'>)[source]

Instantiate and return a Laplace operator

Allows subclasses to change the Laplace operator

iterate2(step, level=0)[source]

Smooths the solution in every multigrid level

class gpaw.xc.functional.XCFunctional(name, type)[source]
calculate(gd, n_sg, v_sg=None, e_g=None)[source]

Calculate energy and potential.

gd: GridDescriptor

Descriptor for 3-d grid.

n_sg: rank-4 ndarray

Spin densities.

v_sg: rank-4 ndarray

Array for potential. The XC potential is added to the values already there.

e_g: rank-3 ndarray

Energy density. Values must be written directly, not added.

The total XC energy is returned.

get_description()[source]

Get long description of functional as a string, or None.

summary(fd)[source]

Write summary of last calculation to file.

todict()[source]

Get dictionary representation of XC functional.

This representation works for libxc kernels; other classes should likely override this function and should probably not rely on this implementation.

tostring()[source]

Get string representation of XC functional.

This will give the name for libxc functionals but other data for hybrids.

class gpaw.xc.gga.GGA(kernel, stencil=2, xp=<module 'numpy' from '/home/docs/checkouts/readthedocs.org/user_builds/gpaw/envs/latest/lib/python3.12/site-packages/numpy/__init__.py'>)[source]
get_description()[source]

Get long description of functional as a string, or None.

todict()[source]

Get dictionary representation of XC functional.

This representation works for libxc kernels; other classes should likely override this function and should probably not rely on this implementation.

gpaw.forces.calculate_forces(wfs, dens, ham, log=None)[source]

Return the atomic forces.

class gpaw.grid_descriptor.GridDescriptor(N_c, cell_cv=[1, 1, 1], pbc_c=True, comm=None, parsize_c=None, allow_empty_domains=False)[source]

Descriptor-class for uniform 3D grid

A GridDescriptor object holds information on how functions, such as wave functions and electron densities, are discreticed in a certain domain in space. The main information here is how many grid points are used in each direction of the unit cell.

There are methods for tasks such as allocating arrays, performing symmetry operations and integrating functions over space. All methods work correctly also when the domain is parallelized via domain decomposition.

This is how a 2x2x2 3D array is laid out in memory:

3-----7
|\    |\
| \   | \
|  1-----5      z
2--|--6  |   y  |
 \ |   \ |    \ |
  \|    \|     \|
   0-----4      +-----x

Example

>>> a = np.zeros((2, 2, 2))
>>> a.ravel()[:] = range(8)
>>> a
array([[[0., 1.],
        [2., 3.]],

       [[4., 5.],
        [6., 7.]]])

Construct grid-descriptor object.

parameters:

N_c: 3 ints

Number of grid points along axes.

cell_cv: 3 float’s or 3x3 floats

Unit cell.

pbc_c: one or three bools

Periodic boundary conditions flag(s).

comm: MPI-communicator

Communicator for domain-decomposition.

parsize_c: tuple of 3 ints, a single int or None

Number of domains.

allow_empty_domains: bool

Allow parallelization that would generate empty domains.

Note that if pbc_c[c] is False, then the actual number of gridpoints along axis c is one less than N_c[c].

Attributes:

dv

Volume per grid point.

h_cv

Array of the grid spacing along the three axes.

N_c

Array of the number of grid points along the three axes.

n_c

Number of grid points on this CPU.

beg_c

Beginning of grid-point indices (inclusive).

end_c

End of grid-point indices (exclusive).

comm

MPI-communicator for domain decomposition.

The length unit is Bohr.

bytecount(dtype=<class 'float'>)[source]

Get the number of bytes used by a grid of specified dtype.

calculate_dipole_moment(rho_g, center=False, origin_c=None)[source]

Calculate dipole moment of density.

coarsen()[source]

Return coarsened \(GridDescriptor\) object.

Reurned descriptor has 2x2x2 fewer grid points.

collect(a_xg, out=None, broadcast=False)[source]

Collect distributed array to master-CPU or all CPU’s.

coords(c, pad=True)[source]

Return coordinates along one of the three axes.

Useful for plotting:

import matplotlib.pyplot as plt
plt.plot(gd.coords(0), data[:, 0, 0])
plt.show()
dipole_moment(rho_R, center_v=None)[source]

Calculate dipole moment of density.

Integration region will be centered on center_v. Default center is center of unit cell.

distribute(B_xg, out=None)[source]

Distribute full array B_xg to subdomains, result in b_xg.

B_xg is not used by the slaves (i.e. it should be None on all slaves) b_xg must be allocated on all nodes and will be overwritten.

empty(n=(), dtype=<class 'float'>, global_array=False, pad=False, xp=<module 'numpy' from '/home/docs/checkouts/readthedocs.org/user_builds/gpaw/envs/latest/lib/python3.12/site-packages/numpy/__init__.py'>)[source]

Return new uninitialized 3D array for this domain.

The type can be set with the dtype keyword (default: float). Extra dimensions can be added with n=dim. A global array spanning all domains can be allocated with global_array=True.

find_center(a_R)[source]

Calculate center of positive function.

get_boxes(spos_c, rcut, cut=True)[source]

Find boxes enclosing sphere.

get_grid_point_coordinates(dtype=<class 'float'>, global_array=False)[source]

Construct cartesian coordinates of grid points in the domain.

get_grid_point_distance_vectors(r_v, mic=True, dtype=<class 'float'>)[source]

Return distances to a given vector in the domain.

mic: if true adopts the mininimum image convention procedure by W. Smith in ‘The Minimum image convention in Non-Cubic MD cells’ March 29, 1989

get_nearest_grid_point(spos_c, force_to_this_domain=False)[source]

Return index of nearest grid point.

The nearest grid point can be on a different CPU than the one the nucleus belongs to (i.e. return can be negative, or larger than gd.end_c), in which case something clever should be done. The point can be forced to the grid descriptors domain to be consistent with self.get_rank_from_position(spos_c).

integrate(a_xg, b_yg=None, global_integral=True, hermitian=False)[source]

Integrate function(s) over domain.

a_xg: ndarray

Function(s) to be integrated.

b_yg: ndarray

If present, integrate a_xg.conj() * b_yg.

global_integral: bool

If the array(s) are distributed over several domains, then the total sum will be returned. To get the local contribution only, use global_integral=False.

hermitian: bool

Result is hermitian.

interpolate_grid_points(spos_nc, vt_g)[source]

Return interpolated values.

Calculate interpolated values from array vt_g based on the scaled coordinates on spos_c.

This doesn’t work in parallel, since it would require communication between neighbouring grids.

is_my_grid_point(R_c)[source]

Check if grid point belongs to this domain.

new_descriptor(N_c=None, cell_cv=None, pbc_c=None, comm=None, parsize_c=None, allow_empty_domains=False)[source]

Create new descriptor based on this one.

The new descriptor will use the same class (possibly a subclass) and all arguments will be equal to those of this descriptor unless new arguments are provided.

plane_wave(k_c)[source]

Evaluate plane wave on grid.

Returns:

  _ _
 ik.r
e    ,

where the wave vector is given by k_c (in units of reciprocal lattice vectors).

refine()[source]

Return refined \(GridDescriptor\) object.

Returned descriptor has 2x2x2 more grid points.

wannier_matrix(psit_nG, psit_nG1, G_c, nbands=None)[source]

Wannier localization integrals

The soft part of Z is given by (Eq. 27 ref1):

~       ~     -i G.r   ~
Z   = <psi | e      |psi >
 nm       n             m

psit_nG and psit_nG1 are the set of wave functions for the two different spin/kpoints in question.

ref1: Thygesen et al, Phys. Rev. B 72, 125119 (2005)

zero_pad(a_xg, global_array=True)[source]

Pad array with zeros as first element along non-periodic directions.

Array may either be local or in standard decomposition.

zeros(n=(), dtype=<class 'float'>, global_array=False, pad=False, xp=<module 'numpy' from '/home/docs/checkouts/readthedocs.org/user_builds/gpaw/envs/latest/lib/python3.12/site-packages/numpy/__init__.py'>)[source]

Return new zeroed 3D array for this domain.

The type can be set with the dtype keyword (default: float). Extra dimensions can be added with n=dim. A global array spanning all domains can be allocated with global_array=True.

class gpaw.scf.SCFLoop(criteria, maxiter=100, niter_fixdensity=None)[source]

Self-consistent field loop.

log(log, converged_items, entries, context)[source]

Output from each iteration.

class gpaw.band_descriptor.BandDescriptor(nbands, comm=None, strided=False)[source]

Descriptor-class for ordered lists of bands

A BandDescriptor object holds information on how functions, such as wave functions and corresponding occupation numbers, are divided into groups according to band indices. The main information here is how many bands are stored on each processor and who gets what.

This is how a 12 band array is laid out in memory on 3 cpu’s:

 a) Blocked groups       b) Strided groups

      3    7   11             9   10   11
myn   2 \  6 \ 10       myn   6    7    8
 |    1  \ 5  \ 9        |    3    4    5
 |    0    4    8        |    0    1    2
 |                       |
 +----- band_rank        +----- band_rank

Example

>>> a = np.zeros((3, 4))
>>> a.ravel()[:] = range(12)
>>> a
array([[ 0.,  1.,  2.,  3.],
       [ 4.,  5.,  6.,  7.],
       [ 8.,  9., 10., 11.]])
>>> b = np.zeros((4, 3))
>>> b.ravel()[:] = range(12)
>>> b.T
array([[ 0.,  3.,  6.,  9.],
       [ 1.,  4.,  7., 10.],
       [ 2.,  5.,  8., 11.]])

Construct band-descriptor object.

Parameters:

nbands: int

Global number of bands.

comm: MPI-communicator

Communicator for band-groups.

strided: bool

Enable strided band distribution for better load balancing with many unoccupied bands.

Note that if comm.size is 1, then all bands are contained on a single CPU and blocked/strided grouping loses its meaning.

Attributes:

nbands

Number of bands in total.

mynbands

Number of bands on this CPU.

beg

Beginning of band indices in group (inclusive).

end

End of band indices in group (exclusive).

step

Stride for band indices between beg and end.

comm

MPI-communicator for band distribution.

collect(a_nx, broadcast=False)[source]

Collect distributed array to master-CPU or all CPU’s.

distribute(B_nx, b_nx)[source]

distribute full array B_nx to band groups, result in b_nx. b_nx must be allocated.

empty(n=(), dtype=<class 'float'>, global_array=False)[source]

Return new uninitialized 3D array for this domain.

The type can be set with the dtype keyword (default: float). Extra dimensions can be added with n=dim. A global array spanning all domains can be allocated with global_array=True.

get_band_indices(band_rank=None)[source]

Return the global band indices which belong to a given rank.

get_band_ranks()[source]

Return array of ranks as a function of global band indices.

get_slice(band_rank=None)[source]

Return the slice of global bands which belong to a given rank.

global_index(myn, band_rank=None)[source]

Convert rank information and local index to global index.

who_has(n)[source]

Convert global band index to rank information and local index.

zeros(n=(), dtype=<class 'float'>, global_array=False)[source]

Return new zeroed 3D array for this domain.

The type can be set with the dtype keyword (default: float). Extra dimensions can be added with n=dim. A global array spanning all domains can be allocated with global_array=True.

class gpaw.spinorbit.BZWaveFunctions(kd, wfs, occ, nelectrons, n_aj, l_aj)[source]

Container for eigenvalues and PAW projections (all of BZ).

calculate_band_energy()[source]

Calculate sum over occupied eigenvalues.

eigenvalues(broadcast=True)[source]

Eigenvalues in eV for the whole BZ.

eigenvectors(broadcast=True)[source]

Eigenvectors for the whole BZ.

get_atomic_density_matrices()[source]

Returns atomic density matrices for each atom.

get_orbital_magnetic_moments()[source]

Returns the orbital magnetic moment vector for each atom a.

occupation_numbers(broadcast=True)[source]

Occupation numbers for the whole BZ.

pdos_weights(a, indices, broadcast=True)[source]

Projections for PDOS.

Returns (nbzkpts, nbands, 2)-shaped ndarray of the square of absolute value of the projections.

spin_projections(broadcast=True)[source]

Spin projections for the whole BZ.

class gpaw.spinorbit.WaveFunction(eigenvalues, projections, bz_index=None)[source]
add_soc(dVL_avii, s_vss, C_ss)[source]

Evaluate H in a basis of S_z eigenstates.

pdos_weights(a, indices)[source]

PDOS weights.

transform(IBZ2BZMap, K)[source]

Transforms PAW projections from IBZ to BZ k-point.

class gpaw.kpt_descriptor.KPointDescriptor(kpts, nspins=1)[source]

Descriptor-class for k-points.

Construct descriptor object for kpoint/spin combinations (ks-pair).

Parameters:

kpts: None, sequence of 3 ints, or (n,3)-shaped array

Specification of the k-point grid. None=Gamma, list of ints=Monkhorst-Pack, ndarray=user specified.

nspins: int

Number of spins.

Attributes =================== ================================================= N_c Number of k-points in the different directions. nspins Number of spins in total. mynspins Number of spins on this CPU. nibzkpts Number of irreducible kpoints in 1st BZ. mynks Number of k-point/spin combinations on this CPU. gamma Boolean indicator for gamma point calculation. comm MPI-communicator for kpoint distribution. weight_k Weights of each k-point ibzk_kc Unknown ibzk_qc Unknown sym_k Unknown time_reversal_k Unknown bz2ibz_k Unknown ibz2bz_k Unknown bz2bz_ks Unknown symmetry Object representing symmetries =================== =================================================

collect(a_ux, broadcast)[source]

Collect distributed data to all.

copy(comm=SerialCommunicator())[source]

Create a copy with shared symmetry object.

create_k_points(sdisp_cd, collinear)[source]

Return a list of KPoints.

find_ibzkpt(symrel, ibzk_kc, bzk_c)[source]

Find index in IBZ and related symmetry operations.

find_k_plus_q(q_c, kpts_k=None)[source]

Find the indices of k+q for all kpoints in the Brillouin zone.

In case that k+q is outside the BZ, the k-point inside the BZ corresponding to k+q is given.

Parameters:
  • q_c (np.ndarray) – Coordinates for the q-vector in units of the reciprocal lattice vectors.

  • kpts_k (Sequence[int]) – Restrict search to specified k-points.

get_bz_q_points(first=False)[source]

Return the q=k1-k2. q-mesh is always Gamma-centered.

get_count(rank=None)[source]

Return the number of ks-pairs which belong to a given rank.

get_ibz_q_points(bzq_qc, op_scc)[source]

Return ibz q points and the corresponding symmetry operations that work for k-mesh as well.

get_indices(rank=None)[source]

Return the global ks-pair indices which belong to a given rank.

get_offset(rank=None)[source]

Return the offset of the first ks-pair on a given rank.

get_rank_and_index(k)[source]

Find rank and local index of k-point/spin combination.

get_transform_wavefunction_index(nG, k)[source]

Get the “wavefunction transform index”.

This is a permutation of the numbers 1, 2, .. N which associates k + q to some k, and where N is the total number of grid points as specified by nG which is a 3D tuple.

Returns index_G and phase_G which are one-dimensional arrays on the grid.

set_communicator(comm)[source]

Set k-point communicator.

set_symmetry(atoms, symmetry, comm=None)[source]

Create symmetry object and construct irreducible Brillouin zone.

atoms: Atoms object

Defines atom positions and types and also unit cell and boundary conditions.

symmetry: Symmetry object

Symmetry object.

transform_wave_function(psit_G, k, index_G=None, phase_G=None)[source]

Transform wave function from IBZ to BZ.

k is the index of the desired k-point in the full BZ.

where_is_q(q_c, bzq_qc)[source]

Find the index of q points in BZ.

who_has(k)[source]

Convert global index to rank information and local index.

class gpaw.projections.Projections(nbands, nproj_a, atom_partition, bcomm=None, collinear=True, spin=0, dtype=None, data=None, bdist=None)[source]
collect()[source]

Collect all bands and atoms to master.

redist(atom_partition)[source]

Redistribute atoms.