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

A nuclear correlation factor class. More...

#include <correlationfactor.h>

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

Public Member Functions

 Slater (World &world, const Molecule &mol, const double a)
 ctor More...
 
corrfactype type () const
 
- Public Member Functions inherited from madness::NuclearCorrelationFactor
 NuclearCorrelationFactor (World &world, const Molecule &mol)
 ctor More...
 
virtual ~NuclearCorrelationFactor ()
 virtual destructor More...
 
virtual real_function_3d apply_U (const real_function_3d &rhs) const
 apply the regularized potential U_nuc on a given function rhs More...
 
coord_3d dsmoothed_unitvec (const coord_3d &xyz, const int axis, double smoothing=0.0) const
 derivative of smoothed unit vector wrt the electronic coordinate More...
 
virtual real_function_3d function () const
 return the nuclear correlation factor More...
 
void initialize (const double vtol1)
 initialize the regularized potentials U1 and U2 More...
 
virtual real_function_3d inverse () const
 return the inverse nuclear correlation factor More...
 
coord_3d smoothed_unitvec (const coord_3d &xyz, double smoothing=0.0) const
 smoothed unit vector for the computation of the U1 potential More...
 
virtual real_function_3d square () const
 return the square of the nuclear correlation factor More...
 
virtual real_function_3d square_times_V_derivative (const int iatom, const int axis) const
 
virtual const real_function_3d U1 (const int axis) const
 return the U1 term of the correlation function More...
 
std::vector< real_function_3dU1vec () const
 return the U1 functions in a vector More...
 
virtual const real_function_3d U2 () const
 return the U2 term of the correlation function More...
 

Private Member Functions

double a_param () const
 
double eprec_param () const
 
double S (const double &r, const double &Z) const
 the nuclear correlation factor More...
 
coord_3d Sp (const coord_3d &vr1A, const double &Z) const
 radial part first derivative of the nuclear correlation factor More...
 
double Spp_div_S (const double &r, const double &Z) const
 second derivative of the nuclear correlation factor More...
 
double Sr_div_S (const double &r, const double &Z) const
 first derivative of the correlation factor wrt (r-R_A) More...
 
double Srr_div_S (const double &r, const double &Z) const
 second derivative of the correlation factor wrt (r-R_A) More...
 
double Srrr_div_S (const double &r, const double &Z) const
 third derivative of the correlation factor wrt (r-R_A) More...
 
double U2X_spherical (const double &r, const double &Z, const double &rcut) const
 derivative of the U2 potential wrt X (scalar part) More...
 

Private Attributes

double a_
 the length scale parameter More...
 
double eprec_
 

Additional Inherited Members

- Public Types inherited from madness::NuclearCorrelationFactor
enum  corrfactype {
  None , GradientalGaussSlater , GaussSlater , LinearSlater ,
  Polynomial , Slater , poly4erfc , Two ,
  Adhoc
}
 
typedef std::shared_ptr< FunctionFunctorInterface< double, 3 > > functorT
 
- Protected Attributes inherited from madness::NuclearCorrelationFactor
std::vector< real_function_3dU1_function
 the three components of the U1 potential More...
 
real_function_3d U2_function
 the purely local U2 potential, having absorbed the nuclear pot V_nuc More...
 

Detailed Description

A nuclear correlation factor class.

Constructor & Destructor Documentation

◆ Slater()

madness::Slater::Slater ( World world,
const Molecule mol,
const double  a 
)
inline

ctor

Parameters
[in]worldthe world
[in]molmolecule with the sites of the nuclei

References a, madness::Molecule::get_eprec(), madness::print(), madness::World::rank(), and madness::NuclearCorrelationFactor::world.

Member Function Documentation

◆ a_param()

double madness::Slater::a_param ( ) const
inlineprivate

◆ eprec_param()

double madness::Slater::eprec_param ( ) const
inlineprivate

◆ S()

double madness::Slater::S ( const double &  r,
const double &  Z 
) const
inlineprivatevirtual

the nuclear correlation factor

Implements madness::NuclearCorrelationFactor.

References a, and Z.

◆ Sp()

coord_3d madness::Slater::Sp ( const coord_3d vr1A,
const double &  Z 
) const
inlineprivatevirtual

radial part first derivative of the nuclear correlation factor

Implements madness::NuclearCorrelationFactor.

References a, madness::Vector< T, N >::normf(), madness::NuclearCorrelationFactor::smoothed_unitvec(), and Z.

◆ Spp_div_S()

double madness::Slater::Spp_div_S ( const double &  r,
const double &  Z 
) const
inlineprivatevirtual

second derivative of the nuclear correlation factor

Implements madness::NuclearCorrelationFactor.

References a, and Z.

◆ Sr_div_S()

double madness::Slater::Sr_div_S ( const double &  r,
const double &  Z 
) const
inlineprivatevirtual

first derivative of the correlation factor wrt (r-R_A)

\[ Sr_div_S = \frac{1}{S(r)}\frac{\partial S(r)}{\partial r} \]

Implements madness::NuclearCorrelationFactor.

References a, and Z.

◆ Srr_div_S()

double madness::Slater::Srr_div_S ( const double &  r,
const double &  Z 
) const
inlineprivatevirtual

second derivative of the correlation factor wrt (r-R_A)

\[ result = \frac{1}{S(r)}\frac{\partial^2 S(r)}{\partial r^2} \]

Implements madness::NuclearCorrelationFactor.

References a, and Z.

◆ Srrr_div_S()

double madness::Slater::Srrr_div_S ( const double &  r,
const double &  Z 
) const
inlineprivatevirtual

third derivative of the correlation factor wrt (r-R_A)

\[ result = \frac{1}{S(r)}\frac{\partial^3 S(r)}{\partial r^3} \]

Implements madness::NuclearCorrelationFactor.

References a, and Z.

◆ type()

corrfactype madness::Slater::type ( ) const
inlinevirtual

◆ U2X_spherical()

double madness::Slater::U2X_spherical ( const double &  r,
const double &  Z,
const double &  rcut 
) const
inlineprivatevirtual

derivative of the U2 potential wrt X (scalar part)

with

\[ \rho = \left| \vec r- \vec R_A \right| \]

returns the term in the parenthesis without the the derivative of rho

\[ \frac{\partial U_2}{\partial X_A} = \frac{\partial \rho}{\partial X} \left(-\frac{1}{2}\frac{S''' S - S'' S'}{S^2} + \frac{1}{\rho^2}\frac{S'}{S} - \frac{1}{\rho} \frac{S''S - S'^2}{S^2} + \frac{Z_A}{\rho^2}\right) \]

Reimplemented from madness::NuclearCorrelationFactor.

References a, a2, madness::r2(), madness::NuclearCorrelationFactor::Sr_div_S(), madness::NuclearCorrelationFactor::Srr_div_S(), madness::NuclearCorrelationFactor::Srrr_div_S(), and Z.

Member Data Documentation

◆ a_

double madness::Slater::a_
private

the length scale parameter

◆ eprec_

double madness::Slater::eprec_
private

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