18.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 Reversed Spin Flip (RSF) calculations with RCCD and RCCSD (see The restricted coupled-cluster module), RLCCD and RLCCSD (see The linearized coupled cluster module), RfpCCD and RfpCCSD (see The restricted frozen-pair coupled-cluster module), and RSFfpLCCD and RSFfpLCCSD (see The restricted frozen-pair coupled-cluster module) reference functions using the Equation of Motion (EOM) formalism. The RSF module is explained in greater detail below.

18.1.1. Supported features

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

18.1.2. How to: RSF

This is a short introduction to the RSF 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 RSF modules, only the spin projection is defined, that is, only the number of excess \(\alpha\) electrons over \(\beta\) electrons is defined.

18.1.2.1. RCCD-type reference function

This version of PyBEST supports RSF-CCD, RSF-LCCD, RSF-RfpCCD and RSF-RfpLCCD calculations for two different states

  • RSF-CCD (RSFCCD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 To Be Announced (TBA)

  • RSF-LCCD (RSFLCCD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

  • RSF-RfpCCD (RSFfpCCD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

  • RSF-RfpLCCD (RSFfpLCCD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

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

18.1.2.2. RCCSD-type reference function

This version of PyBEST supports RSF-CCSD, RSF-LCCSD, RSF-RfpCCSD and RSF-RfpLCCSD calculations for two different states

  • RSF-CCSD (RSFCCSD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

  • RSF-LCCSD (RSFLCCSD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

  • RSF-RfpCCSD (RSFfpCCSD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

  • RSF-RfpLCCSD (RSFfpLCCSD):
    1. quintet state (S_z=2) with Ms=2

    2. triplet state (S_z=1) with Ms=1 (TBA)

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

18.1.2.2.1. RSF-CCD/RSF-CCSD/RSF-LCCD/RSF-LCCSD Triplet and Quintet states

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

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 RSFCCSD calculation. The current version supports the optimization of triplet (TBA) and quintet states. The state can be optimized by setting the keyword argument alpha: 2 for triplet and 4 for quintet when creating an instance of RSFCCSD. The number of targeted states is passed through the nroot keyword argument in the function call.

# Do RSF-CCSD calculation for 4 unpaired electrons (quintet)
rsf = RSFCCSD(lf, occ_model, alpha=4)
rsf_output = rsf(kin, ne, eri, ccsd_output, nroot=3)

The results are returned as a IOData container, while all results are saved to the pybest-results/checkpoint_RSF-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 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.

Note

All remaining flavors (RSFCCD, RSFLCCD, and RSFLCCSD) are accessible in a similar way. To choose a different RSF-CC flavor, simply change the RSFCCSD class to the desired CC model,

  • RSFCCD for RSF-CCD

  • RSFLCCD for RSF-LCCD

  • RSFLCCSD for RSF-LCCSD

18.1.2.2.2. RSF-fpCCD/RSF-fpCCSD/RSF-fpLCCD/RSF-fpLCCSD Triplet and Quintet states

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

We assume that you have performed a restricted fpCC 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 RSFfpCCSD calculation. The current version supports the optimization of triplet (TBA) and quintet states. The state can be optimized by setting the keyword argument alpha: 2 for triplet and 4 for quintet when creating an instance of RSFfpCCSD. The number of targeted states is passed through the nroot keyword argument in the function call.

# Do RSF-fpCCSD calculation for 4 unpaired electrons (quintet)
rsf_fp = RSFfpCCSD(lf, occ_model, alpha=4)
rsf_fp_output = rsf_fp(kin, ne, eri, fpccsd_output, nroot=4)

The results are returned as a IOData container, while all results are saved to the pybest-results/checkpoint_RSF-EOM-fpCCSD.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 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.

Note

All remaining flavors (RSFfpCCD, RSFfpLCCD, and RSFfpLCCSD) are accessible in a similar way. To choose a different RSF-fpCC flavor, simply change the RSFfpCCSD class to the desired CC model,

  • RSFfpCCD for RSF-fpCCD

  • RSFfpLCCD for RSF-fpLCCD

  • RSFfpLCCSD for RSF-fpLCCSD