SIRIUS 7.5.0
Electronic structure library and applications
Classes | Functions
sirius::sf Namespace Reference

Special functions. More...

Classes

class  Spherical_Bessel_functions
 Spherical Bessel functions \( j_{\ell}(q x) \) up to lmax. More...
 

Functions

static void custom_bessel (int lmax, double x, double *result)
 
int lmmax (int lmax)
 Maximum number of \( \ell, m \) combinations for a given \( \ell_{max} \). More...
 
int lm (int l, int m)
 Get composite lm index by angular index l and azimuthal index m. More...
 
int lmax (int lmmax__)
 Get maximum orbital quantum number by the maximum lm index. More...
 
std::vector< int > l_by_lm (int lmax__)
 Get array of orbital quantum numbers for each lm component. More...
 
double hermiteh (int n, double x)
 
std::vector< double > hermiteh_array (int n, double x)
 
double hermiteh_series (int n, double x, const double *a)
 
template<typename T , typename F >
void legendre_plm (int lmax__, double x__, F &&ilm__, T *plm__)
 Generate associated Legendre polynomials. More...
 
template<typename T , typename F >
void legendre_plm_aux (int lmax__, double x__, F &&ilm__, T const *plm__, T *p1lm__, T *p2lm__)
 Generate auxiliary Legendre polynomials which are necessary for spherical harmomic derivatives. More...
 
void spherical_harmonics_ref (int lmax, double theta, double phi, std::complex< double > *ylm)
 Reference implementation of complex spherical harmonics. More...
 
void spherical_harmonics (int lmax, double theta, double phi, std::complex< double > *ylm)
 Optimized implementation of complex spherical harmonics. More...
 
void spherical_harmonics_ref (int lmax, double theta, double phi, double *rlm)
 Reference implementation of real spherical harmonics Rlm. More...
 
void spherical_harmonics (int lmax, double theta, double phi, double *rlm)
 Optimized implementation of real spherical harmonics. More...
 
sddk::mdarray< double, 1 > cosxn (int n__, double x__)
 Generate \( \cos(m x) \) for m in [1, n] using recursion. More...
 
sddk::mdarray< double, 1 > sinxn (int n__, double x__)
 Generate \( \sin(m x) \) for m in [1, n] using recursion. More...
 
void dRlm_dr (int lmax__, r3::vector< double > &r__, sddk::mdarray< double, 2 > &data__, bool divide_by_r__=true)
 Compute the derivatives of real spherical harmonics over the components of cartesian vector. More...
 
void dRlm_dr_numerical (int lmax__, r3::vector< double > &r__, sddk::mdarray< double, 2 > &data__, bool divide_by_r__=true)
 

Detailed Description

Special functions.

Function Documentation

◆ custom_bessel()

static void sirius::sf::custom_bessel ( int  lmax,
double  x,
double *  result 
)
static

Definition at line 14 of file sbessel.cpp.

◆ lmmax()

int sirius::sf::lmmax ( int  lmax)
inline

Maximum number of \( \ell, m \) combinations for a given \( \ell_{max} \).

Definition at line 44 of file specfunc.hpp.

◆ lm()

int sirius::sf::lm ( int  l,
int  m 
)
inline

Get composite lm index by angular index l and azimuthal index m.

Definition at line 50 of file specfunc.hpp.

◆ lmax()

int sirius::sf::lmax ( int  lmmax__)
inline

Get maximum orbital quantum number by the maximum lm index.

Definition at line 56 of file specfunc.hpp.

◆ l_by_lm()

std::vector< int > sirius::sf::l_by_lm ( int  lmax__)
inline

Get array of orbital quantum numbers for each lm component.

Definition at line 69 of file specfunc.hpp.

◆ hermiteh()

double sirius::sf::hermiteh ( int  n,
double  x 
)
inline

Definition at line 80 of file specfunc.hpp.

◆ hermiteh_array()

std::vector< double > sirius::sf::hermiteh_array ( int  n,
double  x 
)
inline

Definition at line 87 of file specfunc.hpp.

◆ hermiteh_series()

double sirius::sf::hermiteh_series ( int  n,
double  x,
const double *  a 
)
inline

Definition at line 94 of file specfunc.hpp.

◆ legendre_plm()

template<typename T , typename F >
void sirius::sf::legendre_plm ( int  lmax__,
double  x__,
F &&  ilm__,
T *  plm__ 
)
inline

Generate associated Legendre polynomials.

Normalised associated Legendre polynomials obey the following recursive relations:

\[ P_{m}^{m}(x) = -\sqrt{1 + \frac{1}{2m}} y P_{m-1}^{m-1}(x) \]

\[ P_{m+1}^{m}(x) = \sqrt{2 m + 3} x P_{m}^{m}(x) \]

\[ P_{\ell}^{m}(x) = a_{\ell}^{m}\big(xP_{\ell-1}^{m}(x) - b_{\ell}^{m}P_{\ell - 2}^{m}(x)\big) \]

where

\begin{eqnarray*} a_{\ell}^{m} &=& \sqrt{\frac{4 \ell^2 - 1}{\ell^2 - m^2}} \\ b_{\ell}^{m} &=& \sqrt{\frac{(\ell-1)^2-m^2}{4(\ell-1)^2-1}} \\ x &=& \cos \theta \\ y &=& \sin \theta \end{eqnarray*}

and

\[ P_{0}^{0} = \sqrt{\frac{1}{4\pi}} \]

Definition at line 123 of file specfunc.hpp.

◆ legendre_plm_aux()

template<typename T , typename F >
void sirius::sf::legendre_plm_aux ( int  lmax__,
double  x__,
F &&  ilm__,
T const *  plm__,
T *  p1lm__,
T *  p2lm__ 
)
inline

Generate auxiliary Legendre polynomials which are necessary for spherical harmomic derivatives.

Generate the following functions:

\begin{eqnarray*} P_{\ell m}^1(x) &=& y P_{\ell}^{m}(x)' \\ P_{\ell m}^2(x) &=& \frac{P_{\ell}^{m}(x)}{y} \end{eqnarray*}

where \( x = \cos \theta \), \( y = \sin \theta \) and \( P_{\ell}^{m}(x) \) are normalized Legendre polynomials.

Both functions obey the recursive relations similar to \( P_{\ell}^{m}(x) \) (see sirius::legendre_plm() for details). For \( P_{\ell m}^1(x) \) this is:

\[ P_{m m}^1(x) = \sqrt{1 + \frac{1}{2m}} \Big( -y P_{m-1,m-1}^{1}(x)' + x P_{m-1}^{m-1}(x) \Big) \]

\[ P_{m+1, m}^1(x) = \sqrt{2 m + 3} \Big( x P_{m,m}^{1}(x) + y P_{m}^{m}(x) \Big) \]

\[ P_{\ell, m}^{1}(x) = a_{\ell}^{m}\Big(x P_{\ell-1,m}^{1}(x) + yP_{\ell-1}^{m}(x) - b_{\ell}^{m} P_{\ell - 2,m}^{1}(x) \Big) \]

where

\[ y = \sqrt{1 - x^2} \]

and

\[ P_{0,0}^1 = 0 \]

And for \( P_{\ell m}^2(x) \) the recursion is:

\[ P_{m,m}^2(x) = -\sqrt{1 + \frac{1}{2m}} P_{m-1}^{m-1}(x) \]

\[ P_{m+1, m}^2(x) = \sqrt{2 m + 3} x P_{m,m}^{2}(x) \]

\[ P_{\ell, m}^2(x) = a_{\ell}^{m}\Big(x P_{\ell-1,m}^{2}(x) - b_{\ell}^{m}P_{\ell - 2,m}^{2}(x) \Big) \]

where

\[ P_{0,0}^2 = 0 \]

See sirius::legendre_plm() for basic definitions.

Definition at line 211 of file specfunc.hpp.

◆ spherical_harmonics_ref() [1/2]

void sirius::sf::spherical_harmonics_ref ( int  lmax,
double  theta,
double  phi,
std::complex< double > *  ylm 
)
inline

Reference implementation of complex spherical harmonics.

Complex spherical harmonics are defined as:

\[ Y_{\ell m}(\theta,\phi) = P_{\ell}^{m}(\cos \theta) e^{im\phi} \]

where \(P_{\ell}^m(x) \) are associated Legendre polynomials.

Mathematica code:

norm[l_, m_] := 4*Pi*Integrate[LegendreP[l, m, x]*LegendreP[l, m, x], {x, 0, 1}]
Ylm[l_, m_, t_, p_] := LegendreP[l, m, Cos[t]]*E^(I*m*p)/Sqrt[norm[l, m]]
Do[Print[ComplexExpand[
 FullSimplify[SphericalHarmonicY[l, m, t, p] - Ylm[l, m, t, p],
  Assumptions -> {0 <= t <= Pi}]]], {l, 0, 5}, {m, -l, l}]

Complex spherical harmonics obey the following symmetry:

\[ Y_{\ell -m}(\theta,\phi) = (-1)^m Y_{\ell m}^{*}(\theta,\phi) \]

Mathematica code:

Do[Print[ComplexExpand[
 FullSimplify[
  SphericalHarmonicY[l, -m, t, p] - (-1)^m*
   Conjugate[SphericalHarmonicY[l, m, t, p]],
    Assumptions -> {0 <= t <= Pi}]]], {l, 0, 4}, {m, 0, l}]

Definition at line 273 of file specfunc.hpp.

◆ spherical_harmonics() [1/2]

void sirius::sf::spherical_harmonics ( int  lmax,
double  theta,
double  phi,
std::complex< double > *  ylm 
)
inline

Optimized implementation of complex spherical harmonics.

Definition at line 298 of file specfunc.hpp.

◆ spherical_harmonics_ref() [2/2]

void sirius::sf::spherical_harmonics_ref ( int  lmax,
double  theta,
double  phi,
double *  rlm 
)
inline

Reference implementation of real spherical harmonics Rlm.

Real spherical harminics are defined as:

\[ R_{\ell m}(\theta,\phi) = \left\{ \begin{array}{lll} \sqrt{2} \Re Y_{\ell m}(\theta,\phi) = \sqrt{2} P_{\ell}^{m}(\cos \theta) \cos m\phi & m > 0 \\ P_{\ell}^{0}(\cos \theta) & m = 0 \\ \sqrt{2} \Im Y_{\ell m}(\theta,\phi) = \sqrt{2} (-1)^{|m|} P_{\ell}^{|m|}(\cos \theta) (-\sin |m|\phi) & m < 0 \end{array} \right. \]

Mathematica code:

(* definition of real spherical harmonics, use Plm(l,m) for m\
\[GreaterEqual]0 only *)

norm[l_, m_] :=
 4*Pi*Integrate[
   LegendreP[l, Abs[m], x]*LegendreP[l, Abs[m], x], {x, 0, 1}]
legendre[l_, m_, x_] := LegendreP[l, Abs[m], x]/Sqrt[norm[l, m]]

(* reference definition *)

RRlm[l_, m_, th_, ph_] :=
 If[m > 0, Sqrt[2]*ComplexExpand[Re[SphericalHarmonicY[l, m, th, ph]]
    ], If[m < 0,
   Sqrt[2]*ComplexExpand[Im[SphericalHarmonicY[l, m, th, ph]]],
   If[m == 0, ComplexExpand[Re[SphericalHarmonicY[l, 0, th, ph]]]]]]

(* definition without ComplexExpand *)

Rlm[l_, m_, th_, ph_] :=
 If[m > 0, legendre[l, m, Cos[th]]*Sqrt[2]*Cos[m*ph],
  If[m < 0, (-1)^m*legendre[l, m, Cos[th]]*Sqrt[2]*(-Sin[Abs[m]*ph]),
   If[m == 0, legendre[l, 0, Cos[th]]]]]

(* check that both definitions are identical *)
Do[
 Print[FullSimplify[Rlm[l, m, a, b] - RRlm[l, m, a, b],
   Assumptions -> {0 <= a <= Pi, 0 <= b <= 2*Pi}]], {l, 0, 5}, {m, -l,
   l}]

Definition at line 374 of file specfunc.hpp.

◆ spherical_harmonics() [2/2]

void sirius::sf::spherical_harmonics ( int  lmax,
double  theta,
double  phi,
double *  rlm 
)
inline

Optimized implementation of real spherical harmonics.

Definition at line 396 of file specfunc.hpp.

◆ cosxn()

sddk::mdarray< double, 1 > sirius::sf::cosxn ( int  n__,
double  x__ 
)
inline

Generate \( \cos(m x) \) for m in [1, n] using recursion.

Definition at line 429 of file specfunc.hpp.

◆ sinxn()

sddk::mdarray< double, 1 > sirius::sf::sinxn ( int  n__,
double  x__ 
)
inline

Generate \( \sin(m x) \) for m in [1, n] using recursion.

Definition at line 446 of file specfunc.hpp.

◆ dRlm_dr()

void sirius::sf::dRlm_dr ( int  lmax__,
r3::vector< double > &  r__,
sddk::mdarray< double, 2 > &  data__,
bool  divide_by_r__ = true 
)
inline

Compute the derivatives of real spherical harmonics over the components of cartesian vector.

The following derivative is computed:

\[ \frac{\partial R_{\ell m}(\theta_r, \phi_r)}{\partial r_{\mu}} = \frac{\partial R_{\ell m}(\theta_r, \phi_r)}{\partial \theta_r} \frac{\partial \theta_r}{\partial r_{\mu}} + \frac{\partial R_{\ell m}(\theta_r, \phi_r)}{\partial \phi_r} \frac{\partial \phi_r}{\partial r_{\mu}} \]

The derivatives of angles are:

\[ \frac{\partial \theta_r}{\partial r_{x}} = \frac{\cos(\phi_r) \cos(\theta_r)}{r} \\ \frac{\partial \theta_r}{\partial r_{y}} = \frac{\cos(\theta_r) \sin(\phi_r)}{r} \\ \frac{\partial \theta_r}{\partial r_{z}} = -\frac{\sin(\theta_r)}{r} \]

and

\[ \frac{\partial \phi_r}{\partial r_{x}} = -\frac{\sin(\phi_r)}{\sin(\theta_r) r} \\ \frac{\partial \phi_r}{\partial r_{y}} = \frac{\cos(\phi_r)}{\sin(\theta_r) r} \\ \frac{\partial \phi_r}{\partial r_{z}} = 0 \]

The derivative of \( \phi \) has discontinuities at \( \theta = 0, \theta=\pi \). This, however, is not a problem, because multiplication by the the derivative of \( R_{\ell m} \) removes it. The following functions have to be hardcoded:

\[ \frac{\partial R_{\ell m}(\theta, \phi)}{\partial \theta} \\ \frac{\partial R_{\ell m}(\theta, \phi)}{\partial \phi} \frac{1}{\sin(\theta)} \]

Spherical harmonics have a separable form:

\[ R_{\ell m}(\theta, \phi) = P_{\ell}^{m}(\cos \theta) f(\phi) \]

The derivative over \( \theta \) is then:

\[ \frac{\partial R_{\ell m}(\theta, \phi)}{\partial \theta} = \frac{\partial P_{\ell}^{m}(x)}{\partial x} \frac{\partial x}{\partial \theta} f(\phi) = -\sin \theta \frac{\partial P_{\ell}^{m}(x)}{\partial x} f(\phi) \]

where \( x = \cos \theta \)

Mathematica script for spherical harmonic derivatives:

Rlm[l_, m_, th_, ph_] :=
 If[m > 0, Sqrt[2]*ComplexExpand[Re[SphericalHarmonicY[l, m, th, ph]]],
   If[m < 0, Sqrt[2]*ComplexExpand[Im[SphericalHarmonicY[l, m, th, ph]]],
     If[m == 0, ComplexExpand[Re[SphericalHarmonicY[l, 0, th, ph]]]]
   ]
 ]
Do[Print[FullSimplify[D[Rlm[l, m, theta, phi], theta]]], {l, 0, 4}, {m, -l, l}]
Do[Print[FullSimplify[TrigExpand[D[Rlm[l, m, theta, phi], phi]/Sin[theta]]]], {l, 0, 4}, {m, -l, l}]

Definition at line 512 of file specfunc.hpp.

◆ dRlm_dr_numerical()

void sirius::sf::dRlm_dr_numerical ( int  lmax__,
r3::vector< double > &  r__,
sddk::mdarray< double, 2 > &  data__,
bool  divide_by_r__ = true 
)
inline

Definition at line 597 of file specfunc.hpp.