DFT-EFE
 
Loading...
Searching...
No Matches
SphericalHarmonicFunctions.h
Go to the documentation of this file.
1#include <vector>
2
3namespace dftefe
4{
5 namespace atoms
6 {
7 void
8 convertCartesianToSpherical(const std::vector<double> &x,
9 double & r,
10 double & theta,
11 double & phi,
12 double polarAngleTolerance);
13
17
18 //
19 // We use the real form of spherical harmonics without the Condon-Shortley
20 // phase (i.e., the (-1)^m prefactor) (see
21 // https://en.wikipedia.org/wiki/Spherical_harmonics) NOTE: 1) The wikipedia
22 // definition has the Condon-Shortley phase.
23 // 2) The definition of the associated Legendre polynomial (P_lm) in
24 // Boost library also contains a Condon-Shortley phase.
25 // Thus, if you're using Boost library, multiply the P_lm
26 // evaluation with (-1)^m to remove the Condon-Shortley phase. Most
27 // Quantum Chemistry codes (e.g., QChem) do not include the
28 // Condon-Shortley phase. So to make it consistent, it is prefered
29 // to remove the Condon-Shortley phase, if there are any to begin
30 // with.
31 // 3) From C++17 onwards, the <cmath> has the associated Legendre
32 // polynomial (see
33 // https://en.cppreference.com/w/cpp/numeric/special_functions/assoc_legendre)
34 // Thus, if you're using C++17 or beyond, you can use the C++
35 // standard's definition of associated Legendre polynomial instead
36 // of Boost. Note that, the C++ standard does not have the
37 // Condon-Shortley phase while Boost has it. So, we do not have to
38 // do anything special to remove it while using the C++ standard.
39 //
40
41 //
42 // Y_lm(theta, phi) = Clm(l,m) * Dm(m) * P_lm(l,m,cos(theta)) * Qm(m,phi),
43 // where theta = polar angle,
44 // phi = azimuthal angle
45 // P_lm is the associated Legendre polynomial of degree l and order m,
46 // Qm is the real form of exp(i*m*phi),
47 // C_lm is the normalization constant for P_lm,
48 // D_m is the normalization constant for Q_m
49 //
50
51 //
52 // For the definition of the associated Legendre polynomials i.e. P_lm and
53 // their derivatives (as used for evaluating the real form of spherical
54 // harmonics and their derivatives) refer:
55 // @article{bosch2000computation,
56 // title={On the computation of derivatives of Legendre functions},
57 // author={Bosch, W},
58 // journal={Physics and Chemistry of the Earth, Part A: Solid Earth
59 // and Geodesy}, volume={25}, number={9-11}, pages={655--659},
60 // year={2000},
61 // publisher={Elsevier}
62 // }
63 // We use the derivative definitions from the above reference because
64 // finding the derivatives on the pole (i.e., theta = 0) is tricky. This is
65 // because the azimuthal angles (phi) is undefined for a point on the pole.
66 // However, the derivative is still well defined on the pole via the
67 // L'Hospital's rule. However, one can avoid implementing tedious
68 // L'Hospital's rule on pole and use much simpler expressions given in the
69 // above reference.
70 //
71
72 double
73 Dm(const int m);
74
75 double
76 Clm(const int l, const int m);
77
78 double
79 Qm(const int m, const double phi);
80
81 double
82 dQmDPhi(const int m, const double phi);
83
84 double
85 Plm(const int l, const int m, const double x);
86
87 double
88 dPlmDTheta(const int l, const int m, const double theta);
89
90 double
91 d2PlmDTheta2(const int l, const int m, const double theta);
95 } // namespace atoms
96} // namespace dftefe
void convertCartesianToSpherical(const std::vector< double > &x, double &r, double &theta, double &phi, double polarAngleTolerance)
Definition: SphericalHarmonicFunctions.cpp:13
double d2PlmDTheta2(const int l, const int m, const double theta)
Definition: SphericalHarmonicFunctions.cpp:276
double Dm(const int m)
Definition: SphericalHarmonicFunctions.cpp:104
double dQmDPhi(const int m, const double phi)
Definition: SphericalHarmonicFunctions.cpp:157
double Qm(const int m, const double phi)
Definition: SphericalHarmonicFunctions.cpp:143
double Clm(const int l, const int m)
Definition: SphericalHarmonicFunctions.cpp:135
double dPlmDTheta(const int l, const int m, const double theta)
Definition: SphericalHarmonicFunctions.cpp:239
double Plm(const int l, const int m, const double x)
Definition: SphericalHarmonicFunctions.cpp:177
dealii includes
Definition: AtomFieldDataSpherical.cpp:31