ksdft_calculation

Wraps the KSDFT calculation and saves results to a .chk file.

This module wraps the KSDFT-calculation as it is implemented in the pyscf library. A Restricted Kohn Sham - class is build and the molecule along with parameters of the calculation are specified. To obtain the intermediate results after each iteration of the self-consistent-field method, which are needed further down the line as additional datapoints for training the model, a callback function was implemented. For the DIIS-coefficients it is further necessary to patch the extrapolate function of the CDIIS - class in pyscf. All values are stored at runtime in a .chk file.

The file format is as follows:

In the “Results” folder the initialization parameters of the calculation are stored:

  • “converged” : bool , if the scf-method converged on an energy.

  • “total_energy” : float , the total energy of the molecule after the last iteration.

  • “occupation_numbers_orbitals” : np.ndarray , the occupation numbers of each orbital after the last iteration.

  • “molecular_coeffs_orbitals” : np.ndarray , the coefficients of the orbitals after the last iteration.

  • “max_cycle” : int , the maximal number of iterations the calculation was initialized with.

  • “name_xc_functional” : string , the specified exchange correlation functional.

  • “init_guess” : string , the specified initial guess method.

  • “convergence_tolerance” : float , the specified convergence tolerance.

  • “diis_start_cycle” : int , the first cycle when diis is being used.

  • “diis_space” : int , the maximal number of Fock matrices to average in the diis scheme.

  • “diis_method” : string , the specified diis method.

  • “grid_level” : int , A parameter specifying the density of the grid used in the xc functional.

  • “prune_method” : string , A parameter specifying the pruning scheme of the grid used in the xc functional. Additionally, in the “Initialization” folder the initial density matrix and the total energy before the first iteration are stored, which are needed for the kinetic energy gradient.

  • “first_density_matrix” : np.ndarray , the density matrix before the first iteration.

  • “first_total_energy” : float , the total energy before the first iteration. for each cycle the intermediate results are saved into a “KS-iteration/{cycle}” folder:

  • “diis_coefficients” : np.ndarray , the diis-coefficients used to construct the current Fock matrix.

  • “occupation_numbers_orbitals” : np.ndarray , the occupations number of the orbitals.

  • “molecular_coeffs_orbitals” : np.ndarray , the coefficients of the orbitals in the specified basis.

  • “total_energy” : float , the total energy of the molecule in this iteration.

  • “coulomb_energy” : float , the total coulomb energy after this iteration (equal to the hartree energy).

  • “exchange_correlation_energy” : float , the energy calculated from the exchange correlation functional.

exception ConvergenceError[source]

Raised when the calculation did not converge.

_save_scf_iteration_callback(envs: dict) None[source]

Save the data of each iteration in the chkfile.

Parameters:

envs – Dictionary with the local environment variables of the iteration.

Returns:

None

ksdft(mol: ~pyscf.gto.mole.Mole, savefile: ~pathlib.Path, xc_functional: str = 'PBE', init_guess: str = 'minao', max_cycle: int = 50, diis_space: int = 8, diis_method: str = 'CDIIS', grid_level: int = 3, prune_method: str | ~typing.Callable | None = <function nwchem_prune>, density_fit_basis: str = 'def2-universal-jfit', density_fit_threshold: int = 30, convergence_tolerance: float = 1e-09, extra_callback: ~typing.Callable = None, use_perturbation: bool = False, perturbation_cfg: ~omegaconf.dictconfig.DictConfig | None = None) None[source]

Calculates the non-relativistic restricted spin calculation as in [M-OFDFT]. XC functional should be PBE. Basis set should be 6-31G(2df,p). DIIS is enabled by default. MINAO initialization is active by default.

Parameters:
  • mol – Molecule object

  • savefile – Path to the savefile

  • xc_functional – The xc functional. Default is PBE.

  • init_guess – the initial guess method (default is MINAO)

  • max_cycle – the maximal number of iterations

  • diis_space – the maximal number of Fock matrices used in the DIIS method

  • diis_method – Either CDIIS(default), EDIIS or ADIIS

  • grid_level – The grid level to use for the integration grid used in the xc-functional.

  • prune_method – The method to prune the integration grid. If None, no pruning is performed.

  • density_fit_basis – The basis set to use for the density fitting.

  • density_fit_threshold – The threshold for the number of atoms in the molecule to use density fitting.

  • convergence_tolerance – The convergence tolerance after which the SCF iteration stops. An alternative value can be 1meV 0.0000367493, see Appendix C.2 in [M-OFDFT].

  • extra_callback – Additional callback function to be called after the original callback is called each iteration.

  • use_perturbation – If True, the Fock matrix is perturbed each iteration.

  • perturbation_cfg – Configuration for the perturbation of the Fock matrix.

Returns:

None

Raises:

ConvergenceError – If the calculation did not converge.

patched_extrapolate(self, n_d: int | None = None, use_perturbation: bool = False, pertubation_function: Callable | None = None, start_cycle: int | None = None, end_cycle: int | None = None) ndarray[source]

Monkeypatch for pyscf.lib.diis.extrapolate to extract the DIIS coefficients. Extrapolate the next Fock matrix based on the previous Fock matrices and errors.

Parameters:

n_d – Number of vectors to be used in the extrapolation. Default is all vectors.

Returns:

A new Fock matrix.

Return type:

np.ndarray

perturb_fock_matrix(fock_matrix: ~numpy.ndarray, mf: <module 'pyscf.dft.rks' from '/home/runner/work/structures25/structures25/.venv/lib/python3.12/site-packages/pyscf/dft/rks.py'>, cycle: int, start_std: float, end_std: float, start_cycle: int, end_cycle: int, of_basis_set: str) tuple[ndarray, ndarray][source]

Perturb the Fock matrix by adding a random perturbation which is sampled in the density basis.

Parameters:
  • fock_matrix – The Fock matrix to perturb.

  • mf – The pyscf mf object.

  • cycle – The current cycle of the SCF iteration.

  • start_std – The standard deviation of the perturbation at the first cycle.

  • end_std – The standard deviation of the perturbation at the last cycle.

  • start_cycle – The cycle at which the perturbation starts.

  • end_cycle – The cycle at which the perturbation ends.

  • of_basis_set – String of the density basis set.