14.1. Quick Guide

Before proceeding with this section, please ensure you have read the general introduction provided in General Remarks concerning Post-Hartree-Fock Calculations. This part of the Documentation builds upon it. The current version of PyBEST offers excited-state calculations with a CC (see The Restricted Coupled Cluster Module) and RpCCD (see The pCCD module) reference function using the Equation of Motion (EOM) formalism for electronic excitations (EE). The EOM module is explained in greater detail below. Specfically, the following EOM flavors are supported (for more details, see below)

14.1.1. Supported features

The EOM module supports spin-restricted orbitals and the DenseLinalgFactory and CholeskyLinalgFactory.

14.1.2. How to: REOM

This is a short introduction to the REOM module. More information on the input and output structure can be found in the following sections. Similar to the previous modules, we will assume the following names for all PyBEST objects

lf

A LinalgFactory instance (see Preliminaries).

occ_model

An Aufbau occupation model of the AufbauOccModel class

kin

The kinetic energy integrals

ne

The nucleus-electron attraction integrals

eri

The two-electron repulsion integrals

14.1.3. Various CC reference functions

PyBEST supports conventional CC flavors, like CCD or CCSD. For unconventional CC methods see the following sections.

14.1.3.1. EOM-CCS

EOM-CCS calculations are equivalent to CIS (configuration interaction singles) calculations (for excited states; see also RCI on top of RHF). Complete examples can be found in the following subsections (Example Python scripts).

We assume that you have performed a restricted CCS calculation, and the results are stored in the IOData container ccs_output. The code snippet below shows how to perform an REOMCCS calculation where only single excitations are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMCCS(lf, occ_model)
eom_output = eom(kin, ne, eri, ccs_output, nroot=3)

Note

The nroot keyword indicates the number of excited states (beyond the ground state).

Note

The order of arguments does not matter because PyBEST assigns them in an automatic fashion.

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-CCS.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (eigenvectors) for each state (column) with the ground state as the first column vector

t_1

The singles amplitudes of CCS

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.3.2. EOM-CCD

We assume that you have performed a restricted CCD calculation, and the results are stored in the IOData container ccd_output. The code snippet below shows how to perform an REOMCCD calculation where only double excitations are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, ccd_output, nroot=3)

Note

The nroot keyword indicates the number of excited states (beyond the ground state).

Note

The order of arguments does not matter because PyBEST assigns them in an automatic fashion.

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-CCD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (eigenvectors) for each state (column) with the ground state as the first column vector

t_2

The doubles amplitudes of CCD

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.3.3. EOM-CCSD

We assume that you have performed a restricted CCSD calculation, and the results are stored in the IOData container ccsd_output. The code snippet below shows how to perform an REOMCCSD calculation where only double excitations are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMCCSD(lf, occ_model)
eom_output = eom(kin, ne, eri, ccsd_output, nroot=3)

Note

The nroot keyword indicates the number of excited states (beyond the ground state).

Note

The order of arguments does not matter because PyBEST assigns them in an automatic fashion.

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-CCSD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (eigenvectors) for each state (column) with the ground state as the first column vector

t_1

The singles amplitudes of CCSD

t_2

The doubles amplitudes of CCSD

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.4. RpCCD reference function

This version of PyBEST supports EOM-pCCD calculations with two different excitation operators

  1. electron-pair excited states only (REOMpCCD)

  2. singly and electron-pair excited states (REOMpCCDS and REOMpCCDCCS)

Complete examples can be found in the following subsections (Example Python scripts).

14.1.4.1. EOM-pCCD: Electron-pair excitations

If you use this module, please cite [boguslawski2016a] and [boguslawski2017c].

We assume that you have performed a restricted pCCD calculation (either with or without orbital optimization), whose results are stored in the IOData container pccd_output (see The pCCD module). The code snippet below shows how to perform an REOMpCCD calculation where only electron-pair excitations are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMpCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, pccd_output, nroot=3)

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-pCCD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (that is, the eigenvectors) for each state (column) (including the ground states as a first column vector)

t_p

The pair amplitudes of pCCD

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.4.2. EOM-pCCD+S and EOM-pCCD-CCS: Single and electron-pair excitations

If you use this module, please cite [boguslawski2016a] and [boguslawski2017c].

To include single excitations as well, you can use the REOMpCCDS class, which performs an EOM-pCCD+S calculation. The overall input structure is equivalent to REOMpCCD,

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMpCCDS(lf, occ_model)
eom_output = eom(kin, ne, eri, pccd_output, nroot=3)

The corresponding IOData container (return value) will also contain all singly-excited states in the civ_ee attribute. All results are saved to the pybest-results/checkpoint_EOM-pCCD+S.h5 checkpoint file.

Note

The EOM-pCCD+S method is not size-intensive. Thus, the ground state energy eigenvalue will not be 0.0 au. All excitation energies (and hence also total energies) have to be adjusted wrt to the ground state energy shift.

The (rows of the) eigenvectors are sorted as follows: reference (pCCD) state, all singly-excited states, all doubly-excited states.

To recover size-intensivity, the REOMpCCDCCS flavor can be used, which exploits a pCCD-CCS reference function. To perform a REOMpCCDCCS, a pCCD-CCS reference function is required. Here, we assume that the corresponding output results are stored in the IOData container pccdccs_output. The code snippet below shows how to perform the corresponding EOM calculation for the 3 energetically lowest excited states.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMpCCDCCS(lf, occ_model)
eom_output = eom(kin, ne, eri, pccdccs_output, nroot=3)

14.1.5. LCC reference function

PyBEST supports EOM-LCC calculations with two different reference functions

  1. linearized coupled-cluster doubles (REOMLCCD)

  2. linearized coupled-cluster singles doubles (REOMLCCSD)

Complete examples can be found in the following subsections (Example Python scripts).

14.1.5.1. EOM-LCCD

We assume that you have performed a restricted LCCD calculation, whose results are stored in the IOData container lcc_output (see RHF reference function). The code snippet below shows how to perform an REOMLCCD calculation where only double excitations are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMLCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, lcc_output, nroot=3)

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-LCCD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (that is, the eigenvectors) for each state (column) (including the ground states as a first column vector)

t_2

The doubles amplitudes of LCCD

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.5.2. EOM-LCCSD

To perform an REOMLCCSD calculation where in addition singles excitations are described within the EOM formalism, you have to use the REOMLCCSD class,

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMLCCSD(lf, occ_model)
eom_output = eom(kin, ne, eri, lcc_output, nroot=3)

Again, we assume that the results of the LCCSD calculations are stored in the IOData container lcc_output. In addition to the EOM-LCCD output, the output container eom_output also contains information on the singles amplitudes,

t_1

The singles amplitudes of LCCSD

14.1.6. pCCD-LCC/fpLCC reference function

If you use these modules, please cite [boguslawski2019].

PyBEST also supports EOM-pCCD-LCC calculations with two different reference functions

  1. pCCD with a linearized coupled-cluster doubles correction (REOMpCCDLCCD)

  2. pCCD with a linearized coupled-cluster singles doubles correction (REOMpCCDLCCSD)

Complete examples can be found in the following subsections (Example Python scripts).

14.1.6.1. EOM-pCCD-LCCD

We assume that you have performed a restricted pCCD-LCCD calculation, whose results are stored in the IOData container pccdlcc_output (see RpCCD reference function). The code snippet below shows how to perform an REOMpCCDLCCD calculation where only double excitations (including electron pairs) are described within the EOM formalism. The number of target states is passed through the nroot keyword argument.

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMpCCDLCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, pccdlcc_output, nroot=3)

The results are returned as an IOData container, while all results are saved to the pybest-results/checkpoint_EOM-pCCD-LCCD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ee

The excitation energies (including the ground state, which equals 0.0)

civ_ee

The CI vector (that is, the eigenvectors) for each state (column) (including the ground states as a first column vector)

t_2

The doubles amplitudes of pCCD-LCCD

t_p

The electron-pair amplitudes of pCCD (fixed in pCCD-LCCD)

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory.

14.1.6.2. EOM-pCCD-LCCSD

To perform an REOMpCCDLCCSD calculation where, in addition, singles excitations are described within the EOM formalism, you have to use the REOMpCCDLCCSD class,

# Calculate 3 lowest-lying roots (excluding the ground state)
eom = REOMpCCDLCCSD(lf, occ_model)
eom_output = eom(kin, ne, eri, pccdlcc_output, nroot=3)

Again, we assume that the results of the pCCD-LCCSD calculations are stored in the IOData container pccdlcc_output. In addition to the EOM-pCCD-LCCD results, the output container eom_output also contains information on the singles amplitudes,

t_1

The singles amplitudes of pCCD-LCCSD

14.1.7. Defining a frozen core

By default, all core orbitals are frozen. If a frozen core has been selected in the CC reference calculation, the same orbitals have to be frozen in the chosen EOM flavor. To freeze some specific (occupied) orbitals, the number of frozen cores has to be specified when an instance of some REOMCC class is created. For instance, if you wish to freeze the first 4 (occupied) orbitals in an EOM-pCCD+S calculation, specify the ncore argument during the initialization of some occupation model class,

# Select 4 frozen core orbital
#-----------------------------
occ_model = AufbauOccModel(basis, ncore=4)
# Perform EOM-CC calculation, ncore is stored in occ_model
#---------------------------------------------------------
eom = REOMpCCDS(lf, occ_model)
eom_output = eom(kin, ne, eri, pccd_output)

This syntax is working for all EOM modules mentioned above.

14.1.8. Restart options

All EOM-CC calculations can be restarted from a previous calculation. To invoke a restart, use the restart keyword to specify the filename that contains the PyBEST-readable restart data (like PyBEST’s checkpoint files).

# Restart EOM-CCSD from a previous, unconverged calculation
# The restart data is assumed to be in PyBEST's default result directory with the default name
#-----------------------------------------------------------------------------------------
eom = REOMCCSD(lf, occ_model, alpha=1)
# Either perform a CCSD calculation or read the corresponding data from a restart file, e.g.,
# cc_output = IOData.from_file("pybest-results/checkpoint_RCCSD.h5")
eom_output = eom(kin, ne, eri, cc_output, restart="pybest-results/checkpoint_EOM-CCSD.h5")

Note

Restarts are only supported for the same molecular geometry. PyBEST does not check if the restart file contains the same geometry information.