lsdalton

Initialize LSDalton interface module.

class daltonproject.lsdalton.OutputParser(filename: str)

Parse LSDalton output files.

property dipole_gradients: ndarray

Extract Cartesian dipole gradients from OpenRSP tensor file.

property electronic_energy: float

Extract the electronic energy from LSDalton output file.

property energy: float

Extract the final energy from LSDalton output file.

property excitation_energies: ndarray

Extract one-photon excitation energies from LSDalton output file.

property filename: str

Name of the LSDalton output files without the extension.

property final_geometry: ndarray

Extract final geometry from LSDalton output file.

property gradients: ndarray

Extract molecular gradient from LSDalton output file.

property hessian: ndarray

Extract molecular hessian from OpenRSP tensor file.

property nuclear_repulsion_energy: float

Extract the nuclear repulsion energy from LSDalton output file.

property oscillator_strengths: ndarray

Extract one-photon oscillator strengths from LSDalton output file.

property polarizability_gradients: PolarizabilityGradients

Extract Cartesian polarizability gradients from OpenRSP tensor file.

daltonproject.lsdalton.compute(molecule: Molecule, basis: Basis, qc_method: QCMethod, properties: Property, environment: Environment | None = None, compute_settings: ComputeSettings | None = None, filename: str | None = None, force_recompute: bool = False) OutputParser

Run a calculation using the LSDalton program.

Parameters
  • molecule – Molecule on which a calculations is performed. This can also be an atom, a fragment, or any collection of atoms.

  • basis – Basis set to use in the calculation.

  • qc_method – Quantum chemistry method, e.g., HF, DFT, or CC, and associated settings.

  • properties – Properties of molecule to be calculated, geometry optimization, excitation energies, etc.

  • environment – Environment description missing.

  • compute_settings – Settings for the calculation, e.g., number of MPI processes and OpenMP threads, work and scratch directory, etc.

  • filename – Optional user-specified filename that will be used for input and output files. If not specified a name will be generated as a hash of the input.

  • force_recompute – Recompute even if the output files already exist.

Returns

OutputParser instance that contains the filename of the output produced in the calculation and can be used to extract results from the output.

daltonproject.lsdalton.compute_farm(molecules: Sequence[Molecule], basis: Basis, qc_method: QCMethod, properties: Property, compute_settings: ComputeSettings | None = None, filenames: Sequence[str] | None = None, force_recompute: bool = False) list[OutputParser]

Run a series of calculations using the LSDalton program.

Parameters
  • molecules – List of molecules on which calculations are performed.

  • basis – Basis set to use in the calculations.

  • qc_method – Quantum chemistry method, e.g., HF, DFT, or CC, and associated settings used for all calculations.

  • properties – Properties of molecule to be calculated, geometry optimization, excitation energies, etc., computed for all molecules.

  • compute_settings – Settings for the calculations, e.g., number of MPI processes and OpenMP threads, work and scratch directory, and more.

  • filenames – Optional list of user-specified filenames that will be used for input and output files. If not specified names will be generated as a hash of the individual inputs.

  • force_recompute – Recompute even if the output files already exist.

Returns

List of OutputParser instances each of which contains the filename of the output produced in the calculation and can be used to extract results from the output.

daltonproject.lsdalton.coulomb_matrix(density_matrix: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the Coulomb matrix.

Parameters
  • density_matrix – Density matrix

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Coulomb matrix

daltonproject.lsdalton.diagonal_density(hamiltonian: ndarray, metric: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod) ndarray

Form an AO density matrix from occupied MOs through diagonalization.

Parameters
  • hamiltonian – Hamiltonian matrix.

  • metric – Metric matrix (i.e. overlap matrix)

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Density matrix

daltonproject.lsdalton.electronic_electrostatic_potential(density_matrix: ndarray, points: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, ep_derivative_order: tuple[int, int] = (0, 0), geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the electronic electrostatic potential at a set of points.

Derivatives of the electrostatic potential can be calculated using the ep_derivative_order argument, e.g., ep_derivative_order = (0, 1) will calculate both the zeroth- and first-order derivatives while ep_derivative_order = (1, 1) will calculate first-order derivatives only.

Parameters
  • density_matrix – Density matrix.

  • points – Set of points where the electrostatic potential will be evaluated.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • ep_derivative_order – Range of orders of the derivative of the electrostatic potential.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Electronic electrostatic potential

daltonproject.lsdalton.electrostatic_potential(density_matrix: ndarray, points: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, ep_derivative_order: tuple[int, int] = (0, 0), geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the molecular electrostatic potential at a set of points.

Derivatives of the electrostatic potential can be calculated using the ep_derivative_order argument, e.g., ep_derivative_order = (0, 1) will calculate both the zeroth- and first-order derivatives while ep_derivative_order = (1, 1) will calculate first-order derivatives only.

Parameters
  • density_matrix – Density matrix.

  • points – Set of points where the electrostatic potential will be evaluated.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • ep_derivative_order – Range of orders of the derivative of the electrostatic potential.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Point-wise electronic multipol-moment potential

daltonproject.lsdalton.electrostatic_potential_integrals(points: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, ep_derivative_order: tuple[int, int] = (0, 0), geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate electrostatic-potential integrals.

Derivatives of the electrostatic potential can be calculated using the ep_derivative_order argument, e.g., ep_derivative_order = (0, 1) will calculate both the zeroth- and first-order derivatives while ep_derivative_order = (1, 1) will calculate first-order derivatives only.

Let a and b be AO basis functions, C a point, and \(\xi\) and \(\zeta\) Cartesian components, then order 0: \((ab|C) = \int a(r) b(r) V(r,C) dr\), with the potential: \(V(r,C) = 1/|r-C|\) order 1: \(\int a(r) b(r) V_\xi(r,C) dr\), with first-order potential derivative: \(V_\xi(r,C) = \partial{V(r,C)}{\partial{C_\xi}}= C_\xi/|r-C|^3\) order 2: \(\int a(r) b(r) V_{\xi,\zeta}(r,C) dr\), with second-order potential derivative: \(V_{\xi,\zeta}(r,C) = \partial^2{V(r,C)}{\partial{C_\xi}\partial{C_\zeta}}\)

Parameters
  • points – Set of points where the electrostatic-potential integrals will be evaluated.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • ep_derivative_order – Range of orders of the derivative of the electrostatic potential.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Electrostatic-potential integrals

daltonproject.lsdalton.eri(specs: str, molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate electron-repulsion integrals for different operator choices.

\[g_{ijkl} = \int\mathrm{d}\mathbf{r}\int\mathrm{d}\mathbf{r}^{\prime} \chi_{i}(\mathbf{r})\chi_{j}(\mathbf{r}^\prime) g(\mathbf{r},\mathbf{r}^\prime) \chi_{k}(\mathbf{r})\chi_{l}(\mathbf{r}^\prime)\]

The first character of the specs string specifies the operator \(g(\mathbf{r}, \mathbf{r}^\prime)\). Valid values are:

  • C for Coulomb (\(g(\mathbf{r}, \mathbf{r}^\prime) = \frac{1}{|\mathbf{r} - \mathbf{r}^\prime|}\))

  • G for Gaussian geminal (\(g\))

  • F geminal divided by the Coulomb operator (\(\frac{g}{|\mathbf{r} - \mathbf{r}^\prime|}\))

  • D double commutator (\([[T,g],g]\))

  • 2 Gaussian geminal operator squared (\(g^2\))

The last four characters of the specs string specify the AO basis to use for the four basis functions \(\chi_{i}\), \(\chi_{j}\), \(\chi_{k}\), and \(\chi_{l}\). Valid values are:

  • R for regular AO basis

  • D for auxiliary (RI) AO basis

  • E for empty AO basis (i.e. for 2- and 3-center ERIs)

  • N for nuclei

Parameters
  • specs – A 5-character-long string with the specification for AO basis and operator to use (see description above).

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Electron-repulsion integral tensor.

daltonproject.lsdalton.eri4(molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate four-center electron-repulsion integrals.

\[(ab|cd) = g_{acbd}\]

with

\[g(\mathbf{r},\mathbf{r}^\prime) = \frac{1}{|\mathbf{r}-\mathbf{r}^\prime|}\]
Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Four-center electron-repulsion integral tensor.

daltonproject.lsdalton.exchange_correlation(density_matrix: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) tuple[float, numpy.ndarray]

Calculate the exchange-correlation energy and matrix.

Parameters
  • density_matrix – Density matrix.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Exchange-correlation energy and matrix

daltonproject.lsdalton.exchange_matrix(density_matrix: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the exchange matrix.

Parameters
  • density_matrix – Density matrix.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Exchange matrix

daltonproject.lsdalton.fock_matrix(density_matrix: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the Fock/KS matrix.

Parameters
  • density_matrix – Density matrix

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Fock/KS matrix

daltonproject.lsdalton.kinetic_energy_matrix(molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the kinetic energy matrix.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Kinetic energy matrix

daltonproject.lsdalton.multipole_interaction_matrix(moments: ndarray, points: ndarray, molecule: Molecule, basis: Basis, qc_method: QCMethod, multipole_orders: tuple[int, int] = (0, 0), geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the electron-multipole electrostatic interaction matrix in atomic-orbital (AO) basis.

Minimum and maximum orders of the multipole moments are provided by the multipole_orders argument, e.g., multipole_orders = (0, 1) includes charges (order 0) and dipoles (order 1) while multipole_orders = (1, 1) includes only dipoles.

Parameters
  • moments – Multipole moments.

  • points – Positions of the multipole moments.

  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • multipole_orders – Minimum and maximum orders of the multipole moments, where 0 = charge, 1 = dipole, etc.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Electron-multipole electrostatic interaction matrix

daltonproject.lsdalton.nuclear_electron_attraction_matrix(molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the nuclear-electron attraction matrix.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Nuclear-electron attraction matrix

daltonproject.lsdalton.nuclear_electrostatic_potential(points: ndarray, molecule: Molecule, ep_derivative_order: tuple[int, int] = (0, 0), geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the nuclear electrostatic potential at a set of points.

Derivatives of the electrostatic potential can be calculated using the ep_derivative_order argument, e.g., ep_derivative_order = (0, 1) will calculate both the zeroth- and first-order derivatives while ep_derivative_order = (1, 1) will calculate first-order derivatives only.

Parameters
  • points – Set of points where the electrostatic potential will be evaluated.

  • molecule – Molecule object.

  • ep_derivative_order – Range of orders of the derivative of the electrostatic potential.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Nuclear electrostatic potential

daltonproject.lsdalton.nuclear_energy(molecule: Molecule) float

Calculate the nuclear-repulsion energy.

Parameters

molecule – Molecule object.

Returns

Nuclear-repulsion energy

daltonproject.lsdalton.num_atoms(molecule: Molecule, basis: Basis, qc_method: QCMethod) int

Return the number of atoms.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Number of atoms

daltonproject.lsdalton.num_basis_functions(molecule: Molecule, basis: Basis, qc_method: QCMethod) int

Return the number of basis functions.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Number of basis functions

daltonproject.lsdalton.num_electrons(molecule: Molecule, basis: Basis, qc_method: QCMethod) int

Return the number of electrons.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Number of electrons

daltonproject.lsdalton.num_ri_basis_functions(molecule: Molecule, basis: Basis, qc_method: QCMethod) int

Return the number of auxiliary (RI) basis functions.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Number of auxiliary (RI) basis functions

daltonproject.lsdalton.overlap_matrix(molecule: Molecule, basis: Basis, qc_method: QCMethod, geometric_derivative_order: int = 0, magnetic_derivative_order: int = 0) ndarray

Calculate the overlap matrix.

Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

  • geometric_derivative_order – Order of the derivative with respect to nuclear displacements.

  • magnetic_derivative_order – Order of the derivative with respect to magnetic field.

Returns

Overlap matrix

daltonproject.lsdalton.ri2(molecule: Molecule, basis: Basis, qc_method: QCMethod) ndarray

Calculate two-center electron-repulsion integrals.

\[(I|J) = g_{I00J}\]

with

\[g(\mathbf{r},\mathbf{r}^\prime) = \frac{1}{|\mathbf{r}-\mathbf{r}^\prime|}\]
Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Two-center RI integrals

daltonproject.lsdalton.ri3(molecule: Molecule, basis: Basis, qc_method: QCMethod) ndarray

Calculate three-center electron-repulsion integrals.

\[(ab|I) = g_{a0bI}\]

with

\[g(\mathbf{r},\mathbf{r}^\prime) = \frac{1}{|\mathbf{r}-\mathbf{r}^\prime|}\]
Parameters
  • molecule – Molecule object.

  • basis – Basis object.

  • qc_method – QCMethod object.

Returns

Three-center RI integrals