Namespace Scine::Sparrow

namespace Sparrow

Enums

enum IntegralMethod

List of the available calculation methods for calculating the density matrix elements.

Values:

ObaraSaika
ClosedForm

Functions

std::vector<std::shared_ptr<Core::Module>> moduleFactory()
class DFTB0MethodWrapper : public Utils::CloneInterface<DFTB0MethodWrapper, DFTBMethodWrapper>
#include <DFTB0MethodWrapper.h>

A method wrapper handling DFTB0 calculations.

class DFTB0Settings : public Settings
#include <DFTB0Settings.h>

The Settings specific to the DFTB0 method, a non SCF method.

Please note that since DFTB0 is not an SCF method, it will not contain SCF options. Furthermore, it is not capable of performing calculation in unrestricted formalism.

class DFTB2MethodWrapper : public Utils::CloneInterface<DFTB2MethodWrapper, DFTBMethodWrapper>
#include <DFTB2MethodWrapper.h>

A method wrapper handling DFTB2 calculations, also known as SCC-DFTB, self-consistent charge DFTB.

class DFTB2Settings : public Settings
#include <DFTB2Settings.h>

The Settings specific to the DFTB2 (SCC-DFTB) method.

class DFTB3MethodWrapper : public Utils::CloneInterface<DFTB3MethodWrapper, DFTBMethodWrapper>
#include <DFTB3MethodWrapper.h>

A method wrapper handling DFTB3 calculations.

class DFTB3Settings : public Settings
#include <DFTB3Settings.h>

The Settings specific to the DFTB3 method.

template<class DFTBMethod>
class DFTBDipoleMatrixCalculator : public Scine::Sparrow::DipoleMatrixCalculator
#include <DFTBDipoleMatrixCalculator.h>

Class responsible for the calculation of the dipole matrix in MO basis for DFTB methods.

This class calculates the dipole matrix in MO basis according to R. Rüger, E. van Lenthe, T. Heine and L. Visscher, Tight-Binding Approximations to Time-Dependent Density Functional Theory

It works with the DFTBO0, DFTB2 and DFTB3 methods in restricted and unrestricted formalism.

template<class DFTBMethod>
class DFTBDipoleMomentCalculator : public Scine::Sparrow::DipoleMomentCalculator
#include <DFTBDipoleMomentCalculator.h>

This class calculates the electrical dipole moment for the DFTB methods.

Right now it calculates the dipole through a Mulliken population analysis. New methods can then be implemented if deemed necessary. In particular, with this approximation transition dipoles from orbitals on the same atom can be vastly underestimated, which is of importance for TD-DFTB.

Template Parameters
  • DFTBMethod: One of the DFTB methods type, i.e. DFTB0, DFTB2, DFTB3.

class DipoleMatrixCalculator
#include <DipoleMatrixCalculator.h>

Interface for the calculation of the dipole matrix in semiempirical methods.

Subclassed by Scine::Sparrow::DFTBDipoleMatrixCalculator< DFTBMethod >, Scine::Sparrow::NDDODipoleMatrixCalculator< NDDOMethod >

class DipoleMomentCalculator
#include <DipoleMomentCalculator.h>

Interface for the calculation of the electrical dipole moment in a semiempirical method.

Subclassed by Scine::Sparrow::DFTBDipoleMomentCalculator< DFTBMethod >, Scine::Sparrow::NDDODipoleMomentCalculator< NDDOMethod >

class GenericMethodWrapper : public Utils::CloneInterface<Utils::Abstract<GenericMethodWrapper>, Core::Calculator>, public WavefunctionOutputGenerator
#include <GenericMethodWrapper.h>

A MethodWrapper running Generic calculations.

class MoldenFileGenerator

Class to create the wavefunction information needed for outputting densities,…

Note that NDDO methods have their own STO-6G expansion, fine tuned according to their parameters, while for DFTB methods the STO-6G expansion of the PM6 method was used. Since STO-6G expansions are very similar one-another, the implementation ease of this was deemed a satisfactory compromise. If the underlying calculator has not been initialized, then care must be taken that outside of this function a calculation is performed. This class handles only s, p and d orbitals.

class AM1TypeSettings : public Settings

The Settings specific to all the AM1-type methods, i.e.

AM1, RM1, PM3. This class uses the curiously recurrent template pattern to have compile-time resolution of the parameters position and settings descriptions, which both are method-dependant. Having independent classes would lead to unwanted code duplication.

Subclassed by Scine::Sparrow::AM1Settings, Scine::Sparrow::PM3Settings, Scine::Sparrow::RM1Settings

template<class AM1Type>
class AM1TypeMethodWrapper : public Utils::CloneInterface<Utils::Abstract<AM1TypeMethodWrapper<AM1Type>>, NDDOMethodWrapper>
#include <AM1TypeMethodWrapper.h>

A parent class for all AM1-type methods using the same formalism but different parameters (AM1, RM1, PM3).

This class uses the curiously recurrent template pattern to have compile-time resolution of the model names, as well as the use of the right constructor. The constructors are method-dependant, having independent classes would lead to unwanted code duplication.

class MNDOMethodWrapper : public Utils::CloneInterface<MNDOMethodWrapper, NDDOMethodWrapper>
#include <MNDOMethodWrapper.h>

A method wrapper handling MNDO calculations.

class NDDOMethodWrapper : public Utils::CloneInterface<Utils::Abstract<NDDOMethodWrapper>, GenericMethodWrapper>
#include <NDDOMethodWrapper.h>

Abstract class acting as a generic Wrapper for NDDO methods.

class PM6MethodWrapper : public Utils::CloneInterface<PM6MethodWrapper, NDDOMethodWrapper>
#include <PM6MethodWrapper.h>

A method wrapper running PM6 calculations.

class PM6Settings : public Settings
#include <PM6Settings.h>

The Settings specific to the PM6 method.

class AtomPairDipole
#include <AtomPairDipole.h>

Class responsible for calculating a block of the dipole matrix.

class GTODipoleMatrixBlock
#include <GTODipoleMatrixBlock.h>

Class responsible for calculating a block of the dipole matrix in GTO ao basis.

The basis must be one of the STO-nG type.

template<class NDDOMethod>
class NDDODipoleMatrixCalculator : public Scine::Sparrow::DipoleMatrixCalculator
#include <NDDODipoleMatrixCalculator.h>

Class responsible for the calculation of the dipole matrix.

Template Parameters
  • NDDOMethod: An NDDO method, must derive from Utils::ScfMethod and must have the method getInitializer(). getInitializer must return a valid NDDOInitializer instance. This class is responsible for the calculation of the dipole matrix in both atomic and molecular orbital basis. The dipole matrix in atomic orbital basis is calculated by explicitly integrating in closed form or using the highly efficient Obara-Saika scheme the dipole integral <|r|>. The dipole matrix in atomic orbital basis is then transformed in molecular orbital basis by D_{MO} = C^T * D_{AO} * C in restricted formalism, D_{MO} = C_^T * D_{AO} * C_ + C_^T * D_{AO} * C_ in unrestricted formalism.

template<class NDDOMethod>
class NDDODipoleMomentCalculator : public Scine::Sparrow::DipoleMomentCalculator
#include <NDDODipoleMomentCalculator.h>

Class resposible for the calculation of the dipole in the NDDO methods.

It must be able to calculate the dipole both with the NDDO approximation and with the use of the dipole matrix.

Template Parameters
  • NDDOMethod: An NDDO method, must derive from Utils::ScfMethod and must have the method getInitializer. getInitializer must return a valid NDDOInitializer instance.

class SparrowModule : public Module
#include <SparrowModule.h>

The SCINE Module implementation for Sparrow.

struct SparrowState : public State
#include <SparrowState.h>

Definition of a calculation state for methods implemented in Sparrow.

The calculation state is defined as the density matrix and the coefficient matrix at some point. If the density matrix is empty, an exception is thrown at loading time.

namespace dftb

Functions

template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>() const
struct RepulsionParameters
#include <RepulsionParameters.h>

DFTB parameters for the repulsion of an element pair.

class ScfFock : public ElectronicContributionCalculator
#include <ScfFock.h>

Implementation of FockMatrixCalculator for SCF-type DFTB methods.

Subclassed by Scine::Sparrow::dftb::SecondOrderFock, Scine::Sparrow::dftb::ThirdOrderFock

class SecondOrderFock : public Scine::Sparrow::dftb::ScfFock
#include <SecondOrderFock.h>

Implementation of FockMatrixCalculator for DFTB2, the SCC-DFTB.

It calculates the electronic contributions to the energy and their derivatives with respect to the nuclear cartesian coordinates.

class ThirdOrderFock : public Scine::Sparrow::dftb::ScfFock
#include <ThirdOrderFock.h>

Implementation of FockMatrixCalculator for DFTB3.

It calculates the electronic contributions to the energy and their derivatives with respect to the nuclear cartesian coordinates.

class ZeroOrderFock : public ElectronicContributionCalculator
#include <ZeroOrderFock.h>

Implementation of FockMatrixCalculator for DFTB0.

class ZeroOrderMatricesCalculator
#include <ZeroOrderMatricesCalculator.h>

This class calculates the matrices resulting from the zeroth order expansion of the DFT energy for the DFTB methods.

namespace nddo

Enums

enum BasisFunctions

Values:

sp
spd
enum sc_t

Values:

F0ss
F0pp
F0dd
F0sp
F0sd
F0pd
F2pp
F2dd
F2pd
F4dd
G1sp
G1pd
G2sd
G3pd
R1sppd
R2sdpp
R2sddd

Functions

template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::zero> integral<Utils::derivOrder::zero>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::one> integral<Utils::derivOrder::one>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::two> integral<Utils::derivOrder::two>(double R) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>() const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::zero> integral<Utils::derivOrder::zero>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::one> integral<Utils::derivOrder::one>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::two> integral<Utils::derivOrder::two>(double R) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>() const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::zero> integral<Utils::derivOrder::zero>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::one> integral<Utils::derivOrder::one>(double R) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::two> integral<Utils::derivOrder::two>(double R) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>() const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>() const
class AM1PairwiseRepulsion
#include <AM1PairwiseRepulsion.h>

This class calculates the core-core repulsion between two atoms.

class AM1RepulsionEnergy : public RepulsionCalculator
#include <AM1RepulsionEnergy.h>

This class sums up the core-core repulsion energies and the corresponding derivatives with respect to the nuclear cartesian coordinate between all pairs of cores.

It inherits from Utils::RepulsionCalculator in order for it to work with the LCAO/ScfMethod polymorphic system.

class MNDOPairwiseRepulsion
#include <MNDOPairwiseRepulsion.h>

This class calculates the core-core repulsion between two atoms.

class MNDORepulsionEnergy : public RepulsionCalculator
#include <MNDORepulsionEnergy.h>

This class sums up the core-core repulsion energies and the corresponding derivatives with respect to the nuclear cartesian coordinate between all pairs of cores.

It inherits from Utils::RepulsionCalculator in order for it to work with the LCAO/ScfMethod polymorphic system.

class PM6PairwiseRepulsion
#include <PM6PairwiseRepulsion.h>

This class calculates the core-core repulsion between two atoms.

class PM6RepulsionEnergy : public RepulsionCalculator
#include <PM6RepulsionEnergy.h>

This class sums up the core-core repulsion energies and the corresponding derivatives with respect to the nuclear cartesian coordinate between all pairs of cores.

It inherits from Utils::RepulsionCalculator in order for it to work with the LCAO/ScfMethod polymorphic system.

template<Utils::derivOrder O>
class AtomPairOverlap
#include <AtomPairOverlap.h>

This class computes a block of the overlap matrix for two atoms.

The actual calculation is done by the GTOOverlapMatrixBlock class, here the blocks that need calculation are identified and scheduled for calculation.

template<Utils::derivOrder O>
class GTOOverlapMatrixBlock
#include <GTOOverlapMatrixBlock.h>

This class calculates the overlap matrix block and its derivatives for two groups of orbitals on different atoms, each group of which shares the same angular momentum.

F.i. s-s, s-p, p-d, d-d, … according to the Obara-Saika method.

class OneCenterIntegralContainer
#include <oneCenterIntegralContainer.h>

This class contains smart pointers to one-center two-electron matrices for the atoms in the AtomCollection. Since they are equal for an element type, they are calculated only once per element type.

class OneCenterSlaterIntegral
#include <oneCenterSlaterIntegral.h>

Calculation of the one-centre integrals (related to Slater-Condon parameters) according to Kumar, Mishra, J. Phys., 1987, 29, 385-390. The calculation is not optimized for performance as it does not need to be fast. NB: this returns U^l(a,b,c,d) = R^k(a,c,b,d) as defined in other papers.

class OneCenterTwoElectronCalculator
#include <OneCenterTwoElectronCalculator.h>

This class implements formulas for the calculation of the one-center two-electron integrals from 14 Slater-Condon parameters and three radial parameters. The main reference is Pelikan, P; Turi Nagy L., Chemical Papers, 1974, 28, 594-598. The indexes used here are the same as in the reference minus one. In formula 17 of the reference, sqrt(12) is used instead of 12. In formula 54, R1sppd would be correct instead of R1spdd In formulas 51, 53, 56, 57, R2sppd is R2sdpp

class OneCenterTwoElectronIntegralExpression
#include <OneCenterTwoElectronIntegralExpression.h>

This class allows the creation of instances containing, so to say, the analytical expression for the calculation of a one-center two-electron integral based on Slater-Type parameters.

class OneCenterTwoElectronIntegrals
#include <oneCenterTwoElectronIntegrals.h>

This class calculates all the unique one-center two-electron integrals for a given element.

class OverlapMatrix : public OverlapCalculator
#include <OverlapMatrix.h>

This class computes the whole overlap matrix and returns it in lower diagonal form.

The basis function overlap, as well as its first and second order derivatives with respect to the nuclear cartesian coordinates is calculated. It inherits from OverlapCalculator in order to make this class compatible with its polymorphic useage.

class TwoCenterIntegralContainer
#include <TwoCenterIntegralContainer.h>

This class contains smart pointers to two-center two-electron matrices for different atoms.

This class stores for each atom pair a shared pointer to a Global2c2eMatrix, containing the ERIs between the two centers.

class TwoElectronIntegralIndexes
#include <TwoElectronIntegralIndexes.h>

This class initializes and stores as a static array the indices of the charge distributions on one center playing a role in the multipole expansion.

The generation of the indices array happens only one time being it declared static. The generated charge distributions are located on one single center, as in the NDDO approximation charge distributions centered on two centers are neglected: \( \langle \phi_\mu\phi_\nu | \phi_\\lambda\phi_\sigma \rangle = \vardelta_{IJ}\vardelta_{KL}\langle \chi_\mu^I\chi_\nu^J | \chi_\\lambda^K\chi_\sigma^L\rangle \) In total, 40 possible charge distributions are possible: out of a symmetric 9x9 matrix, the diagonal and the terms that give rise to a charge distribution are retained. In the end 40 indices are stored. 5 charge distributions do not give rise to a multipole:

  • \( p_x d_{yz} \)

  • \( p_y d_{xz} \)

  • \( p_z d_{xy} \)

  • \( p_z d_{x^2 - y^2} \)

  • \( p_xy d_{x^2 - y^2} \)

class NDDODensityGuess : public DensityMatrixGuessCalculator
#include <NDDODensityGuess.h>

Implementation of DensityMatrixGuessCalculator for the NDDO methods.

class NDDOInitializer : public StructureDependentInitializer
#include <NDDOInitializer.h>

Settings for generic NDDO methods.

Reads the parameters and applies them to the system of interest. Depending on the method, the basis functions used can either be BasisFunction::sp (MNDO, AM1) of BasisFunction::spd (i.e. PM6, AM1*, MNDO/d) and it can have just atomic parameters (i.e. AM1, MNDO) or also diatomic parameters (i.e. PM6)

class OneElectronMatrix
#include <OneElectronMatrix.h>

This class generates the one-electron matrix H for semi-empirical methods.

class AtomicParameters
#include <AtomicParameters.h>

Class for the storage of atomic parameters in the semiempirical methods. (only those needed at runtime).

class DiatomicParameters
#include <DiatomicParameters.h>

Class for the storage of pairwise parameters in the PM6 method. (only those needed at runtime)

Subclassed by Scine::Sparrow::nddo::PM6DiatomicParameters

class ElementPairParameters
#include <ElementPairParameters.h>

This class holds the all the pointers to the parameters for element pairs.

class ElementParameters
#include <ElementParameters.h>

This class holds the pointers to the element-specific runtime parameters

class PM6DiatomicParameters : public Scine::Sparrow::nddo::DiatomicParameters
#include <PM6DiatomicParameters.h>

Class for the storage of pairwise parameters in the PM6 method. (only those needed at runtime)

class RawParameterProcessor
#include <RawParameterProcessor.h>

This class implements functions for the conversion between raw parameters published for PM6 and parameters useful at runtime.

struct gaussianRepulsionParameter
#include <RawParameters.h>

Structure to store the atomic information present in the parameter file for PM6.

class RawDiatomicParameters
#include <RawParameters.h>

Structure to store the diatomic information present in the parameter file for PM6.

class RawParametersContainer
#include <RawParametersContainer.h>

Class containing the PM6 parameters in a raw form.

class SlaterCondonParameters
#include <SlaterCondonParameters.h>

This class is the container for the Slater-Condon parameters.

class TwoElectronMatrix
#include <TwoElectronMatrix.h>

Class to generate the two-electron matrix G for semi-empirical methods. This class is parallelized with OpenMP.

namespace GeneralTypes

Enums

enum orb_t

enum indicating the possible orbitals

orb_t::s A s orbital orb_t::x A p_x orbital orb_t::y A p_y orbital orb_t::z A p_z orbital orb_t::x2y2 A d_{x^2-y^2} orbital orb_t::xz A d_{xz} orbital orb_t::z2 A d_{z^2} orbital orb_t::yz A d_{yz} orbital orb_t::xy A d_{xy} orbital

Values:

s
x
y
z
x2y2
xz
z2
yz
xy
enum twoElIntegral_t

enum listing all of the orbital pairs giving rise to a valid charge distribution

Values:

s_s
s_x
x_x
s_y
x_y
y_y
s_z
x_z
y_z
z_z
s_z2
s_xz
s_yz
s_x2y2
s_xy
x_z2
x_xz
x_x2y2
x_xy
y_z2
y_yz
y_x2y2
y_xy
z_z2
z_xz
z_yz
z2_z2
z2_xz
z2_yz
z2_x2y2
z2_xy
xz_xz
xz_yz
xz_x2y2
xz_xy
yz_yz
yz_x2y2
yz_xy
x2y2_x2y2
xy_xy
enum rotationOrbitalPair

Values:

s_s
x_x
x_y
x_z
y_x
y_y
y_z
z_x
z_y
z_z
x2y2_x2y2
x2y2_xz
x2y2_z2
x2y2_yz
x2y2_xy
xz_x2y2
xz_xz
xz_z2
xz_yz
xz_xy
z2_x2y2
z2_xz
z2_z2
z2_yz
z2_xy
yz_x2y2
yz_xz
yz_z2
yz_yz
yz_xy
xy_x2y2
xy_xz
xy_z2
xy_yz
xy_xy

Functions

int orbitalQN(orb_t o)

gets the orbital quantum number (“l”) of the argument, i.e.

0 for s, 1 for px,py,pz, 2 for all d orbitals.

Return

Returns 0 if the orbital type is invalid.

std::pair<orb_t, orb_t> separatePair(twoElIntegral_t t)

separates an orbital pair into its orbital components, throws InvalidOrbitalPairException if pair not valid

Return

a std::pair with orbital type elements, i.e. for s_x a pair with an orb_t::s and a orb_t::x

Parameters
  • an: orbital pair, i.e. s_x, corresponding to the s and the p_x orbitals

namespace multipole

Enums

enum ChargeDistance

Values:

d0
dp1
dm1
dp2
dm2
dps2
dms2
enum ChargeDistanceSeparation

Values:

d00
d01
d10
d02
d20
d0s2
ds20
p11
m11
p12
p21
m12
m21
p1s2
ps21
m1s2
ms21
p22
m22
p2s2
ps22
m2s2
ms22
ps2s2
ms2s2
enum multipolePair_t

enum listing the possible charge configurations of a multipole.

It is possible i.e. to infer the charge separation D_{sp1} from it. They are separated in monopole (l = 0), dipole (l = 1) and quadrupole (l = 2).

Values:

sp1
pd1
pp2
sd2
dd2
ss0
pp0
dd0
enum multipole_t

Multipole types used in the calculation of the ERI with the multipole expansion approximation.

Confusion might arise by the use of two merged formalisms: the one for just s and p orbitals and the one for s, p and d orbitals.

l: orbital quantum number m: magnetic quantum number \( M00 = M_{0,0} = q^I\) a monopole with l = 0, m = 0 \( M1m1 = M_{1,-1} = \mu_y \) a dipole in y direction with l = 1, m = -1 \( M10 = M_{1,0} = \mu_z \) a dipole in z direction with l = 1, m = 0 \( M11 = M_{1,1} = \mu_x \) a dipole in x direction with l = 1, m = 1 \( Qxx = Q_{x,x} \) a linear quadrupole in x direction with l = 2, m = 0 \( Qyy = Q_{y,y} \) a linear quadrupole in y direction with l = 2, m = 0 \( Qzz = Q_{z,z} \) a linear quadrupole in z direction with l = 2, m = 0 \( M2m2 = M_{2,-2} = Q_{x,y} \) a x,y square quadrupole with l = 2, m = -2 \( M2m1 = M_{2,-1} = Q_{y,z} \) a y,z square quadrupole with l = 2, m = -1 \( M20 = M_{2,0} = -\~{Q}_{x,z} - \frac{1}{2}\~{Q}_{x,y} \) a quadrupole with charges along each axis at \( \sqrt{2} \) distance from the origin with l = 2, m = 0 \( M21 = M_{2,1} = Q_{x,z}\) a x,z square quadrupole with l = 2, m = 1 \( M22 = M_{2,2} = \~{Q}_{x,y} \) a square quadrupole with charges along the x,y axes at \( \sqrt{2} \) distance from the origin with l = 2, m = 2 \( Qzx = \~{Q}_{z,x} = -\~{Q}_{x,z} \) a square quadrupole with charges along the x,z axes at \( \sqrt{2} \) distance from the origin with l = 2, m = 2

Values:

M00
Qxx
Qyy
Qzz
M1m1
M10
M11
M2m2
M2m1
M20
M21
M22
Qzx

Functions

template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::first> getDerivative<Utils::derivativeType::first>(orbPair_index_t op1, orbPair_index_t op2) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_atomic> getDerivative<Utils::derivativeType::second_atomic>(orbPair_index_t op1, orbPair_index_t op2) const
template<>
Utils::AutomaticDifferentiation::DerivativeType<Utils::derivativeType::second_full> getDerivative<Utils::derivativeType::second_full>(orbPair_index_t op1, orbPair_index_t op2) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::zero> expr<Utils::derivOrder::zero>(double f, double, double invsqrt) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::one> expr<Utils::derivOrder::one>(double f, double dz, double invsqrt) const
template<>
Utils::AutomaticDifferentiation::Value1DType<Utils::derivOrder::two> expr<Utils::derivOrder::two>(double f, double dz, double invsqrt) const
GeneralTypes::rotationOrbitalPair getRotPairType(GeneralTypes::orb_t o1, GeneralTypes::orb_t o2)

Given 2 orbitals, gives the corresponding orbital pair.

Throws InvalidOrbitalPairException() if the orbital types given are invalid. Order matters.

Return

a GeneralTypes::rotationOrbitalPair corresponding to the input orbitals

Parameters
  • o1: first orbital

  • o2: second orbital

int MQuantumNumber(multipole_t m)

Returns the magnetic quantum number m of a multipole, i.e.

-1 for a dipole in y direction, 0 for a dipole in z direction and 1 for a dipole in x direction. Throws InvalidMultipoleException() if the multipole is not valid.

int LQuantumNumber(multipole_t m)

Returns the orbital quantum number l of an orbital, i.e.

0 for s, 1 for p and 2 for d type orbitals. Throws InvalidMultipoleException() if the multipole is not a valid one.

multipolePair_t pairType(int l1, int l2, int l)

Function to infer the charge configuration of a multipole.

Return

throws InvalidQuantumNumbersException() if the quantum number is invalid. Otherwise returns a multipolePair_t.

Parameters
  • l1: the orbital quantum number of the first orbital

  • l2: the orbital quantum number of the second orbital

  • l: the multipole orbital quantum number

class ChargesInMultipoles
#include <ChargesInMultipoles.h>

This class defines the point charges of the different multipole.

class Global2c2eMatrix
#include <Global2c2eMatrix.h>

This class calculates the two-center two-electron integrals in the global coordinate system.

class Local2c2eIntegralCalculator
#include <Local2c2eIntegralCalculator.h>

This class is responsible for the calculation of the 2-center-2-electron integrals in the local coordinate system.

template<Utils::derivOrder O>
class Local2c2eMatrix
#include <Local2c2eMatrix.h>

This class creates the local two-center two-electron matrix for an atom pair.

class MultipoleCharge
#include <MultipoleCharge.h>

This class defines an object containing the position and charge of a point charge.

class MultipoleChargePair
#include <MultipoleChargePair.h>

This class stores the information about the distance between two point charges and about the product of their charges

class MultipoleMultipoleInteraction
#include <MultipoleMultipoleInteraction.h>

This header-only class performs the actual calculation of the multipole-multipole interaction.

First all the charge-charge configurations between two multipoles (MultipoleMultipoleTerm) are inferred and stored in a list, then they are calculated by calling MultipoleMultipoleTerm::calculate(…)

class MultipoleMultipoleInteractionContainer
#include <MultipoleMultipoleInteractionContainer.h>

This class keeps a list of terms of charge-charge-interactions for a pair of multipoles.

All the possible interactions can be stored in a 13x13 matrix, as there are 13 multipole types: 1 monopole, 3 dipoles, 3 linear quadrupoles, 3 square quadrupoles with charges between the axis, 3 square quadrupoles with charges along the axis.

class MultipoleMultipoleTerm
#include <MultipoleMultipoleTerm.h>

This header-only class defines an object for the calculation of an interaction between two charges in a multipole.

The total interaction between two electrons is approximated by a classical multipole expansion. This class defines an object that calculates the interaction between two charges of said multipoles. The charge interaction is calculated within the Klopman approximation to retrieve the correct 1-center, 2-electrons interaction in the limit of vanishing distance between the charges. The template functions allow for the analytical calculation of all the values up to the second derivative of the repulsion energy.

class VuvB
#include <VuvB.h>

This class returns the \(V_{\mu\nu,B}\) terms needed in semi-empirical methods.

class ZeroLocal2c2eIntegrals
#include <zeroLocal2c2eIntegrals.h>

Class that specifies which local two-center two-electron integrals are equal to zero in the semi-empirical approximation.

class ChargeSeparationParameter
#include <ChargeSeparationParameter.h>

Charge separation D of semi-empirical models. It describes the separation between two charges of opposite sign in a multipole.

The calculation of the charge separation is described in Thiel, Voityuk, Theor Chim Acta, 1992, 81, 391. The charge separation are stored ad a static c-array of size 5. The charge separations are ordered as in the nddo::multipole::multipolePair_t enum, i.e. sp1, pd1, pp2, sd2, dd2. ss0, pp0, dd0 are not present as there is no charge separation for them (they are monopoles).

class KlopmanParameter
#include <KlopmanParameter.h>

This class is the container for the Klopman-Ohno parameters used for the evaluation of the multipoles.

The Klopman-Ohno are used in the calculation of the two-center ERIs in the NDDO formalism. \( U\left(\Theta^{\mu\nu}_t, \Theta^{\lambda\sigma}_s \right) = \sum^{C_t}_{c=1}\sum^{C_s}_{d=1} \frac{q_c^Iq_d^J}{\sqrt{|r_c^I - r_d^J|^2 + \left( \theta_c^I(\chi_\mu^I\chi_\nu^I) + \theta_d^J(\chi_\lambda^J\chi_\sigma^J) \right)}} \)

namespace PM6Elements

Functions

unsigned int getQuantumNumberForSOrbital(Utils::ElementType e)
unsigned int getQuantumNumberForPOrbital(Utils::ElementType e)
unsigned int getQuantumNumberForDOrbital(Utils::ElementType e)
unsigned int getNumberOfAOs(Utils::ElementType e, BasisFunctions basisFunctions)
unsigned int getNumberOneCenterTwoElectronIntegrals(Utils::ElementType e, BasisFunctions basisFunctions)
double getCoreCharge(Utils::ElementType elementType)