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