Class which handles computing the CAP matrix.
More...
#include <CAP.h>
|
| CAP (System &my_sys, std::map< std::string, std::string > params) |
|
| CAP (System my_sys, py::dict dict, size_t num_states, const std::function< std::vector< double >(std::vector< double > &, std::vector< double > &, std::vector< double > &, std::vector< double > &)> &cap_func={}) |
|
void | compute_ao_cap (py::dict dict, const std::function< std::vector< double >(std::vector< double > &, std::vector< double > &, std::vector< double > &, std::vector< double > &)> &cap_func={}) |
|
void | integrate_cap () |
|
void | compute_projected_cap () |
|
void | verify_method (std::map< std::string, std::string > params) |
|
void | run () |
|
Eigen::MatrixXd | get_ao_cap (std::string ordering="molden", std::string basis_file="") |
|
Eigen::MatrixXd | get_projected_cap () |
|
Eigen::MatrixXd | get_density (size_t row_idx, size_t col_idx, bool beta=false) |
|
Eigen::MatrixXd | get_H () |
|
void | add_tdms (Eigen::MatrixXd &alpha_density, Eigen::MatrixXd &beta_density, size_t row_idx, size_t col_idx, std::string ordering, std::string basis_file="") |
|
void | add_tdm (Eigen::MatrixXd tdm, size_t row_idx, size_t col_idx, std::string ordering, std::string basis_file="") |
|
void | read_electronic_structure_data (py::dict dict) |
|
void | verify_data () |
|
void | renormalize_cap (Eigen::MatrixXd smat, std::string ordering, std::string basis_file="") |
|
void | renormalize () |
|
void | compute_cap_on_grid (py::array_t< double > &x, py::array_t< double > &y, py::array_t< double > &z, py::array_t< double > &grid_w) |
|
Class which handles computing the CAP matrix.
◆ CAP() [1/2]
CAP::CAP |
( |
System & |
my_sys, |
|
|
std::map< std::string, std::string > |
params |
|
) |
| |
Constructs CAP object from System object.
- Parameters
-
my_sys | System object |
params | Map of parameters |
◆ CAP() [2/2]
CAP::CAP |
( |
System |
my_sys, |
|
|
py::dict |
dict, |
|
|
size_t |
num_states, |
|
|
const std::function< std::vector< double >(std::vector< double > &, std::vector< double > &, std::vector< double > &, std::vector< double > &)> & |
cap_func = {} |
|
) |
| |
Constructor for Python.
- Parameters
-
my_sys | System object |
dict | Python dictionary containing parameters |
num_states | number of states |
◆ add_tdm()
void CAP::add_tdm |
( |
Eigen::MatrixXd |
tdm, |
|
|
size_t |
row_idx, |
|
|
size_t |
col_idx, |
|
|
std::string |
ordering, |
|
|
std::string |
basis_file = "" |
|
) |
| |
Adds spin traced transition density matrix from state row_idx --> col_idx.
- Parameters
-
tdm | Spin traced TDM in AO basis of dimension (NBasis,Nbasis) |
row_idx | initial state index |
col_idx | final state index |
ordering | order of GTOs |
basis_file | File containing basis set specification. Required for OpenMolcas. |
◆ add_tdms()
void CAP::add_tdms |
( |
Eigen::MatrixXd & |
alpha_density, |
|
|
Eigen::MatrixXd & |
beta_density, |
|
|
size_t |
row_idx, |
|
|
size_t |
col_idx, |
|
|
std::string |
ordering, |
|
|
std::string |
basis_file = "" |
|
) |
| |
Adds transition density matrices (alpha and beta) from state row_idx --> col_idx.
- Parameters
-
alpha_density | TDM in AO basis of dimension (NBasis,Nbasis) |
beta_density | TDM in AO basis of dimension (NBasis,Nbasis) |
row_idx | initial state index |
col_idx | final state index |
ordering | order of GTOs |
basis_file | File containing basis set specification. Required for OpenMolcas. |
◆ compute_ao_cap()
void CAP::compute_ao_cap |
( |
py::dict |
dict, |
|
|
const std::function< std::vector< double >(std::vector< double > &, std::vector< double > &, std::vector< double > &, std::vector< double > &)> & |
cap_func = {} |
|
) |
| |
Computes CAP in AO basis using parameters in given python dictionary
◆ compute_cap_on_grid()
void CAP::compute_cap_on_grid |
( |
py::array_t< double > & |
x, |
|
|
py::array_t< double > & |
y, |
|
|
py::array_t< double > & |
z, |
|
|
py::array_t< double > & |
grid_w |
|
) |
| |
◆ compute_projected_cap()
void CAP::compute_projected_cap |
( |
| ) |
|
Computes CAP in state basis
◆ get_ao_cap()
Eigen::MatrixXd CAP::get_ao_cap |
( |
std::string |
ordering = "molden" , |
|
|
std::string |
basis_file = "" |
|
) |
| |
Returns CAP matrix in AO basis.
- Parameters
-
ordering | order of GTOs |
basis_file | File containing basis set specification. Required for OpenMolcas. |
◆ get_density()
Eigen::MatrixXd CAP::get_density |
( |
size_t |
row_idx, |
|
|
size_t |
col_idx, |
|
|
bool |
beta = false |
|
) |
| |
Gets TDM for state row_idx --> col_idx.
- Parameters
-
row_idx | initial state index |
col_idx | final state index |
beta | Beta density, false means alpha density |
◆ get_H()
Eigen::MatrixXd CAP::get_H |
( |
| ) |
|
Returns zeroth order Hamiltonian.
◆ get_projected_cap()
Eigen::MatrixXd CAP::get_projected_cap |
( |
| ) |
|
Returns CAP matrix in state basis.
◆ integrate_cap()
void CAP::integrate_cap |
( |
| ) |
|
◆ read_electronic_structure_data()
void CAP::read_electronic_structure_data |
( |
py::dict |
dict | ) |
|
Reads in electronic structure data from file, from python. Valid keywords: method,qc_output,h0_file,rassi_h5, fchk_file,molcas_output.
- Parameters
-
dict | Python dictionary containing keywords for reading in data from disk. |
◆ renormalize()
void CAP::renormalize |
( |
| ) |
|
Renormalizes the AO CAP matrix using the previously parsed overlap matrix.
◆ renormalize_cap()
void CAP::renormalize_cap |
( |
Eigen::MatrixXd |
smat, |
|
|
std::string |
ordering, |
|
|
std::string |
basis_file = "" |
|
) |
| |
Renormalizes the AO CAP matrix using the supplied overlap matrix.
◆ run()
Executes Perturbative CAP method.
◆ verify_data()
void CAP::verify_data |
( |
| ) |
|
Verifies that required electronic structure data is present to perform calculation.
◆ verify_method()
void CAP::verify_method |
( |
std::map< std::string, std::string > |
params | ) |
|
Checks that electronic structure method and package is supported, and that necessary keywords are present.
◆ alpha_dms
std::vector<std::vector<Eigen::MatrixXd> > CAP::alpha_dms |
Transition density matrices in AO basis, alpha densities
◆ AO_CAP_MAT
Eigen::MatrixXd CAP::AO_CAP_MAT |
◆ beta_dms
std::vector<std::vector<Eigen::MatrixXd> > CAP::beta_dms |
Transition density matrices in AO basis, beta densities
◆ cap_integrator
AOCAP CAP::cap_integrator |
◆ CAP_MAT
Eigen::MatrixXd CAP::CAP_MAT |
CAP matrix in correlated many electron basis
◆ nstates
◆ OVERLAP_MAT
Eigen::MatrixXd CAP::OVERLAP_MAT |
Overlap matrix parsed from file, can be empty when constructed from python
◆ parameters
std::map<std::string, std::string> CAP::parameters |
Map containing the parameters defined in the input
◆ python
Set to true when constructed from the python interpreter, important for printing
◆ rotation_matrix
Eigen::MatrixXd CAP::rotation_matrix |
Rotation matrix, required for XMS-CASPT2 variants, otherwise set to I
◆ system
◆ ZERO_ORDER_H
Eigen::MatrixXd CAP::ZERO_ORDER_H |
Zeroth order Hamiltonian. Dimension is (nstates,nstates)
The documentation for this class was generated from the following file:
- /home/runner/work/opencap/opencap/opencap/include/CAP.h