main.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340
  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. def scattcoeffs_(x, m, nmax=-1, pl=-1, mp=False):
  34. if mp:
  35. from scattnlay_mp import mie_mp as mie_
  36. else:
  37. # from scattnlay_dp import mie_dp as mie_
  38. from scattnlay_mp import mie_mp as mie_
  39. mie = mie_()
  40. mie.SetLayersSize(x)
  41. mie.SetLayersIndex(m)
  42. mie.SetPECLayer(pl)
  43. mie.SetMaxTerms(nmax)
  44. mie.calcScattCoeffs()
  45. terms = mie.GetMaxTerms()
  46. a = mie.GetAn()
  47. b = mie.GetBn()
  48. return terms, a, b
  49. def scattcoeffs(x, m, nmax=-1, pl=-1, mp=False):
  50. """
  51. scattcoeffs(x, m[, nmax, pl, mp])
  52. Calculate the scattering coefficients required to calculate both the
  53. near- and far-field parameters.
  54. x: Size parameters (1D or 2D ndarray)
  55. m: Relative refractive indices (1D or 2D ndarray)
  56. nmax: Maximum number of multipolar expansion terms to be used for the
  57. calculations. Only use it if you know what you are doing, otherwise
  58. set this parameter to -1 and the function will calculate it.
  59. pl: Index of PEC layer. If there is none just send -1.
  60. mp: Use multiple (True) or double (False) precision.
  61. Returns: (terms, an, bn)
  62. with
  63. terms: Number of multipolar expansion terms used for the calculations
  64. an, bn: Complex scattering coefficients
  65. """
  66. if len(m.shape) != 1 and len(m.shape) != 2:
  67. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  68. if len(x.shape) == 1:
  69. if len(m.shape) == 1:
  70. return scattcoeffs_(x, m, nmax=nmax, pl=pl, mp=mp)
  71. else:
  72. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  73. elif len(x.shape) != 2:
  74. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  75. # Repeat the same m for all wavelengths
  76. if len(m.shape) == 1:
  77. m = np.repeat(m[np.newaxis, :], x.shape[0], axis=0)
  78. if nmax == -1:
  79. nstore = 0
  80. else:
  81. nstore = nmax
  82. terms = np.zeros((x.shape[0]), dtype=int)
  83. an = np.zeros((0, nstore), dtype=complex)
  84. bn = np.zeros((0, nstore), dtype=complex)
  85. for i, xi in enumerate(x):
  86. terms[i], a, b = scattcoeffs_(xi, m[i], nmax=nmax, pl=pl, mp=mp)
  87. if terms[i] > nstore:
  88. nstore = terms[i]
  89. an.resize((an.shape[0], nstore))
  90. bn.resize((bn.shape[0], nstore))
  91. an = np.vstack((an, a))
  92. bn = np.vstack((bn, b))
  93. return terms, an, bn
  94. #scattcoeffs()
  95. def expancoeffs(x, m, nmax=-1, pl=-1, mp=False):
  96. """
  97. expancoeffs(x, m[, nmax, pl, mp])
  98. Calculate the scattering coefficients required to calculate both the
  99. near- and far-field parameters.
  100. x: Size parameters (1D or 2D ndarray)
  101. m: Relative refractive indices (1D or 2D ndarray)
  102. nmax: Maximum number of multipolar expansion terms to be used for the
  103. calculations. Only use it if you know what you are doing, otherwise
  104. set this parameter to -1 and the function will calculate it.
  105. pl: Index of PEC layer. If there is none just send -1.
  106. mp: Use multiple (True) or double (False) precision.
  107. Returns: (terms, an, bn, cn, dn)
  108. with
  109. terms: Number of multipolar expansion terms used for the calculations
  110. an, bn, cn, dn: Complex expansion coefficients of each layer
  111. """
  112. if mp:
  113. from scattnlay_mp import expancoeffs as expancoeffs_
  114. else:
  115. from scattnlay_dp import expancoeffs as expancoeffs_
  116. if len(m.shape) != 1 and len(m.shape) != 2:
  117. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  118. if len(x.shape) == 1:
  119. if len(m.shape) == 1:
  120. return expancoeffs_(x, m, nmax=nmax, pl=pl)
  121. else:
  122. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  123. elif len(x.shape) != 2:
  124. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  125. # Repeat the same m for all wavelengths
  126. if len(m.shape) == 1:
  127. m = np.repeat(m[np.newaxis, :], x.shape[0], axis=0)
  128. if nmax == -1:
  129. nstore = 0
  130. else:
  131. nstore = nmax
  132. terms = np.zeros((x.shape[0]), dtype=int)
  133. an = np.zeros((0, x.shape[1], nstore), dtype=complex)
  134. bn = np.zeros((0, x.shape[1], nstore), dtype=complex)
  135. cn = np.zeros((0, x.shape[1], nstore), dtype=complex)
  136. dn = np.zeros((0, x.shape[1], nstore), dtype=complex)
  137. for i, xi in enumerate(x):
  138. terms[i], a, b, c, d = expancoeffs_(xi, m[i], nmax=nmax, pl=pl)
  139. if terms[i] > nstore:
  140. nstore = terms[i]
  141. an.resize((an.shape[0], an.shape[1], nstore))
  142. bn.resize((bn.shape[0], bn.shape[1], nstore))
  143. cn.resize((cn.shape[0], cn.shape[1], nstore))
  144. dn.resize((dn.shape[0], dn.shape[1], nstore))
  145. an = np.vstack((an, a))
  146. bn = np.vstack((bn, b))
  147. cn = np.vstack((cn, c))
  148. dn = np.vstack((dn, d))
  149. return terms, an, bn, cn, dn
  150. #expancoeffs()
  151. def scattnlay_(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
  152. if mp:
  153. from scattnlay_mp import mie_mp as mie_
  154. else:
  155. from scattnlay_dp import mie_dp as mie_
  156. mie = mie_()
  157. mie.SetLayersSize(x)
  158. mie.SetLayersIndex(m)
  159. mie.SetAngles(theta)
  160. mie.SetPECLayer(pl)
  161. mie.SetMaxTerms(nmax)
  162. mie.RunMieCalculation()
  163. Qext = mie.GetQext()
  164. Qsca = mie.GetQsca()
  165. Qabs = mie.GetQabs()
  166. Qbk = mie.GetQbk()
  167. Qpr = mie.GetQpr()
  168. g = mie.GetAsymmetryFactor()
  169. Albedo = mie.GetAlbedo()
  170. terms = mie.GetMaxTerms()
  171. S1 = mie.GetS1()
  172. S2 = mie.GetS2()
  173. return terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2
  174. def scattnlay(x, m, theta=np.zeros(0, dtype=float), nmax=-1, pl=-1, mp=False):
  175. """
  176. scattnlay(x, m[, theta, nmax, pl, mp])
  177. Calculate the actual scattering parameters and amplitudes.
  178. x: Size parameters (1D or 2D ndarray)
  179. m: Relative refractive indices (1D or 2D ndarray)
  180. theta: Scattering angles where the scattering amplitudes will be
  181. calculated (optional, 1D ndarray)
  182. nmax: Maximum number of multipolar expansion terms to be used for the
  183. calculations. Only use it if you know what you are doing.
  184. pl: Index of PEC layer. If there is none just send -1.
  185. mp: Use multiple (True) or double (False) precision.
  186. Returns: (terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2)
  187. with
  188. terms: Number of multipolar expansion terms used for the calculations
  189. Qext: Efficiency factor for extinction
  190. Qsca: Efficiency factor for scattering
  191. Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca)
  192. Qbk: Efficiency factor for backscattering
  193. Qpr: Efficiency factor for the radiation pressure
  194. g: Asymmetry factor (g = (Qext-Qpr)/Qsca)
  195. Albedo: Single scattering albedo (Albedo = Qsca/Qext)
  196. S1, S2: Complex scattering amplitudes
  197. """
  198. if len(m.shape) != 1 and len(m.shape) != 2:
  199. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  200. if len(x.shape) == 1:
  201. if len(m.shape) == 1:
  202. return scattnlay_(x, m, theta, nmax=nmax, pl=pl, mp=mp)
  203. else:
  204. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  205. elif len(x.shape) != 2:
  206. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  207. if len(theta.shape) != 1:
  208. raise ValueError('The scattering angles (theta) should be a 1-D NumPy array.')
  209. # Repeat the same m for all wavelengths
  210. if len(m.shape) == 1:
  211. m = np.repeat(m[np.newaxis, :], x.shape[0], axis=0)
  212. terms = np.zeros((x.shape[0]), dtype=int)
  213. Qext = np.zeros((x.shape[0]), dtype=float)
  214. Qsca = np.zeros((x.shape[0]), dtype=float)
  215. Qabs = np.zeros((x.shape[0]), dtype=float)
  216. Qbk = np.zeros((x.shape[0]), dtype=float)
  217. Qpr = np.zeros((x.shape[0]), dtype=float)
  218. g = np.zeros((x.shape[0]), dtype=float)
  219. Albedo = np.zeros((x.shape[0]), dtype=float)
  220. S1 = np.zeros((x.shape[0], theta.shape[0]), dtype=complex)
  221. S2 = np.zeros((x.shape[0], theta.shape[0]), dtype=complex)
  222. for i, xi in enumerate(x):
  223. terms[i], Qext[i], Qsca[i], Qabs[i], Qbk[i], Qpr[i], g[i], Albedo[i], S1[i], S2[i] = scattnlay_(xi, m[i], theta, nmax=nmax, pl=pl, mp=mp)
  224. return terms, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2
  225. #scattnlay()
  226. def fieldnlay_(x, m, xp, yp, zp, nmax=-1, pl=-1, mp=False):
  227. if mp:
  228. from scattnlay_mp import mie_mp as mie_
  229. else:
  230. from scattnlay_dp import mie_dp as mie_
  231. # from scattnlay_mp import mie_mp as mie_
  232. mie = mie_()
  233. mie.SetLayersSize(x)
  234. mie.SetLayersIndex(m)
  235. mie.SetPECLayer(pl)
  236. mie.SetMaxTerms(nmax)
  237. mie.SetFieldCoords(xp, yp, zp)
  238. mie.RunFieldCalculation()
  239. terms = mie.GetMaxTerms()
  240. E = mie.GetFieldE()
  241. H = mie.GetFieldH()
  242. return terms, E, H
  243. def fieldnlay(x, m, xp, yp, zp, nmax=-1, pl=-1, mp=False):
  244. """
  245. fieldnlay(x, m, xp, yp, zp[, nmax, pl, mp])
  246. Calculate the actual scattering parameters and amplitudes.
  247. x: Size parameters (1D or 2D ndarray)
  248. m: Relative refractive indices (1D or 2D ndarray)
  249. xp: Array containing all X coordinates to calculate the complex
  250. electric and magnetic fields (1D* ndarray)
  251. yp: Array containing all Y coordinates to calculate the complex
  252. electric and magnetic fields (1D* ndarray)
  253. zp: Array containing all Z coordinates to calculate the complex
  254. electric and magnetic fields (1D* ndarray)
  255. nmax: Maximum number of multipolar expansion terms to be used for the
  256. calculations. Only use it if you know what you are doing.
  257. pl: Index of PEC layer. If there is none just send -1.
  258. mp: Use multiple (True) or double (False) precision.
  259. Returns: (terms, E, H)
  260. with
  261. terms: Number of multipolar expansion terms used for the calculations
  262. E, H: Complex electric and magnetic field at the provided coordinates
  263. *Note: We assume that the coordinates are referred to the first wavelength
  264. (or structure) and correct it for the following ones
  265. """
  266. if len(m.shape) != 1 and len(m.shape) != 2:
  267. raise ValueError('The relative refractive index (m) should be a 1-D or 2-D NumPy array.')
  268. if len(x.shape) == 1:
  269. if len(m.shape) == 1:
  270. return fieldnlay_(x, m, xp, yp, zp, nmax=nmax, pl=pl, mp=mp)
  271. else:
  272. raise ValueError('The number of of dimensions for the relative refractive index (m) and for the size parameter (x) must be equal.')
  273. elif len(x.shape) != 2:
  274. raise ValueError('The size parameter (x) should be a 1-D or 2-D NumPy array.')
  275. # Repeat the same m for all wavelengths
  276. if len(m.shape) == 1:
  277. m = np.repeat(m[np.newaxis, :], x.shape[0], axis=0)
  278. terms = np.zeros((x.shape[0]), dtype=int)
  279. E = np.zeros((x.shape[0], xp.shape[0], 3), dtype=complex)
  280. H = np.zeros((x.shape[0], xp.shape[0], 3), dtype=complex)
  281. for i, xi in enumerate(x):
  282. # (2020/05/12) We assume that the coordinates are referred to the first wavelength
  283. # (or structure) and correct it for the following ones
  284. terms[i], E[i], H[i] = fieldnlay_(xi, m[i], xp*xi[-1]/x[0, -1], yp*xi[-1]/x[0, -1], zp*xi[-1]/x[0, -1], nmax=nmax, pl=pl, mp=mp)
  285. return terms, E, H
  286. #fieldnlay()