main.py 14 KB

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