MADNESS  0.10.1
Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
madness::MP2 Class Reference

a class for computing the first order wave function and MP2 pair energies More...

#include <mp2.h>

Inheritance diagram for madness::MP2:
Inheritance graph
[legend]
Collaboration diagram for madness::MP2:
Collaboration graph
[legend]

Classes

struct  Pairs
 POD holding all electron pairs with easy access. More...
 
struct  Parameters
 POD for MP2 keywords. More...
 

Public Member Functions

 MP2 (World &world, const commandlineparser &parser)
 ctor More...
 
double asymmetry (const real_function_6d &f, const std::string s) const
 
bool check_core_valence_separation () const
 make sure frozen orbitals don't couple with correlated ones – relocalize if necessary More...
 
double compute_gQf (const int i, const int j, ElectronPair &pair) const
 compute the matrix element <ij | g12 Q12 f12 | phi^0> More...
 
double compute_gQf_cc2interface (const int i, const int j, const real_function_6d &f) const
 compute the matrix element <ij | g12 Q12 f12 | phi^0> More...
 
double coord_chksum () const
 return a checksum for the geometry More...
 
real_function_6d debug_cc2 (const real_function_6d &f, const size_t &i, const size_t &j) const
 
void enforce_core_valence_separation ()
 make sure frozen orbitals don't couple with correlated ones – relocalize if necessary More...
 
HartreeFockget_hf ()
 return the underlying HF reference More...
 
real_function_6d get_residue (const real_function_6d &f, const int i, const int j)
 
void increment (ElectronPair &pair, real_convolution_6d &green)
 compute increments: psi^1 = C + GV C + GVGV C + GVGVGV C + .. More...
 
real_function_6d iterate (const real_function_6d &f) const
 solve the residual equation for electron pair (i,j) More...
 
template<typename T , size_t NDIM>
void load_function (Function< T, NDIM > &f, const std::string name) const
 load a function More...
 
std::string name () const
 
template<typename T >
void print_options (const std::string option, const T val) const
 pretty print the options More...
 
template<typename T , size_t NDIM>
void save_function (const Function< T, NDIM > &f, const std::string name) const
 save a function More...
 
virtual bool selftest ()
 
double solve_coupled_equations (Pairs< ElectronPair > &pairs, const double econv, const double dconv) const
 solve the coupled MP1 equations (e.g. for local orbitals) More...
 
void solve_residual_equations (ElectronPair &pair, const double econv, const double dconv) const
 solve the residual equation for electron pair (i,j) More...
 
double value ()
 return the molecular correlation energy energy (without the HF energy) More...
 
double value (const Tensor< double > &x)
 return the molecular correlation energy as a function of the coordinates More...
 
double zeroth_order_energy (const int i, const int j) const
 return the 0th order energy of pair ij (= sum of orbital energies) More...
 
- Public Member Functions inherited from madness::OptimizationTargetInterface
virtual ~OptimizationTargetInterface ()
 
virtual Tensor< double > gradient (const Tensor< double > &x)
 Should return the derivative of the function. More...
 
virtual bool provides_gradient () const
 Override this to return true if the derivative is implemented. More...
 
double test_gradient (Tensor< double > &x, double value_precision, bool doprint=true)
 Numerical test of the derivative ... optionally prints to stdout, returns max abs error. More...
 
virtual void value_and_gradient (const Tensor< double > &x, double &value, Tensor< double > &gradient)
 Reimplement if more efficient to evaluate both value and gradient in one call. More...
 

Static Public Member Functions

static void help ()
 
static void print_parameters ()
 

Private Member Functions

void add_local_coupling (const Pairs< ElectronPair > &pairs, Pairs< real_function_6d > &coupling) const
 add the coupling terms for local MP2 More...
 
real_function_6d apply_exchange (const real_function_6d &f, const real_function_3d &orbital_ket, const real_function_3d &orbital_bra, const int particle) const
 apply the exchange operator on f More...
 
real_function_6d apply_exchange_vector (const real_function_6d &f, const int particle) const
 
double compute_energy (ElectronPair &pair) const
 compute the singlet and triplet energy for a given electron pair More...
 
void END_TIMER (World &world, const char *msg) const
 
Tensor< double > get_fock_matrix () const
 
void guess_mp1_3 (ElectronPair &pair) const
 compute the first iteration of the residual equations and all intermediates More...
 
real_function_3d J (const real_function_3d &phi) const
 apply the Coulomb operator a on orbital More...
 
real_function_6d JK1phi0_on_demand (const int i, const int j, const bool hc=false) const
 return the function (J(1)-K(1)) |phi0> as on-demand function More...
 
real_function_6d JK2phi0_on_demand (const int i, const int j, const bool hc=false) const
 return the function (J(2)-K(2)) |phi0> as on-demand function More...
 
real_function_3d K (const real_function_3d &phi, const bool hc=false) const
 apply the exchange operator on an orbital More...
 
real_function_6d K (const real_function_6d &phi, const bool is_symmetric=false) const
 apply the exchange operator on a pair function More...
 
std::vector< real_function_3dmake_chi (const real_function_3d &phi, const real_convolution_3d &op, const bool hc=false) const
 make the quantity chi_k More...
 
real_function_6d make_fKphi0 (const int i, const int j) const
 apply the operator K on the reference and multiply with f; fK |phi^0> More...
 
real_function_6d make_KffKphi0 (const ElectronPair &pair) const
 return the function [K,f] phi0; load from disk if available More...
 
ElectronPair make_pair (const int i, const int j) const
 compute some matrix elements that don't change during the SCF More...
 
real_function_6d make_Uphi0 (ElectronPair &pair) const
 return the function Uphi0; load from disk if available More...
 
std::vector< real_function_3dmake_xi (const real_function_3d &phi_i, const real_function_3d &phi_j, const real_convolution_3d &op, const bool hc=false) const
 make the quantity xi_k More...
 
real_function_6d multiply_with_0th_order_Hamiltonian (const real_function_6d &f, const int i, const int j) const
 multiply the given function with the 0th order Hamiltonian, exluding the 0th order energy More...
 
real_function_6d nemo0_on_demand (const int i, const int j) const
 return the function |F1F2> as on-demand function More...
 
real_function_6d phi0_on_demand (const int i, const int j) const
 return the function |phi0> as on-demand function More...
 
void START_TIMER (World &world) const
 

Private Attributes

double coords_sum
 check sum for the geometry More...
 
double correlation_energy
 the correlation energy More...
 
CorrelationFactor corrfac
 correlation factor: Slater More...
 
Tensor< double > fock
 the Fock matrix More...
 
std::shared_ptr< HartreeFockhf
 our reference More...
 
std::shared_ptr< NuclearCorrelationFactornuclear_corrfac
 
Pairs< ElectronPairpairs
 pair functions and energies More...
 
Parameters param
 SCF parameters for MP2. More...
 
std::shared_ptr< real_convolution_3dpoisson
 
StrongOrthogonalityProjector< double, 3 > Q12
 
double sss
 
double ttt
 
Worldworld
 the world More...
 

Detailed Description

a class for computing the first order wave function and MP2 pair energies

Constructor & Destructor Documentation

◆ MP2()

madness::MP2::MP2 ( World world,
const commandlineparser parser 
)

Member Function Documentation

◆ add_local_coupling()

void madness::MP2::add_local_coupling ( const Pairs< ElectronPair > &  pairs,
Pairs< real_function_6d > &  coupling 
) const
private

add the coupling terms for local MP2

\sum_{k\neq i} f_ki |u_kj> + \sum_{l\neq j} f_lj |u_il>

Todo:
Verify this doxygen block.

Referenced by solve_coupled_equations().

◆ apply_exchange()

real_function_6d madness::MP2::apply_exchange ( const real_function_6d f,
const real_function_3d orbital_ket,
const real_function_3d orbital_bra,
const int  particle 
) const
private

apply the exchange operator on f

if the exchange operator is similarity transformed (R-1 K R) the orbital spaces differ for the orbitals underneath the integral sign R-1 K R = \phi-ket(1) * \int \phi-bra(1') * f(1',2)

Parameters
[in]fthe pair function
[in]orbital_brathe orbital underneath the integral sign (typically R2orbitals)
[in]orbital_ketthe orbital to be pre-multiplied with (typically orbitals)
Returns
the pair function, on which the exchange operator has been applied

◆ apply_exchange_vector()

real_function_6d madness::MP2::apply_exchange_vector ( const real_function_6d f,
const int  particle 
) const
private

◆ asymmetry()

double madness::MP2::asymmetry ( const real_function_6d f,
const std::string  s 
) const

◆ check_core_valence_separation()

bool madness::MP2::check_core_valence_separation ( ) const

◆ compute_energy()

double madness::MP2::compute_energy ( ElectronPair pair) const
private

compute the singlet and triplet energy for a given electron pair

Returns
the energy of 1 degenerate triplet and 1 singlet pair

Referenced by increment(), solve_coupled_equations(), and solve_residual_equations().

◆ compute_gQf()

double madness::MP2::compute_gQf ( const int  i,
const int  j,
ElectronPair pair 
) const

compute the matrix element <ij | g12 Q12 f12 | phi^0>

scales quartically. I think I can get this down to cubically by setting Q12 = (1 - O1)(1 - O2) = 1- O1(1 - 0.5 O2) - O2 (1 - 0.5 O1) as for the formulas cf the article mra_molecule

Returns
the energy <ij | g Q f | kl>

References a, bsh_eps, corrfac, madness::CoulombOperator(), e(), madness::f, madness::g, hf, madness::ElectronPair::i, madness::inner(), madness::ElectronPair::j, k, lo, make_xi(), madness::matrix_inner(), madness::OT_F12, madness::OT_G12, param, madness::constants::pi, poisson, madness::print(), Q12, madness::World::rank(), madness::CCTimer::reset(), madness::SlaterF12Operator(), timer(), madness::Tensor< T >::trace(), and world.

Referenced by compute_gQf_cc2interface(), and iterate().

◆ compute_gQf_cc2interface()

double madness::MP2::compute_gQf_cc2interface ( const int  i,
const int  j,
const real_function_6d f 
) const
inline

compute the matrix element <ij | g12 Q12 f12 | phi^0>

scales quartically. I think I can get this down to cubically by setting Q12 = (1 - O1)(1 - O2) = 1- O1(1 - 0.5 O2) - O2 (1 - 0.5 O1) as for the formulas cf the article mra_molecule

Returns
the energy <ij | g Q f | kl>

References compute_gQf(), madness::ElectronPair::constant_term, madness::copy(), madness::f, and madness::ElectronPair::function.

◆ coord_chksum()

double madness::MP2::coord_chksum ( ) const
inline

return a checksum for the geometry

References coords_sum.

Referenced by value().

◆ debug_cc2()

real_function_6d madness::MP2::debug_cc2 ( const real_function_6d f,
const size_t &  i,
const size_t &  j 
) const
inline

◆ END_TIMER()

void madness::MP2::END_TIMER ( World world,
const char *  msg 
) const
inlineprivate

◆ enforce_core_valence_separation()

void madness::MP2::enforce_core_valence_separation ( )

◆ get_fock_matrix()

Tensor<double> madness::MP2::get_fock_matrix ( ) const
inlineprivate

◆ get_hf()

HartreeFock& madness::MP2::get_hf ( )
inline

return the underlying HF reference

References hf.

◆ get_residue()

real_function_6d madness::MP2::get_residue ( const real_function_6d f,
const int  i,
const int  j 
)
inline

◆ guess_mp1_3()

void madness::MP2::guess_mp1_3 ( ElectronPair pair) const
private

compute the first iteration of the residual equations and all intermediates

Referenced by solve_residual_equations().

◆ help()

static void madness::MP2::help ( )
inlinestatic

◆ increment()

void madness::MP2::increment ( ElectronPair pair,
real_convolution_6d green 
)

◆ iterate()

real_function_6d madness::MP2::iterate ( const real_function_6d f) const
inline

solve the residual equation for electron pair (i,j)

Todo:
Parameter documentation. Below are un-doxygenated comments that no longer seem relevant?

References compute_gQf(), madness::ElectronPair::constant_term, madness::copy(), madness::f, madness::ElectronPair::function, madness::ElectronPair::ij_gQf_ij, madness::ElectronPair::ji_gQf_ij, and solve_residual_equations().

◆ J()

real_function_3d madness::MP2::J ( const real_function_3d phi) const
private

apply the Coulomb operator a on orbital

Parameters
[in]phithe orbital
Returns
Jphi

◆ JK1phi0_on_demand()

real_function_6d madness::MP2::JK1phi0_on_demand ( const int  i,
const int  j,
const bool  hc = false 
) const
private

return the function (J(1)-K(1)) |phi0> as on-demand function

Parameters
[in]hccompute hermitian conjugate -> swap bra and ket space

◆ JK2phi0_on_demand()

real_function_6d madness::MP2::JK2phi0_on_demand ( const int  i,
const int  j,
const bool  hc = false 
) const
private

return the function (J(2)-K(2)) |phi0> as on-demand function

◆ K() [1/2]

real_function_3d madness::MP2::K ( const real_function_3d phi,
const bool  hc = false 
) const
private

apply the exchange operator on an orbital

Parameters
[in]phithe orbital
[in]hchermitian conjugate -> swap bra and ket space
Returns
Kphi

◆ K() [2/2]

real_function_6d madness::MP2::K ( const real_function_6d phi,
const bool  is_symmetric = false 
) const
private

apply the exchange operator on a pair function

Parameters
[in]phithe pair function
[in]is_symmetricis the function symmetric wrt particle exchange
Returns
(K1 + K2) |phi >

◆ load_function()

template<typename T , size_t NDIM>
void madness::MP2::load_function ( Function< T, NDIM > &  f,
const std::string  name 
) const

load a function

References madness::f, name(), madness::print(), madness::World::rank(), and world.

Referenced by increment().

◆ make_chi()

std::vector<real_function_3d> madness::MP2::make_chi ( const real_function_3d phi,
const real_convolution_3d op,
const bool  hc = false 
) const
private

make the quantity chi_k

chi is the Poisson kernel applied on an orbital product of the input function and the vector of orbitals

\[ \chi_{k{*} i}(1) = \int dr_2 \frac{k(2) i(2)}{|r_1-r_2|} \]

\[ \chi_{ki{*}}(1) = \int dr_2 \frac{k(2) i(2)}{|r_1-r_2|} \]

if hc

Parameters
[in]phiorbital phi_i
[in]opthe operator in SeparatedConvolution form
[in]hccompute hermitian conjugate -> pass the correct phi!
Returns
a vector of length nocc

◆ make_fKphi0()

real_function_6d madness::MP2::make_fKphi0 ( const int  i,
const int  j 
) const
private

apply the operator K on the reference and multiply with f; fK |phi^0>

Parameters
[in]iindex of orbital i
[in]jindex of orbital j
Returns
the function f12 (K(1)+K(2))|phi^0>

◆ make_KffKphi0()

real_function_6d madness::MP2::make_KffKphi0 ( const ElectronPair pair) const
private

return the function [K,f] phi0; load from disk if available

◆ make_pair()

ElectronPair madness::MP2::make_pair ( const int  i,
const int  j 
) const
private

compute some matrix elements that don't change during the SCF

Referenced by value().

◆ make_Uphi0()

real_function_6d madness::MP2::make_Uphi0 ( ElectronPair pair) const
private

◆ make_xi()

std::vector<real_function_3d> madness::MP2::make_xi ( const real_function_3d phi_i,
const real_function_3d phi_j,
const real_convolution_3d op,
const bool  hc = false 
) const
private

make the quantity xi_k

xi is chi multiplied with an orbital j

\[ \xi_{k{*}i,j}(1) = \chi_{ki}(1) j(1) \]

\[ \xi_{ki{*},j{*}}(1) = \chi_{k{*}i}(1) j(1) \]

if hc

Parameters
[in]phi_iorbital i
[in]phi_jorbital j
[in]opthe operator in SeparatedConvolution form
[in]hccompute hermitian conjugate -> pass the correct phi!
Returns
a vector of length k=0,..,nocc

Referenced by compute_gQf().

◆ multiply_with_0th_order_Hamiltonian()

real_function_6d madness::MP2::multiply_with_0th_order_Hamiltonian ( const real_function_6d f,
const int  i,
const int  j 
) const
private

multiply the given function with the 0th order Hamiltonian, exluding the 0th order energy

Parameters
[in]fthe function we apply H^0 on
Returns
the function g=H^0 f, which is NOT orthogonalized against f

Referenced by debug_cc2(), get_residue(), increment(), solve_coupled_equations(), and solve_residual_equations().

◆ name()

std::string madness::MP2::name ( ) const
inlinevirtual

◆ nemo0_on_demand()

real_function_6d madness::MP2::nemo0_on_demand ( const int  i,
const int  j 
) const
private

return the function |F1F2> as on-demand function

◆ phi0_on_demand()

real_function_6d madness::MP2::phi0_on_demand ( const int  i,
const int  j 
) const
private

return the function |phi0> as on-demand function

◆ print_options()

template<typename T >
void madness::MP2::print_options ( const std::string  option,
const T  val 
) const
inline

pretty print the options

invoke only in (world.rank()==0) !!

◆ print_parameters()

static void madness::MP2::print_parameters ( )
inlinestatic

◆ save_function()

template<typename T , size_t NDIM>
void madness::MP2::save_function ( const Function< T, NDIM > &  f,
const std::string  name 
) const

save a function

References madness::f, name(), madness::print(), madness::World::rank(), and world.

Referenced by solve_residual_equations().

◆ selftest()

virtual bool madness::MP2::selftest ( )
inlinevirtual

Implements QCPropertyInterface.

◆ solve_coupled_equations()

double madness::MP2::solve_coupled_equations ( Pairs< ElectronPair > &  pairs,
const double  econv,
const double  dconv 
) const

solve the coupled MP1 equations (e.g. for local orbitals)

solve the coupled MP1 equations for local orbitals

Parameters
[in]pairsset of (coupled) electron pairs to solve
[in]econvenergy convergence criterion (for all pairs)
[in]dconvdensity convergence criterion (for all pairs)
[in]pairssolved decoupled electron pairs

References std::abs(), add_local_coupling(), bsh_eps, compute_energy(), END_TIMER(), energy, madness::MP2::Parameters::freeze(), hf, lo, madness::MP2::Parameters::maxiter(), multiply_with_0th_order_Hamiltonian(), pairs, param, Q12, madness::World::rank(), residual(), madness::Function< T, NDIM >::scale(), START_TIMER(), madness::stringify(), madness::Function< T, NDIM >::truncate(), madness::truncate(), madness::wall_time(), world, and zeroth_order_energy().

Referenced by value().

◆ solve_residual_equations()

void madness::MP2::solve_residual_equations ( ElectronPair pair,
const double  econv,
const double  dconv 
) const

◆ START_TIMER()

void madness::MP2::START_TIMER ( World world) const
inlineprivate

◆ value() [1/2]

double madness::MP2::value ( )

return the molecular correlation energy energy (without the HF energy)

References check_core_valence_separation(), enforce_core_valence_separation(), and hf.

◆ value() [2/2]

double madness::MP2::value ( const Tensor< double > &  x)
virtual

◆ zeroth_order_energy()

double madness::MP2::zeroth_order_energy ( const int  i,
const int  j 
) const
inline

return the 0th order energy of pair ij (= sum of orbital energies)

References hf.

Referenced by make_Uphi0(), solve_coupled_equations(), and solve_residual_equations().

Member Data Documentation

◆ coords_sum

double madness::MP2::coords_sum
private

check sum for the geometry

Referenced by coord_chksum(), and value().

◆ correlation_energy

double madness::MP2::correlation_energy
private

the correlation energy

Referenced by value().

◆ corrfac

CorrelationFactor madness::MP2::corrfac
private

correlation factor: Slater

Referenced by MP2(), compute_gQf(), and make_Uphi0().

◆ fock

Tensor<double> madness::MP2::fock
mutableprivate

◆ hf

std::shared_ptr<HartreeFock> madness::MP2::hf
private

◆ nuclear_corrfac

std::shared_ptr<NuclearCorrelationFactor> madness::MP2::nuclear_corrfac
private

Referenced by get_residue(), and value().

◆ pairs

Pairs<ElectronPair> madness::MP2::pairs
private

pair functions and energies

Referenced by MP2(), solve_coupled_equations(), and value().

◆ param

Parameters madness::MP2::param
private

◆ poisson

std::shared_ptr<real_convolution_3d> madness::MP2::poisson
private

Referenced by MP2(), and compute_gQf().

◆ Q12

StrongOrthogonalityProjector<double, 3> madness::MP2::Q12
private

◆ sss

double madness::MP2::sss
private

Referenced by END_TIMER(), and START_TIMER().

◆ ttt

double madness::MP2::ttt
mutableprivate

Referenced by END_TIMER(), and START_TIMER().

◆ world

World& madness::MP2::world
private

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