MADNESS
version 0.9
|
#include <SCF.h>
Public Member Functions | |
SCF (World &world, const char *filename) | |
template<std::size_t NDIM> | |
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) |
std::vector< int > | group_orbital_sets (World &world, const tensorT &eps, const tensorT &occ, const int nmo) const |
group orbitals into sets of similar orbital energies for localization More... | |
distmatT | 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=false) const |
compute the unitary localization matrix according to Pipek-Mezey More... | |
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) |
distmatT | kinetic_energy_matrix (World &world, const vecfuncT &v) const |
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) const |
functionT | make_density (World &world, const tensorT &occ, const cvecfuncT &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) const |
apply the HF exchange on a set of orbitals More... | |
functionT | make_lda_potential (World &world, const functionT &arho) |
functionT | make_dft_potential (World &world, const vecfuncT &vf, int ispin, int what) |
double | make_dft_energy (World &world, const vecfuncT &vf, int ispin) |
vecfuncT | apply_potential (World &world, const tensorT &occ, const vecfuncT &amo, const vecfuncT &vf, const vecfuncT &delrho, const functionT &vlocal, double &exc, double &enl, int ispin) |
tensorT | derivatives (World &world) |
tensorT | dipole (World &world) |
void | vector_stats (const std::vector< double > &v, double &rms, double &maxabsval) const |
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) const |
functionT | make_coulomb_potential (const functionT &rho) const |
make the Coulomb potential given the total density More... | |
Tensor< double > | twoint (World &world, const vecfuncT &psi) const |
Compute the two-electron integrals over the provided set of orbitals. More... | |
tensorT | matrix_exponential (const tensorT &A) const |
tensorT | get_fock_transformation (World &world, const tensorT &overlap, tensorT &fock, tensorT &evals, const tensorT &occ, const double thresh_degenerate) const |
compute the unitary transformation that diagonalizes the fock matrix More... | |
tensorT | diag_fock_matrix (World &world, tensorT &fock, vecfuncT &psi, vecfuncT &Vpsi, tensorT &evals, const tensorT &occ, const double thresh) const |
diagonalize the fock matrix, taking care of degenerate states More... | |
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) const |
void | rotate_subspace (World &world, const distmatT &U, subspaceT &subspace, int lo, int nfunc, double trantol) const |
void | update_subspace (World &world, vecfuncT &Vpsia, vecfuncT &Vpsib, tensorT &focka, tensorT &fockb, subspaceT &subspace, tensorT &Q, double &bsh_residual, double &update_residual) |
double | do_step_restriction (World &world, const vecfuncT &mo, vecfuncT &mo_new, std::string spin) const |
perform step restriction following the KAIN solver More... | |
void | orthonormalize (World &world, vecfuncT &amo_new) const |
orthonormalize the vectors More... | |
void | propagate (World &world, double omega, int step0) |
complex_functionT | APPLY (const complex_operatorT *q1d, const complex_functionT &psi) |
void | iterate_trotter (World &world, Convolution1D< double_complex > *G, cvecfuncT &camo, cvecfuncT &cbmo, double t, double time_step, double thresh) |
void | solve (World &world) |
Public Attributes | |
std::shared_ptr< PotentialManager > | potentialmanager |
std::shared_ptr < GTHPseudopotential< double > > | gthpseudopotential |
Molecule | molecule |
CalculationParameters | param |
XCfunctional | xc |
AtomicBasisSet | aobasis |
functionT | mask |
vecfuncT | amo |
alpha and beta molecular orbitals More... | |
vecfuncT | bmo |
std::vector< int > | aset |
std::vector< int > | bset |
vecfuncT | ao |
MRA projection of the minimal basis set. More... | |
std::vector< int > | at_to_bf |
std::vector< int > | at_nbf |
tensorT | aocc |
occupation numbers for alpha and beta orbitals More... | |
tensorT | bocc |
tensorT | aeps |
orbital energies for alpha and beta orbitals More... | |
tensorT | beps |
poperatorT | coulop |
std::vector< std::shared_ptr < real_derivative_3d > > | gradop |
double | vtol |
double | current_energy |
Static Public Attributes | |
static const int | vnucextra = 12 |
SCF::SCF | ( | World & | world, |
const char * | filename | ||
) |
References madness::CalculationParameters::aobasis, aobasis, madness::WorldGopInterface::broadcast_serializable(), madness::CalculationParameters::core_type, madness::CalculationParameters::econv, madness::Molecule::get_atom_number(), madness::get_charge_from_file(), madness::World::gop, gthpseudopotential, madness::Molecule::guess_file(), madness::XCfunctional::initialize(), madness::CalculationParameters::L, molecule, madness::Molecule::n_core_orb_all(), madness::Molecule::natom(), madness::CalculationParameters::no_orient, madness::Molecule::orient(), param, potentialmanager, madness::CalculationParameters::psp_calc, madness::World::rank(), madness::Molecule::read_core_file(), madness::Molecule::read_file(), madness::CalculationParameters::read_file(), madness::AtomicBasisSet::read_file(), madness::Molecule::set_atom_charge(), madness::CalculationParameters::set_molecular_info(), madness::CalculationParameters::spin_restricted, TAU_START, TAU_STOP, xc, and madness::CalculationParameters::xc_data.
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>() |
||
) |
References ao, aobasis, C, madness::END_TIMER(), energy(), madness::gesvp(), madness::inner(), madness::matrix_inner(), molecule, madness::mul_sparse(), madness::AtomicBasisSet::print_anal(), madness::World::rank(), sqrt(), madness::START_TIMER(), TAU_START, TAU_STOP, madness::transpose(), and vtol.
Referenced by solve().
|
inline |
vecfuncT SCF::apply_hf_exchange | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | psi, | ||
const vecfuncT & | f | ||
) | const |
apply the HF exchange on a set of orbitals
[in] | world | the world |
[in] | occ | occupation numbers |
[in] | psi | the orbitals in the exchange operator |
[in] | f | the orbitals |i> that the operator is applied on |
Important this is consistent with Coulomb
References madness::apply(), madness::compress(), coulop, madness::f, madness::mul_sparse(), madness::norm_tree(), madness::reconstruct(), and madness::truncate().
Referenced by SCF::apply_potential(), and apply_potential().
vecfuncT SCF::apply_potential | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | amo, | ||
const vecfuncT & | vf, | ||
const vecfuncT & | delrho, | ||
const functionT & | vlocal, | ||
double & | exc, | ||
double & | enl, | ||
int | ispin | ||
) |
References apply_hf_exchange(), core_projection(), madness::CalculationParameters::core_type, madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, gthpseudopotential, madness::XCfunctional::hf_exchange_coefficient(), madness::inner(), madness::XCfunctional::is_dft(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), make_dft_energy(), make_dft_potential(), madness::mul_sparse(), madness::CalculationParameters::nbeta, param, potentialmanager, madness::print_meminfo(), madness::CalculationParameters::psp_calc, madness::World::rank(), madness::Function< T, NDIM >::scale(), madness::START_TIMER(), TAU_START, TAU_STOP, madness::truncate(), madness::Function< T, NDIM >::truncate(), vtol, and xc.
Referenced by solve().
vecfuncT SCF::compute_residual | ( | World & | world, |
tensorT & | occ, | ||
tensorT & | fock, | ||
const vecfuncT & | psi, | ||
vecfuncT & | Vpsi, | ||
double & | err | ||
) |
References madness::apply(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, make_bsh_operators(), mpfr::min(), madness::norm2s(), madness::print(), madness::World::rank(), madness::scale(), madness::set_thresh(), madness::START_TIMER(), madness::sub(), TAU_START, TAU_STOP, madness::transform(), madness::truncate(), vector_stats(), and vtol.
Referenced by update_subspace().
References madness::Atom::atomic_number, c, madness::Function< T, NDIM >::fence(), madness::Molecule::get_atom(), madness::Molecule::get_core_bc(), madness::Molecule::get_core_l(), madness::inner(), m, molecule, madness::Molecule::n_core_orb(), madness::Molecule::natom(), madness::print(), psi(), and madness::Function< T, NDIM >::scale().
Referenced by SCF::apply_potential(), and apply_potential().
References amo, aocc, bmo, bocc, madness::CalculationParameters::core_type, madness::END_TIMER(), func(), madness::Function< T, NDIM >::gaxpy(), madness::Molecule::get_atom(), madness::inner(), madness::Molecule::is_potential_defined_atom(), make_density(), molecule, madness::Molecule::natom(), madness::CalculationParameters::nbeta, madness::Molecule::nuclear_repulsion_derivative(), param, potentialmanager, madness::print(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), TAU_START, TAU_STOP, madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by madness::HartreeFock< T, NDIM >::gradient(), madness::MolecularEnergy::gradient(), and main().
tensorT SCF::diag_fock_matrix | ( | World & | world, |
tensorT & | fock, | ||
vecfuncT & | psi, | ||
vecfuncT & | Vpsi, | ||
tensorT & | evals, | ||
const tensorT & | occ, | ||
const double | thresh | ||
) | const |
diagonalize the fock matrix, taking care of degenerate states
Vpsi is passed in to make sure orbitals and Vpsi are in phase
[in] | world | the world |
[in,out] | fock | the fock matrix (diagonal upon exit) |
[in,out] | psi | the orbitals |
[in,out] | Vpsi | the orbital times the potential |
[out] | evals | the orbital energies |
[in] | occ | occupation numbers |
[in] | thresh | threshold for rotation and truncation |
References madness::END_TIMER(), get_fock_transformation(), madness::matrix_inner(), mpfr::min(), madness::normalize(), TAU_STOP, madness::transform(), madness::truncate(), and vtol.
Referenced by solve().
|
inline |
Referenced by SCF::localize_boys().
References amo, aocc, bmo, bocc, madness::END_TIMER(), make_density(), molecule, mu, madness::CalculationParameters::nbeta, madness::Molecule::nuclear_dipole(), param, madness::print(), madness::World::rank(), madness::Function< T, NDIM >::scale(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), TAU_START, and TAU_STOP.
void SCF::do_plots | ( | World & | world | ) |
References amo, aocc, madness::apply(), bmo, bocc, madness::copy(), coulop, madness::END_TIMER(), make_density(), madness::CalculationParameters::nalpha, madness::CalculationParameters::nbeta, madness::CalculationParameters::npt_plot, param, madness::CalculationParameters::plot_cell, madness::CalculationParameters::plotcoul, madness::CalculationParameters::plotdens, madness::plotdx(), madness::CalculationParameters::plothi, madness::CalculationParameters::plotlo, potentialmanager, madness::Function< T, NDIM >::reconstruct(), madness::Function< T, NDIM >::scale(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), TAU_START, TAU_STOP, and madness::Function< T, NDIM >::truncate().
Referenced by main().
double SCF::do_step_restriction | ( | World & | world, |
const vecfuncT & | mo, | ||
vecfuncT & | mo_new, | ||
std::string | spin | ||
) | const |
perform step restriction following the KAIN solver
undo the rotation from the KAIN solver if the rotation exceeds the maxrotn parameter
[in] | world | the world |
[in] | mo | vector of orbitals from previous iteration |
[in,out] | new_mo | vector of orbitals from the KAIN solver |
[in] | spin | "alpha" or "beta" for user information |
References madness::WorldGopInterface::fence(), madness::World::gop, madness::CalculationParameters::maxrotn, madness::norm2s(), param, madness::print(), madness::World::rank(), madness::sub(), and vector_stats().
Referenced by update_subspace().
tensorT SCF::get_fock_transformation | ( | World & | world, |
const tensorT & | overlap, | ||
tensorT & | fock, | ||
tensorT & | evals, | ||
const tensorT & | occ, | ||
const double | thresh_degenerate | ||
) | const |
compute the unitary transformation that diagonalizes the fock matrix
[in] | world | the world |
[in] | overlap | the overlap matrix of the orbitals |
[in,out] | fock | the fock matrix; diagonal upon exit |
[out] | evals | the orbital energies |
[in] | occ | the occupation numbers |
[in] | thresh_degenerate | threshold for orbitals being degenerate |
References madness::WorldGopInterface::broadcast(), madness::copy(), madness::END_TIMER(), mpfr::fabs(), madness::World::gop, madness::inner(), matrix_exponential(), max, rot, madness::START_TIMER(), std::swap(), madness::sygvp(), TAU_START, TAU_STOP, and madness::transpose().
Referenced by diag_fock_matrix().
std::vector< int > SCF::group_orbital_sets | ( | World & | world, |
const tensorT & | eps, | ||
const tensorT & | occ, | ||
const int | nmo | ||
) | const |
group orbitals into sets of similar orbital energies for localization
[in] | eps | orbital energies |
[in] | occ | occupation numbers |
[in] | nmo | number of MOs for the given spin |
References madness::print(), and madness::World::rank().
Referenced by initial_guess().
void SCF::initial_guess | ( | World & | world | ) |
References madness::LoadBalanceDeux< NDIM >::add_tree(), aeps, amo, ao, aobasis, aocc, madness::apply(), aset, beps, bmo, bocc, bset, c, madness::Function< T, NDIM >::clear(), madness::compress(), madness::DistributedMatrix< T >::copy_to_replicated(), madness::CalculationParameters::core_type, coulop, madness::END_TIMER(), group_orbital_sets(), gthpseudopotential, kinetic_energy_matrix(), madness::LoadBalanceDeux< NDIM >::load_balance(), load_mos(), make_lda_potential(), madness::matrix_inner(), molecule, madness::mul_sparse(), madness::Molecule::n_core_orb_all(), madness::CalculationParameters::nalpha, madness::CalculationParameters::nbeta, madness::CalculationParameters::nmo_alpha, madness::CalculationParameters::nmo_beta, madness::normalize(), param, potential(), potentialmanager, madness::print(), madness::print_meminfo(), madness::CalculationParameters::psp_calc, madness::World::rank(), madness::reconstruct(), madness::Function< T, NDIM >::reconstruct(), madness::CalculationParameters::restart, save(), madness::Function< T, NDIM >::scale(), madness::World::size(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), madness::sygvp(), TAU_START, TAU_STOP, madness::Molecule::total_nuclear_charge(), madness::Function< T, NDIM >::trace(), madness::transform(), madness::transpose(), madness::truncate(), madness::Function< T, NDIM >::truncate(), vnucextra, and vtol.
Referenced by madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
void SCF::initial_load_bal | ( | World & | world | ) |
References madness::LoadBalanceDeux< NDIM >::add_tree(), gthpseudopotential, madness::LoadBalanceDeux< NDIM >::load_balance(), param, potentialmanager, madness::CalculationParameters::psp_calc, and vnucextra.
Referenced by main(), and propagate().
void SCF::iterate_trotter | ( | World & | world, |
Convolution1D< double_complex > * | G, | ||
cvecfuncT & | camo, | ||
cvecfuncT & | cbmo, | ||
double | t, | ||
double | time_step, | ||
double | thresh | ||
) |
References amo, aocc, madness::apply(), APPLY(), bmo, bocc, coulop, make_density(), madness::mul_sparse(), madness::CalculationParameters::nalpha, madness::CalculationParameters::nbeta, param, potentialmanager, madness::CalculationParameters::spin_restricted, and vtol.
Referenced by propagate().
References madness::apply(), madness::DistributedMatrixDistribution::distribution(), gradop, madness::matrix_inner(), and madness::reconstruct().
Referenced by initial_guess(), SCF::make_fock_matrix(), and make_fock_matrix().
void SCF::load_mos | ( | World & | world | ) |
References aeps, amo, aocc, aset, beps, bmo, bocc, bset, madness::copy(), current_energy, madness::WorldGopInterface::fence(), madness::World::gop, k, molecule, madness::Molecule::n_core_orb_all(), madness::CalculationParameters::nmo_alpha, madness::CalculationParameters::nmo_beta, param, madness::project(), madness::reconstruct(), madness::CalculationParameters::spin_restricted, TAU_START, TAU_STOP, and thresh.
Referenced by initial_guess(), propagate(), madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
void SCF::loadbal | ( | World & | world, |
functionT & | arho, | ||
functionT & | brho, | ||
functionT & | arho_old, | ||
functionT & | brho_old, | ||
subspaceT & | subspace | ||
) |
References madness::LoadBalanceDeux< NDIM >::add_tree(), amo, bmo, madness::WorldGopInterface::fence(), madness::World::gop, gthpseudopotential, madness::LoadBalanceDeux< NDIM >::load_balance(), madness::CalculationParameters::nbeta, param, potentialmanager, madness::CalculationParameters::psp_calc, madness::World::size(), madness::CalculationParameters::spin_restricted, and vnucextra.
Referenced by solve().
distmatT 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 = false |
||
) | const |
compute the unitary localization matrix according to Pipek-Mezey
[in] | world | the world |
[in] | mo | the MOs |
[in] | set | only orbitals within the same set will be mixed |
[in] | thresh | the localization threshold |
[in] | thetamax | ?? |
[in] | randomize | ?? |
References ao, at_nbf, at_to_bf, madness::distributed_localize_PM(), madness::END_TIMER(), madness::START_TIMER(), and TAU_START.
Referenced by solve().
std::vector< poperatorT > SCF::make_bsh_operators | ( | World & | world, |
const tensorT & | evals | ||
) |
References madness::CalculationParameters::lo, param, madness::print(), madness::World::rank(), and sqrt().
Referenced by compute_residual().
make the Coulomb potential given the total density
References madness::apply().
Referenced by madness::HartreeFock< T, NDIM >::get_coulomb_potential().
References madness::compress(), madness::Function< T, NDIM >::compress(), madness::WorldGopInterface::fence(), madness::Function< T, NDIM >::gaxpy(), madness::World::gop, and madness::square().
Referenced by derivatives(), dipole(), do_plots(), madness::HartreeFock< T, NDIM >::get_coulomb_potential(), iterate_trotter(), main(), propagate(), and solve().
References madness::Function< T, NDIM >::trace().
Referenced by SCF::apply_potential(), and apply_potential().
|
inline |
Referenced by SCF::apply_potential(), and apply_potential().
tensorT SCF::make_fock_matrix | ( | World & | world, |
const vecfuncT & | psi, | ||
const vecfuncT & | Vpsi, | ||
const tensorT & | occ, | ||
double & | ekinetic | ||
) | const |
References madness::DistributedMatrix< T >::copy_to_replicated(), madness::END_TIMER(), k, kinetic_energy_matrix(), madness::matrix_inner(), madness::START_TIMER(), TAU_START, TAU_STOP, and madness::transpose().
Referenced by solve().
References madness::copy(), madness::Function< T, NDIM >::reconstruct(), and madness::Function< T, NDIM >::unaryop().
Referenced by initial_guess().
void SCF::make_nuclear_potential | ( | World & | world | ) |
References B, madness::inner(), k, and madness::scale().
Referenced by SCF::diag_fock_matrix(), and get_fock_transformation().
orthonormalize the vectors
[in] | world | the world |
[in,out] | amo_new | the vectors to be orthonormalized |
References amo, madness::END_TIMER(), madness::matrix_inner(), max, mpfr::min(), madness::normalize(), madness::print(), madness::Q2(), madness::World::rank(), madness::START_TIMER(), TAU_START, TAU_STOP, madness::transform(), madness::truncate(), and vtol.
Referenced by update_subspace().
void SCF::project | ( | World & | world | ) |
References amo, bmo, madness::WorldGopInterface::fence(), madness::World::gop, madness::CalculationParameters::nbeta, madness::normalize(), param, madness::project(), madness::reconstruct(), madness::CalculationParameters::spin_restricted, and madness::truncate().
Referenced by madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
void SCF::project_ao_basis | ( | World & | world | ) |
References ao, aobasis, at_nbf, at_to_bf, madness::AtomicBasisSet::atoms_to_bfn(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::AtomicBasisSet::get_atomic_basis_function(), madness::World::gop, molecule, madness::AtomicBasisSet::nbf(), madness::normalize(), madness::print_meminfo(), madness::World::rank(), madness::START_TIMER(), TAU_START, TAU_STOP, and madness::truncate().
Referenced by madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
void SCF::propagate | ( | World & | world, |
double | omega, | ||
int | step0 | ||
) |
References amo, aocc, bmo, c, madness::copy(), mpfr::exp(), G, initial_load_bal(), madness::inner(), iterate_trotter(), madness::CalculationParameters::L, load_mos(), make_density(), make_nuclear_potential(), madness::CalculationParameters::nalpha, madness::CalculationParameters::nbeta, param, madness::constants::pi, madness::qm_1d_free_particle_propagator(), madness::CalculationParameters::spin_restricted, and thresh.
Referenced by main().
void SCF::rotate_subspace | ( | World & | world, |
const tensorT & | U, | ||
subspaceT & | subspace, | ||
int | lo, | ||
int | nfunc, | ||
double | trantol | ||
) | const |
References madness::transform().
Referenced by solve().
void SCF::rotate_subspace | ( | World & | world, |
const distmatT & | U, | ||
subspaceT & | subspace, | ||
int | lo, | ||
int | nfunc, | ||
double | trantol | ||
) | const |
References madness::transform().
void SCF::save_mos | ( | World & | world | ) |
References aeps, amo, aocc, aset, beps, bmo, bocc, bset, current_energy, madness::CalculationParameters::nio, param, madness::CalculationParameters::spin_restricted, TAU_START, and TAU_STOP.
Referenced by madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
|
inline |
References madness::CalculationParameters::dconv, madness::f, madness::Molecule::get_eprec(), madness::FunctionDefaults< NDIM >::get_thresh(), k, madness::CalculationParameters::k, madness::CalculationParameters::L, madness::CalculationParameters::lo, max, mpfr::min(), NDIM, madness::print(), madness::World::rank(), madness::FunctionDefaults< NDIM >::set_apply_randomize(), madness::FunctionDefaults< NDIM >::set_autorefine(), madness::FunctionDefaults< NDIM >::set_cubic_cell(), madness::Molecule::set_eprec(), madness::FunctionDefaults< NDIM >::set_initial_level(), madness::FunctionDefaults< NDIM >::set_k(), madness::FunctionDefaults< NDIM >::set_project_randomize(), madness::FunctionDefaults< NDIM >::set_refine(), madness::FunctionDefaults< NDIM >::set_thresh(), and madness::FunctionDefaults< NDIM >::set_truncate_mode().
Referenced by main(), madness::MP2::MP2(), madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
void SCF::solve | ( | World & | world | ) |
References aeps, amo, analyze_vectors(), aocc, madness::apply(), apply_potential(), aset, beps, bmo, bocc, bset, madness::Function< T, NDIM >::clear(), madness::CalculationParameters::conv_only_dens, coulop, current_energy, madness::DistributedMatrix< T >::data(), madness::CalculationParameters::dconv, diag_fock_matrix(), dipole(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::World::gop, gradop, gthpseudopotential, madness::inner(), madness::XCfunctional::is_dft(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), loadbal(), madness::CalculationParameters::localize, localize_PM(), make_density(), make_fock_matrix(), madness::matrix_inner(), max, madness::CalculationParameters::maxiter, madness::CalculationParameters::maxsub, mpfr::min(), molecule, madness::Molecule::natom(), madness::CalculationParameters::nbeta, madness::norm2(), madness::normalize(), madness::Molecule::nuclear_repulsion_energy(), madness::Molecule::nuclear_repulsion_energy_pseudo(), param, potentialmanager, madness::print(), madness::print_meminfo(), madness::CalculationParameters::psp_calc, madness::World::rank(), madness::reconstruct(), madness::Function< T, NDIM >::reconstruct(), madness::Function< T, NDIM >::refine_to_common_level(), rotate_subspace(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), madness::sygvp(), TAU_START, TAU_STOP, madness::transform(), madness::truncate(), madness::Function< T, NDIM >::truncate(), update_subspace(), vtol, madness::wall_time(), and xc.
Referenced by madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
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(), coulop, madness::WorldGopInterface::fence(), 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 | ||
) |
References amo, aocc, bmo, bocc, madness::WorldGopInterface::broadcast(), madness::WorldGopInterface::broadcast_serializable(), c, madness::compress(), compute_residual(), do_step_restriction(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, madness::KAIN(), m, max, madness::CalculationParameters::maxsub, madness::CalculationParameters::nbeta, madness::CalculationParameters::nmo_alpha, madness::CalculationParameters::nmo_beta, orthonormalize(), param, madness::print(), madness::World::rank(), madness::CalculationParameters::spin_restricted, madness::START_TIMER(), madness::WorldGopInterface::sum(), TAU_START, and TAU_STOP.
Referenced by solve().
void SCF::vector_stats | ( | const std::vector< double > & | v, |
double & | rms, | ||
double & | maxabsval | ||
) | const |
References sqrt().
Referenced by compute_residual(), and do_step_restriction().
tensorT madness::SCF::aeps |
orbital energies for alpha and beta orbitals
Referenced by madness::CIS::CIS(), initial_guess(), load_mos(), main(), madness::HartreeFock< T, NDIM >::orbital_energy(), save_mos(), madness::TDA::setup(), and solve().
vecfuncT madness::SCF::amo |
alpha and beta molecular orbitals
Referenced by madness::TDA_DFT::apply_kernel(), madness::CIS::CIS(), derivatives(), dipole(), do_plots(), initial_guess(), iterate_trotter(), load_mos(), loadbal(), main(), madness::HartreeFock< T, NDIM >::nemo(), madness::HartreeFock< T, NDIM >::nemos(), orthonormalize(), project(), propagate(), save_mos(), madness::CIS::solve(), solve(), madness::CIS::solve_internal_par(), madness::TDA_DFT::TDA_DFT(), and update_subspace().
vecfuncT madness::SCF::ao |
MRA projection of the minimal basis set.
Referenced by analyze_vectors(), initial_guess(), localize_PM(), and project_ao_basis().
AtomicBasisSet madness::SCF::aobasis |
Referenced by analyze_vectors(), initial_guess(), project_ao_basis(), and SCF().
tensorT madness::SCF::aocc |
occupation numbers for alpha and beta orbitals
Referenced by derivatives(), dipole(), do_plots(), madness::HartreeFock< T, NDIM >::get_coulomb_potential(), initial_guess(), iterate_trotter(), load_mos(), main(), propagate(), save_mos(), solve(), and update_subspace().
std::vector<int> madness::SCF::aset |
sets of orbitals grouped by their orbital energies (for localization?) only orbitals within the same set will be mixed to localize
Referenced by initial_guess(), load_mos(), save_mos(), and solve().
std::vector<int> madness::SCF::at_nbf |
Referenced by localize_PM(), SCF::localize_PM_task_kernel(), and project_ao_basis().
std::vector<int> madness::SCF::at_to_bf |
Referenced by localize_PM(), SCF::localize_PM_task_kernel(), and project_ao_basis().
tensorT madness::SCF::beps |
Referenced by initial_guess(), load_mos(), save_mos(), and solve().
vecfuncT madness::SCF::bmo |
Referenced by derivatives(), dipole(), do_plots(), initial_guess(), iterate_trotter(), load_mos(), loadbal(), project(), propagate(), save_mos(), solve(), and update_subspace().
tensorT madness::SCF::bocc |
Referenced by derivatives(), dipole(), do_plots(), initial_guess(), iterate_trotter(), load_mos(), save_mos(), solve(), and update_subspace().
std::vector<int> madness::SCF::bset |
Referenced by initial_guess(), load_mos(), save_mos(), and solve().
poperatorT madness::SCF::coulop |
Referenced by SCF::apply_hf_exchange(), apply_hf_exchange(), do_plots(), initial_guess(), iterate_trotter(), solve(), SCF::twoint(), and twoint().
double madness::SCF::current_energy |
Referenced by load_mos(), save_mos(), solve(), madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
std::vector< std::shared_ptr<real_derivative_3d> > madness::SCF::gradop |
Referenced by SCF::apply_potential(), SCF::kinetic_energy_matrix(), kinetic_energy_matrix(), and solve().
std::shared_ptr<GTHPseudopotential <double> > madness::SCF::gthpseudopotential |
Referenced by apply_potential(), initial_guess(), initial_load_bal(), loadbal(), make_nuclear_potential(), SCF(), and solve().
functionT madness::SCF::mask |
Molecule madness::SCF::molecule |
Referenced by analyze_vectors(), core_projection(), core_projector_derivative(), madness::create_nuclear_correlation_factor(), derivatives(), dipole(), initial_guess(), load_mos(), main(), madness::MP2::MP2(), project_ao_basis(), SCF(), solve(), madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
CalculationParameters madness::SCF::param |
Referenced by SCF::apply_potential(), apply_potential(), madness::CIS::CIS(), madness::create_nuclear_correlation_factor(), derivatives(), dipole(), do_plots(), do_step_restriction(), madness::HartreeFock< T, NDIM >::get_coulomb_potential(), initial_guess(), initial_load_bal(), iterate_trotter(), load_mos(), loadbal(), main(), make_bsh_operators(), make_nuclear_potential(), madness::MP2::MP2(), madness::HartreeFock< T, NDIM >::nemo(), madness::HartreeFock< T, NDIM >::nemos(), madness::HartreeFock< T, NDIM >::nocc(), madness::HartreeFock< T, NDIM >::orbital(), madness::HartreeFock< T, NDIM >::orbital_energy(), madness::HartreeFock< T, NDIM >::orbitals(), project(), propagate(), madness::HartreeFock< T, NDIM >::R2orbital(), madness::HartreeFock< T, NDIM >::R2orbitals(), save_mos(), SCF(), madness::TDA::setup(), solve(), update_subspace(), madness::HartreeFock< T, NDIM >::value(), and madness::MolecularEnergy::value().
std::shared_ptr<PotentialManager> madness::SCF::potentialmanager |
Referenced by apply_potential(), apply_V(), madness::create_nuclear_correlation_factor(), derivatives(), do_plots(), madness::HartreeFock< T, NDIM >::get_nuclear_potential(), initial_guess(), initial_load_bal(), iterate_trotter(), loadbal(), main(), make_nuclear_potential(), SCF(), solve(), and madness::HartreeFock< T, NDIM >::value().
|
static |
Referenced by initial_guess(), initial_load_bal(), and loadbal().
double madness::SCF::vtol |
XCfunctional madness::SCF::xc |
Referenced by SCF::apply_potential(), apply_potential(), SCF::make_dft_energy(), SCF::make_dft_potential(), SCF(), and solve().