ampsci
High-precision calculations for one- and two-valence atomic systems
SphericalBessel Namespace Reference

Detailed Description

Spherical Bessel functions j_L(x) and y_L(x), and related utilities.

Provides an approximate version (good to ~1e-9) and exact GSL-backed versions, plus cell-averaged vectors for integration on coarse grids, and a lookup table.

Classes

class  JL_table
 Lookup table of spherical Bessel functions: j_L(q*r) = J[L][q][r]. More...
 

Functions

double jL (int L, double x)
 Spherical Bessel function j_L(x); uses low-x approximation where appropriate.
 
double exactGSL_JL (int L, double x)
 Spherical Bessel function of first kind via GSL.
 
double exactGSL_YL (int L, double x)
 Spherical Bessel functions of second kind.
 
double yL (int L, double x)
 Spherical Bessel function of second kind, y_L(x)
 
double PhiL (int L, double x, bool tilde=false)
 Phi_L(x) = [(2L+1)!! / x^L] * j_L(x). If tilde=true, returns 1 - Phi_L(x)
 
double PsiL (int L, double x, bool tilde=false)
 Psi_L(x) = [x^{L+1} / (2L-1)!!] * y_L(x). If tilde=true, returns 1 - Psi_L(x)
 
std::vector< double > fillBesselVec (int l, const std::vector< double > &xvec, bool cell_average=true)
 Fills a vector with the spherical Bessel function j_L sampled on the supplied grid; oscillation-resolving cell-average is used where needed.
 
std::vector< double > fillBesselVec_kr (int l, double k, const std::vector< double > &rvec, bool cell_average=true)
 As fillBesselVec, with the argument of j_L being k*r instead of x.
 
void fillBesselVec_kr (int l, double k, const std::vector< double > &r, std::vector< double > *jl, bool cell_average=true)
 In-place variant of fillBesselVec_kr; writes into a caller-owned vector.
 

Function Documentation

◆ jL()

double SphericalBessel::jL ( int  L,
double  x 
)

Spherical Bessel function j_L(x); uses low-x approximation where appropriate.

◆ exactGSL_JL()

double SphericalBessel::exactGSL_JL ( int  L,
double  x 
)

Spherical Bessel function of first kind via GSL.

◆ exactGSL_YL()

double SphericalBessel::exactGSL_YL ( int  L,
double  x 
)

Spherical Bessel functions of second kind.

◆ yL()

double SphericalBessel::yL ( int  L,
double  x 
)

Spherical Bessel function of second kind, y_L(x)

◆ PhiL()

double SphericalBessel::PhiL ( int  L,
double  x,
bool  tilde 
)

Phi_L(x) = [(2L+1)!! / x^L] * j_L(x). If tilde=true, returns 1 - Phi_L(x)

◆ PsiL()

double SphericalBessel::PsiL ( int  L,
double  x,
bool  tilde 
)

Psi_L(x) = [x^{L+1} / (2L-1)!!] * y_L(x). If tilde=true, returns 1 - Psi_L(x)

◆ fillBesselVec()

std::vector< double > SphericalBessel::fillBesselVec ( int  l,
const std::vector< double > &  xvec,
bool  cell_average = true 
)

Fills a vector with the spherical Bessel function j_L sampled on the supplied grid; oscillation-resolving cell-average is used where needed.

Returns a vector v such that v[i] is the value (or cell-average) of j_L on the grid point xvec[i]. The "cell" at index i is the interval bounded by the midpoints to its neighbours (asymmetric for the first/last point).

This vector is intended for integration of the form \( \int j_L(x)\,\psi(x)\,dx \approx \sum_i w_i\,v_i\,\psi(x_i) \) against a smooth wavefunction \(\psi\) sampled on the same grid. When the grid cell at index i resolves the oscillation of j_L ( \(\Delta x_i \le \lambda_j / N\) with \(\lambda_j = 2\pi\) asymptotically and N = 10), v[i] is the point value j_L(xvec[i]). When the cell is too coarse, v[i] is the cell average

\[ v_i = \frac{1}{\Delta x_i}\int_{\text{cell}_i} j_L(x)\,dx , \]

computed by trapezoidal sub-sampling with enough sub-points to resolve j_L on the cell. Using these cell averages in the Riemann-sum quadrature is equivalent to a Filon-style integration: the (smooth) factor \(\psi\) is treated as piecewise-constant per cell while the oscillation of j_L is integrated exactly within each cell. This recovers the correct Riemann-Lebesgue (~1/x) suppression at large x and removes the aliasing seen with pure point-sampling on a coarse grid.

Parameters
lOrder of the spherical Bessel function.
xvecGrid (monotonic) on which j_L is sampled.
cell_averageIf true (default), enables the cell averaging described above. If false, every entry is a pure point sample j_L(xvec[i]); use this for inspection / plotting or to reproduce legacy behaviour.
Returns
Vector v with v[i] = j_L(xvec[i]) or its cell average; same size as xvec.
Warning
The returned values are NOT pure point samples on cells where averaging triggers. They are only meaningful in conjunction with a Riemann-sum-style quadrature using the same grid weights.

◆ fillBesselVec_kr() [1/2]

std::vector< double > SphericalBessel::fillBesselVec_kr ( int  l,
double  k,
const std::vector< double > &  rvec,
bool  cell_average = true 
)

As fillBesselVec, with the argument of j_L being k*r instead of x.

Returns a vector v with v[i] = j_L(k*rvec[i]) or, where the cell at index i is too coarse for the oscillation of j_L(kr) (cell width \(\Delta r_i > 2\pi / (k\,N)\) with N = 10), the cell average

\[ v_i = \frac{1}{\Delta r_i}\int_{\text{cell}_i} j_L(k r)\,dr . \]

See fillBesselVec for the usage and warning.

Parameters
lOrder of the spherical Bessel function.
kMomentum (or other) factor multiplying r in the argument.
rvecRadial grid on which j_L(kr) is sampled.
cell_averageIf true (default), enables cell averaging; if false, point-samples j_L(k*rvec[i]) at every grid point.
Returns
Vector v of same length as rvec.

◆ fillBesselVec_kr() [2/2]

void SphericalBessel::fillBesselVec_kr ( int  l,
double  k,
const std::vector< double > &  r,
std::vector< double > *  jl,
bool  cell_average = true 
)

In-place variant of fillBesselVec_kr; writes into a caller-owned vector.

Identical semantics to the returning overload (including the cell-averaging behaviour); jl is resized to r.size().

Parameters
lOrder of the spherical Bessel function.
kMomentum factor multiplying r in the argument.
rRadial grid.
jlOutput vector (resized internally); must be non-null.
cell_averageIf true (default), enables cell averaging; if false, point-samples j_L(k*r[i]) at every grid point.