m2m模型翻译
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

340 lines
11 KiB

6 months ago
  1. from sympy.core.expr import Expr
  2. from sympy.core.function import Function, ArgumentIndexError
  3. from sympy.core.numbers import I, pi
  4. from sympy.core.singleton import S
  5. from sympy.core.symbol import Dummy
  6. from sympy.core.sympify import sympify
  7. from sympy.functions import assoc_legendre
  8. from sympy.functions.combinatorial.factorials import factorial
  9. from sympy.functions.elementary.complexes import Abs
  10. from sympy.functions.elementary.exponential import exp
  11. from sympy.functions.elementary.miscellaneous import sqrt
  12. from sympy.functions.elementary.trigonometric import sin, cos, cot
  13. _x = Dummy("x")
  14. class Ynm(Function):
  15. r"""
  16. Spherical harmonics defined as
  17. .. math::
  18. Y_n^m(\theta, \varphi) := \sqrt{\frac{(2n+1)(n-m)!}{4\pi(n+m)!}}
  19. \exp(i m \varphi)
  20. \mathrm{P}_n^m\left(\cos(\theta)\right)
  21. Explanation
  22. ===========
  23. ``Ynm()`` gives the spherical harmonic function of order $n$ and $m$
  24. in $\theta$ and $\varphi$, $Y_n^m(\theta, \varphi)$. The four
  25. parameters are as follows: $n \geq 0$ an integer and $m$ an integer
  26. such that $-n \leq m \leq n$ holds. The two angles are real-valued
  27. with $\theta \in [0, \pi]$ and $\varphi \in [0, 2\pi]$.
  28. Examples
  29. ========
  30. >>> from sympy import Ynm, Symbol, simplify
  31. >>> from sympy.abc import n,m
  32. >>> theta = Symbol("theta")
  33. >>> phi = Symbol("phi")
  34. >>> Ynm(n, m, theta, phi)
  35. Ynm(n, m, theta, phi)
  36. Several symmetries are known, for the order:
  37. >>> Ynm(n, -m, theta, phi)
  38. (-1)**m*exp(-2*I*m*phi)*Ynm(n, m, theta, phi)
  39. As well as for the angles:
  40. >>> Ynm(n, m, -theta, phi)
  41. Ynm(n, m, theta, phi)
  42. >>> Ynm(n, m, theta, -phi)
  43. exp(-2*I*m*phi)*Ynm(n, m, theta, phi)
  44. For specific integers $n$ and $m$ we can evaluate the harmonics
  45. to more useful expressions:
  46. >>> simplify(Ynm(0, 0, theta, phi).expand(func=True))
  47. 1/(2*sqrt(pi))
  48. >>> simplify(Ynm(1, -1, theta, phi).expand(func=True))
  49. sqrt(6)*exp(-I*phi)*sin(theta)/(4*sqrt(pi))
  50. >>> simplify(Ynm(1, 0, theta, phi).expand(func=True))
  51. sqrt(3)*cos(theta)/(2*sqrt(pi))
  52. >>> simplify(Ynm(1, 1, theta, phi).expand(func=True))
  53. -sqrt(6)*exp(I*phi)*sin(theta)/(4*sqrt(pi))
  54. >>> simplify(Ynm(2, -2, theta, phi).expand(func=True))
  55. sqrt(30)*exp(-2*I*phi)*sin(theta)**2/(8*sqrt(pi))
  56. >>> simplify(Ynm(2, -1, theta, phi).expand(func=True))
  57. sqrt(30)*exp(-I*phi)*sin(2*theta)/(8*sqrt(pi))
  58. >>> simplify(Ynm(2, 0, theta, phi).expand(func=True))
  59. sqrt(5)*(3*cos(theta)**2 - 1)/(4*sqrt(pi))
  60. >>> simplify(Ynm(2, 1, theta, phi).expand(func=True))
  61. -sqrt(30)*exp(I*phi)*sin(2*theta)/(8*sqrt(pi))
  62. >>> simplify(Ynm(2, 2, theta, phi).expand(func=True))
  63. sqrt(30)*exp(2*I*phi)*sin(theta)**2/(8*sqrt(pi))
  64. We can differentiate the functions with respect
  65. to both angles:
  66. >>> from sympy import Ynm, Symbol, diff
  67. >>> from sympy.abc import n,m
  68. >>> theta = Symbol("theta")
  69. >>> phi = Symbol("phi")
  70. >>> diff(Ynm(n, m, theta, phi), theta)
  71. m*cot(theta)*Ynm(n, m, theta, phi) + sqrt((-m + n)*(m + n + 1))*exp(-I*phi)*Ynm(n, m + 1, theta, phi)
  72. >>> diff(Ynm(n, m, theta, phi), phi)
  73. I*m*Ynm(n, m, theta, phi)
  74. Further we can compute the complex conjugation:
  75. >>> from sympy import Ynm, Symbol, conjugate
  76. >>> from sympy.abc import n,m
  77. >>> theta = Symbol("theta")
  78. >>> phi = Symbol("phi")
  79. >>> conjugate(Ynm(n, m, theta, phi))
  80. (-1)**(2*m)*exp(-2*I*m*phi)*Ynm(n, m, theta, phi)
  81. To get back the well known expressions in spherical
  82. coordinates, we use full expansion:
  83. >>> from sympy import Ynm, Symbol, expand_func
  84. >>> from sympy.abc import n,m
  85. >>> theta = Symbol("theta")
  86. >>> phi = Symbol("phi")
  87. >>> expand_func(Ynm(n, m, theta, phi))
  88. sqrt((2*n + 1)*factorial(-m + n)/factorial(m + n))*exp(I*m*phi)*assoc_legendre(n, m, cos(theta))/(2*sqrt(pi))
  89. See Also
  90. ========
  91. Ynm_c, Znm
  92. References
  93. ==========
  94. .. [1] https://en.wikipedia.org/wiki/Spherical_harmonics
  95. .. [2] http://mathworld.wolfram.com/SphericalHarmonic.html
  96. .. [3] http://functions.wolfram.com/Polynomials/SphericalHarmonicY/
  97. .. [4] http://dlmf.nist.gov/14.30
  98. """
  99. @classmethod
  100. def eval(cls, n, m, theta, phi):
  101. n, m, theta, phi = [sympify(x) for x in (n, m, theta, phi)]
  102. # Handle negative index m and arguments theta, phi
  103. if m.could_extract_minus_sign():
  104. m = -m
  105. return S.NegativeOne**m * exp(-2*I*m*phi) * Ynm(n, m, theta, phi)
  106. if theta.could_extract_minus_sign():
  107. theta = -theta
  108. return Ynm(n, m, theta, phi)
  109. if phi.could_extract_minus_sign():
  110. phi = -phi
  111. return exp(-2*I*m*phi) * Ynm(n, m, theta, phi)
  112. # TODO Add more simplififcation here
  113. def _eval_expand_func(self, **hints):
  114. n, m, theta, phi = self.args
  115. rv = (sqrt((2*n + 1)/(4*pi) * factorial(n - m)/factorial(n + m)) *
  116. exp(I*m*phi) * assoc_legendre(n, m, cos(theta)))
  117. # We can do this because of the range of theta
  118. return rv.subs(sqrt(-cos(theta)**2 + 1), sin(theta))
  119. def fdiff(self, argindex=4):
  120. if argindex == 1:
  121. # Diff wrt n
  122. raise ArgumentIndexError(self, argindex)
  123. elif argindex == 2:
  124. # Diff wrt m
  125. raise ArgumentIndexError(self, argindex)
  126. elif argindex == 3:
  127. # Diff wrt theta
  128. n, m, theta, phi = self.args
  129. return (m * cot(theta) * Ynm(n, m, theta, phi) +
  130. sqrt((n - m)*(n + m + 1)) * exp(-I*phi) * Ynm(n, m + 1, theta, phi))
  131. elif argindex == 4:
  132. # Diff wrt phi
  133. n, m, theta, phi = self.args
  134. return I * m * Ynm(n, m, theta, phi)
  135. else:
  136. raise ArgumentIndexError(self, argindex)
  137. def _eval_rewrite_as_polynomial(self, n, m, theta, phi, **kwargs):
  138. # TODO: Make sure n \in N
  139. # TODO: Assert |m| <= n ortherwise we should return 0
  140. return self.expand(func=True)
  141. def _eval_rewrite_as_sin(self, n, m, theta, phi, **kwargs):
  142. return self.rewrite(cos)
  143. def _eval_rewrite_as_cos(self, n, m, theta, phi, **kwargs):
  144. # This method can be expensive due to extensive use of simplification!
  145. from sympy.simplify import simplify, trigsimp
  146. # TODO: Make sure n \in N
  147. # TODO: Assert |m| <= n ortherwise we should return 0
  148. term = simplify(self.expand(func=True))
  149. # We can do this because of the range of theta
  150. term = term.xreplace({Abs(sin(theta)):sin(theta)})
  151. return simplify(trigsimp(term))
  152. def _eval_conjugate(self):
  153. # TODO: Make sure theta \in R and phi \in R
  154. n, m, theta, phi = self.args
  155. return S.NegativeOne**m * self.func(n, -m, theta, phi)
  156. def as_real_imag(self, deep=True, **hints):
  157. # TODO: Handle deep and hints
  158. n, m, theta, phi = self.args
  159. re = (sqrt((2*n + 1)/(4*pi) * factorial(n - m)/factorial(n + m)) *
  160. cos(m*phi) * assoc_legendre(n, m, cos(theta)))
  161. im = (sqrt((2*n + 1)/(4*pi) * factorial(n - m)/factorial(n + m)) *
  162. sin(m*phi) * assoc_legendre(n, m, cos(theta)))
  163. return (re, im)
  164. def _eval_evalf(self, prec):
  165. # Note: works without this function by just calling
  166. # mpmath for Legendre polynomials. But using
  167. # the dedicated function directly is cleaner.
  168. from mpmath import mp, workprec
  169. n = self.args[0]._to_mpmath(prec)
  170. m = self.args[1]._to_mpmath(prec)
  171. theta = self.args[2]._to_mpmath(prec)
  172. phi = self.args[3]._to_mpmath(prec)
  173. with workprec(prec):
  174. res = mp.spherharm(n, m, theta, phi)
  175. return Expr._from_mpmath(res, prec)
  176. def Ynm_c(n, m, theta, phi):
  177. r"""
  178. Conjugate spherical harmonics defined as
  179. .. math::
  180. \overline{Y_n^m(\theta, \varphi)} := (-1)^m Y_n^{-m}(\theta, \varphi).
  181. Examples
  182. ========
  183. >>> from sympy import Ynm_c, Symbol, simplify
  184. >>> from sympy.abc import n,m
  185. >>> theta = Symbol("theta")
  186. >>> phi = Symbol("phi")
  187. >>> Ynm_c(n, m, theta, phi)
  188. (-1)**(2*m)*exp(-2*I*m*phi)*Ynm(n, m, theta, phi)
  189. >>> Ynm_c(n, m, -theta, phi)
  190. (-1)**(2*m)*exp(-2*I*m*phi)*Ynm(n, m, theta, phi)
  191. For specific integers $n$ and $m$ we can evaluate the harmonics
  192. to more useful expressions:
  193. >>> simplify(Ynm_c(0, 0, theta, phi).expand(func=True))
  194. 1/(2*sqrt(pi))
  195. >>> simplify(Ynm_c(1, -1, theta, phi).expand(func=True))
  196. sqrt(6)*exp(I*(-phi + 2*conjugate(phi)))*sin(theta)/(4*sqrt(pi))
  197. See Also
  198. ========
  199. Ynm, Znm
  200. References
  201. ==========
  202. .. [1] https://en.wikipedia.org/wiki/Spherical_harmonics
  203. .. [2] http://mathworld.wolfram.com/SphericalHarmonic.html
  204. .. [3] http://functions.wolfram.com/Polynomials/SphericalHarmonicY/
  205. """
  206. from sympy.functions.elementary.complexes import conjugate
  207. return conjugate(Ynm(n, m, theta, phi))
  208. class Znm(Function):
  209. r"""
  210. Real spherical harmonics defined as
  211. .. math::
  212. Z_n^m(\theta, \varphi) :=
  213. \begin{cases}
  214. \frac{Y_n^m(\theta, \varphi) + \overline{Y_n^m(\theta, \varphi)}}{\sqrt{2}} &\quad m > 0 \\
  215. Y_n^m(\theta, \varphi) &\quad m = 0 \\
  216. \frac{Y_n^m(\theta, \varphi) - \overline{Y_n^m(\theta, \varphi)}}{i \sqrt{2}} &\quad m < 0 \\
  217. \end{cases}
  218. which gives in simplified form
  219. .. math::
  220. Z_n^m(\theta, \varphi) =
  221. \begin{cases}
  222. \frac{Y_n^m(\theta, \varphi) + (-1)^m Y_n^{-m}(\theta, \varphi)}{\sqrt{2}} &\quad m > 0 \\
  223. Y_n^m(\theta, \varphi) &\quad m = 0 \\
  224. \frac{Y_n^m(\theta, \varphi) - (-1)^m Y_n^{-m}(\theta, \varphi)}{i \sqrt{2}} &\quad m < 0 \\
  225. \end{cases}
  226. Examples
  227. ========
  228. >>> from sympy import Znm, Symbol, simplify
  229. >>> from sympy.abc import n, m
  230. >>> theta = Symbol("theta")
  231. >>> phi = Symbol("phi")
  232. >>> Znm(n, m, theta, phi)
  233. Znm(n, m, theta, phi)
  234. For specific integers n and m we can evaluate the harmonics
  235. to more useful expressions:
  236. >>> simplify(Znm(0, 0, theta, phi).expand(func=True))
  237. 1/(2*sqrt(pi))
  238. >>> simplify(Znm(1, 1, theta, phi).expand(func=True))
  239. -sqrt(3)*sin(theta)*cos(phi)/(2*sqrt(pi))
  240. >>> simplify(Znm(2, 1, theta, phi).expand(func=True))
  241. -sqrt(15)*sin(2*theta)*cos(phi)/(4*sqrt(pi))
  242. See Also
  243. ========
  244. Ynm, Ynm_c
  245. References
  246. ==========
  247. .. [1] https://en.wikipedia.org/wiki/Spherical_harmonics
  248. .. [2] http://mathworld.wolfram.com/SphericalHarmonic.html
  249. .. [3] http://functions.wolfram.com/Polynomials/SphericalHarmonicY/
  250. """
  251. @classmethod
  252. def eval(cls, n, m, theta, phi):
  253. n, m, th, ph = [sympify(x) for x in (n, m, theta, phi)]
  254. if m.is_positive:
  255. zz = (Ynm(n, m, th, ph) + Ynm_c(n, m, th, ph)) / sqrt(2)
  256. return zz
  257. elif m.is_zero:
  258. return Ynm(n, m, th, ph)
  259. elif m.is_negative:
  260. zz = (Ynm(n, m, th, ph) - Ynm_c(n, m, th, ph)) / (sqrt(2)*I)
  261. return zz