main.py 14 KB

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