DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
SphericalFunctionUtil.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17// @author Bikash Kanungo
18//
19
20#ifndef DFTFE_SPHERICALFUNCTIONUTIL_H
21#define DFTFE_SPHERICALFUNCTIONUTIL_H
22
23#include <vector>
24#include <TypeConfig.h>
25namespace dftfe
26{
27 namespace utils
28 {
29 namespace sphUtils
30 {
31 /*
32 * @brief Function to convert cartesian coordinates into spherical coordinates
33 * @param[in] x Vector of size 3 containing the cartesian coordinates of
34 * the point
35 * @param[out] r Stores the computed radius of the point
36 * @param[out]theta Stores the computed polar angle of the point
37 * @param[out] phiStores the computed azimuthal angle of the point
38 * @param[in] rTol Defines a tolerance for the radius,
39 * below which theta and phi are set to zero.
40 * To elaborate, for radius->0, the angles are undefined.
41 * We set them to zero as a convenient choice
42 * @param[in] angleTol Defines a tolerance for the polar angle (theta)
43 * approching the poles, below which the azimuthal angle is
44 * set to zero. To elaborate, for a point on the pole
45 * (theta = 0 or theta = pi), the azimuthal angle is undefined.
46 * We set phi to zero if |theta - 0| < angleTol or |theta - pi|
47 * < angleTol.
48 */
49 void
50 convertCartesianToSpherical(const std::vector<double> &x,
51 double &r,
52 double &theta,
53 double &phi,
54 const double rTol,
55 const double angleTol);
56
57 double
58 Clm(const dftfe::Int l, const dftfe::Int m);
59
60 /*
61 * @brief Function to compute the azimuthal angle (phi) dependent
62 * part of the real spherical harmonics. It does not include
63 * any normalization constant
64 * Q(m,phi) = cos(m*phi) if m > 0
65 * Qm(m,phi) = 1 , if m = 0
66 * Qm(m, phi = sin(|m|phi), if m < 0
67 * @param[in] m order of the spherical harmonic (m-quantum number) for
68 * which Qm is to be evaluated
69 * @param[in] phi Azimuthal angle at which Qm is to be evaluated
70 * @return Value of the Qm function
71 */
72 double
73 Qm(const dftfe::Int m, const double phi);
74
75 /*
76 * @brief Function to compute derivative of the Qm(m,phi) function (defined above) with respect to phi.
77 * @param[in] m order of the spherical harmonic (m-quantum number) for
78 * which the derivative of Qm is to be evaluated
79 * @param[in] phi Azimuthal angle at which the derivative of Qm is to be
80 * evaluated
81 * @return Value of the derivative of Qm with respect to phi function
82 */
83 double
84 dQmDPhi(const dftfe::Int m, const double phi);
85
86 /*
87 * @brief Function to compute double derivative of the Qm(m,phi) function (defined above) with respect to phi.
88 * @param[in] m order of the spherical harmonic (m-quantum number) for
89 * which the double derivative of Qm is to be evaluated
90 * @param[in] phi Azimuthal angle at which the double derivative of Qm is
91 * to be evaluated
92 * @return Value of the double derivative of Qm with respect to phi function
93 */
94 double
95 d2QmDPhi2(const dftfe::Int m, const double phi);
96
97 /*
98 * @brief Function to compute the polar angle (theta) dependent
99 * part of the real spherical harmonics. Given the degree l and
100 * order m (i.e., l and m quantum numbers), this amounts to just the
101 * P_{l,|m|}, that is the associated Legendre function evaluated at
102 * |m|. It does not include any normalization constant or Condon-Shockley
103 * factor.
104 * @param[in] l degree of the spherical harmonic (l-quantum number)
105 * @param[in] m order of the spherical harmonic (m-quantum number)
106 * @param[in] theta Polar angle
107 * @return Value of the P_{l,|m|}
108 */
109 double
110 Plm(const dftfe::Int l, const dftfe::Int m, const double theta);
111
112 /*
113 * @brief Function to compute the derivative of the P_{l,|m|} function (defined above)
114 * with respect to the polar angle (theta).
115 * @param[in] l degree of the spherical harmonic (l-quantum number)
116 * @param[in] m order of the spherical harmonic (m-quantum number)
117 * @param[in] theta Polar angle
118 * @return Value of the derivative of P_{l,|m|} with respect to theta
119 */
120 double
121 dPlmDTheta(const dftfe::Int l, const dftfe::Int m, const double theta);
122
123 /*
124 * @brief Function to compute the double derivative of the P_{l,|m|} function (defined above)
125 * with respect to the polar angle (theta).
126 * @param[in] l degree of the spherical harmonic (l-quantum number)
127 * @param[in] m order of the spherical harmonic (m-quantum number)
128 * @param[in] theta Polar angle
129 * @return Value of the double derivative of P_{l,|m|} with respect to theta
130 */
131 double
132 d2PlmDTheta2(const dftfe::Int l, const dftfe::Int m, const double theta);
133
134 /*
135 * @brief Function to evaluate the real spherical harmonics YlmReal for a
136 * given degree (l), order (m), polar angle (theta), azimuthal
137 * angle (phi)
138 * @param[in] l degree of the spherical harmonic (l-quantum number)
139 * @param[in] m order of the spherical harmonic (m-quantum number)
140 * @param[in] theta Polar angle
141 * @param[in] phi Azimuthal angle
142 * @return Value of YlmReal
143 */
144 double
146 const dftfe::Int m,
147 const double theta,
148 const double phi);
149
150 /*
151 * @brief Function to evaluate the parial derivatives of the YlmReal function defined above with respect to
152 * polar angle (theta) and azimuthal angle (phi).
153 * @param[in] l degree of the spherical harmonic (l-quantum number)
154 * @param[in] m order of the spherical harmonic (m-quantum number)
155 * @param[in] theta Polar angle
156 * @param[in] phi Azimuthal angle
157 * @return Vector containing the partial derivatives of YlmReal with respect to theta and phi, in that order.
158 */
159 std::vector<double>
161 const dftfe::Int m,
162 const double theta,
163 const double phi);
164
165 /*
166 * @brief Function to evaluate the second-order parial derivatives of the YlmReal function defined above with respect to
167 * polar angle (theta) and azimuthal angle (phi).
168 * @param[in] l degree of the spherical harmonic (l-quantum number)
169 * @param[in] m order of the spherical harmonic (m-quantum number)
170 * @param[in] theta Polar angle
171 * @param[in] phi Azimuthal angle
172 * @return Vector containing the second-order partial derivatives of YlmReal with respect to theta and phi, in that order.
173 */
174 std::vector<double>
176 const dftfe::Int m,
177 const double theta,
178 const double phi);
179
180 /*
181 * @brief Function to evaluate the inverse of the Jacobian for the transform from cartesian to spherical coordinates
182 * @param[in] r Radius of the point
183 * @param[in] theta Polar angle of the point
184 * @param[in] phi Azimuthal angle of the point
185 * @return 2D Vector containing the inverse of the Jacobian
186 */
187 std::vector<std::vector<double>>
188 getJInv(const double r, const double theta, const double phi);
189 } // end of namespace sphUtils
190 } // end of namespace utils
191} // end of namespace dftfe
192#endif // DFTFE_SPHERICALFUNCTIONUTIL_H
Definition SphericalFunctionUtil.h:30
double d2PlmDTheta2(const dftfe::Int l, const dftfe::Int m, const double theta)
std::vector< double > dYlmReal(const dftfe::Int l, const dftfe::Int m, const double theta, const double phi)
double d2QmDPhi2(const dftfe::Int m, const double phi)
double dQmDPhi(const dftfe::Int m, const double phi)
std::vector< double > d2YlmReal(const dftfe::Int l, const dftfe::Int m, const double theta, const double phi)
double YlmReal(const dftfe::Int l, const dftfe::Int m, const double theta, const double phi)
double Plm(const dftfe::Int l, const dftfe::Int m, const double theta)
std::vector< std::vector< double > > getJInv(const double r, const double theta, const double phi)
double Qm(const dftfe::Int m, const double phi)
double dPlmDTheta(const dftfe::Int l, const dftfe::Int m, const double theta)
void convertCartesianToSpherical(const std::vector< double > &x, double &r, double &theta, double &phi, const double rTol, const double angleTol)
double Clm(const dftfe::Int l, const dftfe::Int m)
Definition Cell.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
std::int32_t Int
Definition TypeConfig.h:11