18.3. Example Python scripts
Several complete examples can be found in the directory data/examples/rsf
.
18.3.1. RSF-RCCD calculations on the C atom with 4 unpaired electrons
This is a basic example of how to perform an RSF-RCCD calculation in PyBEST. This script performs an RCCD calculation and one RSF flavor with an RCCD reference function on the C molecule using the cc-pVDZ basis set. Only two unpaired electrons are considered (quintet state).
#!/usr/bin/env python3
from pybest import context
from pybest.cc import RCCD
from pybest.gbasis import (
compute_cholesky_eri,
compute_kinetic,
compute_nuclear,
compute_nuclear_repulsion,
compute_overlap,
get_gobasis,
)
from pybest.linalg import CholeskyLinalgFactory
from pybest.occ_model import AufbauOccModel
from pybest.rsf_eom import RSFCCD
from pybest.wrappers import RHF
#
# Set up molecule, define basis set
#
# get the XYZ file from PyBEST's test data directory
c_xyz = context.get_fn("test/c.xyz")
basis = get_gobasis("cc-pvdz", c_xyz)
#
# Define Occupation model, orbitals and overlap
#
lf = CholeskyLinalgFactory(basis.nbasis)
orb_a = lf.create_orbital(basis.nbasis)
olp = compute_overlap(basis)
occ_model = AufbauOccModel(basis)
#
# Construct Hamiltonian
#
kin = compute_kinetic(basis)
ne = compute_nuclear(basis)
eri = compute_cholesky_eri(basis)
external = compute_nuclear_repulsion(basis)
#
# Do a Hartree-Fock calculation
#
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
#
# Do CCD optimization
#
ccd = RCCD(lf, occ_model)
ccd_output = ccd(kin, ne, eri, hf_output)
#
# Do RSF-CCD calculation with 4 unpaired electrons
#
rsf = RSFCCD(lf, occ_model, alpha=4)
rsf_output = rsf(kin, ne, eri, ccd_output, nroot=3)
18.3.2. RSF-RLCCSD calculations on the C atom with 4 unpaired electrons
This is a basic example of how to perform an RSF-RLCCSD calculation in PyBEST. This script performs an RLCCSD calculation and one RSF flavor with an RLCCSD reference function on the C molecule using the cc-pVDZ basis set. Only two unpaired electrons are considered (quintet state).
from pybest import context
from pybest.cc import RLCCSD
from pybest.gbasis import (
compute_eri,
compute_kinetic,
compute_nuclear,
compute_nuclear_repulsion,
compute_overlap,
get_gobasis,
)
from pybest.linalg import DenseLinalgFactory
from pybest.occ_model import AufbauOccModel
from pybest.rsf_eom import RSFLCCD
from pybest.wrappers import RHF
#
# Set up molecule, define basis set
#
# get the XYZ file from PyBEST's test data directory
c_xyz = context.get_fn("test/c.xyz")
basis = get_gobasis("cc-pvdz", c_xyz)
#
# Define Occupation model, orbitals and overlap
#
lf = DenseLinalgFactory(basis.nbasis)
orb_a = lf.create_orbital(basis.nbasis)
olp = compute_overlap(basis)
occ_model = AufbauOccModel(basis)
#
# Construct Hamiltonian
#
kin = compute_kinetic(basis)
ne = compute_nuclear(basis)
eri = compute_eri(basis)
external = compute_nuclear_repulsion(basis)
#
# Do a Hartree-Fock calculation
#
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
#
# Do LCCSD optimization
#
lccd = RLCCSD(lf, occ_model)
lccd_output = lccd(kin, ne, eri, hf_output)
#
# Do RSF-LCCSD calculation with 4 unpaired electrons
#
rsf_L = RSFLCCD(lf, occ_model, alpha=4)
rsf_L_output = rsf_L(kin, ne, eri, lccd_output, nroot=8)
18.3.3. RSF-RfpCCSD calculations on the C atom with 4 unpaired electrons
This is a basic example of how to perform an RSF-RfpCCSD calculation in PyBEST. This script performs an RpCCD and then RfpCCSD calculation and one RSF flavor with an RfpCCSD reference function on the C molecule using the cc-pVDZ basis set. Only two unpaired electrons are considered (quintet state).
from pybest import context
from pybest.cc import RfpCCSD
from pybest.gbasis import (
compute_eri,
compute_kinetic,
compute_nuclear,
compute_nuclear_repulsion,
compute_overlap,
get_gobasis,
)
from pybest.geminals import RpCCD
from pybest.linalg import DenseLinalgFactory
from pybest.occ_model import AufbauOccModel
from pybest.rsf_eom import RSFfpCCSD
from pybest.wrappers import RHF
#
# Set up molecule, define basis set
#
# get the XYZ file from PyBEST's test data directory
c_xyz = context.get_fn("test/c.xyz")
basis = get_gobasis("cc-pvdz", c_xyz)
#
# Define Occupation model, orbitals and overlap
#
lf = DenseLinalgFactory(basis.nbasis)
orb_a = lf.create_orbital(basis.nbasis)
olp = compute_overlap(basis)
occ_model = AufbauOccModel(basis)
#
# Construct Hamiltonian
#
kin = compute_kinetic(basis)
ne = compute_nuclear(basis)
eri = compute_eri(basis)
external = compute_nuclear_repulsion(basis)
#
# Do a Hartree-Fock calculation
#
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
#
# Do pCCD optimization
#
pccd = RpCCD(lf, occ_model)
pccd_output = pccd(kin, ne, eri, hf_output)
#
# Do fpCCSD optimization
#
fpccsd = RfpCCSD(lf, occ_model)
fpccsd_output = fpccsd(kin, ne, eri, pccd_output)
#
# Do RSF-fpCCSD calculation with 4 unpaired electrons
#
rsf_fp = RSFfpCCSD(lf, occ_model, alpha=4)
rsf_fp_output = rsf_fp(kin, ne, eri, fpccsd_output, nroot=8)
18.3.4. RSF-RfpLCCSD calculations on the C atom with 4 unpaired electrons
This is a basic example of how to perform an RSF-RfpLCCSD calculation in PyBEST. This script performs an RpCCD and then RpCCDLCCSD calculation and one RSF flavor with an RfpLCCSD reference function on the C molecule using the cc-pVDZ basis set. Only two unpaired electrons are considered (quintet state).
from pybest import context
from pybest.cc import RpCCDLCCSD
from pybest.gbasis import (
compute_eri,
compute_kinetic,
compute_nuclear,
compute_nuclear_repulsion,
compute_overlap,
get_gobasis,
)
from pybest.geminals import RpCCD
from pybest.linalg import DenseLinalgFactory
from pybest.occ_model import AufbauOccModel
from pybest.rsf_eom import RSFfpLCCSD
from pybest.wrappers import RHF
#
# Set up molecule, define basis set
#
# get the XYZ file from PyBEST's test data directory
c_xyz = context.get_fn("test/c.xyz")
basis = get_gobasis("cc-pvdz", c_xyz)
#
# Define Occupation model, orbitals and overlap
#
lf = DenseLinalgFactory(basis.nbasis)
orb_a = lf.create_orbital(basis.nbasis)
olp = compute_overlap(basis)
occ_model = AufbauOccModel(basis)
#
# Construct Hamiltonian
#
kin = compute_kinetic(basis)
ne = compute_nuclear(basis)
eri = compute_eri(basis)
external = compute_nuclear_repulsion(basis)
#
# Do a Hartree-Fock calculation
#
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
#
# Do pCCD optimization
#
pccd = RpCCD(lf, occ_model)
pccd_output = pccd(kin, ne, eri, hf_output)
#
# Do fpLCCSD optimization
#
fplccsd = RpCCDLCCSD(lf, occ_model)
fplccsd_output = fplccsd(kin, ne, eri, pccd_output)
#
# Do RSF-fpLCCSD calculation with 4 unpaired electrons
#
rsf_fpl = RSFfpLCCSD(lf, occ_model, alpha=4)
rsf_fpl_output = rsf_fpl(kin, ne, eri, fplccsd_output, nroot=8)