17.1. Quick Guide

Before reading this section, please read the general introduction mentioned in General Remarks concerning Post-Hartree-Fock Calculations. This part of the Documentation builds upon it. The current version of PyBEST offers Electron Attachment (EA) calculations with a RpCCD (see The pCCD module) reference function and various RCC (see The restricted coupled-cluster module) and fpCC (see The restricted frozen-pair coupled-cluster module) models using the Equation of Motion (EOM) formalism. The EA module is explained in greater detail below.

17.1.1. Supported features

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

17.1.2. How to: REA

This is a short introduction to the REA 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

Note

In all EA modules, only the spin projection is defined, that is, only the number of excess \(\alpha\) electrons over \(\beta\) electrons is defined.

17.1.2.1. RHF reference function

This version of PyBEST supports single EA-CCD, EA-CCSD, EA-LCCD, and EA-LCCSD with two different excitation operators

  • EA-CCD (REACCD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-CCSD (REACCSD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-LCCD (REALCCD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-LCCSD (REALCCSD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

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

17.1.2.2. RpCCD reference function

This version of PyBEST supports single EA-pCCD, EA-fpCCD, EA-fpCCSD, EA-fpLCCD, and EA-fpLCCSD and double (D)EA-pCCD calculations with two different excitation operators

  • EA-pCCD (REApCCD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-fpCCD (REAfpCCD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-fpCCSD (REAfpCCSD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-fpLCCD (REAfpLCCD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • EA-fpLCCSD (REAfpLCCSD):
    1. 1 particle (1p) operators

    2. 2 particle, 1 hole (2p1h) operators

  • DEA-pCCD (RDEApCCD):
    1. 2 particle (2p) operators

    2. 3 particle, 1 hole (3p1h) operators

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

17.1.2.2.1. EA-CCD/EA-CCSD/EA-LCCD/EA-LCCSD: Doublet and quartet states

If you use this module, please cite [behjou2025].

We assume that you have performed a restricted CC calculation, whose results are stored in the IOData container cc_output (see The linearized coupled cluster module for LCC flavors and The restricted coupled-cluster module for conventional CC). The code snippet below shows how to perform a REACCSD calculation. The current version supports the optimization of doublet and quartet states. Both states can be optimized by setting the keyword argument alpha to 1 when creating an instance of REACCSD. The number of targeted states is passed through the nroot keyword argument in the function call.

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5)
# Both doublet and quartet states can be targeted
ea = REACCSD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, cc_output, nroot=3)

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

e_ea

The electron affinities in \(E_h\)

civ_ea_alpha

The CI vectors (that is, the eigenvectors) for each state (column) for a given alpha value

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory. Since for one spin projection \(S_z\), several spin multiplets are accessible, PyBEST suggests the multiplicity for each targeted state in the output file.

By default, REACCSD includes up to 2p1h terms. For alpha=1, the 2p1h terms can be neglected and only the 1p terms are considered during the diagonalization. The number of particle/hole operators is specified by the keyword argument nparticle. The following code snippet shows how to perform a REACCSD calculation for 1 unpaired electron, including only 1p terms,

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5) and 1p terms
# Both doublet and quartet states can be targeted
ea = REACCSD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, cc_output, nroot=3, nparticle=1)

Note

All remaining flavors (REACCD, REALCCD, and REALCCSD) are accessible in a similar way. To choose a different EA-CC flavor, simply change the REACCSD class to the desired CC model,

  • REACCD for EA-CCD

  • REALCCD for EA-LCCD

  • REALCCSD for EA-LCCSD

17.1.2.2.2. EA-pCCD: Doublet and quartet states

If you use this module, please cite [galynska2024b].

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 a REApCCD calculation. The current version supports the optimization of doublet and quartet states and just the quartet states. The former can be optimized by setting the keyword argument alpha to 1 when creating an instance of REApCCD, while the latter are accessible by choosing alpha=3. The number of targeted states is passed through the nroot keyword argument in the function call.

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5)
# Both doublet and quartet states can be targeted
ea = REApCCD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, pccd_output, nroot=3)

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

e_ea

The electron affinities in \(E_h\)

civ_ea_alpha

The CI vectors (that is, the eigenvectors) for each state (column) for a given alpha value

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory. Since for one spin projection \(S_z\), several spin multiplets are accessible, PyBEST suggests the multiplicity for each targeted state in the output file.

By default, REApCCD includes up to 2p1h terms. For alpha=1, the 2p1h terms can be neglected and only the 1p terms are considered during the diagonalization. The number of particle/hole operators is specified by the keyword argument nparticle. The following code snippet shows how to perform a REApCCD calculation for 1 unpaired electron including only 1p terms,

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5) and 1p terms
# Only doublet states can be targeted
ea = REApCCD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, pccd_output, nroot=3, nparticle=1)

Note

For three unpaired electrons, only 2p1h terms are supported. Specifying nparticle=1 will raise an error.

Note

A spin-free implementation is also available. It might help to eliminate convergence difficulties in some cases. If set to True, only doublet state can be targeted.

# Calculate 3 electron-attached states for 1 unpaired electron (S_z=0.5)
# Only doublet states can be targeted
ea = REApCCD(lf, occ_model, alpha=1, spinfree=True)
ea_output = ea(kin, ne, eri, cc_output, nroot=3)

17.1.2.2.3. EA-fpCCD/EA-fpCCSD/EA-fpLCCD/EA-fpLCCSD: Doublet and quartet states

If you use this module, please cite [behjou2025].

We assume that you have performed a restricted frozen-pair (fp)CC calculation, whose results are stored in the IOData container cc_output (see The linearized coupled cluster module for LCC flavors on top of pCCD and The restricted frozen-pair coupled-cluster module for supported fpCC models). The code snippet below shows how to perform a REAfpCCSD calculation. The current version supports the optimization of doublet and quartet states. Both states can be optimized by setting the keyword argument alpha to 1 when creating an instance of REAfpCCSD. The number of targeted states is passed through the nroot keyword argument in the function call.

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5)
# Both doublet and quartet states can be targeted
ea = REAfpCCSD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, cc_output, nroot=3)

The results are returned as a IOData container, while all results are saved to the pybest-results/checkpoint_EA-EOM-fpCCSD.h5 file. Specifically, the IOData container contains (amongst others) the following attributes

e_ea

The electron affinities in \(E_h\)

civ_ea_alpha

The CI vectors (that is, the eigenvectors) for each state (column) for a given alpha value

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory. Since for one spin projection \(S_z\), several spin multiplets are accessible, PyBEST suggests the multiplicity for each targeted state in the output file.

By default, REAfpCCSD includes up to 2p1h terms. For alpha=1, the 2p1h terms can be neglected and only the 1p terms are considered during the diagonalization. The number of particle/hole operators is specified by the keyword argument nparticle. The following code snippet shows how to perform a REAfpCCSD calculation for 1 unpaired electron including only 1p terms,

# Calculate 3 lowest electron-attached states for 1 unpaired electron (S_z=0.5) and 1p terms
# Both doublet and quartet states can be targeted
ea = REAfpCCSD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, cc_output, nroot=3, nparticle=1)

Note

All remaining flavors (REAfpCCD, REAfpLCCD, and REAfpLCCSD) are accessible in a similar way. To choose a different EA-CC flavor, simply change the REAfpCCSD class to the desired CC model,

  • REAfpCCD for EA-fpCCD

  • REAfpLCCD for EA-fpLCCD

  • REAfpLCCSD for EA-fpLCCSD

17.1.2.2.4. DEA-pCCD: Singlet, triplet, and quintet states

If you use this module, please cite [galynska2024b].

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 a RDEApCCD calculation. The current version supports the optimization of singlet, triplet, and quintet states, triplet and quintet states, or just quintet states. The first case can be optimized by setting the keyword argument alpha to 0 when creating an instance of RDEApCCD, the second case is accessible by choosing alpha=2, while the last spin multiplicity can be chosen by setting alpha=4. The number of targeted states is passed through the nroot keyword argument in the function call.

# Calculate 3 lowest double electron-attached states for 0 unpaired electron (S_z=0.0)
# Singlet, triplet, and quintet states can be targeted
dea = RDEApCCD(lf, occ_model, alpha=0)
dea_output = dea(kin, ne, eri, pccd_output, nroot=3)

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

e_ea

The ionization energies in \(E_h\)

civ_ea_alpha

The CI vectors (that is, the eigenvectors) for each state (column) for a given alpha value

The eigenvalues and eigenvectors are stored as numpy arrays, not as instances of the LinalgFactory. Since for one spin projection \(S_z\), several spin multiplets are accessible, PyBEST suggests the multiplicity for each targeted state in the output file.

By default, RDEApCCD includes up to 3p1h terms. For alpha=0 and alpha=2, the 3p1h terms can be neglected so that only the 2p terms are considered during the diagonalization. The number of particle/hole operators is specified by the keyword argument nparticle. The following code snippet shows how to perform a RDEApCCD calculation for 0 unpaired electron including only 2p terms,

# Calculate 3 lowest double electron-attached states for 0 unpaired electron (S_z=0.0) and 2p terms
# Singlet and triplet states can be targeted
dea = RDEApCCD(lf, occ_model, alpha=0)
dea_output = dea(kin, ne, eri, pccd_output, nroot=3, nparticle=2)

Note

For four unpaired electrons (qunitet states), only 3p1h terms are supported. Specifying nparticle=2 will raise an error.

17.1.2.3. 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 EA flavor. To freeze some specific (occupied) orbitals, the number of frozen cores has to be specified when an instance of some OccupationModel class is created. For instance, if you wish to freeze the first 4 (occupied) orbitals in an EA-pCCD 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 EA-CC calculation for 1 unpaired electron (S_z=0.5),
# ncore is stored in occ_model
#-------------------------------------------------------------
ea = REApCCD(lf, occ_model, alpha=1)
ea_output = ea(kin, ne, eri, pccd_output)

This syntax is working for all EA modules mentioned above.

17.1.2.4. Restart options

Restart options are not supported yet.