14.3. Example Python scripts
Several complete examples can be found in the directory data/examples/eom
.
14.3.1. EOM-pCCD and EOM-pCCD+S calculations on the water molecule
This is a basic example of how to perform an EOM-pCCD and EOM-pCCD+S calculation in PyBEST. This script performs an RHF calculation followed by a pCCD calculation and two different EOM flavors with a pCCD reference function on the water molecule using the cc-pVDZ basis set.
from pybest import context
from pybest.ee_eom import REOMpCCD, REOMpCCDS
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.wrappers import RHF
###############################################################################
## Set up molecule, define basis set ##########################################
###############################################################################
# get the XYZ file from PyBEST's test data directory
fn_xyz = context.get_fn("test/water.xyz")
obasis = get_gobasis("cc-pvdz", fn_xyz)
###############################################################################
## Define Occupation model, expansion coefficients and overlap ################
###############################################################################
lf = DenseLinalgFactory(obasis.nbasis)
occ_model = AufbauOccModel(obasis)
orb_a = lf.create_orbital(obasis.nbasis)
olp = compute_overlap(obasis)
###############################################################################
## Construct Hamiltonian ######################################################
###############################################################################
kin = compute_kinetic(obasis)
ne = compute_nuclear(obasis)
eri = compute_eri(obasis)
external = compute_nuclear_repulsion(obasis)
###############################################################################
## Do a Hartree-Fock calculation ##############################################
###############################################################################
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
###############################################################################
## Do pCCD optimization #######################################################
###############################################################################
oopccd = RpCCD(lf, occ_model)
oopccd_output = oopccd(kin, ne, eri, hf_output)
###############################################################################
### Do REOM-pCCD calculation ##################################################
###############################################################################
eom = REOMpCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, oopccd_output, nroot=3)
###############################################################################
### Do REOM-pCCD+S calculation ################################################
###############################################################################
eom = REOMpCCDS(lf, occ_model)
eom_output = eom(kin, ne, eri, oopccd_output, nroot=3)
###############################################################################
### Do REOM-pCCD+S calculation using exact diagonalization ####################
###############################################################################
eom = REOMpCCDS(lf, occ_model)
eom_output = eom(kin, ne, eri, oopccd_output, nroot=3, davidson=False)
14.3.2. EOM-CCS and EOM-pCCD-CCS calculations on the water molecule
This is a basic example of how to perform an EOM-CCS and EOM-pCCD-CCS calculation in PyBEST. This script performs an RHF calculation followed by a CCS and pCCD calculation and two different EOM flavors, one with a CCS and one with a pCCD reference function on the water molecule using the cc-pVDZ basis set.
from pybest import context
from pybest.cc import RCCS, RpCCDCCS
from pybest.ee_eom import REOMCCS, REOMpCCDCCS
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.wrappers import RHF
###############################################################################
## Set up molecule, define basis set ##########################################
###############################################################################
# get the XYZ file from PyBEST's test data directory
fn_xyz = context.get_fn("test/water.xyz")
obasis = get_gobasis("cc-pvdz", fn_xyz)
###############################################################################
## Define Occupation model, expansion coefficients and overlap ################
###############################################################################
lf = DenseLinalgFactory(obasis.nbasis)
occ_model = AufbauOccModel(obasis)
orb_a = lf.create_orbital(obasis.nbasis)
olp = compute_overlap(obasis)
###############################################################################
## Construct Hamiltonian ######################################################
###############################################################################
kin = compute_kinetic(obasis)
ne = compute_nuclear(obasis)
eri = compute_eri(obasis)
external = compute_nuclear_repulsion(obasis)
###############################################################################
## Do a Hartree-Fock calculation ##############################################
###############################################################################
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
###############################################################################
## Do CCS optimization ########################################################
###############################################################################
ccs = RCCS(lf, occ_model)
ccs_output = ccs(kin, ne, eri, hf_output)
###############################################################################
### Do REOM-CCS calculation ###################################################
###############################################################################
eom = REOMCCS(lf, occ_model)
eom_output = eom(kin, ne, eri, ccs_output, nroot=3)
###############################################################################
## Do pCCD optimization #######################################################
###############################################################################
oopccd = RpCCD(lf, occ_model)
oopccd_output = oopccd(kin, ne, eri, hf_output)
###############################################################################
## Do pCCD-CCS optimization ###################################################
###############################################################################
pccdccs = RpCCDCCS(lf, occ_model)
pccdccs_output = pccdccs(kin, ne, eri, oopccd_output)
###############################################################################
### Do REOM-pCCD-CCS calculation ##############################################
###############################################################################
eom = REOMpCCDCCS(lf, occ_model)
eom_output = eom(kin, ne, eri, pccdccs_output, nroot=3)
14.3.3. EOM-LCCD and EOM-LCCSD calculations on the water molecule
This is a basic example of how to perform an EOM-LCCD and EOM-LCCSD calculation in PyBEST. This script performs an RHF calculation followed by an LCCD and LCCSD calculation and two different EOM flavors, one with an LCCD and one with an LCCSD reference function on the water molecule using the cc-pVDZ basis set.
from pybest import context
from pybest.cc import RLCCD, RLCCSD
from pybest.ee_eom import REOMLCCD, REOMLCCSD
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.wrappers import RHF
###############################################################################
## Set up molecule, define basis set ##########################################
###############################################################################
# get the XYZ file from PyBEST's test data directory
fn_xyz = context.get_fn("test/water.xyz")
obasis = get_gobasis("cc-pvdz", fn_xyz)
###############################################################################
## Define Occupation model, expansion coefficients and overlap ################
###############################################################################
lf = DenseLinalgFactory(obasis.nbasis)
occ_model = AufbauOccModel(obasis)
orb_a = lf.create_orbital(obasis.nbasis)
olp = compute_overlap(obasis)
###############################################################################
## Construct Hamiltonian ######################################################
###############################################################################
kin = compute_kinetic(obasis)
ne = compute_nuclear(obasis)
eri = compute_eri(obasis)
external = compute_nuclear_repulsion(obasis)
###############################################################################
## Do a Hartree-Fock calculation ##############################################
###############################################################################
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
###############################################################################
## Do LCCD optimization #######################################################
###############################################################################
lccd = RLCCD(lf, occ_model)
lccd_output = lccd(kin, ne, eri, hf_output)
###############################################################################
### Do REOM-LCCD calculation ##################################################
###############################################################################
eom = REOMLCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, lccd_output, nroot=3)
###############################################################################
## Do LCCSD optimization ######################################################
###############################################################################
lccsd = RLCCSD(lf, occ_model)
lccsd_output = lccsd(kin, ne, eri, hf_output)
###############################################################################
### Do REOM-LCCSD calculation #################################################
###############################################################################
eom = REOMLCCSD(lf, occ_model)
eom_output = eom(kin, ne, eri, lccsd_output, nroot=3)
14.3.4. EOM-pCCD-LCCD and EOM-pCCD-LCCSD calculations on the water molecule
This is a basic example of how to perform an EOM-pCCD-LCCD and EOM-pCCD-LCCSD calculation in PyBEST. This script performs an RHF and pCCD calculation followed by a pCCD-LCCD and pCCD-LCCSD calculation and two different EOM flavors, one with a pCCD-LCCD and one with a pCCD-LCCSD reference function on the water molecule using the cc-pVDZ basis set.
from pybest import context
from pybest.cc import RpCCDLCCD, RpCCDLCCSD
from pybest.ee_eom import REOMpCCDLCCD, REOMpCCDLCCSD
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.wrappers import RHF
###############################################################################
## Set up molecule, define basis set ##########################################
###############################################################################
# get the XYZ file from PyBEST's test data directory
fn_xyz = context.get_fn("test/water.xyz")
obasis = get_gobasis("cc-pvdz", fn_xyz)
###############################################################################
## Define Occupation model, expansion coefficients and overlap ################
###############################################################################
lf = DenseLinalgFactory(obasis.nbasis)
occ_model = AufbauOccModel(obasis)
orb_a = lf.create_orbital(obasis.nbasis)
olp = compute_overlap(obasis)
###############################################################################
## Construct Hamiltonian ######################################################
###############################################################################
kin = compute_kinetic(obasis)
ne = compute_nuclear(obasis)
eri = compute_eri(obasis)
external = compute_nuclear_repulsion(obasis)
###############################################################################
## Do a Hartree-Fock calculation ##############################################
###############################################################################
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
###############################################################################
## Do pCCD optimization #######################################################
###############################################################################
oopccd = RpCCD(lf, occ_model)
oopccd_output = oopccd(kin, ne, eri, hf_output)
###############################################################################
## Do pCCD-LCCD optimization ##################################################
###############################################################################
lccd = RpCCDLCCD(lf, occ_model)
lccd_output = lccd(kin, ne, eri, oopccd_output)
###############################################################################
### Do REOM-pCCD-LCCD calculation #############################################
###############################################################################
eom = REOMpCCDLCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, lccd_output, nroot=3)
###############################################################################
## Do pCCD-LCCSD optimization #################################################
###############################################################################
lccsd = RpCCDLCCSD(lf, occ_model)
lccsd_output = lccsd(kin, ne, eri, oopccd_output)
###############################################################################
### Do REOM-pCCD-LCCSD calculation ############################################
###############################################################################
eom = REOMpCCDLCCSD(lf, occ_model)
eom_output = eom(kin, ne, eri, lccsd_output, nroot=3)
14.3.5. EOM-pCCD and EOM-pCCD+S calculations on the water molecule using a frozen core
This is a basic example on how to perform a EOM-pCCD and EOM-pCCD+S calculation in PyBEST. This script performs an RHF calculation followed by a pCCD calculation and two different EOM flavors with a pCCD reference function on the water molecule using the cc-pVDZ basis set and freezing the lowest orbital.
from pybest import context
from pybest.ee_eom import REOMpCCD, REOMpCCDS
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.wrappers import RHF
###############################################################################
## Set up molecule, define basis set ##########################################
###############################################################################
# get the XYZ file from PyBEST's test data directory
fn_xyz = context.get_fn("test/water.xyz")
obasis = get_gobasis("cc-pvdz", fn_xyz)
###############################################################################
## Define Occupation model, expansion coefficients and overlap ################
###############################################################################
lf = DenseLinalgFactory(obasis.nbasis)
occ_model = AufbauOccModel(obasis, ncore=1)
orb_a = lf.create_orbital(obasis.nbasis)
olp = compute_overlap(obasis)
###############################################################################
## Construct Hamiltonian ######################################################
###############################################################################
kin = compute_kinetic(obasis)
ne = compute_nuclear(obasis)
eri = compute_eri(obasis)
external = compute_nuclear_repulsion(obasis)
###############################################################################
## Do a Hartree-Fock calculation ##############################################
###############################################################################
hf = RHF(lf, occ_model)
hf_output = hf(kin, ne, eri, external, olp, orb_a)
###############################################################################
## Do pCCD optimization #######################################################
###############################################################################
oopccd = RpCCD(lf, occ_model)
oopccd_output = oopccd(kin, ne, eri, hf_output)
################################################################################
### Do REOM-pCCD calculation ###################################################
################################################################################
eom = REOMpCCD(lf, occ_model)
eom_output = eom(kin, ne, eri, oopccd_output, nroot=3)
################################################################################
### Do REOM-pCCD+S calculation #################################################
################################################################################
eom = REOMpCCDS(lf, occ_model)
eom_output = eom(kin, ne, eri, oopccd_output, nroot=3)