ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
NumericalRadial Class Reference

A class that represents a numerical radial function. More...

#include <numerical_radial.h>

Collaboration diagram for NumericalRadial:

Public Member Functions

 NumericalRadial ()=default
 
 NumericalRadial (NumericalRadial const &)
 Deep-copy grid & values.
 
NumericalRadialoperator= (NumericalRadial const &)
 Deep-copy grid & values.
 
 ~NumericalRadial ()
 
void build (const int l, const bool for_r_space, const int ngrid, const double *const grid, const double *const value, const int p=0, const int izeta=0, const std::string symbol="", const int itype=0, const bool init_sbt=true)
 Initializes the object by providing the grid & values in one space.
 
void to_numerical_orbital_lm (Numerical_Orbital_Lm &orbital_lm, const int nk_legacy=4005, const double lcao_dk=0.01) const
 Overwrites the content of a Numerical_Orbital_Lm object with the current object.
 
void set_transformer (ModuleBase::SphericalBesselTransformer sbt, int update=0)
 Sets a SphericalBesselTransformer.
 
void set_grid (const bool for_r_space, const int ngrid, const double *const grid, const char mode='i')
 Sets up a grid.
 
void set_uniform_grid (const bool for_r_space, const int ngrid, const double cutoff, const char mode='i', const bool enable_fft=false)
 Sets up a uniform grid.
 
void set_value (const bool for_r_space, const double *const value, const int p)
 Updates values on an existing grid.
 
void wipe (const bool r_space=true, const bool k_space=true)
 Removes the grid & values in r or k space.
 
void radtab (const char op, const NumericalRadial &ket, const int l, double *const table, const int nr_tab, const double rmax_tab, const bool deriv=false) const
 Computes the radial table for two-center integrals.
 
void normalize (bool for_r_space=true)
 Normalizes the radial function.
 
Getters
std::string const & symbol () const
 padded zeros ignored
 
int itype () const
 padded zeros ignored
 
int izeta () const
 padded zeros ignored
 
int l () const
 padded zeros ignored
 
int nr () const
 padded zeros ignored
 
int nk () const
 padded zeros ignored
 
double rcut () const
 padded zeros ignored
 
double kcut () const
 padded zeros ignored
 
double rmax () const
 padded zeros considered
 
double kmax () const
 padded zeros ignored
 
const double * rgrid () const
 padded zeros ignored
 
const double * kgrid () const
 padded zeros ignored
 
const double * rvalue () const
 padded zeros ignored
 
const double * kvalue () const
 padded zeros ignored
 
double pr () const
 padded zeros ignored
 
double pk () const
 padded zeros ignored
 
bool is_fft_compliant () const
 padded zeros ignored
 
ModuleBase::SphericalBesselTransformer sbt () const
 padded zeros ignored
 
double rgrid (int ir) const
 padded zeros ignored
 
double kgrid (int ik) const
 padded zeros ignored
 
double rvalue (int ir) const
 padded zeros ignored
 
double kvalue (int ik) const
 padded zeros ignored
 

Private Member Functions

void transform (const bool forward)
 Transforms the r-space values to get k-space values, or vice versa.
 
void set_icut (const bool for_r_space, const bool for_k_space, const double tol=1e-15)
 Updates ircut_ and/or ikcut_.
 

Static Private Member Functions

static bool is_uniform (const int n, const double *const grid, const double tol=1e-15)
 Checks whether a grid is uniform.
 
static bool is_fft_compliant (const int nr, const double *const rgrid, const int nk, const double *const kgrid, const double tol=1e-15)
 Checks whether the given two grids are FFT-compliant.
 

Private Attributes

std::string symbol_ = ""
 chemical symbol
 
int itype_ = 0
 element index in calculation
 
int l_ = -1
 angular momentum
 
int izeta_ = 0
 further index for NumericalRadial objects with the same itype_and l_
 
int nr_ = 0
 number of r-space grid points
 
int nk_ = 0
 number of k-space grid points
 
double * rgrid_ = nullptr
 r-space grid
 
double * kgrid_ = nullptr
 k-space grid
 
int ircut_ = 0
 Index of the first trailing zero.
 
int ikcut_ = 0
 
double * rvalue_ = nullptr
 r-space value
 
double * kvalue_ = nullptr
 k-space value
 
bool is_fft_compliant_ = false
 A flag that tells whether the r & k grids are FFT-compliant.
 
ModuleBase::SphericalBesselTransformer sbt_
 An object that provides spherical Bessel transforms.
 
Implicit exponents in values

Sometimes a radial function is given in the form of pow(r,p) * F(r) rather than F(r) (same applies to k). For example, the Kleinman-Bylander beta functions are often given as r*beta(r) instead of bare beta(r). Very often using r*beta(r) is adequate; there's no need to get bare beta(r) at all.

This class takes care of this situation. When building the object, one can optionally provide an exponent p so that the input values are interpreted as pow(x[i],p) * F(x[i]), where F(x) is what the object actually represents. pr_ & pk_ keep track of these exponents within r & k values. They are taken taken account automatically during spherical Bessel transforms.

int pr_ = 0
 implicit exponent in rvalues_
 
int pk_ = 0
 implicit exponent in kvalues_
 

Detailed Description

A class that represents a numerical radial function.

This class is designed to be the container for the radial part of numerical atomic orbitals, Kleinman-Bylander beta functions, and all other similar numerical radial functions in three-dimensional space, each of which is associated with some angular momentum l and whose r and k space values are related by an l-th order spherical Bessel transform.

A NumericalRadial object can be initialized by "build", which requires the angular momentum, the number of grid points, the grid and the corresponding values. Grid does not have to be uniform. One can initialize the object in either r or k space. After initialization, one can set the grid in the other space via set_grid or set_uniform_grid. Values in the other space are automatically computed by a spherical Bessel transform.

Usage:

int l = 1;
int itype = 3;
int izeta = 5;
std::string symbol = "Au";

// Prepares the grid & values to initialize the objects
int sz = 2000;
double dr = 0.01;
double* grid = new double[sz];
for (int ir = 0; ir != sz; ++ir) {
    grid[ir] = ir * dr; 
    f[ir] = std::exp(-grid[ir] * grid[ir]);
}
// grid does not necessarily have to be uniform; it just
// has to be positive and strictly increasing.

// The class will interpret the input values as r^p * F(r)
// where F is the underlying radial function that the class object
// actually represents.
int p1 = 0;
int p2 = -2;

NumericalRadial chi1, chi2;
chi1.build(0, true, sz, grid, f, p1);
chi2.build(2, true, sz, grid, f, p2);

// Now chi1 represents exp(-r^2); chi2 actually represents
// r^2*exp(-r^2), even though the values stored is also exp(-r^2).

// Adds the k-space grid.
chi1.set_uniform_grid(false, sz, PI/dr, 't');
chi2.set_uniform_grid(false, sz, PI/dr, 't');
// k-space values are automatically computed above

Constructor & Destructor Documentation

◆ NumericalRadial() [1/2]

NumericalRadial::NumericalRadial ( )
default

◆ NumericalRadial() [2/2]

NumericalRadial::NumericalRadial ( NumericalRadial const &  other)

Deep-copy grid & values.

Here is the call graph for this function:

◆ ~NumericalRadial()

NumericalRadial::~NumericalRadial ( )

Member Function Documentation

◆ build()

void NumericalRadial::build ( const int  l,
const bool  for_r_space,
const int  ngrid,
const double *const  grid,
const double *const  value,
const int  p = 0,
const int  izeta = 0,
const std::string  symbol = "",
const int  itype = 0,
const bool  init_sbt = true 
)

Initializes the object by providing the grid & values in one space.

Parameters
[in]lAngular momentum.
[in]for_r_spaceSpecifies whether the input corresponds to r or k space.
[in]ngridNumber of grid points.
[in]gridGrid points, must be positive & strictly increasing.
[in]valueValues on the grid.
[in]pImplicit exponent in input values, see pr_ & pk_.
[in]izetaMultiplicity index of radial functions of the same itype and l.
[in]symbolChemical symbol.
[in]itypeIndex for the element in calculation.
[in]init_sbtIf true, internal SphericalBesselTransformer will be initialized.
Note
init_sbt is only useful when the internal SphericalBesselTransformer (sbt_) is null-initialized; The function will NOT reset sbt_ if it is already usable.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ is_fft_compliant() [1/2]

bool NumericalRadial::is_fft_compliant ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ is_fft_compliant() [2/2]

bool NumericalRadial::is_fft_compliant ( const int  nr,
const double *const  rgrid,
const int  nk,
const double *const  kgrid,
const double  tol = 1e-15 
)
staticprivate

Checks whether the given two grids are FFT-compliant.

See also
is_fft_compliant_
Here is the call graph for this function:

◆ is_uniform()

bool NumericalRadial::is_uniform ( const int  n,
const double *const  grid,
const double  tol = 1e-15 
)
staticprivate

Checks whether a grid is uniform.

Here is the caller graph for this function:

◆ itype()

int NumericalRadial::itype ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ izeta()

int NumericalRadial::izeta ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ kcut()

double NumericalRadial::kcut ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ kgrid() [1/2]

const double * NumericalRadial::kgrid ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ kgrid() [2/2]

double NumericalRadial::kgrid ( int  ik) const
inline

padded zeros ignored

◆ kmax()

double NumericalRadial::kmax ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ kvalue() [1/2]

const double * NumericalRadial::kvalue ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ kvalue() [2/2]

double NumericalRadial::kvalue ( int  ik) const
inline

padded zeros ignored

◆ l()

int NumericalRadial::l ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ nk()

int NumericalRadial::nk ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ normalize()

void NumericalRadial::normalize ( bool  for_r_space = true)

Normalizes the radial function.

The radial function is normalized such that

 / +inf     2
 |      dx x  f(x) = 1
 /  0

where x is r or k. The integral is evaluated with Simpson's rule. Values in the other space are updated automatically via a spherical Bessel transform.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nr()

int NumericalRadial::nr ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ operator=()

NumericalRadial & NumericalRadial::operator= ( NumericalRadial const &  rhs)

Deep-copy grid & values.

Here is the call graph for this function:

◆ pk()

double NumericalRadial::pk ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ pr()

double NumericalRadial::pr ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ radtab()

void NumericalRadial::radtab ( const char  op,
const NumericalRadial ket,
const int  l,
double *const  table,
const int  nr_tab,
const double  rmax_tab,
const bool  deriv = false 
) const

Computes the radial table for two-center integrals.

TODO This function shall be removed from this class in the future; its functionality should be moved to TwoCenterTable class.

This function requires that both "this" and "ket" have existing kgrid_, and the two kgrid_ be identical.

On finish, table is filled with values on tabgrid. If tabgrid is not given, the rgrid_ of "this" is assumed (it would be an error if rgrid_ does not exist in this case).

op could be:

  • 'S' or 'I': overlap integral
                      / +inf     2
          table[ir] = |      dk k  f(k) g(k) j (k*r[ir])
                      /  0                    l
    
  • 'T': kinetic integral table
                      / +inf     4
          table[ir] = |      dk k  f(k) g(k) j (k*r[ir])
                      /  0                    l
    
  • 'U': Coulomb integral table. This is slightly different from overlap or kinetic integral that in this case the two-center integral is a double integral:
         /        f(r) g(r'-R)
         | dr dr' ------------
         /          |r-r'|
    
    The corresponding table is
                 / +inf
     table[ir] = |      dk  f(k) g(k) j (k*r[ir])
                 /  0                  l
    
Parameters
[in]opoperator, could be:
  • 'S' or 'I': overlap
  • 'T': kinetic
  • 'U': Coulomb
[in]ketthe other NumericalRadial object with which the two-center integral is computed
[in]langular momentum of the table
[out]tableon finish, contain the computed table
[in]nr_tabsize of table grid
[in]rmax_tabcutoff radius of table grid
[in]derivif true, calculates the derivative of the table
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rcut()

double NumericalRadial::rcut ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ rgrid() [1/2]

const double * NumericalRadial::rgrid ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ rgrid() [2/2]

double NumericalRadial::rgrid ( int  ir) const
inline

padded zeros ignored

◆ rmax()

double NumericalRadial::rmax ( ) const
inline

padded zeros considered

Here is the caller graph for this function:

◆ rvalue() [1/2]

const double * NumericalRadial::rvalue ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ rvalue() [2/2]

double NumericalRadial::rvalue ( int  ir) const
inline

padded zeros ignored

◆ sbt()

ModuleBase::SphericalBesselTransformer NumericalRadial::sbt ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ set_grid()

void NumericalRadial::set_grid ( const bool  for_r_space,
const int  ngrid,
const double *const  grid,
const char  mode = 'i' 
)

Sets up a grid.

This function can be used to set up the grid which is absent in "build" (in which case values on the new grid are automatically computed by a spherical Bessel transform) or interpolate on an existing grid to a new grid.

Parameters
[in]for_r_spaceSpecifies whether to set grid for the r or k space.
[in]ngridNumber of grid points.
[in]gridGrid points, must be positive & strictly increasing.
[in]modeSpecifies how values are updated, could be 'i' or 't':
  • 'i': New values are obtained by interpolating and zero-padding the existing values from current space. With this option, it is an error if the designated space does not have a grid;
  • 't': New values are obtained via transform from the other space. With this option, it is an error if the other space does not have a grid.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_icut()

void NumericalRadial::set_icut ( const bool  for_r_space,
const bool  for_k_space,
const double  tol = 1e-15 
)
private

Updates ircut_ and/or ikcut_.

Here is the caller graph for this function:

◆ set_transformer()

void NumericalRadial::set_transformer ( ModuleBase::SphericalBesselTransformer  sbt,
int  update = 0 
)

Sets a SphericalBesselTransformer.

By default the class uses an internal SphericalBesselTransformer, but one can optionally use a shared one. This could be beneficial when there are a lot of NumericalRadial objects whose grids have the same size.

Parameters
[in]sbtAn external transformer.
[in]updateSpecifies whether and how values are recomputed with the new transformer. Accepted values are:
  • 0: does not recompute values;
  • 1: calls a forward transform;
  • -1: calls a backward transform.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_uniform_grid()

void NumericalRadial::set_uniform_grid ( const bool  for_r_space,
const int  ngrid,
const double  cutoff,
const char  mode = 'i',
const bool  enable_fft = false 
)

Sets up a uniform grid.

The functionality of this function is similar to set_grid, except that the new grid is a uniform grid specified by the cutoff and the number of grid points given by

               cutoff
 grid[i] = i * -------
               ngrid-1
See also
set_grid

If enable_fft is true, this function will not only set up the grid & values in the designated space, but also sets the grid in the other space such that the r & k grids are FFT-compliant (and updates values via a FFT-based spherical Bessel transform).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_value()

void NumericalRadial::set_value ( const bool  for_r_space,
const double *const  value,
const int  p 
)

Updates values on an existing grid.

This function does not alter the grid; it merely updates values on the existing grid. The number of values to read from "value" is nr_ or nk_. Values of the other space will be automatically updated if they exist.

Note
This function does not check the index bound; use with care!
Here is the call graph for this function:

◆ symbol()

std::string const & NumericalRadial::symbol ( ) const
inline

padded zeros ignored

Here is the caller graph for this function:

◆ to_numerical_orbital_lm()

void NumericalRadial::to_numerical_orbital_lm ( Numerical_Orbital_Lm orbital_lm,
const int  nk_legacy = 4005,
const double  lcao_dk = 0.01 
) const

Overwrites the content of a Numerical_Orbital_Lm object with the current object.

This function provides an interface to the corresponding object in module_ao. Due to algorithmic difference (FFT vs. Simpson's integration), it is inappropriate to use the k grid of NumericalRadial (which is FFT-compliant with r grid) to initialize the k grid of Numerical_Orbital_Lm.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ transform()

void NumericalRadial::transform ( const bool  forward)
private

Transforms the r-space values to get k-space values, or vice versa.

The grid & values where the transform is initiated must exist; this function does nothing if grid in the destination space does not exist.

forward : r to k backward: k to r

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wipe()

void NumericalRadial::wipe ( const bool  r_space = true,
const bool  k_space = true 
)

Removes the grid & values in r or k space.

Here is the caller graph for this function:

Member Data Documentation

◆ ikcut_

int NumericalRadial::ikcut_ = 0
private

◆ ircut_

int NumericalRadial::ircut_ = 0
private

Index of the first trailing zero.

A numerical radial function might be zero-padded for the sake of FFT-based spherical Bessel transform algorithms. The following two variables keep track of the actual cutoff radius. Specifically, if there are no trailing zeros in rvalues_, then ircut_ = nr_; if there are trailing zeros, then ircut_ is the index of the first trailing zero. For example, rvalues_ = {1, 2, 3, 0, 0, 0} -> ircut_ = 3 rvalues_ = {1, 2, 3, 4, 5, 6} -> ircut_ = 6 rvalues_ = {0, 0, 0, 0, 0, 0} -> ircut_ = 0

◆ is_fft_compliant_

bool NumericalRadial::is_fft_compliant_ = false
private

A flag that tells whether the r & k grids are FFT-compliant.

r & k grids are considered FFT-compliant if they

  1. have the same number of grid points;
  2. are both uniform;
  3. both starts from 0;
  4. satisfy dr*dk = pi/(N-1) where N >= 2 is the number of each grid points.

If the grids are FFT-compliant, spherical Bessel transforms are performed with an FFT-based algorithm. Otherwise, the transforms are performed with numerical integration (Simpson's rule).

◆ itype_

int NumericalRadial::itype_ = 0
private

element index in calculation

◆ izeta_

int NumericalRadial::izeta_ = 0
private

further index for NumericalRadial objects with the same itype_and l_

◆ kgrid_

double* NumericalRadial::kgrid_ = nullptr
private

k-space grid

◆ kvalue_

double* NumericalRadial::kvalue_ = nullptr
private

k-space value

◆ l_

int NumericalRadial::l_ = -1
private

angular momentum

◆ nk_

int NumericalRadial::nk_ = 0
private

number of k-space grid points

◆ nr_

int NumericalRadial::nr_ = 0
private

number of r-space grid points

◆ pk_

int NumericalRadial::pk_ = 0
private

implicit exponent in kvalues_

This parameter affects how this class interprets kvalues_. Specifically, kvalues_[ik] is interpreted as pow(kgrid_[ik], pk_) * F(kgrid_[ik]) during spherical Bessel transforms.

◆ pr_

int NumericalRadial::pr_ = 0
private

implicit exponent in rvalues_

This parameter affects how this class interprets rvalues_. Specifically, rvalues_[ir] is interpreted as pow(rgrid_[ir], pr_) * F(rgrid_[ir]) during spherical Bessel transforms.

◆ rgrid_

double* NumericalRadial::rgrid_ = nullptr
private

r-space grid

◆ rvalue_

double* NumericalRadial::rvalue_ = nullptr
private

r-space value

◆ sbt_

ModuleBase::SphericalBesselTransformer NumericalRadial::sbt_
private

An object that provides spherical Bessel transforms.

◆ symbol_

std::string NumericalRadial::symbol_ = ""
private

chemical symbol


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