MADNESS
0.10.1
|
#include <SCF.h>
Public Member Functions | |
SCF (World &world, const commandlineparser &parser) | |
collective constructor for SCF uses contents of file filename and broadcasts to all nodes More... | |
complex_functionT | APPLY (const complex_operatorT *q1d, const complex_functionT &psi) |
vecfuncT | apply_potential (World &world, const tensorT &occ, const vecfuncT &amo, const functionT &vlocal, double &exc, double &enl, int ispin) |
vecfuncT | compute_residual (World &world, tensorT &occ, tensorT &fock, const vecfuncT &psi, vecfuncT &Vpsi, double &err) |
void | copy_data (World &world, const SCF &other) |
tensorT | derivatives (World &world, const functionT &rho) const |
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... | |
tensorT | dipole (World &world, const functionT &rho) const |
compute the total dipole moment of the molecule More... | |
void | do_plots (World &world) |
double | do_step_restriction (World &world, const vecfuncT &mo, vecfuncT &mo_new, std::string spin) const |
perform step restriction following the KAIN solver More... | |
const vecfuncT & | get_amo () const |
getter for the molecular orbitals, alpha spin More... | |
const tensorT & | get_aocc () const |
getter for the occupation numbers, alpha spin More... | |
const vecfuncT & | get_bmo () const |
getter for the molecular orbitals, beta spin More... | |
const tensorT & | get_bocc () const |
getter for the occupation numbers, alpha spin More... | |
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... | |
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... | |
void | initial_guess (World &world) |
void | initial_load_bal (World &world) |
bool | is_spin_restricted () const |
distmatT | kinetic_energy_matrix (World &world, const vecfuncT &v) const |
void | load_mos (World &world) |
void | loadbal (World &world, functionT &arho, functionT &brho, functionT &arho_old, functionT &brho_old, subspaceT &subspace) |
std::vector< poperatorT > | make_bsh_operators (World &world, const tensorT &evals) const |
functionT | make_coulomb_potential (const functionT &rho) const |
make the Coulomb potential given the total density More... | |
functionT | make_density (World &world, const tensorT &occ, const cvecfuncT &v) |
functionT | make_density (World &world, const tensorT &occ, const vecfuncT &v) const |
double | make_dft_energy (World &world, const vecfuncT &vf, int ispin) |
tensorT | make_fock_matrix (World &world, const vecfuncT &psi, const vecfuncT &Vpsi, const tensorT &occ, double &ekinetic) const |
void | make_nuclear_potential (World &world) |
void | orthonormalize (World &world, vecfuncT &amo_new) const |
orthonormalize the vectors More... | |
void | orthonormalize (World &world, vecfuncT &amo_new, int nocc) const |
orthonormalize the vectors (symmetric in occupied spaced, gramm-schmidt for virt to occ) More... | |
void | output_calc_info_schema () const |
void | project (World &world) |
vecfuncT | project_ao_basis (World &world, const AtomicBasisSet &aobasis) |
void | reset_aobasis (const std::string &aobasisname) |
bool | restart_aos (World &world) |
void | rotate_subspace (World &world, const distmatT &U, subspaceT &subspace, int lo, int nfunc, double trantol) const |
void | rotate_subspace (World &world, const tensorT &U, subspaceT &subspace, int lo, int nfunc, double trantol) const |
void | save_mos (World &world) |
void | set_print_timings (const bool value) |
template<std::size_t NDIM> | |
void | set_protocol (World &world, double thresh) |
void | solve (World &world) |
Tensor< double > | twoint (World &world, const vecfuncT &psi) const |
Compute the two-electron integrals over the provided set of orbitals. More... | |
void | update_subspace (World &world, vecfuncT &Vpsia, vecfuncT &Vpsib, tensorT &focka, tensorT &fockb, subspaceT &subspace, tensorT &Q, double &bsh_residual, double &update_residual) |
void | vector_stats (const std::vector< double > &v, double &rms, double &maxabsval) const |
Static Public Member Functions | |
static void | analyze_vectors (World &world, const vecfuncT &mo, const vecfuncT &ao, double vtol, const Molecule &molecule, const int print_level, const AtomicBasisSet &aobasis, const tensorT &occ=tensorT(), const tensorT &energy=tensorT(), const std::vector< int > &set=std::vector< int >()) |
static void | help () |
static functionT | make_lda_potential (World &world, const functionT &arho) |
static void | print_parameters () |
static vecfuncT | project_ao_basis_only (World &world, const AtomicBasisSet &aobasis, const Molecule &molecule) |
Public Attributes | |
tensorT | aeps |
orbital energies for alpha and beta orbitals More... | |
vecfuncT | amo |
alpha and beta molecular orbitals More... | |
vecfuncT | ao |
MRA projection of the minimal basis set. More... | |
AtomicBasisSet | aobasis |
tensorT | aocc |
occupation numbers for alpha and beta orbitals More... | |
std::vector< int > | aset |
std::vector< int > | at_nbf |
std::vector< int > | at_to_bf |
tensorT | beps |
vecfuncT | bmo |
tensorT | bocc |
std::vector< int > | bset |
double | converged_for_thresh =1.e10 |
mos are converged for this threshold More... | |
poperatorT | coulop |
double | current_energy |
scf_data | e_data |
std::vector< std::shared_ptr< real_derivative_3d > > | gradop |
std::shared_ptr< GTHPseudopotential< double > > | gthpseudopotential |
functionT | mask |
Molecule | molecule |
CalculationParameters | param |
PCM | pcm |
std::shared_ptr< PotentialManager > | potentialmanager |
double | vtol |
XCfunctional | xc |
madness::SCF::SCF | ( | World & | world, |
const commandlineparser & | parser | ||
) |
collective constructor for SCF uses contents of file filename
and broadcasts to all nodes
collective constructor, reads input
on rank 0, broadcasts to all
collective constructor for SCF uses contents of stream input
and broadcasts to all nodes
References madness::CalculationParameters::aobasis(), aobasis, madness::WorldGopInterface::broadcast_serializable(), madness::Molecule::GeometryParameters::core_type(), madness::World::gop, madness::XCfunctional::initialize(), madness::CalculationParameters::L(), molecule, madness::CalculationParameters::nwfile(), param, madness::Molecule::parameters, madness::CalculationParameters::print_level(), madness::print_timings, PROFILE_MEMBER_FUNC, madness::World::rank(), madness::Molecule::read_core_file(), madness::AtomicBasisSet::read_nw_file(), reset_aobasis(), madness::FunctionDefaults< NDIM >::set_cubic_cell(), madness::CalculationParameters::set_derived_values(), madness::FunctionDefaults< NDIM >::set_truncate_mode(), madness::CalculationParameters::spin_restricted(), madness::CalculationParameters::xc(), and xc.
|
static |
References madness::_(), ao, aobasis, axis, C, madness::END_TIMER(), energy, madness::gesvp(), madness::inner(), madness::matrix_inner(), molecule, madness::mul_sparse(), madness::AtomicBasisSet::print_anal(), PROFILE_MEMBER_FUNC, madness::World::rank(), madness::rsquared(), madness::BaseTensor::size(), madness::START_TIMER(), madness::transpose(), and vtol.
Referenced by madness::Znemo::analyze(), madness::OEP::iterate(), and solve().
|
inline |
vecfuncT madness::SCF::apply_potential | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | amo, | ||
const functionT & | vlocal, | ||
double & | exc, | ||
double & | enl, | ||
int | ispin | ||
) |
References amo, madness::XCOperator< T, NDIM >::compute_xc_energy(), madness::copy(), madness::CalculationParameters::dft_deriv(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, gthpseudopotential, madness::XCfunctional::hf_exchange_coefficient(), madness::CalculationParameters::hfexalg(), madness::inner(), madness::XCfunctional::is_dft(), madness::XCfunctional::is_spin_polarized(), K(), madness::XCOperator< T, NDIM >::make_xc_potential(), molecule, madness::mul_sparse(), param, madness::Molecule::parameters, madness::CalculationParameters::print_level(), madness::print_meminfo(), PROFILE_MEMBER_FUNC, madness::Molecule::GeometryParameters::pure_ae(), madness::World::rank(), madness::START_TIMER(), madness::Function< T, NDIM >::truncate(), madness::truncate(), vtol, and xc.
Referenced by solve().
vecfuncT madness::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(), madness::norm2s(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, psi(), madness::World::rank(), madness::scale(), madness::set_thresh(), madness::START_TIMER(), madness::sub(), madness::transform(), madness::truncate(), vector_stats(), and vtol.
Referenced by update_subspace().
References amo, aocc, axis, bmo, bocc, madness::Molecule::GeometryParameters::core_type(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::func(), madness::Molecule::get_atom(), madness::World::gop, madness::inner(), madness::Molecule::is_potential_defined_atom(), molecule, madness::Molecule::natom(), madness::CalculationParameters::nbeta(), madness::Molecule::nuclear_repulsion_derivative(), param, madness::Molecule::parameters, potentialmanager, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, madness::World::rank(), madness::BaseTensor::size(), madness::CalculationParameters::spin_restricted(), madness::START_TIMER(), truncate_mode, madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by madness::MolecularEnergy::energy_and_gradient(), madness::MolecularEnergy::gradient(), and main().
tensorT madness::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 |
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::BaseTensor::dim(), madness::END_TIMER(), get_fock_transformation(), madness::matrix_inner(), madness::CalculationParameters::nalpha(), madness::normalize(), param, PROFILE_MEMBER_FUNC, psi(), madness::START_TIMER(), thresh, madness::transform(), madness::truncate(), and vtol.
Referenced by solve().
compute the total dipole moment of the molecule
[in] | rho | the total (alpha + beta) density |
References axis, madness::END_TIMER(), molecule, mu, madness::Molecule::nuclear_dipole(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, madness::World::rank(), and madness::START_TIMER().
void madness::SCF::do_plots | ( | World & | world | ) |
References amo, aocc, madness::apply(), bmo, bocc, bufsize, madness::copy(), coulop, madness::END_TIMER(), madness::QCCalculationParametersBase::get(), make_density(), madness::CalculationParameters::nalpha(), madness::CalculationParameters::nbeta(), param, madness::CalculationParameters::plot_cell(), madness::plotdx(), potentialmanager, PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), madness::Function< T, NDIM >::scale(), madness::BaseTensor::size(), madness::CalculationParameters::spin_restricted(), madness::START_TIMER(), madness::Function< T, NDIM >::truncate(), and vnuc().
Referenced by main().
double madness::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] | mo_new | vector of orbitals from the KAIN solver |
[in] | spin | "alpha" or "beta" for user information |
Limit maximum step size to make convergence more robust
[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::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, madness::World::rank(), madness::sub(), and vector_stats().
Referenced by update_subspace().
|
inline |
getter for the molecular orbitals, alpha spin
References amo.
Referenced by madness::Coulomb< T, NDIM >::compute_density().
|
inline |
getter for the occupation numbers, alpha spin
References aocc.
Referenced by madness::Coulomb< T, NDIM >::compute_density().
|
inline |
getter for the molecular orbitals, beta spin
References bmo.
Referenced by madness::Coulomb< T, NDIM >::compute_density().
|
inline |
getter for the occupation numbers, alpha spin
References bocc.
Referenced by madness::Coulomb< T, NDIM >::compute_density().
tensorT madness::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 |
[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::World::gop, PROFILE_MEMBER_FUNC, madness::Tensor< T >::ptr(), madness::BaseTensor::size(), madness::sygvp(), madness::Localizer::undo_degenerate_rotations(), and madness::Localizer::undo_reordering().
Referenced by diag_fock_matrix().
std::vector< int > madness::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 lo, madness::CalculationParameters::localize_pm(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, and madness::World::rank().
Referenced by initial_guess().
|
inlinestatic |
References madness::print(), and madness::print_header2().
Referenced by main().
void madness::SCF::initial_guess | ( | World & | world | ) |
References madness::_(), madness::LoadBalanceDeux< NDIM >::add_tree(), aeps, amo, ao, aobasis, aocc, madness::apply(), aset, madness::Atom::atomic_number, madness::atomic_number_to_symbol(), slymer::ES_Interface::atoms, slymer::Properties::Basis, slymer::ES_Interface::basis_set, beps, slymer::ES_Interface::beta_energies, slymer::ES_Interface::beta_MOs, slymer::ES_Interface::beta_occupancies, bmo, bocc, bset, c, madness::Function< T, NDIM >::clear(), madness::compress(), madness::copy(), madness::DistributedMatrix< T >::copy_to_replicated(), madness::Molecule::GeometryParameters::core_type(), coulop, e(), madness::END_TIMER(), slymer::Properties::Energies, slymer::ES_Interface::energies, slymer::ES_Interface::err, madness::QCCalculationParametersBase::get(), madness::Molecule::get_atom(), madness::Molecule::get_atom_charge(), madness::Molecule::get_atomic_number(), madness::Molecule::get_pseudo_atom(), group_orbital_sets(), gthpseudopotential, madness::inner(), kinetic_energy_matrix(), madness::LoadBalanceDeux< NDIM >::load_balance(), load_mos(), MADNESS_ASSERT, make_lda_potential(), make_nuclear_potential(), madness::matrix_inner(), madness::AtomicBasisSet::modify_dmat_psp(), molecule, slymer::Properties::MOs, slymer::ES_Interface::MOs, madness::mul_sparse(), madness::Molecule::n_core_orb_all(), madness::CalculationParameters::nalpha(), madness::Molecule::natom(), madness::CalculationParameters::nbeta(), madness::CalculationParameters::nmo_alpha(), madness::CalculationParameters::nmo_beta(), madness::normalize(), madness::CalculationParameters::nwfile(), slymer::Properties::Occupancies, slymer::ES_Interface::occupancies, param, madness::Molecule::parameters, potential(), potentialmanager, madness::print(), madness::CalculationParameters::print_level(), madness::print_meminfo(), PROFILE_MEMBER_FUNC, madness::Molecule::GeometryParameters::psp_calc(), madness::Molecule::GeometryParameters::pure_ae(), q(), madness::World::rank(), slymer::NWChem_Interface::read(), madness::Function< T, NDIM >::reconstruct(), madness::reconstruct(), madness::FunctionDefaults< NDIM >::redistribute(), madness::CalculationParameters::restart(), madness::Molecule::rotate(), madness::Function< T, NDIM >::scale(), sigma, madness::World::size(), madness::CalculationParameters::spin_restricted(), madness::START_TIMER(), madness::svd(), madness::sygvp(), madness::symbol_to_atomic_number(), madness::Function< T, NDIM >::trace(), madness::transform(), madness::Molecule::translate(), madness::transpose(), madness::Function< T, NDIM >::truncate(), madness::truncate(), vnuc(), madness::CalculationParameters::vnucextra(), vtol, madness::Atom::x, madness::Atom::y, and madness::Atom::z.
Referenced by madness::MolecularEnergy::value().
void madness::SCF::initial_load_bal | ( | World & | world | ) |
References madness::LoadBalanceDeux< NDIM >::add_tree(), gthpseudopotential, madness::LoadBalanceDeux< NDIM >::load_balance(), madness::CalculationParameters::loadbalparts(), molecule, param, madness::Molecule::parameters, potentialmanager, PROFILE_MEMBER_FUNC, madness::Molecule::GeometryParameters::psp_calc(), madness::Molecule::GeometryParameters::pure_ae(), madness::FunctionDefaults< NDIM >::redistribute(), vnuc(), and madness::CalculationParameters::vnucextra().
Referenced by main().
|
inline |
References madness::QCCalculationParametersBase::get(), and param.
Referenced by madness::Coulomb< T, NDIM >::compute_density().
References madness::apply(), madness::compress(), madness::DistributedMatrixDistribution::distribution(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::World::gop, gradop, madness::matrix_inner(), PROFILE_MEMBER_FUNC, madness::reconstruct(), madness::START_TIMER(), and v.
Referenced by initial_guess(), and make_fock_matrix().
void madness::SCF::load_mos | ( | World & | world | ) |
References aeps, amo, aocc, aset, beps, bmo, bocc, bset, converged_for_thresh, madness::copy(), current_energy, madness::WorldGopInterface::fence(), madness::FunctionDefaults< NDIM >::get_k(), madness::FunctionDefaults< NDIM >::get_thresh(), madness::World::gop, k, k1, L, madness::CalculationParameters::L(), madness::CalculationParameters::localize_method(), MADNESS_ASSERT, molecule, madness::Molecule::n_core_orb_all(), madness::CalculationParameters::nmo_alpha(), madness::CalculationParameters::nmo_beta(), param, madness::CalculationParameters::prefix(), madness::Molecule::print(), madness::print(), PROFILE_MEMBER_FUNC, madness::project(), madness::World::rank(), madness::reconstruct(), madness::set_thresh(), madness::CalculationParameters::spin_restricted(), thresh, madness::info::version(), and madness::CalculationParameters::xc().
Referenced by initial_guess(), and madness::MolecularEnergy::value().
void madness::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::loadbalparts(), molecule, madness::CalculationParameters::nbeta(), param, madness::Molecule::parameters, potentialmanager, madness::Molecule::GeometryParameters::psp_calc(), madness::Molecule::GeometryParameters::pure_ae(), madness::FunctionDefaults< NDIM >::redistribute(), madness::World::size(), madness::CalculationParameters::spin_restricted(), vnuc(), and madness::CalculationParameters::vnucextra().
Referenced by solve().
std::vector< poperatorT > madness::SCF::make_bsh_operators | ( | World & | world, |
const tensorT & | evals | ||
) | const |
References madness::BSHOperatorPtr3D(), madness::BaseTensor::dim(), madness::FunctionDefaults< NDIM >::get_thresh(), madness::CalculationParameters::lo(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, and madness::World::rank().
Referenced by compute_residual().
make the Coulomb potential given the total density
References madness::apply(), and coulop.
functionT madness::SCF::make_density | ( | World & | world, |
const tensorT & | occ, | ||
const vecfuncT & | v | ||
) | const |
References madness::Function< T, NDIM >::compress(), madness::compress(), madness::WorldGopInterface::fence(), madness::Function< T, NDIM >::gaxpy(), madness::World::gop, PROFILE_MEMBER_FUNC, madness::square(), and v.
Referenced by madness::Coulomb< double, 3 >::Coulomb(), madness::Coulomb< T, NDIM >::compute_density(), do_plots(), madness::MolecularEnergy::energy_and_gradient(), madness::MolecularEnergy::gradient(), main(), and solve().
References madness::Function< T, NDIM >::trace(), and xc.
tensorT madness::SCF::make_fock_matrix | ( | World & | world, |
const vecfuncT & | psi, | ||
const vecfuncT & | Vpsi, | ||
const tensorT & | occ, | ||
double & | ekinetic | ||
) | const |
References madness::LoadBalanceDeux< NDIM >::add_tree(), madness::copy(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::Tensor< T >::gaxpy(), madness::FunctionDefaults< NDIM >::get_pmap(), madness::World::gop, k, kinetic_energy_matrix(), madness::LoadBalanceDeux< NDIM >::load_balance(), madness::CalculationParameters::loadbalparts(), madness::matrix_inner(), param, PROFILE_MEMBER_FUNC, psi(), madness::FunctionDefaults< NDIM >::set_pmap(), madness::BaseTensor::size(), madness::World::size(), madness::START_TIMER(), and madness::transpose().
Referenced by solve().
References madness::copy(), PROFILE_MEMBER_FUNC, madness::Function< T, NDIM >::reconstruct(), and madness::Function< T, NDIM >::unaryop().
Referenced by initial_guess().
void madness::SCF::make_nuclear_potential | ( | World & | world | ) |
References madness::Molecule::GeometryParameters::core_type(), madness::END_TIMER(), gthpseudopotential, molecule, madness::Molecule::parameters, potentialmanager, PROFILE_MEMBER_FUNC, madness::Molecule::GeometryParameters::psp_calc(), madness::Molecule::GeometryParameters::pure_ae(), and madness::START_TIMER().
Referenced by dnuclear_anchor_test(), initial_guess(), main(), nuclear_anchor_test(), and madness::MolecularEnergy::value().
orthonormalize the vectors
orthonormalize the vectors ignoring occupied/virtual distinctions
[in] | world | the world |
[in,out] | amo_new | the vectors to be orthonormalized |
[in] | world | the world |
[in,out] | amo_new | the vectors to be orthonormalized |
References std::abs(), amo, madness::END_TIMER(), madness::matrix_inner(), max, madness::normalize(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, Q(), madness::Q2(), madness::World::rank(), madness::START_TIMER(), madness::transform(), madness::truncate(), and vtol.
Referenced by restart_aos(), and update_subspace().
orthonormalize the vectors (symmetric in occupied spaced, gramm-schmidt for virt to occ)
[in] | world | the world |
[in,out] | amo_new | the vectors to be orthonormalized |
References std::abs(), madness::END_TIMER(), madness::matrix_inner(), max, madness::normalize(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, Q(), madness::Q2(), madness::World::rank(), madness::START_TIMER(), madness::transform(), madness::truncate(), and vtol.
void madness::SCF::output_calc_info_schema | ( | ) | const |
References molecule, param, InputParameters::prefix, madness::World::rank(), and madness::to_json().
|
inlinestatic |
References param, madness::QCCalculationParametersBase::print(), madness::print(), and madness::Molecule::print_parameters().
Referenced by main().
void madness::SCF::project | ( | World & | world | ) |
References amo, bmo, madness::WorldGopInterface::fence(), madness::World::gop, madness::CalculationParameters::nbeta(), madness::normalize(), param, PROFILE_MEMBER_FUNC, madness::project(), madness::reconstruct(), madness::CalculationParameters::spin_restricted(), and madness::truncate().
Referenced by madness::MolecularEnergy::value().
vecfuncT madness::SCF::project_ao_basis | ( | World & | world, |
const AtomicBasisSet & | aobasis | ||
) |
References aobasis, at_nbf, at_to_bf, madness::AtomicBasisSet::atoms_to_bfn(), molecule, PROFILE_MEMBER_FUNC, and project_ao_basis_only().
Referenced by madness::MolecularEnergy::value().
|
static |
References ao, aobasis, madness::WorldGopInterface::fence(), madness::AtomicBasisSet::get_atomic_basis_function(), madness::World::gop, molecule, madness::AtomicBasisSet::nbf(), madness::normalize(), and madness::truncate().
Referenced by madness::Znemo::analyze(), madness::Znemo::hcore_guess(), main(), project_ao_basis(), madness::MolecularOrbitals< T, NDIM >::read_restartaodata(), madness::MolecularOrbitals< T, NDIM >::save_restartaodata(), and test_read_restartaodata().
|
inline |
References aobasis, and madness::AtomicBasisSet::read_file().
Referenced by SCF(), and madness::MolecularEnergy::value().
bool madness::SCF::restart_aos | ( | World & | world | ) |
References aeps, amo, ao, aocc, aset, beps, bmo, bocc, madness::WorldGopInterface::broadcast(), madness::WorldGopInterface::broadcast_serializable(), bset, c, madness::BaseTensor::dim(), fred(), madness::gesvp(), madness::World::gop, madness::matrix_inner(), madness::CalculationParameters::nalpha(), madness::CalculationParameters::nbeta(), madness::CalculationParameters::nmo_alpha(), madness::CalculationParameters::nmo_beta(), orthonormalize(), param, madness::CalculationParameters::prefix(), madness::print(), madness::World::rank(), madness::CalculationParameters::spin_restricted(), madness::transform(), madness::truncate(), and vtol.
Referenced by madness::MolecularEnergy::value().
void madness::SCF::rotate_subspace | ( | World & | world, |
const distmatT & | U, | ||
subspaceT & | subspace, | ||
int | lo, | ||
int | nfunc, | ||
double | trantol | ||
) | const |
References madness::WorldGopInterface::fence(), madness::World::gop, lo, PROFILE_MEMBER_FUNC, subspace, madness::transform(), and v.
void madness::SCF::rotate_subspace | ( | World & | world, |
const tensorT & | U, | ||
subspaceT & | subspace, | ||
int | lo, | ||
int | nfunc, | ||
double | trantol | ||
) | const |
References madness::WorldGopInterface::fence(), madness::World::gop, lo, PROFILE_MEMBER_FUNC, subspace, madness::transform(), and v.
void madness::SCF::save_mos | ( | World & | world | ) |
References aeps, amo, ao, aocc, aset, beps, bmo, bocc, bset, converged_for_thresh, current_energy, madness::QCCalculationParametersBase::get(), madness::FunctionDefaults< NDIM >::get_k(), madness::CalculationParameters::L(), madness::CalculationParameters::localize_method(), madness::matrix_inner(), molecule, madness::CalculationParameters::nwfile(), param, madness::CalculationParameters::prefix(), PROFILE_MEMBER_FUNC, madness::World::rank(), madness::CalculationParameters::spin_restricted(), madness::info::version(), and madness::CalculationParameters::xc().
Referenced by madness::MolecularEnergy::value().
void madness::SCF::set_print_timings | ( | const bool | value | ) |
References madness::print_timings.
|
inline |
References madness::CoulombOperatorPtr(), coulop, madness::CalculationParameters::dconv(), madness::CalculationParameters::deriv(), madness::f, madness::FunctionDefaults< NDIM >::get_thresh(), gradop, k, madness::CalculationParameters::k(), madness::CalculationParameters::L(), madness::CalculationParameters::lo(), mask, madness::mask3(), max, NDIM, param, madness::print(), madness::CalculationParameters::print_level(), madness::World::rank(), madness::FunctionDefaults< NDIM >::set_apply_randomize(), madness::FunctionDefaults< NDIM >::set_autorefine(), madness::FunctionDefaults< NDIM >::set_cubic_cell(), 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(), thresh, and vtol.
Referenced by main(), run_all_calculations(), test_read_restartaodata(), test_read_restartdata(), and madness::MolecularEnergy::value().
void madness::SCF::solve | ( | World & | world | ) |
References madness::scf_data::add_data(), aeps, amo, analyze_vectors(), ao, aobasis, aocc, madness::apply(), apply_potential(), aset, beps, bmo, bocc, bset, madness::Function< T, NDIM >::clear(), madness::Localizer::compute_localization_matrix(), madness::PCM::compute_pcm_energy(), madness::PCM::compute_pcm_potential(), converged_for_thresh, coulop, current_energy, madness::CalculationParameters::dconv(), diag_fock_matrix(), dipole(), madness::CalculationParameters::do_localize(), e_data, madness::CalculationParameters::econv(), madness::END_TIMER(), madness::f, madness::QCCalculationParametersBase::get(), gthpseudopotential, madness::inner(), loadbal(), madness::CalculationParameters::localize_method(), make_density(), make_fock_matrix(), madness::matrix_inner(), max, madness::CalculationParameters::maxiter(), molecule, madness::name(), madness::Molecule::natom(), madness::CalculationParameters::nbeta(), madness::norm2(), madness::normalize(), madness::Molecule::nuclear_repulsion_energy(), madness::CalculationParameters::nwfile(), param, madness::Molecule::parameters, pcm, madness::CalculationParameters::pcm_data(), potentialmanager, madness::CalculationParameters::prefix(), madness::print(), madness::CalculationParameters::print_level(), madness::print_meminfo(), PROFILE_MEMBER_FUNC, madness::Molecule::GeometryParameters::psp_calc(), madness::Molecule::GeometryParameters::pure_ae(), Q(), madness::World::rank(), madness::Function< T, NDIM >::scale(), madness::Localizer::set_method(), madness::CalculationParameters::spin_restricted(), madness::START_TIMER(), subspace, madness::sygvp(), madness::transform(), madness::transpose(), madness::Function< T, NDIM >::truncate(), madness::truncate(), update_subspace(), vnuc(), vtol, and madness::wall_time().
Referenced by 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::FunctionDefaults< NDIM >::get_thresh(), madness::World::gop, madness::matrix_inner(), madness::mul_sparse(), madness::norm_tree(), PROFILE_MEMBER_FUNC, psi(), madness::reconstruct(), and madness::truncate().
void madness::SCF::update_subspace | ( | World & | world, |
vecfuncT & | Vpsia, | ||
vecfuncT & | Vpsib, | ||
tensorT & | focka, | ||
tensorT & | fockb, | ||
subspaceT & | subspace, | ||
tensorT & | Q, | ||
double & | bsh_residual, | ||
double & | update_residual | ||
) |
References madness::_(), amo, aocc, bmo, bocc, madness::WorldGopInterface::broadcast(), madness::WorldGopInterface::broadcast_serializable(), c, madness::compress(), compute_residual(), do_step_restriction(), e(), madness::END_TIMER(), madness::WorldGopInterface::fence(), madness::gaxpy(), madness::World::gop, madness::KAIN(), m, max, madness::CalculationParameters::maxsub(), madness::CalculationParameters::nalpha(), madness::CalculationParameters::nbeta(), madness::CalculationParameters::nmo_alpha(), madness::CalculationParameters::nmo_beta(), orthonormalize(), param, madness::print(), madness::CalculationParameters::print_level(), PROFILE_MEMBER_FUNC, madness::Tensor< T >::ptr(), Q(), madness::World::rank(), madness::CalculationParameters::spin_restricted(), madness::START_TIMER(), subspace, and madness::WorldGopInterface::sum().
Referenced by solve().
void madness::SCF::vector_stats | ( | const std::vector< double > & | v, |
double & | rms, | ||
double & | maxabsval | ||
) | const |
References std::abs(), PROFILE_MEMBER_FUNC, and v.
Referenced by compute_residual(), and do_step_restriction().
tensorT madness::SCF::aeps |
orbital energies for alpha and beta orbitals
Referenced by compare_calc_and_mos(), copy_data(), initial_guess(), load_mos(), madness::MolecularEnergy::output_calc_info_schema(), restart_aos(), save_mos(), and solve().
vecfuncT madness::SCF::amo |
alpha and beta molecular orbitals
Referenced by madness::Exchange< T, NDIM >::ExchangeImpl< T, NDIM >::ExchangeImpl(), apply_potential(), compare_calc_and_mos(), copy_data(), derivatives(), do_plots(), madness::MolecularEnergy::energy_and_gradient(), get_amo(), madness::MolecularEnergy::gradient(), initial_guess(), load_mos(), loadbal(), main(), orthonormalize(), project(), restart_aos(), save_mos(), solve(), update_subspace(), and madness::MolecularEnergy::value().
vecfuncT madness::SCF::ao |
MRA projection of the minimal basis set.
Referenced by analyze_vectors(), copy_data(), initial_guess(), project_ao_basis_only(), restart_aos(), save_mos(), solve(), and madness::MolecularEnergy::value().
AtomicBasisSet madness::SCF::aobasis |
tensorT madness::SCF::aocc |
occupation numbers for alpha and beta orbitals
Referenced by compare_calc_and_mos(), copy_data(), derivatives(), do_plots(), madness::MolecularEnergy::energy_and_gradient(), get_aocc(), madness::MolecularEnergy::gradient(), initial_guess(), load_mos(), main(), restart_aos(), save_mos(), solve(), update_subspace(), and madness::MolecularEnergy::value().
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 compare_calc_and_mos(), copy_data(), initial_guess(), load_mos(), restart_aos(), save_mos(), and solve().
std::vector<int> madness::SCF::at_nbf |
Referenced by copy_data(), and project_ao_basis().
std::vector<int> madness::SCF::at_to_bf |
Referenced by copy_data(), and project_ao_basis().
tensorT madness::SCF::beps |
Referenced by copy_data(), initial_guess(), load_mos(), madness::MolecularEnergy::output_calc_info_schema(), restart_aos(), save_mos(), and solve().
vecfuncT madness::SCF::bmo |
Referenced by madness::Coulomb< double, 3 >::Coulomb(), madness::Exchange< T, NDIM >::ExchangeImpl< T, NDIM >::ExchangeImpl(), copy_data(), derivatives(), do_plots(), madness::MolecularEnergy::energy_and_gradient(), get_bmo(), madness::MolecularEnergy::gradient(), initial_guess(), load_mos(), loadbal(), main(), project(), restart_aos(), save_mos(), solve(), update_subspace(), and madness::MolecularEnergy::value().
tensorT madness::SCF::bocc |
Referenced by madness::Coulomb< double, 3 >::Coulomb(), copy_data(), derivatives(), do_plots(), madness::MolecularEnergy::energy_and_gradient(), get_bocc(), madness::MolecularEnergy::gradient(), initial_guess(), load_mos(), main(), restart_aos(), save_mos(), solve(), update_subspace(), and madness::MolecularEnergy::value().
std::vector<int> madness::SCF::bset |
Referenced by copy_data(), initial_guess(), load_mos(), restart_aos(), save_mos(), and solve().
double madness::SCF::converged_for_thresh =1.e10 |
mos are converged for this threshold
Referenced by load_mos(), save_mos(), and solve().
poperatorT madness::SCF::coulop |
Referenced by do_plots(), initial_guess(), make_coulomb_potential(), set_protocol(), solve(), and twoint().
double madness::SCF::current_energy |
Referenced by madness::MolecularEnergy::energy_and_gradient(), load_mos(), save_mos(), solve(), and madness::MolecularEnergy::value().
scf_data madness::SCF::e_data |
Referenced by madness::MolecularEnergy::output_calc_info_schema(), and solve().
std::vector<std::shared_ptr<real_derivative_3d> > madness::SCF::gradop |
Referenced by kinetic_energy_matrix(), and set_protocol().
std::shared_ptr<GTHPseudopotential<double> > madness::SCF::gthpseudopotential |
Referenced by apply_potential(), initial_guess(), initial_load_bal(), loadbal(), make_nuclear_potential(), and solve().
functionT madness::SCF::mask |
Referenced by set_protocol().
Molecule madness::SCF::molecule |
Referenced by madness::DNuclear< T, NDIM >::DNuclear(), madness::Nuclear< T, NDIM >::Nuclear(), SCF(), analyze_vectors(), apply_potential(), derivatives(), dipole(), dnuclear_anchor_test(), initial_guess(), initial_load_bal(), load_mos(), loadbal(), main(), make_nuclear_potential(), nuclear_anchor_test(), madness::MolecularEnergy::output_calc_info_schema(), project_ao_basis(), project_ao_basis_only(), run_all_calculations(), save_mos(), solve(), test_read_restartaodata(), test_read_restartdata(), and madness::MolecularEnergy::value().
CalculationParameters madness::SCF::param |
Referenced by madness::Coulomb< T, NDIM >::Coulomb(), SCF(), madness::XCOperator< T, NDIM >::XCOperator(), apply_potential(), compute_residual(), derivatives(), diag_fock_matrix(), dipole(), dnuclear_anchor_test(), do_plots(), do_step_restriction(), madness::MolecularEnergy::energy_and_gradient(), madness::MolecularEnergy::gradient(), group_orbital_sets(), initial_guess(), initial_load_bal(), is_spin_restricted(), load_mos(), loadbal(), main(), make_bsh_operators(), make_fock_matrix(), nuclear_anchor_test(), orthonormalize(), madness::MolecularEnergy::output_calc_info_schema(), print_parameters(), project(), restart_aos(), save_mos(), set_protocol(), solve(), test_read_restartaodata(), test_read_restartdata(), update_subspace(), and madness::MolecularEnergy::value().
PCM madness::SCF::pcm |
Referenced by solve(), and madness::MolecularEnergy::value().
std::shared_ptr<PotentialManager> madness::SCF::potentialmanager |
double madness::SCF::vtol |
Referenced by analyze_vectors(), apply_potential(), compute_residual(), diag_fock_matrix(), initial_guess(), orthonormalize(), restart_aos(), set_protocol(), and solve().
XCfunctional madness::SCF::xc |
Referenced by SCF(), apply_potential(), and make_dft_energy().