Hex  1.0
Hydrogen-electron collision solver
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions
RadialIntegrals Class Reference

#include <radial.h>

Public Member Functions

 RadialIntegrals (Bspline const &bspline)
 
void setupOneElectronIntegrals ()
 
void setupTwoElectronIntegrals (Parallel const &par, CommandLine const &cmd, Array< bool > const &lambdas)
 
Complex computeD_iknot (int i, int j, int iknot) const
 
Complex computeD (int i, int j, int maxknot=-1) const
 
Complex computeM_iknot (int a, int i, int j, int iknot, Complex R) const
 
Complex computeM (int a, int i, int j, int maxknot=0) const
 
cArray computeMi (int a, int iknotmax=0) const
 
rArray computeScale (int a, int iknotmax=0) const
 
void M_integrand (int n, Complex *in, Complex *out, int i, int j, int a, int iknot, int iknotmax, double &logscale) const
 
Complex computeR (int lambda, int i, int j, int k, int l, cArray const &Mtr_L, cArray const &Mtr_mLm1) const
 
Complex computeRdiag (int L, int a, int b, int c, int d, int iknot, int iknotmax) const
 
Complex computeRtri (int L, int k, int l, int m, int n, int iknot, int iknotmax) const
 
void R_inner_integrand (int n, Complex *in, Complex *out, int i, int j, int L, int iknot, int iknotmax, Complex x) const
 
void R_outer_integrand (int n, Complex *in, Complex *out, int i, int j, int k, int l, int L, int iknot, int iknotmax) const
 
void allSymmetries (int i, int j, int k, int l, Complex Rijkl_tr, NumberArray< long > &R_tr_i, NumberArray< long > &R_tr_j, NumberArray< Complex > &R_tr_v) const
 
template<class Functor >
cArray overlapP (int n, int l, Functor weightf) const
 
template<class Functor >
cArray overlapj (int maxell, const rArrayView vk, Functor weightf) const
 Compute j-overlaps. More...
 
Bspline const & bspline () const
 
SymDiaMatrix const & D () const
 
SymDiaMatrix const & S () const
 
SymDiaMatrix const & Mm1 () const
 
SymDiaMatrix const & Mm1_tr () const
 
SymDiaMatrix const & Mm2 () const
 
SymDiaMatrix const & R_tr_dia (unsigned i) const
 
size_t maxlambda () const
 

Constructor & Destructor Documentation

RadialIntegrals::RadialIntegrals ( Bspline const &  bspline)
inline

Member Function Documentation

void RadialIntegrals::allSymmetries ( int  i,
int  j,
int  k,
int  l,
Complex  Rijkl_tr,
NumberArray< long > &  R_tr_i,
NumberArray< long > &  R_tr_j,
NumberArray< Complex > &  R_tr_v 
) const
Bspline const& RadialIntegrals::bspline ( ) const
inline
Complex RadialIntegrals::computeD ( int  i,
int  j,
int  maxknot = -1 
) const

Compute derivative overlap for B-splines \( B_i \) and \( B_j \).

Parameters
iB-spline index.
jB-spline index.
maxknotRight-most knot of any integration.
Complex RadialIntegrals::computeD_iknot ( int  i,
int  j,
int  iknot 
) const

Compute derivative overlap of B-splines \( B_i \) and \( B_j \) over the knot "iknot", using Gauss-Legendre integration.

Parameters
iB-spline index.
jB-spline index.
iknotInterval index.
Complex RadialIntegrals::computeM ( int  a,
int  i,
int  j,
int  maxknot = 0 
) const

Compute integral moment of coordinate power between the B-splines \( B_i \) and \( B_j \)

Parameters
aExponent.
iB-spline index.
jB-spline index.
maxknotRight-most knot of any integration.
Complex RadialIntegrals::computeM_iknot ( int  a,
int  i,
int  j,
int  iknot,
Complex  R 
) const

Compute integral moment of coordinate power between the B-splines \( B_i \) and \( B_j \) over the knot "iknot", using Gauss-Legendre integration.

Parameters
aExponent.
iB-spline index.
jB-spline index.
iknotInterval index.
RPotential truncation point.
cArray RadialIntegrals::computeMi ( int  a,
int  iknotmax = 0 
) const

Compute logarithms of integral moment of degree "a" for every B-spline pair and every interknot sub-interval. Store in 1-D array of shape

* [ Nspline × Nspline × Nintval ]
*
Parameters
aMoment degree.
iknotmaxIndex of knot that terminates the integration range.
Complex RadialIntegrals::computeR ( int  lambda,
int  i,
int  j,
int  k,
int  l,
cArray const &  Mtr_L,
cArray const &  Mtr_mLm1 
) const

Compute the two-electron (Slater-type) four-B-spline integral.

Parameters
lambdaMultipole degree.
iFirst (x-dependent) B-spline index.
jSecond (y-dependent) B-spline index.
kThird (x-dependent) B-spline index.
lFourth (y-dependent) B-spline index.
Mtr_LLogarithms of R₀-truncated partial moments.
Mtr_mLm1Logarithms of R₀-truncated partial moments.

Given the R-type integral symmetry, following calls will produce identical results:

* computeR(lambda, i, j, k, l);
* computeR(lambda, j, i, l, k);
* computeR(lambda, k, j, i, l);
* computeR(lambda, i, l, k, j);
* computeR(lambda, k, l, i, j);
*
Returns
Value of Rtr.
Complex RadialIntegrals::computeRdiag ( int  L,
int  a,
int  b,
int  c,
int  d,
int  iknot,
int  iknotmax 
) const

Diagonal part of Slater integral

Complex RadialIntegrals::computeRtri ( int  L,
int  k,
int  l,
int  m,
int  n,
int  iknot,
int  iknotmax 
) const

Triangle integral

rArray RadialIntegrals::computeScale ( int  a,
int  iknotmax = 0 
) const
SymDiaMatrix const& RadialIntegrals::D ( ) const
inline
void RadialIntegrals::M_integrand ( int  n,
Complex in,
Complex out,
int  i,
int  j,
int  a,
int  iknot,
int  iknotmax,
double &  logscale 
) const
size_t RadialIntegrals::maxlambda ( ) const
inline
SymDiaMatrix const& RadialIntegrals::Mm1 ( ) const
inline
SymDiaMatrix const& RadialIntegrals::Mm1_tr ( ) const
inline
SymDiaMatrix const& RadialIntegrals::Mm2 ( ) const
inline
template<class Functor >
cArray RadialIntegrals::overlapj ( int  maxell,
const rArrayView  vk,
Functor  weightf 
) const
inline

Compute B-spline overlap integrals for Riccati-Bessel function.

Parameters
maxellMaximal degree of the Riccati-Bessel function.
vkArray containing linear momenta.
weightfWeight function to multiply every value of the Bessel function. It is expected to have the "double operator() (Complex z)" interface, where the sent value is the complex coordinate.
Returns
Array of shape [vk.size() × (maxell + 1) × Nspline] in column-major format.
template<class Functor >
cArray RadialIntegrals::overlapP ( int  n,
int  l,
Functor  weightf 
) const
inline

Compute P-overlaps Compute overlap vector of B-splines vs. hydrogen Pnl function.

Parameters
nPrincipal quantum number.
lOrbital quantum number.
weightfWeight function to multiply every value of the hydrogenic function. It is expected to have the "double operator() (Complex z)" interface, where the sent value is the complex coordinate.
void RadialIntegrals::R_inner_integrand ( int  n,
Complex in,
Complex out,
int  i,
int  j,
int  L,
int  iknot,
int  iknotmax,
Complex  x 
) const
void RadialIntegrals::R_outer_integrand ( int  n,
Complex in,
Complex out,
int  i,
int  j,
int  k,
int  l,
int  L,
int  iknot,
int  iknotmax 
) const
SymDiaMatrix const& RadialIntegrals::R_tr_dia ( unsigned  i) const
inline
SymDiaMatrix const& RadialIntegrals::S ( ) const
inline
void RadialIntegrals::setupOneElectronIntegrals ( )
void RadialIntegrals::setupTwoElectronIntegrals ( Parallel const &  par,
CommandLine const &  cmd,
Array< bool > const &  lambdas 
)

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