API

DFT

class pydft.DFT(mol: Molecule, basis: str | list[cgf] = 'sto3g', functional: str = 'svwn5', nshells: int = 32, nangpts: int = 110, lmax: int = 8, fdpts: int = 7, normalize: bool = True, parallel: bool = False)
__init__(mol: Molecule, basis: str | list[cgf] = 'sto3g', functional: str = 'svwn5', nshells: int = 32, nangpts: int = 110, lmax: int = 8, fdpts: int = 7, normalize: bool = True, parallel: bool = False)

Constructs the DFT class

Parameters:
  • mol (Molecule) – molecule

  • basis (str, optional) – basis set, can be either a string or a list of cgf objects, by default ‘sto3g’

  • functional (str, optional) – exchange-correlation function, valid options are swn5 and pbe, by default ‘svwn5’

  • nshells (int, optional) – number of radial shells, by default 32

  • nangpts (int, optional) – number of angular sampling points, by default 110

  • lmax (int, optional) – maximum value of l in the spherical harmonic expansion, by default 8

  • fdpts (int, optional) – number of grid point in finite difference scheme, by default 7

  • normalize (whether to perform intermediary normalization of the electron) – density

  • verbose (bool, optional) – whether to provide verbose output, by default False

  • parallel (bool, optional) – whether to use multiprocessing features (only for Linux), by default False

get_construction_times() dict

Return construct times dictionary from MolecularGrid to get insights into the time to construct particular properties

Returns:

dictionary of construction times for the MolecularGrid objects

Return type:

dict

get_data() dict

Get relevant memory objects, only valid after successful SCF calculation

Returns:

Dictionary containing relevant memory objects (see below)

Return type:

dict

Notes

The dictionary contains the following elements:

  • S: overlap matrix

  • T: kinetic energy matrix

  • V: nuclear attraction matrix

  • C: coefficient matrix

  • J: Hartree matrix

  • P: density matrix

  • XC: exchange-correlation matrix

  • F: Fock-matrix

  • Exc: exchange-correlation energy

  • Ex: exchange energy

  • Ec: correlation energy

  • energies: all energies

  • energy: last energy

  • orbc: coeffient matrix (duplicate)

  • orbe: molecular orbital eigenvalues

  • enucrep: electrostatic repulsion of the nuclei

  • timedata: computation time for various parts of the calculation

get_density_at_points(spoints: ndarray) ndarray

Get the electron density at provided points

Parameters:

spoints (ndarray) – Set of sampling points (\(X \times 3\) array)

Returns:

Electron density scalar field (\(N \times 1\) array)

Return type:

ndarray

get_gradient_at_points(spoints: ndarray) ndarray

Get the electron density gradient at provided points

Parameters:

spoints (ndarray) – Set of sampling points (\(X \times 3\) array)

Returns:

Electron density gradient vector field (\(N \times 3\) array)

Return type:

ndarray

get_molgrid_copy() MolecularGrid

Get a copy of the underlying molgrid object

Returns:

copy of the MolecularGrid object

Return type:

MolecularGrid

scf(tol: float = 1e-05, verbose: bool = False) float

Perform the self-consistent field procedure

Parameters:

tol (float, optional) – electronic convergence criterion, by default 1e-5

Returns:

total electronic energy (in Hartrees)

Return type:

float

MolecularGrid

class pydft.MolecularGrid(atoms: list, cgfs: list[cgf], nshells: int = 32, nangpts: int = 110, lmax: int = 8, fdpts: int = 7, functional: str = 'svwn5', parallel: bool = False)
__init__(atoms: list, cgfs: list[cgf], nshells: int = 32, nangpts: int = 110, lmax: int = 8, fdpts: int = 7, functional: str = 'svwn5', parallel: bool = False)

Construct MolecularGrid

Parameters:
  • atoms (list) – list of atoms and their charges

  • cgfs (list[pq.cgf]) – list of contracted Gaussian functions (basis set)

  • nshells (int, optional) – number of radial shells, by default 32

  • nangpts (int, optional) – number of angular sampling points per shell, by default 110

  • lmax (int, optional) – maximum value for l for projection of spherical harmonics, by default 8

  • fdpts (int, optional) – number of grid point in finite difference scheme, by default 7

  • functional (str, optional) – exchange-correlation functional, by default ‘svwn5’

  • parallel (bool, optional) – whether to use multiprocessing features, by default False

build_density(P: ndarray, normalize: bool = True)

Build the electron density from the density matrix

Parameters:
  • P (np.ndarray) – density matrix

  • normalize (bool, optional) – whether to normalize the density matrix based on the total number of electrons, by default True

static build_ylmgpts_atom(Omega, lmax, npts)

Get the spherical harmonic parameters for a single atom

This function is ‘pickable’ and can be used in a multiprocessing pool

calculate_correlation() -> (<class 'numpy.ndarray'>, <class 'float'>)

Build the correlation matrix and give the correlation energy for the molecule

Returns:

Correlation matrix and exchange energy

Return type:

np.ndarray,float

calculate_coulomb_potential_at_points(pts: ndarray) ndarray

Collect the spherical harmonic expansion coefficients of the hartree potential for all the atoms at the grid points

Parameters:

pts (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

Returns:

Coulomb potential at grid points \((\left[l_{\textrm{max}}+1\right]^{2} \times N)\) with \(N\) the number of grid points

Return type:

np.array

calculate_coulombic_matrix() array

Build the coulomb matrix \(\mathbf{J}\)

Returns:

Coulomb matrix \(\mathbf{J}\)

Return type:

np.array

calculate_dfa_coulomb() float

Get density functional approximation for the coulombic repulsion between electrons

Returns:

Total coulombic repulsion (Hartree)

Return type:

float

calculate_dfa_coulomb_no_interpolation() float

Calculate the Coulombic self-repulsion in the limit of the sum of atomic centers, without any interpolation of the contribution of the other atoms

Returns:

Approximate coulombic repulsion (Hartree)

Return type:

float

calculate_dfa_exchange()

Get density functional approximation for the exchange energy

Returns:

Total exchange energy (Hartree)

Return type:

float

calculate_dfa_kinetic() float

Get density functional approximation for the kinetic energy

Returns:

Total kinetic energy (Hartree)

Return type:

float

calculate_dfa_nuclear_attraction_full() array

Get density functional approximation for the nuclear attraction energy

Returns:

Total nuclear attraction energy (Hartree)

Return type:

float

calculate_dfa_nuclear_attraction_local() float

Get density functional approximation for the nuclear attraction energy using only the local potential for each atomic grid and ignoring the influence of atomic attraction due to other atoms

Returns:

Approximate nuclear attraction energy (Hartree)

Return type:

float

Notes

This function will of course not give an accurate result and is only meant for educational purposes

calculate_exchange() -> (<class 'numpy.ndarray'>, <class 'float'>)

Build the exchange matrix and give the exchange energy for the molecule

Returns:

Exchange matrix \(\mathbf{X}\) and exchange energy

Return type:

np.ndarray,float

calculate_weights_at_points(points: ndarray, k: int = 3) ndarray

Custom function that (re-)calculates weights at given points

Parameters:
  • points (np.ndarray) – \(N \times 3\) array of grid points

  • k (int, optional) – smoothing value, by default 3

Returns:

\(N_{a} \times N\) array with \(N_{a}\) the number of atoms and \(N\) the number of grid points

Return type:

np.ndarray

count_electrons() float

Count number of electrons in the system

Returns:

Total number of electrons

Return type:

float

count_electrons_from_rho_lm() float

Count number of electrons from the spherical harmonic expansion; this function mainly acts to verify that the spherical harmonic expansion works correctly

Returns:

Total number of electrons

Return type:

float

Notes

The Becke weights are already used when generating the density coefficients and hence for the back-transformation, we should not multiply with the Becke weights

get_amplitude_at_points(spoints: ndarray, c: ndarray) ndarray

Calculate the wave function amplitude at the points as given by spoints using solution vector c

Parameters:
  • spoints (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

  • c (np.ndarray) – Coefficient vector \(\vec{c}\) of a single molecular orbital

Returns:

Electron density specified at grid points \((N \times 1)\) with \(N\) the number of grid points

Return type:

np.array

get_atomic_grid(i: int) AtomicGrid

Get an atomic grid of atom \(i\)

Parameters:

i (int) – atom index

Returns:

pydft.AtomicGrid of atom \(i\)

Return type:

AtomicGrid

get_becke_weights() ndarray

Get the Becke coefficients

Returns:

\((N \times G \times 3)\) array with \(N\) number of atoms and \(G\) number of grid points

Return type:

np.array

get_correlation_potential_at_points(pts: ndarray, P: ndarray) ndarray

Get the correlation potential at points pts

Parameters:
  • pts (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

  • P (np.ndarray) – Density matrix \((K \times K)\) with \(K\) the number of basis functions

Returns:

Correlation potential \(\nu_{c}\) at grid points \((N \times 1)\) with \(N\) the number of grid points

Return type:

np.ndarray

get_densities() list[array]

Get the densities at the grid points

Returns:

list over atoms with for each atom a \((G \times 1)\) array of electron densities where \(G\) is the number of grid points per atom

Return type:

list[np.array]

get_density_at_points(spoints: ndarray, P: ndarray) array

Calculate the electron density at the points as given by spoints using density matrix P

Parameters:
  • spoints (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

  • P (np.ndarray) – Density matrix \((K \times K)\) with \(K\) the number of basis functions

Returns:

Electron density specified at grid points \((N \times 1)\) with \(N\) the number of grid points

Return type:

np.array

get_exchange_potential_at_points(pts: ndarray, P: ndarray) ndarray

Get the exchange potential at points pts

Parameters:
  • pts (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

  • P (np.ndarray) – Density matrix \((K \times K)\) with \(K\) the number of basis functions

Returns:

Exchange potential \(\nu_{x}\) at grid points \((N \times 1)\) with \(N\) the number of grid points

Return type:

np.ndarray

get_gradient_at_points(spoints: ndarray, P: ndarray)

Calculate the gradient of the electron density at the points as given by spoints using density matrix \(\mathbf{P}\)

Parameters:
  • spoints (np.ndarray) – Array of grid points \((N \times 3)\) with \(N\) the number of grid points

  • P (np.ndarray) – Density matrix \((K \times K)\) with \(K\) the number of basis functions

Returns:

Electron density gradients at grid points \((N \times 3)\) with \(N\) the number of grid points

Return type:

np.array

get_gradients() list[array]

Get the gradient magnitude of the electron density field

Returns:

list over atoms with for each atom a \((G \times 1)\) array of electron gradient magnitudes where \(G\) is the number of grid points per atom

Return type:

list[np.array]

get_grid_coordinates() list[array]

Get the grid coordinates

Returns:

list over atoms with for each atom a \((G \times 1)\) array of electron densities where \(G\) is the number of grid points per atom

Return type:

list[np.array]

get_hartree_potential() array

Get the Hartree potential at each grid cell

Returns:

\((G \times 1)\) array with \(G\) number of molecular grid points

Return type:

np.array

get_rho_lm_atoms() array

Get the electron density projected onto spherical harmonics per atom

Returns:

\(\rho_{n,lm}\) where \(n\) loops over the atoms and the pairs \(lm\) over the spherical harmonics

Return type:

np.array

get_spherical_harmonic_expansion_of_amplitude(c: ndarray, radial_factor=False) list[ndarray]

Get the spherical harmonic expansion representation of a wavefunction amplitude using solution vector c

Parameters:
  • c (np.ndarray) – Wave function expansion coefficient

  • radial_factor (bool, optional) – whether to multiply result by \(r^{2}\), by default False

Returns:

List of spherical harmonic expansions per atom stored as \((N_{r} \times N_{lm})\) arrays where \(N_{r}\) is the number of radial points and \(N_{lm}\) the set of spherical harmonics used in the expansion.

Return type:

list[np.ndarray]

initialize()

Initialize the molecular grid

This is a relatively expensive procedure and thus the user can delay the initialization of the molecular grid until.