|
SIRIUS 7.5.0
Electronic structure library and applications
|
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... | |
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
|
inline |
| void sirius::SHT::check | ( | ) | const |
| 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) \]
| [in] | ld | Size of leading dimension of flm. |
| [in] | flm | Raw pointer to \( f_{\ell m}(r) \). |
| [in] | nr | Number of radial points. |
| [in] | lmmax | Maximum number of lm- harmonics to take into sum. |
| [out] | ftp | Raw pointer to \( f(\theta, \phi, r) \). |
| 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 \]
| [in] | ftp | Raw pointer to \( f(\theta, \phi, r) \). |
| [in] | nr | Number of radial points. |
| [in] | lmmax | Maximum number of lm- coefficients to generate. |
| [in] | ld | Size of leading dimension of flm. |
| [out] | flm | Raw pointer to \( f_{\ell m}(r) \). |
|
inlinestatic |
|
inlinestatic |
|
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}]
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inlinestatic |
|
inline |
|
inline |
| void sirius::SHT::backward_transform | ( | int | ld, |
| double const * | flm, | ||
| int | nr, | ||
| int | lmmax, | ||
| double * | ftp | ||
| ) | const |
| void sirius::SHT::forward_transform | ( | double const * | ftp, |
| int | nr, | ||
| int | lmmax, | ||
| int | ld, | ||
| double * | flm | ||
| ) | const |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |