MADNESS  0.10.1
Classes | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
madness::XCOperator< T, NDIM > Class Template Reference

operator class for the handling of DFT exchange-correlation functionals More...

#include <SCFOperators.h>

Inheritance diagram for madness::XCOperator< T, NDIM >:
Inheritance graph
[legend]
Collaboration diagram for madness::XCOperator< T, NDIM >:
Collaboration graph
[legend]

Classes

struct  expme
 simple structure to take the pointwise exponential of a function, shifted by +14 More...
 
struct  logme
 simple structure to take the pointwise logarithm of a function, shifted by +14 More...
 

Public Member Functions

 XCOperator (World &world)
 default ctor without information about the XC functional More...
 
 XCOperator (World &world, const Nemo *nemo, int ispin=0)
 ctor with a Nemo calculation, will initialize the necessary intermediates More...
 
 XCOperator (World &world, const Nemo *scf, const real_function_3d &arho, const real_function_3d &brho, int ispin=0)
 ctor with an Nemo calculation, will initialize the necessary intermediates More...
 
 XCOperator (World &world, const SCF *scf, const real_function_3d &arho, const real_function_3d &brho, int ispin=0, std::string deriv="abgv")
 ctor with an SCF calculation, will initialize the necessary intermediates More...
 
 XCOperator (World &world, const SCF *scf, int ispin=0, std::string deriv="abgv")
 ctor with an SCF calculation, will initialize the necessary intermediates More...
 
 XCOperator (World &world, std::string xc_data, const bool spin_polarized, const real_function_3d &arho, const real_function_3d &brho, std::string deriv="abgv")
 custom ctor with information about the XC functional More...
 
real_function_3d apply_xc_kernel (const real_function_3d &density, const vecfuncT grad_dens_pt=vecfuncT()) const
 construct the xc kernel and apply it directly on the (response) density More...
 
double compute_xc_energy () const
 compute the xc energy using the precomputed intermediates vf and delrho More...
 
std::string info () const
 print some information about this operator More...
 
real_function_3d make_xc_potential () const
 return the local xc potential More...
 
T operator() (const Function< T, NDIM > &bra, const Function< T, NDIM > &ket) const
 compute the matrix element <bra | op | ket> More...
 
Function< T, NDIMoperator() (const Function< T, NDIM > &ket) const
 apply the xc potential on an orbitals More...
 
std::vector< Function< T, NDIM > > operator() (const std::vector< Function< T, NDIM > > &vket) const
 apply the xc potential on a set of orbitals More...
 
Tensor< Toperator() (const std::vector< Function< T, NDIM >> &vbra, const std::vector< Function< T, NDIM >> &vket) const
 
XCOperatorset_extra_truncation (const double &fac)
 
void set_ispin (const int i) const
 set the spin state this operator is acting on More...
 
- Public Member Functions inherited from madness::SCFOperatorBase< T, NDIM >
 SCFOperatorBase ()=default
 
 SCFOperatorBase (std::shared_ptr< MacroTaskQ > taskq)
 
virtual ~SCFOperatorBase ()
 
virtual tensorT operator() (const vecfuncT &vbra, const vecfuncT &vket) const =0
 compute the matrix <vbra | op | vket> More...
 
virtual vecfuncT operator() (const vecfuncT &vket) const =0
 apply this operator on the argument vector of functions More...
 

Public Attributes

std::shared_ptr< XCfunctionalxc
 interface to the actual XC functionals More...
 
- Public Attributes inherited from madness::SCFOperatorBase< T, NDIM >
std::shared_ptr< MacroTaskQtaskq =0
 

Private Member Functions

bool is_initialized () const
 check if the intermediates are initialized More...
 
vecfuncT prep_xc_args (const real_function_3d &arho, const real_function_3d &brho) const
 compute the intermediates for the XC functionals More...
 
void prep_xc_args_response (const real_function_3d &dens_pt, vecfuncT &xc_args, vecfuncT &ddens_pt) const
 compute the intermediates for the XC functionals More...
 

Private Attributes

std::string dft_deriv
 which derivative operator to use More...
 
double extra_truncation
 additional truncation for the densities in the XC kernel More...
 
int ispin
 the XC functionals depend on the spin of the orbitals they act on More...
 
int nbeta
 number of beta orbitals More...
 
std::shared_ptr< NuclearCorrelationFactorncf
 the nuclear correlation factor, if it exists, for computing derivatives for GGA More...
 
Worldworld
 the world More...
 
vecfuncT xc_args
 functions that are need for the computation of the XC operator More...
 

Additional Inherited Members

- Public Types inherited from madness::SCFOperatorBase< T, NDIM >
typedef Function< T, NDIMfunctionT
 
typedef Tensor< TtensorT
 
typedef std::vector< functionTvecfuncT
 

Detailed Description

template<typename T, std::size_t NDIM>
class madness::XCOperator< T, NDIM >

operator class for the handling of DFT exchange-correlation functionals

Constructor & Destructor Documentation

◆ XCOperator() [1/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world)
inline

default ctor without information about the XC functional

◆ XCOperator() [2/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world,
std::string  xc_data,
const bool  spin_polarized,
const real_function_3d arho,
const real_function_3d brho,
std::string  deriv = "abgv" 
)

◆ XCOperator() [3/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world,
const SCF scf,
int  ispin = 0,
std::string  deriv = "abgv" 
)

◆ XCOperator() [4/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world,
const Nemo nemo,
int  ispin = 0 
)

ctor with a Nemo calculation, will initialize the necessary intermediates

◆ XCOperator() [5/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world,
const SCF scf,
const real_function_3d arho,
const real_function_3d brho,
int  ispin = 0,
std::string  deriv = "abgv" 
)

◆ XCOperator() [6/6]

template<typename T , std::size_t NDIM>
madness::XCOperator< T, NDIM >::XCOperator ( World world,
const Nemo scf,
const real_function_3d arho,
const real_function_3d brho,
int  ispin = 0 
)

ctor with an Nemo calculation, will initialize the necessary intermediates

Member Function Documentation

◆ apply_xc_kernel()

template<typename T , std::size_t NDIM>
real_function_3d madness::XCOperator< T, NDIM >::apply_xc_kernel ( const real_function_3d dens_pt,
const vecfuncT  grad_dens_pt = vecfuncT() 
) const

construct the xc kernel and apply it directly on the (response) density

apply the xc kernel on a perturbed density

the xc kernel is the second derivative of the xc functions wrt the density

Parameters
[in]densitythe (response) density on which the kernel is applied
Returns
kernel * density

cf Eq. (13) of T. Yanai, R. J. Harrison, and N. Handy, “Multiresolution quantum chemistry in multiwavelet bases: time-dependent density functional theory with asymptotically corrected potentials in local density and generalized gradient approximations,” Mol. Phys., vol. 103, no. 2, pp. 413–424, 2005.

the application of the xc kernel is (RHF only)

\[ \frac{\partial^2E_{xc}}{\partial \rho_\alpha^2}\circ\tilde\rho = second_{local} + second_{semilocal} + first_{semilocal} \]

where the second partial derivatives are

\[ second_{local} = \frac{\partial^2 f_{xc}}{\partial \rho_\alpha^2}\tilde \rho + 2\frac{\partial^2 f_{xc}}{\partial \rho_\alpha\sigma_{\alpha\alpha}} \left(\vec\nabla \rho_a\cdot \vec \nabla\tilde\rho\right) \]

the second partial derivatives that need to be multiplied with the density gradients

\[ second_{semilocal} = -\vec\nabla\cdot\left((\vec\nabla\rho) \left[2\frac{\partial^2 f_{xc}}{\partial\rho_\alpha\partial\sigma_{\alpha\alpha}}\tilde\rho + 4\frac{\partial^2 f_{xc}}{\partial\sigma_{\alpha\alpha}^2} \left(\vec\nabla\rho_\alpha\cdot\vec\nabla\tilde\rho\right)\right]\right) \]

and the first derivatives that need to be multiplied with the density gradients

\[ first_{semilocal} = -\vec\nabla\cdot\left(2\frac{\partial f_{xc}}{\partial\sigma_{\alpha\alpha}}\vec\nabla\tilde\rho\right) \]

References div(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), MADNESS_ASSERT, MADNESS_EXCEPTION, madness::multi_to_multi_op_values(), op(), madness::refine_to_common_level(), madness::Function< T, NDIM >::truncate(), truncate(), and xc.

◆ compute_xc_energy()

template<typename T , std::size_t NDIM>
double madness::XCOperator< T, NDIM >::compute_xc_energy

◆ info()

template<typename T , std::size_t NDIM>
std::string madness::XCOperator< T, NDIM >::info ( ) const
inlinevirtual

print some information about this operator

Implements madness::SCFOperatorBase< T, NDIM >.

◆ is_initialized()

template<typename T , std::size_t NDIM>
bool madness::XCOperator< T, NDIM >::is_initialized ( ) const
inlineprivate

check if the intermediates are initialized

References madness::XCOperator< T, NDIM >::xc_args.

◆ make_xc_potential()

template<typename T , std::size_t NDIM>
real_function_3d madness::XCOperator< T, NDIM >::make_xc_potential

◆ operator()() [1/4]

template<typename T , std::size_t NDIM>
T madness::XCOperator< T, NDIM >::operator() ( const Function< T, NDIM > &  bra,
const Function< T, NDIM > &  ket 
) const
inlinevirtual

compute the matrix element <bra | op | ket>

Parameters
brabra state
ketket state
Returns
the matrix element <bra | op | ket>

Implements madness::SCFOperatorBase< T, NDIM >.

References MADNESS_EXCEPTION.

◆ operator()() [2/4]

template<typename T , std::size_t NDIM>
Function<T,NDIM> madness::XCOperator< T, NDIM >::operator() ( const Function< T, NDIM > &  ket) const
inlinevirtual

apply the xc potential on an orbitals

Implements madness::SCFOperatorBase< T, NDIM >.

References madness::XCOperator< T, NDIM >::operator()().

◆ operator()() [3/4]

template<typename T , std::size_t NDIM>
std::vector< Function< T, NDIM > > madness::XCOperator< T, NDIM >::operator() ( const std::vector< Function< T, NDIM > > &  vket) const

apply the xc potential on a set of orbitals

Referenced by madness::XCOperator< T, NDIM >::operator()().

◆ operator()() [4/4]

template<typename T , std::size_t NDIM>
Tensor<T> madness::XCOperator< T, NDIM >::operator() ( const std::vector< Function< T, NDIM >> &  vbra,
const std::vector< Function< T, NDIM >> &  vket 
) const
inline

References MADNESS_EXCEPTION.

◆ prep_xc_args()

template<typename T , std::size_t NDIM>
vecfuncT madness::XCOperator< T, NDIM >::prep_xc_args ( const real_function_3d arho,
const real_function_3d brho 
) const
private

compute the intermediates for the XC functionals

prepare xc args

Parameters
[in]arhodensity of the alpha orbitals
[in]brhodensity of the beta orbitals (necessary only if spin-polarized)
Returns
xc_args vector of intermediates as described above

References copy(), madness::dot(), madness::WorldGopInterface::fence(), madness::World::gop, madness::grad(), madness::grad_ble_one(), madness::grad_bspline_one(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), madness::Function< T, NDIM >::reconstruct(), truncate(), madness::unary_op(), madness::Function< T, NDIM >::world(), and xc.

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

◆ prep_xc_args_response()

template<typename T , std::size_t NDIM>
void madness::XCOperator< T, NDIM >::prep_xc_args_response ( const real_function_3d dens_pt,
vecfuncT xc_args,
vecfuncT ddens_pt 
) const
private

compute the intermediates for the XC functionals

add intermediates for the response kernels to xc_args

Parameters
[in]dens_ptperturbed densities from CPHF or TDDFT equations
[in,out]xc_argsvector of intermediates as described above
[out]ddens_ptxyz-derivatives of dens_pt

References madness::dot(), madness::WorldGopInterface::fence(), madness::World::gop, madness::grad(), madness::XCfunctional::is_gga(), madness::XCfunctional::is_spin_polarized(), print(), truncate(), madness::Function< T, NDIM >::world(), xc, and zeta.

◆ set_extra_truncation()

template<typename T , std::size_t NDIM>
XCOperator& madness::XCOperator< T, NDIM >::set_extra_truncation ( const double &  fac)
inline

◆ set_ispin()

template<typename T , std::size_t NDIM>
void madness::XCOperator< T, NDIM >::set_ispin ( const int  i) const
inline

set the spin state this operator is acting on

References madness::XCOperator< T, NDIM >::ispin.

Referenced by test_XCOperator().

Member Data Documentation

◆ dft_deriv

template<typename T , std::size_t NDIM>
std::string madness::XCOperator< T, NDIM >::dft_deriv
private

which derivative operator to use

◆ extra_truncation

template<typename T , std::size_t NDIM>
double madness::XCOperator< T, NDIM >::extra_truncation
private

additional truncation for the densities in the XC kernel

the densities in the DFT kernal are processed as their inverses, so noise in the small density regions might amplify and lead to inaccurate results. Extra truncation will tighten the truncation threshold by a specified factor, default is 0.01.

Referenced by madness::XCOperator< T, NDIM >::set_extra_truncation().

◆ ispin

template<typename T , std::size_t NDIM>
int madness::XCOperator< T, NDIM >::ispin
mutableprivate

the XC functionals depend on the spin of the orbitals they act on

Referenced by madness::XCOperator< T, NDIM >::set_ispin().

◆ nbeta

template<typename T , std::size_t NDIM>
int madness::XCOperator< T, NDIM >::nbeta
private

number of beta orbitals

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

◆ ncf

template<typename T , std::size_t NDIM>
std::shared_ptr<NuclearCorrelationFactor> madness::XCOperator< T, NDIM >::ncf
private

the nuclear correlation factor, if it exists, for computing derivatives for GGA

◆ world

template<typename T , std::size_t NDIM>
World& madness::XCOperator< T, NDIM >::world
private

◆ xc

template<typename T , std::size_t NDIM>
std::shared_ptr<XCfunctional> madness::XCOperator< T, NDIM >::xc

interface to the actual XC functionals

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

◆ xc_args

template<typename T , std::size_t NDIM>
vecfuncT madness::XCOperator< T, NDIM >::xc_args
mutableprivate

functions that are need for the computation of the XC operator

the ordering of the intermediates is fixed, but the code can handle non-initialized functions, so if e.g. no GGA is requested, all the corresponding vector components may be left empty. For the ordering of the intermediates see xcfunctional::xc_arg

Referenced by madness::XCOperator< T, NDIM >::XCOperator(), and madness::XCOperator< T, NDIM >::is_initialized().


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