main.py 8.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219
  1. #!/usr/bin/env python
  2. # -*- coding: UTF-8 -*-
  3. #
  4. # Copyright (C) 2009-2019 Ovidio Peña Rodríguez <ovidio@bytesfall.com>
  5. # Copyright (C) 2013-2019 Konstantin Ladutenko <kostyfisik@gmail.com>
  6. #
  7. # This file is part of scattnlay
  8. #
  9. # This program is free software: you can redistribute it and/or modify
  10. # it under the terms of the GNU General Public License as published by
  11. # the Free Software Foundation, either version 3 of the License, or
  12. # (at your option) any later version.
  13. #
  14. # This program is distributed in the hope that it will be useful,
  15. # but WITHOUT ANY WARRANTY; without even the implied warranty of
  16. # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  17. # GNU General Public License for more details.
  18. #
  19. # The only additional remark is that we expect that all publications
  20. # describing work using this software, or all commercial products
  21. # using it, cite at least one of the following references:
  22. # [1] O. Peña and U. Pal, "Scattering of electromagnetic radiation by
  23. # a multilayered sphere," Computer Physics Communications,
  24. # vol. 180, Nov. 2009, pp. 2348-2354.
  25. # [2] K. Ladutenko, U. Pal, A. Rivera, and O. Peña-Rodríguez, "Mie
  26. # calculation of electromagnetic near-field for a multilayered
  27. # sphere," Computer Physics Communications, vol. 214, May 2017,
  28. # pp. 225-230.
  29. #
  30. # You should have received a copy of the GNU General Public License
  31. # along with this program. If not, see <http://www.gnu.org/licenses/>.
  32. import numpy as np
  33. import sys
  34. def scattcoeffs(x, m, nmax=-1, pl=-1, mp=False):
  35. """
  36. scattcoeffs(x, m[, nmax, pl])
  37. Calculate the scattering coefficients required to calculate both the
  38. near- and far-field parameters.
  39. x: Size parameters (1D or 2D ndarray)
  40. m: Relative refractive indices (1D or 2D ndarray)
  41. nmax: Maximum number of multipolar expansion terms to be used for the
  42. calculations. Only use it if you know what you are doing, otherwise
  43. set this parameter to -1 and the function will calculate it.
  44. pl: Index of PEC layer. If there is none just send -1.
  45. mp: Use multiple (True) or double (False) precision.
  46. Returns: (terms, an, bn)
  47. with
  48. terms: Number of multipolar expansion terms used for the calculations
  49. an, bn: Complex scattering coefficients
  50. """
  51. if mp:
  52. from scattnlay_mp import scattcoeffs, scattnlay, fieldnlay
  53. else:
  54. from scattnlay_dp import scattcoeffs, scattnlay, fieldnlay
  55. if len(m.shape) != 1 and len(m.shape) != 2:
  56. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  57. if len(x.shape) == 1:
  58. if len(m.shape) == 1:
  59. return scattcoeffs(x, m, nmax=nmax, pl=pl)
  60. else:
  61. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  62. elif len(x.shape) != 2:
  63. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  64. if nmax == -1:
  65. nstore = 0
  66. else:
  67. nstore = nmax
  68. terms = np.zeros((x.shape[0]), dtype=int)
  69. an = np.zeros((0, nstore), dtype=complex)
  70. bn = np.zeros((0, nstore), dtype=complex)
  71. for i, xi in enumerate(x):
  72. if len(m.shape) == 1:
  73. mi = m
  74. else:
  75. mi = m[i]
  76. terms[i], a, b = scattcoeffs(xi, mi, nmax=nmax, pl=pl)
  77. if terms[i] > nstore:
  78. nstore = terms[i]
  79. an.resize((an.shape[0], nstore))
  80. bn.resize((bn.shape[0], nstore))
  81. an = np.vstack((an, a))
  82. bn = np.vstack((bn, b))
  83. return terms, an, bn
  84. #scattcoeffs()
  85. def scattnlay(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
  86. """
  87. scattnlay(x, m[, theta, nmax, pl])
  88. Calculate the actual scattering parameters and amplitudes.
  89. x: Size parameters (1D or 2D ndarray)
  90. m: Relative refractive indices (1D or 2D ndarray)
  91. theta: Scattering angles where the scattering amplitudes will be
  92. calculated (optional, 1D ndarray)
  93. nmax: Maximum number of multipolar expansion terms to be used for the
  94. calculations. Only use it if you know what you are doing.
  95. pl: Index of PEC layer. If there is none just send -1.
  96. mp: Use multiple (True) or double (False) precision.
  97. Returns: (terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2)
  98. with
  99. terms: Number of multipolar expansion terms used for the calculations
  100. Qext: Efficiency factor for extinction
  101. Qsca: Efficiency factor for scattering
  102. Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca)
  103. Qbk: Efficiency factor for backscattering
  104. Qpr: Efficiency factor for the radiation pressure
  105. g: Asymmetry factor (g = (Qext-Qpr)/Qsca)
  106. Albedo: Single scattering albedo (Albedo = Qsca/Qext)
  107. S1, S2: Complex scattering amplitudes
  108. """
  109. if mp:
  110. from scattnlay_mp import scattcoeffs, scattnlay, fieldnlay
  111. else:
  112. from scattnlay_dp import scattcoeffs, scattnlay, fieldnlay
  113. if len(m.shape) != 1 and len(m.shape) != 2:
  114. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  115. if len(x.shape) == 1:
  116. if len(m.shape) == 1:
  117. return scattnlay(x, m, theta, nmax=nmax, pl=pl)
  118. else:
  119. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  120. elif len(x.shape) != 2:
  121. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  122. if len(theta.shape) != 1:
  123. raise ValueError('The scattering angles (theta) should be a 1-D NumPy array.')
  124. terms = np.zeros((x.shape[0]), dtype=int)
  125. Qext = np.zeros((x.shape[0]), dtype=float)
  126. Qsca = np.zeros((x.shape[0]), dtype=float)
  127. Qabs = np.zeros((x.shape[0]), dtype=float)
  128. Qbk = np.zeros((x.shape[0]), dtype=float)
  129. Qpr = np.zeros((x.shape[0]), dtype=float)
  130. g = np.zeros((x.shape[0]), dtype=float)
  131. Albedo = np.zeros((x.shape[0]), dtype=float)
  132. S1 = np.zeros((x.shape[0], theta.shape[0]), dtype=complex)
  133. S2 = np.zeros((x.shape[0], theta.shape[0]), dtype=complex)
  134. for i, xi in enumerate(x):
  135. if len(m.shape) == 1:
  136. mi = m
  137. else:
  138. mi = m[i]
  139. terms[i], Qext[i], Qsca[i], Qabs[i], Qbk[i], Qpr[i], g[i], Albedo[i], S1[i], S2[i] = scattnlay(xi, mi, theta, nmax=nmax, pl=pl)
  140. return terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2
  141. #scattnlay()
  142. def fieldnlay(x, m, xp, yp, zp, nmax=-1, pl=-1, mp=False):
  143. """
  144. fieldnlay(x, m, xp, yp, zp[, theta, nmax, pl])
  145. Calculate the actual scattering parameters and amplitudes.
  146. x: Size parameters (1D or 2D ndarray)
  147. m: Relative refractive indices (1D or 2D ndarray)
  148. xp: Array containing all X coordinates to calculate the complex
  149. electric and magnetic fields (1D ndarray)
  150. nmax: Maximum number of multipolar expansion terms to be used for the
  151. calculations. Only use it if you know what you are doing.
  152. pl: Index of PEC layer. If there is none just send -1.
  153. mp: Use multiple (True) or double (False) precision.
  154. Returns: (terms, E, H)
  155. with
  156. terms: Number of multipolar expansion terms used for the calculations
  157. E, H: Complex electric and magnetic field at the provided coordinates
  158. """
  159. if mp:
  160. from scattnlay_mp import scattcoeffs, scattnlay, fieldnlay
  161. else:
  162. from scattnlay_dp import scattcoeffs, scattnlay, fieldnlay
  163. if len(m.shape) != 1 and len(m.shape) != 2:
  164. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  165. if len(x.shape) == 1:
  166. if len(m.shape) == 1:
  167. return fieldnlay(x, m, xp, yp, zp, nmax=nmax, pl=pl)
  168. else:
  169. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  170. elif len(x.shape) != 2:
  171. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  172. terms = np.zeros((x.shape[0]), dtype=int)
  173. E = np.zeros((x.shape[0], xp.shape[0], 3), dtype=complex)
  174. H = np.zeros((x.shape[0], xp.shape[0], 3), dtype=complex)
  175. for i, xi in enumerate(x):
  176. if len(m.shape) == 1:
  177. mi = m
  178. else:
  179. mi = m[i]
  180. terms[i], E[i], H[i] = fieldnlay(xi, mi, xp, yp, zp, nmax=nmax, pl=pl)
  181. return terms, E, H
  182. #fieldnlay()