CCWavefunction

class psi4.core.CCWavefunction

Bases: psi4.core.Wavefunction

docstring

Methods Summary

get_amplitudes(self) Get dict of converged T amplitudes

Methods Documentation

get_amplitudes(self: psi4.core.CCWavefunction) → Dict[str, psi4.core.Matrix]

Get dict of converged T amplitudes

amps : dict (spacestr, SharedMatrix)

spacestr is a description of the amplitude set using the following conventions.

I,J,K -> alpha occupied i,j,k -> beta occupied A,B,C -> alpha virtual a,b,c -> beta virtual

The following entries are stored in the amps, depending on the reference type

RHF: “tIA”, “tIjAb” UHF: tIA, tia, tIjAb, tIJAB, tijab ROHF: tIA, tia, tIjAb, tIJAB, tijab

Examples

RHF T1 diagnostic = sqrt(sum_ia (T_ia * T_ia)/nelec) >>> mol = “”” … 0 1 … Ne 0.0 0.0 0.0 … symmetry c1”“” >>> e, wfn = psi4.energy(“CCSD/cc-pvdz”, return_wfn=True) >>> t1 = wfn.get_amplitudes()[‘tia’].to_array() >>> t1_diagnostic = np.sqrt(np.dot(t1.ravel(),t1.ravel())/ (2 * wfn.nalpha()) >>> t1_diagnostic == psi4.variable(“CC T1 DIAGNOSTIC”) True

Warning

Symmetry free calculations only (nirreps > 1 will cause error)

Warning

No checks that the amplitudes will fit in core. Do not use for proteins

Ca(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Alpha Orbitals.

Ca_subset(self: psi4.core.Wavefunction, arg0: str, arg1: str) → psi4.core.Matrix

Returns the requested Alpha Orbital subset.

Cb(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Beta Orbitals.

Cb_subset(self: psi4.core.Wavefunction, arg0: str, arg1: str) → psi4.core.Matrix

Returns the requested Beta Orbital subset.

Da(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Alpha Density Matrix.

Da_subset(self: psi4.core.Wavefunction, arg0: str) → psi4.core.Matrix

Returns the requested Alpha Density subset.

Db(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Beta Density Matrix.

Db_subset(self: psi4.core.Wavefunction, arg0: str) → psi4.core.Matrix

Returns the requested Beta Density subset.

Fa(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Alpha Fock Matrix.

Fa_subset(self: psi4.core.Wavefunction, arg0: str) → psi4.core.Matrix

Returns the Alpha Fock Matrix in the requested basis (AO,SO).

Fb(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Beta Fock Matrix.

Fb_subset(self: psi4.core.Wavefunction, arg0: str) → psi4.core.Matrix

Returns the Beta Fock Matrix in the requested basis (AO,SO).

H(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the ‘Core’ Matrix (Potential + Kinetic) Integrals.

PCM_enabled(self: psi4.core.Wavefunction) → bool

Whether running a PCM calculation

S(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the One-electron Overlap Matrix.

X(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Lagrangian Matrix.

alpha_orbital_space(self: psi4.core.Wavefunction, arg0: str, arg1: str, arg2: str) → psi4.core.OrbitalSpace

docstring

aotoso(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Atomic Orbital to Symmetry Orbital transformer.

array_variable(self: psi4.core.Wavefunction, arg0: str) → psi4.core.Matrix

Returns copy of the requested (case-insensitive) Matrix QC variable.

array_variables(self: psi4.core.Wavefunction) → Dict[str, psi4.core.Matrix]

Returns the dictionary of all Matrix QC variables.

arrays()
atomic_point_charges(self: psi4.core.Wavefunction) → psi4.core.Vector

Returns the set atomic point charges.

basis_projection(self: psi4.core.Wavefunction, arg0: psi4.core.Matrix, arg1: psi4.core.Dimension, arg2: psi4.core.BasisSet, arg3: psi4.core.BasisSet) → psi4.core.Matrix

Projects a orbital matrix from one basis to another.

basisset(self: psi4.core.Wavefunction) → psi4.core.BasisSet

Returns the current orbital basis.

beta_orbital_space(self: psi4.core.Wavefunction, arg0: str, arg1: str, arg2: str) → psi4.core.OrbitalSpace

docstring

static build(mol, basis=None)
c1_deep_copy(self: psi4.core.Wavefunction, basis: psi4.core.BasisSet) → psi4.core.Wavefunction

Returns a new wavefunction with internal data converted to C_1 symmetry, using pre-c1-constructed BasisSet basis

compute_energy(self: psi4.core.Wavefunction) → float

Computes the energy of the Wavefunction.

compute_gradient(self: psi4.core.Wavefunction) → psi4.core.Matrix

Computes the gradient of the Wavefunction

compute_hessian(self: psi4.core.Wavefunction) → psi4.core.Matrix

Computes the Hessian of the Wavefunction.

deep_copy(self: psi4.core.Wavefunction, arg0: psi4.core.Wavefunction) → None

Deep copies the internal data.

del_array_variable(self: psi4.core.Wavefunction, arg0: str) → int

Removes the requested (case-insensitive) Matrix QC variable.

del_scalar_variable(self: psi4.core.Wavefunction, arg0: str) → int

Removes the requested (case-insensitive) double QC variable.

del_variable(key)
density_fitted(self: psi4.core.Wavefunction) → bool

Returns whether this wavefunction was obtained using density fitting or not.

doccpi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of doubly occupied orbitals per irrep.

efzc(self: psi4.core.Wavefunction) → float

Returns the frozen-core energy

energy(self: psi4.core.Wavefunction) → float

Returns the Wavefunction’s energy.

epsilon_a(self: psi4.core.Wavefunction) → psi4.core.Vector

Returns the Alpha Eigenvalues.

epsilon_a_subset(self: psi4.core.Wavefunction, arg0: str, arg1: str) → psi4.core.Vector

Returns the requested Alpha Eigenvalues subset.

epsilon_b(self: psi4.core.Wavefunction) → psi4.core.Vector

Returns the Beta Eigenvalues.

epsilon_b_subset(self: psi4.core.Wavefunction, arg0: str, arg1: str) → psi4.core.Vector

Returns the requested Beta Eigenvalues subset.

esp_at_nuclei(self: psi4.core.Wavefunction) → psi4.core.Vector

returns electrostatic potentials at nuclei

force_doccpi(self: psi4.core.Wavefunction, arg0: psi4.core.Dimension) → None

Specialized expert use only. Sets the number of doubly occupied oribtals per irrep. Note that this results in inconsistent Wavefunction objects for SCF, so caution is advised.

force_soccpi(self: psi4.core.Wavefunction, arg0: psi4.core.Dimension) → None

Specialized expert use only. Sets the number of singly occupied oribtals per irrep. Note that this results in inconsistent Wavefunction objects for SCF, so caution is advised.

frequencies()
static from_file(wfn_data)

Summary

Parameters:wfn_data (str or dict) – If a str reads a Wavefunction from a disk otherwise, assumes the data is passed in.
Returns:A deserialized Wavefunction object
Return type:Wavefunction
frzcpi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of frozen core orbitals per irrep.

frzvpi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of frozen virtual orbitals per irrep.

get_amplitudes(self: psi4.core.CCWavefunction) → Dict[str, psi4.core.Matrix]

Get dict of converged T amplitudes

amps : dict (spacestr, SharedMatrix)

spacestr is a description of the amplitude set using the following conventions.

I,J,K -> alpha occupied i,j,k -> beta occupied A,B,C -> alpha virtual a,b,c -> beta virtual

The following entries are stored in the amps, depending on the reference type

RHF: “tIA”, “tIjAb” UHF: tIA, tia, tIjAb, tIJAB, tijab ROHF: tIA, tia, tIjAb, tIJAB, tijab

Examples

RHF T1 diagnostic = sqrt(sum_ia (T_ia * T_ia)/nelec) >>> mol = “”” … 0 1 … Ne 0.0 0.0 0.0 … symmetry c1”“” >>> e, wfn = psi4.energy(“CCSD/cc-pvdz”, return_wfn=True) >>> t1 = wfn.get_amplitudes()[‘tia’].to_array() >>> t1_diagnostic = np.sqrt(np.dot(t1.ravel(),t1.ravel())/ (2 * wfn.nalpha()) >>> t1_diagnostic == psi4.variable(“CC T1 DIAGNOSTIC”) True

Warning

Symmetry free calculations only (nirreps > 1 will cause error)

Warning

No checks that the amplitudes will fit in core. Do not use for proteins

get_array(key)
get_basisset(self: psi4.core.Wavefunction, arg0: str) → psi4.core.BasisSet

Returns the requested auxiliary basis.

get_dipole_field_strength(self: psi4.core.Wavefunction) → List[float[3]]

Returns a vector of length 3, containing the x, y, and z dipole field strengths.

get_print(self: psi4.core.Wavefunction) → int

Get the print level of the Wavefunction.

get_scratch_filename(filenumber)

Given a wavefunction and a scratch file number, canonicalizes the name so that files can be consistently written and read

get_variable(key)
gradient(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Wavefunction’s gradient.

has_array_variable(self: psi4.core.Wavefunction, arg0: str) → bool

Is the Matrix QC variable (case-insensitive) set?

has_scalar_variable(self: psi4.core.Wavefunction, arg0: str) → bool

Is the double QC variable (case-insensitive) set?

has_variable(key)
hessian(self: psi4.core.Wavefunction) → psi4.core.Matrix

Returns the Wavefunction’s Hessian.

legacy_frequencies()
mo_extents(self: psi4.core.Wavefunction) → List[psi4.core.Vector]

returns the wavefunction’s electronic orbital extents.

molecule(self: psi4.core.Wavefunction) → psi4.core.Molecule

Returns the Wavefunction’s molecule.

nalpha(self: psi4.core.Wavefunction) → int

Number of Alpha electrons.

nalphapi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of alpha orbitals per irrep.

name(self: psi4.core.Wavefunction) → str

The level of theory this wavefunction corresponds to.

nbeta(self: psi4.core.Wavefunction) → int

Number of Beta electrons.

nbetapi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of beta orbitals per irrep.

nfrzc(self: psi4.core.Wavefunction) → int

Number of frozen core electrons.

nirrep(self: psi4.core.Wavefunction) → int

Number of irreps in the system.

nmo(self: psi4.core.Wavefunction) → int

Number of molecule orbitals.

nmopi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of molecular orbitals per irrep.

no_occupations(self: psi4.core.Wavefunction) → List[List[Tuple[float, int, int]]]

returns the natural orbital occupations on the wavefunction.

nso(self: psi4.core.Wavefunction) → int

Number of symmetry orbitals.

nsopi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of symmetry orbitals per irrep.

reference_wavefunction(self: psi4.core.Wavefunction) → psi4.core.Wavefunction

Returns the reference wavefunction.

same_a_b_dens(self: psi4.core.Wavefunction) → bool

Returns true if the alpha and beta densities are the same.

same_a_b_orbs(self: psi4.core.Wavefunction) → bool

Returns true if the alpha and beta orbitals are the same.

scalar_variable(self: psi4.core.Wavefunction, arg0: str) → float

Returns the requested (case-insensitive) double QC variable.

scalar_variables(self: psi4.core.Wavefunction) → Dict[str, float]

Returns the dictionary of all double QC variables.

set_array(key, val)
set_array_variable(self: psi4.core.Wavefunction, arg0: str, arg1: psi4.core.Matrix) → None

Sets the requested (case-insensitive) Matrix QC variable.

set_basisset(self: psi4.core.Wavefunction, arg0: str, arg1: psi4.core.BasisSet) → None

Sets the requested auxiliary basis.

set_energy(self: psi4.core.Wavefunction, arg0: float) → None

Sets the Wavefunction’s energy.

set_external_potential(self: psi4.core.Wavefunction, arg0: psi4.core.ExternalPotential) → None

Sets the requested external potential.

set_frequencies(val)
set_gradient(self: psi4.core.Wavefunction, arg0: psi4.core.Matrix) → None

Sets the Wavefunction’s gradient.

set_hessian(self: psi4.core.Wavefunction, arg0: psi4.core.Matrix) → None

Sets the Wavefunction’s Hessian.

set_legacy_frequencies(self: psi4.core.Wavefunction, arg0: psi4.core.Vector) → None

Sets the frequencies of the Hessian.

set_name(self: psi4.core.Wavefunction, arg0: str) → None

Sets the level of theory this wavefunction corresponds to.

set_print(self: psi4.core.Wavefunction, arg0: int) → None

Sets the print level of the Wavefunction.

set_reference_wavefunction(self: psi4.core.Wavefunction, arg0: psi4.core.Wavefunction) → None

docstring

set_scalar_variable(self: psi4.core.Wavefunction, arg0: str, arg1: float) → None

Sets the requested (case-insensitive) double QC variable.

set_variable(key, val)
shallow_copy(self: psi4.core.Wavefunction, arg0: psi4.core.Wavefunction) → None

Copies the pointers to the internal data.

sobasisset(self: psi4.core.Wavefunction) → psi4.core.SOBasisSet

Returns the symmetry orbitals basis.

soccpi(self: psi4.core.Wavefunction) → psi4.core.Dimension

Returns the number of singly occupied orbitals per irrep.

to_file(filename=None)

Converts a Wavefunction object to a base class

Parameters:
  • wfn (Wavefunction) – A Wavefunction or inherited class
  • filename (None, optional) – An optional filename to write the data to
Returns:

A dictionary and NumPy representation of the Wavefunction.

Return type:

dict

variable(key)
variables()