Processing math: 100%

qml.qchem.electron_integrals

electron_integrals(mol, core=None, active=None)[source]

Return a function that computes the one- and two-electron integrals in the molecular orbital basis.

The one- and two-electron integrals are required to construct a molecular Hamiltonian in the second-quantized form

H=pqhpqcpcq+12pqrshpqrscpcqcrcs,

where c and c are the creation and annihilation operators, respectively, and hpq and hpqrs are the one- and two-electron integrals. These integrals can be computed by integrating over molecular orbitals ϕ as

hpq=ϕp(r)(2r2iZi|rRi|)ϕq(r)dr,

and

hpqrs=ϕp(r1)ϕq(r2)ϕr(r2)ϕs(r1)|r1r2|dr1dr2.

The molecular orbitals are constructed as a linear combination of atomic orbitals as

ϕi=νciνχν.

The one- and two-electron integrals can be written in the molecular orbital basis as

hpq=μνCpμhμνCνq,

and

hpqrs=μνρσCpμCqνhμνρσCρrCσs.

The hμν and hμνρσ terms refer to the elements of the core matrix and the electron repulsion tensor, respectively, and C is the molecular orbital expansion coefficient matrix.

Parameters
  • mol (Molecule) – the molecule object

  • core (list[int]) – indices of the core orbitals

  • active (list[int]) – indices of the active orbitals

Returns

function that computes the core constant and the one- and two-electron integrals

Return type

function

Example

>>> symbols  = ['H', 'H']
>>> geometry = np.array([[0.0, 0.0, 0.0], [0.0, 0.0, 1.0]], requires_grad = False)
>>> alpha = np.array([[3.42525091, 0.62391373, 0.1688554],
>>>                   [3.42525091, 0.62391373, 0.1688554]], requires_grad=True)
>>> mol = qml.qchem.Molecule(symbols, geometry, alpha=alpha)
>>> args = [alpha]
>>> electron_integrals(mol)(*args)
(1.0,
 array([[-1.3902192695e+00,  0.0000000000e+00],
        [-4.4408920985e-16, -2.9165331336e-01]]),
 array([[[[ 7.1443907755e-01, -2.7755575616e-17],
          [ 5.5511151231e-17,  1.7024144301e-01]],
         [[ 5.5511151231e-17,  1.7024144301e-01],
          [ 7.0185315353e-01,  6.6613381478e-16]]],
        [[[-1.3877787808e-16,  7.0185315353e-01],
          [ 1.7024144301e-01,  2.2204460493e-16]],
         [[ 1.7024144301e-01, -4.4408920985e-16],
          [ 6.6613381478e-16,  7.3883668974e-01]]]]))