MADNESS  version 0.9
Public Member Functions | Public Attributes | Static Public Attributes | List of all members
madness::SCF Class Reference

#include <SCF.h>

Collaboration diagram for madness::SCF:
Collaboration graph
[legend]

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< poperatorTmake_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< PotentialManagerpotentialmanager
 
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
 

Constructor & Destructor Documentation

SCF::SCF ( World world,
const char *  filename 
)

Member Function Documentation

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>() 
)
complex_functionT madness::SCF::APPLY ( const complex_operatorT q1d,
const complex_functionT psi 
)
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

Parameters
[in]worldthe world
[in]occoccupation numbers
[in]psithe orbitals in the exchange operator
[in]fthe orbitals |i> that the operator is applied on
Returns
a vector of orbitals K| i>

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 
)
vecfuncT SCF::compute_residual ( World world,
tensorT occ,
tensorT fock,
const vecfuncT psi,
vecfuncT Vpsi,
double &  err 
)
vecfuncT SCF::core_projection ( World world,
const vecfuncT psi,
const bool  include_Bc = true 
)
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,
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

Parameters
[in]worldthe world
[in,out]fockthe fock matrix (diagonal upon exit)
[in,out]psithe orbitals
[in,out]Vpsithe orbital times the potential
[out]evalsthe orbital energies
[in]occoccupation numbers
[in]threshthreshold for rotation and truncation
Returns
the unitary matrix U: U^T F U = evals

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().

double SCF::DIP ( const tensorT dip,
int  i,
int  j,
int  k,
int  l 
)
inline

Referenced by SCF::localize_boys().

tensorT SCF::dipole ( World world)
void SCF::do_plots ( World world)
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

Parameters
[in]worldthe world
[in]movector of orbitals from previous iteration
[in,out]new_movector of orbitals from the KAIN solver
[in]spin"alpha" or "beta" for user information
Returns
max residual

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

Parameters
[in]worldthe world
[in]overlapthe overlap matrix of the orbitals
[in,out]fockthe fock matrix; diagonal upon exit
[out]evalsthe orbital energies
[in]occthe occupation numbers
[in]thresh_degeneratethreshold for orbitals being degenerate
Returns
the unitary matrix U: U^T F U = evals

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

Parameters
[in]epsorbital energies
[in]occoccupation numbers
[in]nmonumber of MOs for the given spin
Returns
vector of length nmo with the set index for each MO

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)
void SCF::iterate_trotter ( World world,
Convolution1D< double_complex > *  G,
cvecfuncT camo,
cvecfuncT cbmo,
double  t,
double  time_step,
double  thresh 
)
distmatT SCF::kinetic_energy_matrix ( World world,
const vecfuncT v 
) const
void SCF::load_mos ( World world)
void SCF::loadbal ( World world,
functionT arho,
functionT brho,
functionT arho_old,
functionT brho_old,
subspaceT subspace 
)
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

Parameters
[in]worldthe world
[in]mothe MOs
[in]setonly orbitals within the same set will be mixed
[in]threshthe 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 
)
functionT madness::SCF::make_coulomb_potential ( const functionT rho) const
inline

make the Coulomb potential given the total density

References madness::apply().

Referenced by madness::HartreeFock< T, NDIM >::get_coulomb_potential().

functionT SCF::make_density ( World world,
const tensorT occ,
const vecfuncT v 
) const
functionT SCF::make_density ( World world,
const tensorT occ,
const cvecfuncT v 
)
double madness::SCF::make_dft_energy ( World world,
const vecfuncT vf,
int  ispin 
)
inline
functionT madness::SCF::make_dft_potential ( World world,
const vecfuncT vf,
int  ispin,
int  what 
)
inline
tensorT SCF::make_fock_matrix ( World world,
const vecfuncT psi,
const vecfuncT Vpsi,
const tensorT occ,
double &  ekinetic 
) const
functionT SCF::make_lda_potential ( World world,
const functionT arho 
)
void SCF::make_nuclear_potential ( World world)
tensorT SCF::matrix_exponential ( const tensorT A) const
void SCF::orthonormalize ( World world,
vecfuncT amo_new 
) const

orthonormalize the vectors

Parameters
[in]worldthe world
[in,out]amo_newthe 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)
void SCF::project_ao_basis ( World world)
void SCF::propagate ( World world,
double  omega,
int  step0 
)
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)
template<std::size_t NDIM>
void madness::SCF::set_protocol ( World world,
double  thresh 
)
inline
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().

Tensor< double > SCF::twoint ( World world,
const vecfuncT psi 
) const

Compute the two-electron integrals over the provided set of orbitals.

Returned is a replicated tensor of $(ij|kl)$ with $i>=j$ and $k>=l$. The symmetry $(ij|kl)=(kl|ij)$ 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 
)
void SCF::vector_stats ( const std::vector< double > &  v,
double &  rms,
double &  maxabsval 
) const

References sqrt().

Referenced by compute_residual(), and do_step_restriction().

Member Data Documentation

tensorT madness::SCF::aeps
vecfuncT madness::SCF::amo
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
tensorT madness::SCF::aocc
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
std::vector<int> madness::SCF::at_to_bf
tensorT madness::SCF::beps
vecfuncT madness::SCF::bmo
tensorT madness::SCF::bocc
std::vector<int> madness::SCF::bset
poperatorT madness::SCF::coulop
double madness::SCF::current_energy
std::vector< std::shared_ptr<real_derivative_3d> > madness::SCF::gradop
std::shared_ptr<GTHPseudopotential <double> > madness::SCF::gthpseudopotential
functionT madness::SCF::mask
Molecule madness::SCF::molecule
CalculationParameters madness::SCF::param
std::shared_ptr<PotentialManager> madness::SCF::potentialmanager
const int madness::SCF::vnucextra = 12
static
double madness::SCF::vtol
XCfunctional madness::SCF::xc

The documentation for this class was generated from the following files: