SIRIUS 7.5.0
Electronic structure library and applications
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
sirius::SHT Class Reference

Spherical harmonics transformations and related oprtations. More...

#include <sht.hpp>

Public Member Functions

 SHT (sddk::device_t pu__, int lmax__, int mesh_type__=0)
 Default constructor. More...
 
void check () const
 Check the transformations. More...
 
template<typename T >
void backward_transform (int ld, T const *flm, int nr, int lmmax, T *ftp) const
 Perform a backward transformation from spherical harmonics to spherical coordinates. More...
 
template<typename T >
void forward_transform (T const *ftp, int nr, int lmmax, int ld, T *flm) const
 Perform a forward transformation from spherical coordinates to spherical harmonics. More...
 
void uniform_coverage ()
 
auto ylm_backward (int lm, int itp) const
 
auto rlm_backward (int lm, int itp) const
 
auto coord (int x, int itp) const
 
auto coord (int idx__) const
 
auto theta (int idx__) const
 
auto phi (int idx__) const
 
auto weight (int idx__) const
 
auto num_points () const
 
auto lmax () const
 
auto lmmax () const
 
template<>
void backward_transform (int ld, double const *flm, int nr, int lmmax, double *ftp) const
 
template<>
void forward_transform (double const *ftp, int nr, int lmmax, int ld, double *flm) const
 

Static Public Member Functions

static void convert (int lmax__, double const *f_rlm__, std::complex< double > *f_ylm__)
 Convert form Rlm to Ylm representation. More...
 
static void convert (int lmax__, std::complex< double > const *f_ylm__, double *f_rlm__)
 Convert from Ylm to Rlm representation. More...
 
static std::complex< double > ylm_dot_rlm (int l, int m1, int m2)
 Compute element of the transformation matrix from complex to real spherical harmonics. More...
 
static std::complex< double > rlm_dot_ylm (int l, int m1, int m2)
 
static double gaunt_yyy (int l1, int l2, int l3, int m1, int m2, int m3)
 Gaunt coefficent of three complex spherical harmonics. More...
 
static double gaunt_rrr (int l1, int l2, int l3, int m1, int m2, int m3)
 Gaunt coefficent of three real spherical harmonics. More...
 
static double gaunt_rlm_ylm_rlm (int l1, int l2, int l3, int m1, int m2, int m3)
 Gaunt coefficent of two real spherical harmonics with a complex one. More...
 
static std::complex< double > gaunt_hybrid (int l1, int l2, int l3, int m1, int m2, int m3)
 Gaunt coefficent of two complex and one real spherical harmonics. More...
 
static double clebsch_gordan (int l1, int l2, int l3, int m1, int m2, int m3)
 Return Clebsch-Gordan coefficient. More...
 

Private Attributes

sddk::device_t pu_
 Type of processing unit. More...
 
int lmax_
 Maximum \( \ell \) of spherical harmonics. More...
 
int lmmax_
 Maximum number of \( \ell, m \) components. More...
 
int num_points_
 Number of real-space \( (\theta, \phi) \) points on the sphere. More...
 
sddk::mdarray< double, 2 > coord_
 Cartesian coordinates of points (normalized to 1). More...
 
sddk::mdarray< double, 2 > tp_
 \( (\theta, \phi) \) angles of points. More...
 
std::vector< double > w_
 Point weights. More...
 
sddk::mdarray< std::complex< double >, 2 > ylm_backward_
 Backward transformation from Ylm to spherical coordinates. More...
 
sddk::mdarray< std::complex< double >, 2 > ylm_forward_
 Forward transformation from spherical coordinates to Ylm. More...
 
sddk::mdarray< double, 2 > rlm_backward_
 Backward transformation from Rlm to spherical coordinates. More...
 
sddk::mdarray< double, 2 > rlm_forward_
 Forward transformation from spherical coordinates to Rlm. More...
 
int mesh_type_ {0}
 Type of spherical grid (0: Lebedev-Laikov, 1: uniform). More...
 

Detailed Description

Spherical harmonics transformations and related oprtations.

This class is responsible for the generation of complex and real spherical harmonics, generation of transformation matrices, transformation between spectral and real-space representations, generation of Gaunt and Clebsch-Gordan coefficients and calculation of spherical harmonic derivatives

Definition at line 80 of file sht.hpp.

Constructor & Destructor Documentation

◆ SHT()

sirius::SHT::SHT ( sddk::device_t  pu__,
int  lmax__,
int  mesh_type__ = 0 
)
inline

Default constructor.

Definition at line 121 of file sht.hpp.

Member Function Documentation

◆ check()

void sirius::SHT::check ( ) const

Check the transformations.

Definition at line 283 of file sht.cpp.

◆ backward_transform() [1/2]

template<typename T >
void sirius::SHT::backward_transform ( int  ld,
T const *  flm,
int  nr,
int  lmmax,
T *  ftp 
) const

Perform a backward transformation from spherical harmonics to spherical coordinates.

\[ f(\theta, \phi, r) = \sum_{\ell m} f_{\ell m}(r) Y_{\ell m}(\theta, \phi) \]

Parameters
[in]ldSize of leading dimension of flm.
[in]flmRaw pointer to \( f_{\ell m}(r) \).
[in]nrNumber of radial points.
[in]lmmaxMaximum number of lm- harmonics to take into sum.
[out]ftpRaw pointer to \( f(\theta, \phi, r) \).

◆ forward_transform() [1/2]

template<typename T >
void sirius::SHT::forward_transform ( T const *  ftp,
int  nr,
int  lmmax,
int  ld,
T *  flm 
) const

Perform a forward transformation from spherical coordinates to spherical harmonics.

\[ f_{\ell m}(r) = \iint f(\theta, \phi, r) Y_{\ell m}^{*}(\theta, \phi) \sin \theta d\phi d\theta = \sum_{i} f(\theta_i, \phi_i, r) Y_{\ell m}^{*}(\theta_i, \phi_i) w_i \]

Parameters
[in]ftpRaw pointer to \( f(\theta, \phi, r) \).
[in]nrNumber of radial points.
[in]lmmaxMaximum number of lm- coefficients to generate.
[in]ldSize of leading dimension of flm.
[out]flmRaw pointer to \( f_{\ell m}(r) \).

◆ convert() [1/2]

static void sirius::SHT::convert ( int  lmax__,
double const *  f_rlm__,
std::complex< double > *  f_ylm__ 
)
inlinestatic

Convert form Rlm to Ylm representation.

Definition at line 264 of file sht.hpp.

◆ convert() [2/2]

static void sirius::SHT::convert ( int  lmax__,
std::complex< double > const *  f_ylm__,
double *  f_rlm__ 
)
inlinestatic

Convert from Ylm to Rlm representation.

Definition at line 281 of file sht.hpp.

◆ ylm_dot_rlm()

static std::complex< double > sirius::SHT::ylm_dot_rlm ( int  l,
int  m1,
int  m2 
)
inlinestatic

Compute element of the transformation matrix from complex to real spherical harmonics.

Real spherical harmonic can be written as a linear combination of complex harmonics:

\[ R_{\ell m}(\theta, \phi) = \sum_{m'} a^{\ell}_{m' m}Y_{\ell m'}(\theta, \phi) \]

where

\[ a^{\ell}_{m' m} = \langle Y_{\ell m'} | R_{\ell m} \rangle \]

which gives the name to this function.

Transformation from real to complex spherical harmonics is conjugate transpose:

\[ Y_{\ell m}(\theta, \phi) = \sum_{m'} a^{\ell*}_{m m'}R_{\ell m'}(\theta, \phi) \]

Mathematica code:

b[m1_, m2_] :=
 If[m1 == 0, 1,
 If[m1 < 0 && m2 < 0, -I/Sqrt[2],
 If[m1 > 0 && m2 < 0, (-1)^m1*I/Sqrt[2],
 If[m1 < 0 && m2 > 0, (-1)^m2/Sqrt[2],
 If[m1 > 0 && m2 > 0, 1/Sqrt[2]]]]]]

a[m1_, m2_] := If[Abs[m1] == Abs[m2], b[m1, m2], 0]

Rlm[l_, m_, t_, p_] := Sum[a[m1, m]*SphericalHarmonicY[l, m1, t, p], {m1, -l, l}]

Definition at line 363 of file sht.hpp.

◆ rlm_dot_ylm()

static std::complex< double > sirius::SHT::rlm_dot_ylm ( int  l,
int  m1,
int  m2 
)
inlinestatic

Definition at line 392 of file sht.hpp.

◆ gaunt_yyy()

static double sirius::SHT::gaunt_yyy ( int  l1,
int  l2,
int  l3,
int  m1,
int  m2,
int  m3 
)
inlinestatic

Gaunt coefficent of three complex spherical harmonics.

\[ \langle Y_{\ell_1 m_1} | Y_{\ell_2 m_2} | Y_{\ell_3 m_3} \rangle \]

Definition at line 403 of file sht.hpp.

◆ gaunt_rrr()

static double sirius::SHT::gaunt_rrr ( int  l1,
int  l2,
int  l3,
int  m1,
int  m2,
int  m3 
)
inlinestatic

Gaunt coefficent of three real spherical harmonics.

\[ \langle R_{\ell_1 m_1} | R_{\ell_2 m_2} | R_{\ell_3 m_3} \rangle \]

Definition at line 423 of file sht.hpp.

◆ gaunt_rlm_ylm_rlm()

static double sirius::SHT::gaunt_rlm_ylm_rlm ( int  l1,
int  l2,
int  l3,
int  m1,
int  m2,
int  m3 
)
inlinestatic

Gaunt coefficent of two real spherical harmonics with a complex one.

\[ \langle R_{\ell_1 m_1} | Y_{\ell_2 m_2} | R_{\ell_3 m_3} \rangle \]

Definition at line 452 of file sht.hpp.

◆ gaunt_hybrid()

static std::complex< double > sirius::SHT::gaunt_hybrid ( int  l1,
int  l2,
int  l3,
int  m1,
int  m2,
int  m3 
)
inlinestatic

Gaunt coefficent of two complex and one real spherical harmonics.

\[ \langle Y_{\ell_1 m_1} | R_{\ell_2 m_2} | Y_{\ell_3 m_3} \rangle \]

Definition at line 478 of file sht.hpp.

◆ uniform_coverage()

void sirius::SHT::uniform_coverage ( )
inline

Definition at line 495 of file sht.hpp.

◆ clebsch_gordan()

static double sirius::SHT::clebsch_gordan ( int  l1,
int  l2,
int  l3,
int  m1,
int  m2,
int  m3 
)
inlinestatic

Return Clebsch-Gordan coefficient.

Clebsch-Gordan coefficients arise when two angular momenta are combined into a total angular momentum.

Definition at line 515 of file sht.hpp.

◆ ylm_backward()

auto sirius::SHT::ylm_backward ( int  lm,
int  itp 
) const
inline

Definition at line 528 of file sht.hpp.

◆ rlm_backward()

auto sirius::SHT::rlm_backward ( int  lm,
int  itp 
) const
inline

Definition at line 533 of file sht.hpp.

◆ coord() [1/2]

auto sirius::SHT::coord ( int  x,
int  itp 
) const
inline

Definition at line 538 of file sht.hpp.

◆ coord() [2/2]

auto sirius::SHT::coord ( int  idx__) const
inline

Definition at line 543 of file sht.hpp.

◆ theta()

auto sirius::SHT::theta ( int  idx__) const
inline

Definition at line 548 of file sht.hpp.

◆ phi()

auto sirius::SHT::phi ( int  idx__) const
inline

Definition at line 553 of file sht.hpp.

◆ weight()

auto sirius::SHT::weight ( int  idx__) const
inline

Definition at line 558 of file sht.hpp.

◆ num_points()

auto sirius::SHT::num_points ( ) const
inline

Definition at line 563 of file sht.hpp.

◆ lmax()

auto sirius::SHT::lmax ( ) const
inline

Definition at line 568 of file sht.hpp.

◆ lmmax()

auto sirius::SHT::lmmax ( ) const
inline

Definition at line 573 of file sht.hpp.

◆ backward_transform() [2/2]

template<>
void sirius::SHT::backward_transform ( int  ld,
double const *  flm,
int  nr,
int  lmmax,
double *  ftp 
) const

Definition at line 245 of file sht.cpp.

◆ forward_transform() [2/2]

template<>
void sirius::SHT::forward_transform ( double const *  ftp,
int  nr,
int  lmmax,
int  ld,
double *  flm 
) const

Definition at line 265 of file sht.cpp.

Member Data Documentation

◆ pu_

sddk::device_t sirius::SHT::pu_
private

Type of processing unit.

Definition at line 84 of file sht.hpp.

◆ lmax_

int sirius::SHT::lmax_
private

Maximum \( \ell \) of spherical harmonics.

Definition at line 87 of file sht.hpp.

◆ lmmax_

int sirius::SHT::lmmax_
private

Maximum number of \( \ell, m \) components.

Definition at line 90 of file sht.hpp.

◆ num_points_

int sirius::SHT::num_points_
private

Number of real-space \( (\theta, \phi) \) points on the sphere.

Definition at line 93 of file sht.hpp.

◆ coord_

sddk::mdarray<double, 2> sirius::SHT::coord_
private

Cartesian coordinates of points (normalized to 1).

Definition at line 96 of file sht.hpp.

◆ tp_

sddk::mdarray<double, 2> sirius::SHT::tp_
private

\( (\theta, \phi) \) angles of points.

Definition at line 99 of file sht.hpp.

◆ w_

std::vector<double> sirius::SHT::w_
private

Point weights.

Definition at line 102 of file sht.hpp.

◆ ylm_backward_

sddk::mdarray<std::complex<double>, 2> sirius::SHT::ylm_backward_
private

Backward transformation from Ylm to spherical coordinates.

Definition at line 105 of file sht.hpp.

◆ ylm_forward_

sddk::mdarray<std::complex<double>, 2> sirius::SHT::ylm_forward_
private

Forward transformation from spherical coordinates to Ylm.

Definition at line 108 of file sht.hpp.

◆ rlm_backward_

sddk::mdarray<double, 2> sirius::SHT::rlm_backward_
private

Backward transformation from Rlm to spherical coordinates.

Definition at line 111 of file sht.hpp.

◆ rlm_forward_

sddk::mdarray<double, 2> sirius::SHT::rlm_forward_
private

Forward transformation from spherical coordinates to Rlm.

Definition at line 114 of file sht.hpp.

◆ mesh_type_

int sirius::SHT::mesh_type_ {0}
private

Type of spherical grid (0: Lebedev-Laikov, 1: uniform).

Definition at line 117 of file sht.hpp.


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