pyopencap.CAP

The CAP class is used to compute the CAP matrix first in AO basis, and then in wave function basis using the one-particle densities which are passed in. It is also capable of parsing OpenMolcas output files to obtain the zeroth order Hamiltonian and return it to the user.

class pyopencap.CAP[source]
__init__(*args, **kwargs)

Overloaded function.

  1. __init__(self: pyopencap.pyopencap_cpp.CAP, system: pyopencap.pyopencap_cpp.System, cap_dict: dict, nstates: int) -> None

Constructs CAP object from system, cap dictionary, and number of states.

  1. __init__(self: pyopencap.pyopencap_cpp.CAP, system: pyopencap.pyopencap_cpp.System, cap_dict: dict, nstates: int, cap_func: Callable[[List[float], List[float], List[float], List[float]], List[float]]) -> None

Constructs CAP object from system, cap dictionary, number of states, and cap function.

add_tdm(self: pyopencap.pyopencap_cpp.CAP, tdm: numpy.ndarray[numpy.float64[m, n]], initial_idx: int, final_idx: int, ordering: str, basis_file: str = '') None

Adds spin-traced tdm to CAP object at specified indices. The optional argument basis_file is required when using the OpenMolcas interface, and it must point to the path to the rassi.5 file. Supported orderings: pyscf, openmolcas, qchem, psi4, molden.

add_tdms(self: pyopencap.pyopencap_cpp.CAP, alpha_density: numpy.ndarray[numpy.float64[m, n]], beta_density: numpy.ndarray[numpy.float64[m, n]], initial_idx: int, final_idx: int, ordering: str, basis_file: str = '') None

Adds alpha/beta tdms to CAP object at specified indices. The optional argument basis_file is required when using the OpenMolcas interface, and it must point to the path to the rassi.5 file. Supported orderings: pyscf, openmolcas, qchem, psi4, molden.

compute_ao_cap(self: pyopencap.pyopencap_cpp.CAP, cap_dict: dict, cap_func: Callable[[List[float], List[float], List[float], List[float]], List[float]] = None) None

Computes CAP matrix in AO basis.

compute_projected_cap(self: pyopencap.pyopencap_cpp.CAP) None

Computes CAP matrix in state basis using transition density matrices.

get_H(self: pyopencap.pyopencap_cpp.CAP) numpy.ndarray[numpy.float64[m, n]]

Returns zeroth order Hamiltonian read from file.

get_ao_cap(self: pyopencap.pyopencap_cpp.CAP, ordering: str = 'molden', basis_file: str = '') numpy.ndarray[numpy.float64[m, n]]

Returns CAP matrix in AO basis. Supported orderings: pyscf, openmolcas, qchem, psi4, molden.

get_projected_cap(self: pyopencap.pyopencap_cpp.CAP) numpy.ndarray[numpy.float64[m, n]]

Returns CAP matrix in state basis.

read_data(self: pyopencap.pyopencap_cpp.CAP, es_dict: dict) None

Reads electronic structure data specified in dictionary.

renormalize(self: pyopencap.pyopencap_cpp.CAP) None

Re-normalizes AO CAP using electronic structure data.

renormalize_cap(self: pyopencap.pyopencap_cpp.CAP, smat: numpy.ndarray[numpy.float64[m, n]], ordering: str, basis_file: str = '') None

Re-normalizes AO CAP matrix using input overlap matrix.

compute_cap_on_grid(self: pyopencap.pyopencap_cpp.CAP, x: numpy.ndarray[numpy.float64], y: numpy.ndarray[numpy.float64], z: numpy.ndarray[numpy.float64], w: numpy.ndarray[numpy.float64]) None

Computes CAP matrix on supplied grid. Sum will cumulated for each successive grid until compute_projected_cap is called.