MADNESS
version 0.9
|
#include <scf.h>
Public Types | |
typedef madness::Function < double, 3 > | functionT |
typedef std::vector< functionT > | vecfuncT |
typedef std::pair< vecfuncT, vecfuncT > | pairvecfuncT |
typedef std::vector< pairvecfuncT > | subspaceT |
typedef madness::Tensor< double > | tensorT |
typedef madness::SeparatedConvolution < double, 3 > | operatorT |
typedef std::shared_ptr < operatorT > | poperatorT |
Public Member Functions | |
SCF (World &world, const char *filename) | |
void | set_protocol (World &world, double thresh) |
void | save_mos (World &world) |
void | load_mos (World &world) |
void | do_plots (World &world) |
void | project (World &world) |
void | make_nuclear_potential (World &world) |
void | project_ao_basis (World &world) |
double | PM_q (const tensorT &S, const tensorT &C, int i, int j, int lo, int nbf) |
void | localize_PM_task_kernel (tensorT &Q, std::vector< tensorT > &Svec, tensorT &C, const bool &doprint, const std::vector< int > &set, const double thetamax, tensorT &U, const double thresh) |
tensorT | localize_PM (World &world, const vecfuncT &mo, const std::vector< int > &set, const double thresh=1e-9, const double thetamax=0.5, const bool randomize=true, const bool doprint=true) |
void | analyze_vectors (World &world, const vecfuncT &mo, const tensorT &occ=tensorT(), const tensorT &energy=tensorT(), const std::vector< int > &set=std::vector< int >()) |
double | DIP (const tensorT &dip, int i, int j, int k, int l) |
tensorT | localize_boys (World &world, const vecfuncT &mo, const std::vector< int > &set, const double thresh=1e-9, const double thetamax=0.5, const bool randomize=true) |
tensorT | kinetic_energy_matrix (World &world, const vecfuncT &v) |
vecfuncT | core_projection (World &world, const vecfuncT &psi, const bool include_Bc=true) |
double | core_projector_derivative (World &world, const vecfuncT &mo, const tensorT &occ, int atom, int axis) |
void | initial_guess (World &world) |
void | initial_load_bal (World &world) |
functionT | make_density (World &world, const tensorT &occ, const vecfuncT &v) |
std::vector< poperatorT > | make_bsh_operators (World &world, const tensorT &evals) |
vecfuncT | apply_hf_exchange (World &world, const tensorT &occ, const vecfuncT &psi, const vecfuncT &f) |
functionT | make_lda_potential (World &world, const functionT &arho) |
functionT | make_dft_potential (World &world, const vecfuncT &vf, int what) |
double | make_dft_energy (World &world, const vecfuncT &vf) |
vecfuncT | apply_potential (World &world, const tensorT &occ, const vecfuncT &amo, const vecfuncT &vf, const vecfuncT &delrho, const functionT &vlocal, double &exc, int ispin) |
tensorT | derivatives (World &world) |
tensorT | dipole (World &world) |
void | vector_stats (const std::vector< double > &v, double &rms, double &maxabsval) |
vecfuncT | compute_residual (World &world, tensorT &occ, tensorT &fock, const vecfuncT &psi, vecfuncT &Vpsi, double &err) |
tensorT | make_fock_matrix (World &world, const vecfuncT &psi, const vecfuncT &Vpsi, const tensorT &occ, double &ekinetic) |
Tensor< double > | twoint (World &world, const vecfuncT &psi) |
Compute the two-electron integrals over the provided set of orbitals. More... | |
tensorT | matrix_exponential (const tensorT &A) |
tensorT | diag_fock_matrix (World &world, tensorT &fock, vecfuncT &psi, vecfuncT &Vpsi, tensorT &evals, const tensorT &occ, double thresh) |
void | loadbal (World &world, functionT &arho, functionT &brho, functionT &arho_old, functionT &brho_old, subspaceT &subspace) |
void | rotate_subspace (World &world, const tensorT &U, subspaceT &subspace, int lo, int nfunc, double trantol) |
void | update_subspace (World &world, vecfuncT &Vpsia, vecfuncT &Vpsib, tensorT &focka, tensorT &fockb, subspaceT &subspace, tensorT &Q, double &bsh_residual, double &update_residual) |
void | solve (World &world) |
Public Attributes | |
MolecularSystem | molsys |
Molecule | molecule |
SCFParameters | param |
XCfunctional | xc |
AtomicBasisSet | aobasis |
functionT | vnuc |
functionT | mask |
vecfuncT | ao |
std::vector< int > | at_to_bf |
std::vector< int > | at_nbf |
poperatorT | coulop |
std::vector< std::shared_ptr < real_derivative_3d > > | gradop |
double | vtol |
double | current_energy |
typedef madness::Function<double,3> SCF::functionT |
typedef madness::SeparatedConvolution<double,3> SCF::operatorT |
typedef std::pair<vecfuncT,vecfuncT> SCF::pairvecfuncT |
typedef std::shared_ptr<operatorT> SCF::poperatorT |
typedef std::vector<pairvecfuncT> SCF::subspaceT |
typedef madness::Tensor<double> SCF::tensorT |
typedef std::vector<functionT> SCF::vecfuncT |
SCF::SCF | ( | World & | world, |
const char * | filename | ||
) |
void SCF::analyze_vectors | ( | World & | world, |
const vecfuncT & | mo, | ||
const tensorT & | occ = tensorT() , |
||
const tensorT & | energy = tensorT() , |
||
const std::vector< int > & | set = std::vector< int >() |
||
) |
vecfuncT SCF::apply_hf_exchange | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | psi, | ||
const vecfuncT & | f | ||
) |
Important this is consistent with Coulomb
References madness::apply(), madness::compress(), madness::SCF::coulop, madness::f, madness::FunctionDefaults< NDIM >::get_thresh(), madness::mul_sparse(), madness::norm_tree(), madness::reconstruct(), and madness::truncate().
vecfuncT SCF::apply_potential | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | amo, | ||
const vecfuncT & | vf, | ||
const vecfuncT & | delrho, | ||
const functionT & | vlocal, | ||
double & | exc, | ||
int | ispin | ||
) |
References madness::SCF::apply_hf_exchange(), madness::SCF::core_projection(), madness::CalculationParameters::core_type, madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, madness::SCF::gradop, madness::XCfunctional::hf_exchange_coefficient(), madness::inner(), madness::XCfunctional::is_dft(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), madness::SCF::make_dft_energy(), madness::SCF::make_dft_potential(), madness::mul_sparse(), madness::SCF::param, madness::START_TIMER(), madness::truncate(), madness::SCF::vtol, and madness::SCF::xc.
vecfuncT SCF::compute_residual | ( | World & | world, |
tensorT & | occ, | ||
tensorT & | fock, | ||
const vecfuncT & | psi, | ||
vecfuncT & | Vpsi, | ||
double & | err | ||
) |
double SCF::core_projector_derivative | ( | World & | world, |
const vecfuncT & | mo, | ||
const tensorT & | occ, | ||
int | atom, | ||
int | axis | ||
) |
tensorT SCF::derivatives | ( | World & | world | ) |
tensorT SCF::diag_fock_matrix | ( | World & | world, |
tensorT & | fock, | ||
vecfuncT & | psi, | ||
vecfuncT & | Vpsi, | ||
tensorT & | evals, | ||
const tensorT & | occ, | ||
double | thresh | ||
) |
References madness::WorldGopInterface::broadcast(), madness::copy(), madness::END_TIMER(), mpfr::fabs(), madness::World::gop, madness::inner(), madness::SCF::matrix_exponential(), madness::matrix_inner(), max, mpfr::min(), madness::normalize(), rot, madness::START_TIMER(), std::swap(), madness::sygv(), madness::transform(), madness::transpose(), madness::truncate(), and madness::SCF::vtol.
|
inline |
tensorT SCF::dipole | ( | World & | world | ) |
void SCF::do_plots | ( | World & | world | ) |
void SCF::initial_guess | ( | World & | world | ) |
void SCF::initial_load_bal | ( | World & | world | ) |
References madness::apply(), madness::SCF::gradop, madness::matrix_inner(), and madness::reconstruct().
void SCF::load_mos | ( | World & | world | ) |
void SCF::loadbal | ( | World & | world, |
functionT & | arho, | ||
functionT & | brho, | ||
functionT & | arho_old, | ||
functionT & | brho_old, | ||
subspaceT & | subspace | ||
) |
tensorT SCF::localize_boys | ( | World & | world, |
const vecfuncT & | mo, | ||
const std::vector< int > & | set, | ||
const double | thresh = 1e-9 , |
||
const double | thetamax = 0.5 , |
||
const bool | randomize = true |
||
) |
References madness::WorldGopInterface::broadcast(), c, madness::copy(), mpfr::cos(), madness::SCF::DIP(), doit(), madness::drot(), madness::drot3(), madness::END_TIMER(), mpfr::fabs(), madness::g, madness::World::gop, madness::matrix_inner(), max, madness::mul_sparse(), madness::print(), madness::RandomValue< double >(), madness::World::rank(), mpfr::sin(), madness::START_TIMER(), mpfr::sum(), madness::transpose(), and madness::SCF::vtol.
tensorT SCF::localize_PM | ( | World & | world, |
const vecfuncT & | mo, | ||
const std::vector< int > & | set, | ||
const double | thresh = 1e-9 , |
||
const double | thetamax = 0.5 , |
||
const bool | randomize = true , |
||
const bool | doprint = true |
||
) |
void SCF::localize_PM_task_kernel | ( | tensorT & | Q, |
std::vector< tensorT > & | Svec, | ||
tensorT & | C, | ||
const bool & | doprint, | ||
const std::vector< int > & | set, | ||
const double | thetamax, | ||
tensorT & | U, | ||
const double | thresh | ||
) |
References a(), mpfr::acos(), madness::SCF::at_nbf, madness::SCF::at_to_bf, C, c, mpfr::cos(), madness::drot(), mpfr::fabs(), max, mpfr::min(), madness::Molecule::natom(), madness::print(), mpfr::sin(), sqrt(), and mpfr::sum().
std::vector<poperatorT> SCF::make_bsh_operators | ( | World & | world, |
const tensorT & | evals | ||
) |
References madness::Function< T, NDIM >::trace(), and madness::SCF::xc.
References madness::SCF::xc.
void SCF::make_nuclear_potential | ( | World & | world | ) |
References B, madness::inner(), and madness::scale().
void SCF::project | ( | World & | world | ) |
void SCF::project_ao_basis | ( | World & | world | ) |
void SCF::rotate_subspace | ( | World & | world, |
const tensorT & | U, | ||
subspaceT & | subspace, | ||
int | lo, | ||
int | nfunc, | ||
double | trantol | ||
) |
References madness::transform().
void SCF::save_mos | ( | World & | world | ) |
void SCF::set_protocol | ( | World & | world, |
double | thresh | ||
) |
void SCF::solve | ( | World & | world | ) |
Compute the two-electron integrals over the provided set of orbitals.
Returned is a replicated tensor of with and . The symmetry is enforced.
Important this is consistent with Coulomb
References madness::apply(), madness::SCF::coulop, madness::WorldGopInterface::fence(), madness::FunctionDefaults< NDIM >::get_thresh(), madness::World::gop, madness::matrix_inner(), madness::mul_sparse(), madness::norm_tree(), madness::reconstruct(), and madness::truncate().
void SCF::update_subspace | ( | World & | world, |
vecfuncT & | Vpsia, | ||
vecfuncT & | Vpsib, | ||
tensorT & | focka, | ||
tensorT & | fockb, | ||
subspaceT & | subspace, | ||
tensorT & | Q, | ||
double & | bsh_residual, | ||
double & | update_residual | ||
) |
void SCF::vector_stats | ( | const std::vector< double > & | v, |
double & | rms, | ||
double & | maxabsval | ||
) |
References sqrt().
vecfuncT SCF::ao |
AtomicBasisSet SCF::aobasis |
std::vector<int> SCF::at_nbf |
std::vector<int> SCF::at_to_bf |
poperatorT SCF::coulop |
double SCF::current_energy |
std::vector< std::shared_ptr<real_derivative_3d> > SCF::gradop |
functionT SCF::mask |
Molecule SCF::molecule |
MolecularSystem SCF::molsys |
SCFParameters SCF::param |
Referenced by madness::Nemo::Nemo().
functionT SCF::vnuc |
double SCF::vtol |
XCfunctional SCF::xc |