nmie.cc 47 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005
  1. //**********************************************************************************//
  2. // Copyright (C) 2009-2015 Ovidio Pena <ovidio@bytesfall.com> //
  3. // //
  4. // This file is part of scattnlay //
  5. // //
  6. // This program is free software: you can redistribute it and/or modify //
  7. // it under the terms of the GNU General Public License as published by //
  8. // the Free Software Foundation, either version 3 of the License, or //
  9. // (at your option) any later version. //
  10. // //
  11. // This program is distributed in the hope that it will be useful, //
  12. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  13. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  14. // GNU General Public License for more details. //
  15. // //
  16. // The only additional remark is that we expect that all publications //
  17. // describing work using this software, or all commercial products //
  18. // using it, cite the following reference: //
  19. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  20. // a multilayered sphere," Computer Physics Communications, //
  21. // vol. 180, Nov. 2009, pp. 2348-2354. //
  22. // //
  23. // You should have received a copy of the GNU General Public License //
  24. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  25. //**********************************************************************************//
  26. //**********************************************************************************//
  27. // This library implements the algorithm for a multilayered sphere described by: //
  28. // [1] W. Yang, "Improved recursive algorithm for light scattering by a //
  29. // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. //
  30. // //
  31. // You can find the description of all the used equations in: //
  32. // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  33. // a multilayered sphere," Computer Physics Communications, //
  34. // vol. 180, Nov. 2009, pp. 2348-2354. //
  35. // //
  36. // Hereinafter all equations numbers refer to [2] //
  37. //**********************************************************************************//
  38. #include <math.h>
  39. #include <stdlib.h>
  40. #include <stdio.h>
  41. #include "nmie.h"
  42. #define round(x) ((x) >= 0 ? (int)((x) + 0.5):(int)((x) - 0.5))
  43. const double PI=3.14159265358979323846;
  44. // light speed [m s-1]
  45. double const cc = 2.99792458e8;
  46. // assume non-magnetic (MU=MU0=const) [N A-2]
  47. double const mu = 4.0*PI*1.0e-7;
  48. // Calculate Nstop - equation (17)
  49. int Nstop(double xL) {
  50. int result;
  51. if (xL <= 8) {
  52. result = round(xL + 4*pow(xL, 1.0/3.0) + 1);
  53. } else if (xL <= 4200) {
  54. result = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  55. } else {
  56. result = round(xL + 4*pow(xL, 1.0/3.0) + 2);
  57. }
  58. return result;
  59. }
  60. //**********************************************************************************//
  61. int Nmax(int L, int fl, int pl,
  62. std::vector<double> x,
  63. std::vector<std::complex<double> > m) {
  64. int i, result, ri, riM1;
  65. result = Nstop(x[L - 1]);
  66. for (i = fl; i < L; i++) {
  67. if (i > pl) {
  68. ri = round(std::abs(x[i]*m[i]));
  69. } else {
  70. ri = 0;
  71. }
  72. if (result < ri) {
  73. result = ri;
  74. }
  75. if ((i > fl) && ((i - 1) > pl)) {
  76. riM1 = round(std::abs(x[i - 1]* m[i]));
  77. } else {
  78. riM1 = 0;
  79. }
  80. if (result < riM1) {
  81. result = riM1;
  82. }
  83. }
  84. return result + 15;
  85. }
  86. //**********************************************************************************//
  87. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  88. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  89. // //
  90. // Input parameters: //
  91. // z: Real argument to evaluate jn and h1n //
  92. // nmax: Maximum number of terms to calculate jn and h1n //
  93. // //
  94. // Output parameters: //
  95. // jn, h1n: Spherical Bessel and Hankel functions //
  96. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  97. // //
  98. // The implementation follows the algorithm by I.J. Thompson and A.R. Barnett, //
  99. // Comp. Phys. Comm. 47 (1987) 245-257. //
  100. // //
  101. // Complex spherical Bessel functions from n=0..nmax-1 for z in the upper half //
  102. // plane (Im(z) > -3). //
  103. // //
  104. // j[n] = j/n(z) Regular solution: j[0]=sin(z)/z //
  105. // j'[n] = d[j/n(z)]/dz //
  106. // h1[n] = h[0]/n(z) Irregular Hankel function: //
  107. // h1'[n] = d[h[0]/n(z)]/dz h1[0] = j0(z) + i*y0(z) //
  108. // = (sin(z)-i*cos(z))/z //
  109. // = -i*exp(i*z)/z //
  110. // Using complex CF1, and trigonometric forms for n=0 solutions. //
  111. //**********************************************************************************//
  112. int sbesjh(std::complex<double> z, int nmax, std::vector<std::complex<double> >& jn, std::vector<std::complex<double> >& jnp, std::vector<std::complex<double> >& h1n, std::vector<std::complex<double> >& h1np) {
  113. const int limit = 20000;
  114. double const accur = 1.0e-12;
  115. double const tm30 = 1e-30;
  116. int n;
  117. double absc;
  118. std::complex<double> zi, w;
  119. std::complex<double> pl, f, b, d, c, del, jn0, jndb, h1nldb, h1nbdb;
  120. absc = std::abs(std::real(z)) + std::abs(std::imag(z));
  121. if ((absc < accur) || (std::imag(z) < -3.0)) {
  122. return -1;
  123. }
  124. zi = 1.0/z;
  125. w = zi + zi;
  126. pl = double(nmax)*zi;
  127. f = pl + zi;
  128. b = f + f + zi;
  129. d = 0.0;
  130. c = f;
  131. for (n = 0; n < limit; n++) {
  132. d = b - d;
  133. c = b - 1.0/c;
  134. absc = std::abs(std::real(d)) + std::abs(std::imag(d));
  135. if (absc < tm30) {
  136. d = tm30;
  137. }
  138. absc = std::abs(std::real(c)) + std::abs(std::imag(c));
  139. if (absc < tm30) {
  140. c = tm30;
  141. }
  142. d = 1.0/d;
  143. del = d*c;
  144. f = f*del;
  145. b += w;
  146. absc = std::abs(std::real(del - 1.0)) + std::abs(std::imag(del - 1.0));
  147. if (absc < accur) {
  148. // We have obtained the desired accuracy
  149. break;
  150. }
  151. }
  152. if (absc > accur) {
  153. // We were not able to obtain the desired accuracy
  154. return -2;
  155. }
  156. jn[nmax - 1] = tm30;
  157. jnp[nmax - 1] = f*jn[nmax - 1];
  158. // Downward recursion to n=0 (N.B. Coulomb Functions)
  159. for (n = nmax - 2; n >= 0; n--) {
  160. jn[n] = pl*jn[n + 1] + jnp[n + 1];
  161. jnp[n] = pl*jn[n] - jn[n + 1];
  162. pl = pl - zi;
  163. }
  164. // Calculate the n=0 Bessel Functions
  165. jn0 = zi*std::sin(z);
  166. h1n[0] = std::exp(std::complex<double>(0.0, 1.0)*z)*zi*(-std::complex<double>(0.0, 1.0));
  167. h1np[0] = h1n[0]*(std::complex<double>(0.0, 1.0) - zi);
  168. // Rescale j[n], j'[n], converting to spherical Bessel functions.
  169. // Recur h1[n], h1'[n] as spherical Bessel functions.
  170. w = 1.0/jn[0];
  171. pl = zi;
  172. for (n = 0; n < nmax; n++) {
  173. jn[n] = jn0*(w*jn[n]);
  174. jnp[n] = jn0*(w*jnp[n]) - zi*jn[n];
  175. if (n != 0) {
  176. h1n[n] = (pl - zi)*h1n[n - 1] - h1np[n - 1];
  177. // check if hankel is increasing (upward stable)
  178. if (std::abs(h1n[n]) < std::abs(h1n[n - 1])) {
  179. jndb = z;
  180. h1nldb = h1n[n];
  181. h1nbdb = h1n[n - 1];
  182. }
  183. pl += zi;
  184. h1np[n] = -(pl*h1n[n]) + h1n[n - 1];
  185. }
  186. }
  187. // success
  188. return 0;
  189. }
  190. //**********************************************************************************//
  191. // This function calculates the spherical Bessel functions (bj and by) and the //
  192. // logarithmic derivative (bd) for a given complex value z. See pag. 87 B&H. //
  193. // //
  194. // Input parameters: //
  195. // z: Complex argument to evaluate bj, by and bd //
  196. // nmax: Maximum number of terms to calculate bj, by and bd //
  197. // //
  198. // Output parameters: //
  199. // bj, by: Spherical Bessel functions //
  200. // bd: Logarithmic derivative //
  201. //**********************************************************************************//
  202. void sphericalBessel(std::complex<double> z, int nmax, std::vector<std::complex<double> >& bj, std::vector<std::complex<double> >& by, std::vector<std::complex<double> >& bd) {
  203. std::vector<std::complex<double> > jn, jnp, h1n, h1np;
  204. jn.resize(nmax);
  205. jnp.resize(nmax);
  206. h1n.resize(nmax);
  207. h1np.resize(nmax);
  208. // TODO verify that the function succeeds
  209. int ifail = sbesjh(z, nmax, jn, jnp, h1n, h1np);
  210. for (int n = 0; n < nmax; n++) {
  211. bj[n] = jn[n];
  212. by[n] = (h1n[n] - jn[n])/std::complex<double>(0.0, 1.0);
  213. bd[n] = jnp[n]/jn[n] + 1.0/z;
  214. }
  215. }
  216. // external scattering field = incident + scattered
  217. // BH p.92 (4.37), 94 (4.45), 95 (4.50)
  218. // assume: medium is non-absorbing; refim = 0; Uabs = 0
  219. void fieldExt(int nmax, double Rho, double Phi, double Theta, std::vector<double> Pi, std::vector<double> Tau,
  220. std::vector<std::complex<double> > an, std::vector<std::complex<double> > bn,
  221. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  222. int i, n, n1;
  223. double rn;
  224. std::complex<double> ci, zn, xxip, encap;
  225. std::vector<std::complex<double> > vm3o1n, vm3e1n, vn3o1n, vn3e1n;
  226. vm3o1n.resize(3);
  227. vm3e1n.resize(3);
  228. vn3o1n.resize(3);
  229. vn3e1n.resize(3);
  230. std::vector<std::complex<double> > Ei, Hi, Es, Hs;
  231. Ei.resize(3);
  232. Hi.resize(3);
  233. Es.resize(3);
  234. Hs.resize(3);
  235. for (i = 0; i < 3; i++) {
  236. Ei[i] = std::complex<double>(0.0, 0.0);
  237. Hi[i] = std::complex<double>(0.0, 0.0);
  238. Es[i] = std::complex<double>(0.0, 0.0);
  239. Hs[i] = std::complex<double>(0.0, 0.0);
  240. }
  241. std::vector<std::complex<double> > bj, by, bd;
  242. bj.resize(nmax+1);
  243. by.resize(nmax+1);
  244. bd.resize(nmax+1);
  245. // Calculate spherical Bessel and Hankel functions
  246. sphericalBessel(Rho, nmax, bj, by, bd);
  247. ci = std::complex<double>(0.0, 1.0);
  248. for (n = 0; n < nmax; n++) {
  249. n1 = n + 1;
  250. rn = double(n + 1);
  251. zn = bj[n1] + ci*by[n1];
  252. xxip = Rho*(bj[n] + ci*by[n]) - rn*zn;
  253. vm3o1n[0] = std::complex<double>(0.0, 0.0);
  254. vm3o1n[1] = std::cos(Phi)*Pi[n]*zn;
  255. vm3o1n[2] = -std::sin(Phi)*Tau[n]*zn;
  256. vm3e1n[0] = std::complex<double>(0.0, 0.0);
  257. vm3e1n[1] = -std::sin(Phi)*Pi[n]*zn;
  258. vm3e1n[2] = -std::cos(Phi)*Tau[n]*zn;
  259. vn3o1n[0] = std::sin(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  260. vn3o1n[1] = std::sin(Phi)*Tau[n]*xxip/Rho;
  261. vn3o1n[2] = std::cos(Phi)*Pi[n]*xxip/Rho;
  262. vn3e1n[0] = std::cos(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  263. vn3e1n[1] = std::cos(Phi)*Tau[n]*xxip/Rho;
  264. vn3e1n[2] = -std::sin(Phi)*Pi[n]*xxip/Rho;
  265. // scattered field: BH p.94 (4.45)
  266. encap = std::pow(ci, rn)*(2.0*rn + 1.0)/(rn*rn + rn);
  267. for (i = 0; i < 3; i++) {
  268. Es[i] = Es[i] + encap*(ci*an[n]*vn3e1n[i] - bn[n]*vm3o1n[i]);
  269. Hs[i] = Hs[i] + encap*(ci*bn[n]*vn3o1n[i] + an[n]*vm3e1n[i]);
  270. }
  271. }
  272. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  273. // basis unit vectors = er, etheta, ephi
  274. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  275. Ei[0] = eifac*std::sin(Theta)*std::cos(Phi);
  276. Ei[1] = eifac*std::cos(Theta)*std::cos(Phi);
  277. Ei[2] = -eifac*std::sin(Phi);
  278. // magnetic field
  279. double hffact = 1.0/(cc*mu);
  280. for (i = 0; i < 3; i++) {
  281. Hs[i] = hffact*Hs[i];
  282. }
  283. // incident H field: BH p.26 (2.43), p.89 (4.21)
  284. std::complex<double> hffacta = hffact;
  285. std::complex<double> hifac = eifac*hffacta;
  286. Hi[0] = hifac*std::sin(Theta)*std::sin(Phi);
  287. Hi[1] = hifac*std::cos(Theta)*std::sin(Phi);
  288. Hi[2] = hifac*std::cos(Phi);
  289. for (i = 0; i < 3; i++) {
  290. // electric field E [V m-1] = EF*E0
  291. E[i] = Ei[i] + Es[i];
  292. H[i] = Hi[i] + Hs[i];
  293. }
  294. }
  295. // Calculate an - equation (5)
  296. std::complex<double> calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  297. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  298. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  299. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  300. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  301. return Num/Denom;
  302. }
  303. // Calculate bn - equation (6)
  304. std::complex<double> calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  305. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  306. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  307. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  308. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  309. return Num/Denom;
  310. }
  311. // Calculates S1 - equation (25a)
  312. std::complex<double> calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  313. double Pi, double Tau) {
  314. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  315. }
  316. // Calculates S2 - equation (25b) (it's the same as (25a), just switches Pi and Tau)
  317. std::complex<double> calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  318. double Pi, double Tau) {
  319. return calc_S1(n, an, bn, Tau, Pi);
  320. }
  321. //**********************************************************************************//
  322. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  323. // real argument (x). //
  324. // Equations (20a) - (21b) //
  325. // //
  326. // Input parameters: //
  327. // x: Real argument to evaluate Psi and Zeta //
  328. // nmax: Maximum number of terms to calculate Psi and Zeta //
  329. // //
  330. // Output parameters: //
  331. // Psi, Zeta: Riccati-Bessel functions //
  332. //**********************************************************************************//
  333. void calcPsiZeta(double x, int nmax,
  334. std::vector<std::complex<double> > D1,
  335. std::vector<std::complex<double> > D3,
  336. std::vector<std::complex<double> >& Psi,
  337. std::vector<std::complex<double> >& Zeta) {
  338. int n;
  339. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  340. Psi[0] = std::complex<double>(sin(x), 0);
  341. Zeta[0] = std::complex<double>(sin(x), -cos(x));
  342. for (n = 1; n <= nmax; n++) {
  343. Psi[n] = Psi[n - 1]*(n/x - D1[n - 1]);
  344. Zeta[n] = Zeta[n - 1]*(n/x - D3[n - 1]);
  345. }
  346. }
  347. //**********************************************************************************//
  348. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  349. // functions (D1 and D3) for a complex argument (z). //
  350. // Equations (16a), (16b) and (18a) - (18d) //
  351. // //
  352. // Input parameters: //
  353. // z: Complex argument to evaluate D1 and D3 //
  354. // nmax: Maximum number of terms to calculate D1 and D3 //
  355. // //
  356. // Output parameters: //
  357. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  358. //**********************************************************************************//
  359. void calcD1D3(std::complex<double> z, int nmax,
  360. std::vector<std::complex<double> >& D1,
  361. std::vector<std::complex<double> >& D3) {
  362. int n;
  363. std::complex<double> nz, PsiZeta;
  364. // Downward recurrence for D1 - equations (16a) and (16b)
  365. D1[nmax] = std::complex<double>(0.0, 0.0);
  366. for (n = nmax; n > 0; n--) {
  367. nz = double(n)/z;
  368. D1[n - 1] = nz - 1.0/(D1[n] + nz);
  369. }
  370. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  371. PsiZeta = 0.5*(1.0 - std::complex<double>(cos(2.0*z.real()), sin(2.0*z.real()))*exp(-2.0*z.imag()));
  372. D3[0] = std::complex<double>(0.0, 1.0);
  373. for (n = 1; n <= nmax; n++) {
  374. nz = double(n)/z;
  375. PsiZeta = PsiZeta*(nz - D1[n - 1])*(nz - D3[n - 1]);
  376. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta;
  377. }
  378. }
  379. //**********************************************************************************//
  380. // This function calculates Pi and Tau for all values of Theta. //
  381. // Equations (26a) - (26c) //
  382. // //
  383. // Input parameters: //
  384. // nmax: Maximum number of terms to calculate Pi and Tau //
  385. // nTheta: Number of scattering angles //
  386. // Theta: Array containing all the scattering angles where the scattering //
  387. // amplitudes will be calculated //
  388. // //
  389. // Output parameters: //
  390. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  391. //**********************************************************************************//
  392. void calcPiTau(int nmax, double Theta, std::vector<double>& Pi, std::vector<double>& Tau) {
  393. int n;
  394. //****************************************************//
  395. // Equations (26a) - (26c) //
  396. //****************************************************//
  397. // Initialize Pi and Tau
  398. Pi[0] = 1.0;
  399. Tau[0] = cos(Theta);
  400. // Calculate the actual values
  401. if (nmax > 1) {
  402. Pi[1] = 3*Tau[0]*Pi[0];
  403. Tau[1] = 2*Tau[0]*Pi[1] - 3*Pi[0];
  404. for (n = 2; n < nmax; n++) {
  405. Pi[n] = ((n + n + 1)*Tau[0]*Pi[n - 1] - (n + 1)*Pi[n - 2])/n;
  406. Tau[n] = (n + 1)*Tau[0]*Pi[n] - (n + 2)*Pi[n - 1];
  407. }
  408. }
  409. }
  410. //**********************************************************************************//
  411. // This function calculates the scattering coefficients required to calculate //
  412. // both the near- and far-field parameters. //
  413. // //
  414. // Input parameters: //
  415. // L: Number of layers //
  416. // pl: Index of PEC layer. If there is none just send -1 //
  417. // x: Array containing the size parameters of the layers [0..L-1] //
  418. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  419. // nmax: Maximum number of multipolar expansion terms to be used for the //
  420. // calculations. Only use it if you know what you are doing, otherwise //
  421. // set this parameter to -1 and the function will calculate it. //
  422. // //
  423. // Output parameters: //
  424. // an, bn: Complex scattering amplitudes //
  425. // //
  426. // Return value: //
  427. // Number of multipolar expansion terms used for the calculations //
  428. //**********************************************************************************//
  429. int ScattCoeffs(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int nmax,
  430. std::vector<std::complex<double> >& an, std::vector<std::complex<double> >& bn) {
  431. //************************************************************************//
  432. // Calculate the index of the first layer. It can be either 0 (default) //
  433. // or the index of the outermost PEC layer. In the latter case all layers //
  434. // below the PEC are discarded. //
  435. //************************************************************************//
  436. int fl = (pl > 0) ? pl : 0;
  437. if (nmax <= 0) {
  438. nmax = Nmax(L, fl, pl, x, m);
  439. }
  440. std::complex<double> z1, z2;
  441. std::complex<double> Num, Denom;
  442. std::complex<double> G1, G2;
  443. std::complex<double> Temp;
  444. int n, l;
  445. //**************************************************************************//
  446. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  447. // means that index = layer number - 1 or index = n - 1. The only exception //
  448. // are the arrays for representing D1, D3 and Q because they need a value //
  449. // for the index 0 (zero), hence it is important to consider this shift //
  450. // between different arrays. The change was done to optimize memory usage. //
  451. //**************************************************************************//
  452. // Allocate memory to the arrays
  453. std::vector<std::vector<std::complex<double> > > D1_mlxl, D1_mlxlM1;
  454. D1_mlxl.resize(L);
  455. D1_mlxlM1.resize(L);
  456. std::vector<std::vector<std::complex<double> > > D3_mlxl, D3_mlxlM1;
  457. D3_mlxl.resize(L);
  458. D3_mlxlM1.resize(L);
  459. std::vector<std::vector<std::complex<double> > > Q;
  460. Q.resize(L);
  461. std::vector<std::vector<std::complex<double> > > Ha, Hb;
  462. Ha.resize(L);
  463. Hb.resize(L);
  464. for (l = 0; l < L; l++) {
  465. D1_mlxl[l].resize(nmax + 1);
  466. D1_mlxlM1[l].resize(nmax + 1);
  467. D3_mlxl[l].resize(nmax + 1);
  468. D3_mlxlM1[l].resize(nmax + 1);
  469. Q[l].resize(nmax + 1);
  470. Ha[l].resize(nmax);
  471. Hb[l].resize(nmax);
  472. }
  473. an.resize(nmax);
  474. bn.resize(nmax);
  475. std::vector<std::complex<double> > D1XL, D3XL;
  476. D1XL.resize(nmax + 1);
  477. D3XL.resize(nmax + 1);
  478. std::vector<std::complex<double> > PsiXL, ZetaXL;
  479. PsiXL.resize(nmax + 1);
  480. ZetaXL.resize(nmax + 1);
  481. //*************************************************//
  482. // Calculate D1 and D3 for z1 in the first layer //
  483. //*************************************************//
  484. if (fl == pl) { // PEC layer
  485. for (n = 0; n <= nmax; n++) {
  486. D1_mlxl[fl][n] = std::complex<double>(0.0, -1.0);
  487. D3_mlxl[fl][n] = std::complex<double>(0.0, 1.0);
  488. }
  489. } else { // Regular layer
  490. z1 = x[fl]* m[fl];
  491. // Calculate D1 and D3
  492. calcD1D3(z1, nmax, D1_mlxl[fl], D3_mlxl[fl]);
  493. }
  494. //******************************************************************//
  495. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  496. //******************************************************************//
  497. for (n = 0; n < nmax; n++) {
  498. Ha[fl][n] = D1_mlxl[fl][n + 1];
  499. Hb[fl][n] = D1_mlxl[fl][n + 1];
  500. }
  501. //*****************************************************//
  502. // Iteration from the second layer to the last one (L) //
  503. //*****************************************************//
  504. for (l = fl + 1; l < L; l++) {
  505. //************************************************************//
  506. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  507. //************************************************************//
  508. z1 = x[l]*m[l];
  509. z2 = x[l - 1]*m[l];
  510. //Calculate D1 and D3 for z1
  511. calcD1D3(z1, nmax, D1_mlxl[l], D3_mlxl[l]);
  512. //Calculate D1 and D3 for z2
  513. calcD1D3(z2, nmax, D1_mlxlM1[l], D3_mlxlM1[l]);
  514. //*********************************************//
  515. //Calculate Q, Ha and Hb in the layers fl+1..L //
  516. //*********************************************//
  517. // Upward recurrence for Q - equations (19a) and (19b)
  518. Num = exp(-2.0*(z1.imag() - z2.imag()))*std::complex<double>(cos(-2.0*z2.real()) - exp(-2.0*z2.imag()), sin(-2.0*z2.real()));
  519. Denom = std::complex<double>(cos(-2.0*z1.real()) - exp(-2.0*z1.imag()), sin(-2.0*z1.real()));
  520. Q[l][0] = Num/Denom;
  521. for (n = 1; n <= nmax; n++) {
  522. Num = (z1*D1_mlxl[l][n] + double(n))*(double(n) - z1*D3_mlxl[l][n - 1]);
  523. Denom = (z2*D1_mlxlM1[l][n] + double(n))*(double(n) - z2*D3_mlxlM1[l][n - 1]);
  524. Q[l][n] = (((x[l - 1]*x[l - 1])/(x[l]*x[l])* Q[l][n - 1])*Num)/Denom;
  525. }
  526. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  527. for (n = 1; n <= nmax; n++) {
  528. //Ha
  529. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  530. G1 = -D1_mlxlM1[l][n];
  531. G2 = -D3_mlxlM1[l][n];
  532. } else {
  533. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[l][n]);
  534. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[l][n]);
  535. }
  536. Temp = Q[l][n]*G1;
  537. Num = (G2*D1_mlxl[l][n]) - (Temp*D3_mlxl[l][n]);
  538. Denom = G2 - Temp;
  539. Ha[l][n - 1] = Num/Denom;
  540. //Hb
  541. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  542. G1 = Hb[l - 1][n - 1];
  543. G2 = Hb[l - 1][n - 1];
  544. } else {
  545. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[l][n]);
  546. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[l][n]);
  547. }
  548. Temp = Q[l][n]*G1;
  549. Num = (G2*D1_mlxl[l][n]) - (Temp* D3_mlxl[l][n]);
  550. Denom = (G2- Temp);
  551. Hb[l][n - 1] = (Num/ Denom);
  552. }
  553. }
  554. //**************************************//
  555. //Calculate D1, D3, Psi and Zeta for XL //
  556. //**************************************//
  557. // Calculate D1XL and D3XL
  558. calcD1D3(x[L - 1], nmax, D1XL, D3XL);
  559. // Calculate PsiXL and ZetaXL
  560. calcPsiZeta(x[L - 1], nmax, D1XL, D3XL, PsiXL, ZetaXL);
  561. //*********************************************************************//
  562. // Finally, we calculate the scattering coefficients (an and bn) and //
  563. // the angular functions (Pi and Tau). Note that for these arrays the //
  564. // first layer is 0 (zero), in future versions all arrays will follow //
  565. // this convention to save memory. (13 Nov, 2014) //
  566. //*********************************************************************//
  567. for (n = 0; n < nmax; n++) {
  568. //********************************************************************//
  569. //Expressions for calculating an and bn coefficients are not valid if //
  570. //there is only one PEC layer (ie, for a simple PEC sphere). //
  571. //********************************************************************//
  572. if (pl < (L - 1)) {
  573. an[n] = calc_an(n + 1, x[L - 1], Ha[L - 1][n], m[L - 1], PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  574. bn[n] = calc_bn(n + 1, x[L - 1], Hb[L - 1][n], m[L - 1], PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  575. } else {
  576. an[n] = calc_an(n + 1, x[L - 1], std::complex<double>(0.0, 0.0), std::complex<double>(1.0, 0.0), PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  577. bn[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  578. }
  579. }
  580. return nmax;
  581. }
  582. //**********************************************************************************//
  583. // This function calculates the actual scattering parameters and amplitudes //
  584. // //
  585. // Input parameters: //
  586. // L: Number of layers //
  587. // pl: Index of PEC layer. If there is none just send -1 //
  588. // x: Array containing the size parameters of the layers [0..L-1] //
  589. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  590. // nTheta: Number of scattering angles //
  591. // Theta: Array containing all the scattering angles where the scattering //
  592. // amplitudes will be calculated //
  593. // nmax: Maximum number of multipolar expansion terms to be used for the //
  594. // calculations. Only use it if you know what you are doing, otherwise //
  595. // set this parameter to -1 and the function will calculate it //
  596. // //
  597. // Output parameters: //
  598. // Qext: Efficiency factor for extinction //
  599. // Qsca: Efficiency factor for scattering //
  600. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  601. // Qbk: Efficiency factor for backscattering //
  602. // Qpr: Efficiency factor for the radiation pressure //
  603. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  604. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  605. // S1, S2: Complex scattering amplitudes //
  606. // //
  607. // Return value: //
  608. // Number of multipolar expansion terms used for the calculations //
  609. //**********************************************************************************//
  610. int nMie(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m,
  611. int nTheta, std::vector<double> Theta, int nmax,
  612. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  613. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  614. int i, n, t;
  615. std::vector<std::complex<double> > an, bn;
  616. std::complex<double> Qbktmp;
  617. // Calculate scattering coefficients
  618. nmax = ScattCoeffs(L, pl, x, m, nmax, an, bn);
  619. std::vector<double> Pi, Tau;
  620. Pi.resize(nmax);
  621. Tau.resize(nmax);
  622. double x2 = x[L - 1]*x[L - 1];
  623. // Initialize the scattering parameters
  624. *Qext = 0;
  625. *Qsca = 0;
  626. *Qabs = 0;
  627. *Qbk = 0;
  628. Qbktmp = std::complex<double>(0.0, 0.0);
  629. *Qpr = 0;
  630. *g = 0;
  631. *Albedo = 0;
  632. // Initialize the scattering amplitudes
  633. for (t = 0; t < nTheta; t++) {
  634. S1[t] = std::complex<double>(0.0, 0.0);
  635. S2[t] = std::complex<double>(0.0, 0.0);
  636. }
  637. // By using downward recurrence we avoid loss of precision due to float rounding errors
  638. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  639. // http://en.wikipedia.org/wiki/Loss_of_significance
  640. for (i = nmax - 2; i >= 0; i--) {
  641. n = i + 1;
  642. // Equation (27)
  643. *Qext += (n + n + 1)*(an[i].real() + bn[i].real());
  644. // Equation (28)
  645. *Qsca += (n + n + 1)*(an[i].real()*an[i].real() + an[i].imag()*an[i].imag() + bn[i].real()*bn[i].real() + bn[i].imag()*bn[i].imag());
  646. // Equation (29) TODO We must check carefully this equation. If we
  647. // remove the typecast to double then the result changes. Which is
  648. // the correct one??? Ovidio (2014/12/10) With cast ratio will
  649. // give double, without cast (n + n + 1)/(n*(n + 1)) will be
  650. // rounded to integer. Tig (2015/02/24)
  651. *Qpr += ((n*(n + 2)/(n + 1))*((an[i]*std::conj(an[n]) + bn[i]*std::conj(bn[n])).real()) + ((double)(n + n + 1)/(n*(n + 1)))*(an[i]*std::conj(bn[i])).real());
  652. // Equation (33)
  653. Qbktmp = Qbktmp + (double)(n + n + 1)*(1 - 2*(n % 2))*(an[i]- bn[i]);
  654. //****************************************************//
  655. // Calculate the scattering amplitudes (S1 and S2) //
  656. // Equations (25a) - (25b) //
  657. //****************************************************//
  658. for (t = 0; t < nTheta; t++) {
  659. calcPiTau(nmax, Theta[t], Pi, Tau);
  660. S1[t] += calc_S1(n, an[i], bn[i], Pi[i], Tau[i]);
  661. S2[t] += calc_S2(n, an[i], bn[i], Pi[i], Tau[i]);
  662. }
  663. }
  664. *Qext = 2*(*Qext)/x2; // Equation (27)
  665. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  666. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  667. *Qabs = *Qext - *Qsca; // Equation (30)
  668. *Albedo = *Qsca / *Qext; // Equation (31)
  669. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  670. *Qbk = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  671. return nmax;
  672. }
  673. //**********************************************************************************//
  674. // This function is just a wrapper to call the full 'nMie' function with fewer //
  675. // parameters, it is here mainly for compatibility with older versions of the //
  676. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  677. // any limit for the maximum number of terms. //
  678. // //
  679. // Input parameters: //
  680. // L: Number of layers //
  681. // x: Array containing the size parameters of the layers [0..L-1] //
  682. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  683. // nTheta: Number of scattering angles //
  684. // Theta: Array containing all the scattering angles where the scattering //
  685. // amplitudes will be calculated //
  686. // //
  687. // Output parameters: //
  688. // Qext: Efficiency factor for extinction //
  689. // Qsca: Efficiency factor for scattering //
  690. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  691. // Qbk: Efficiency factor for backscattering //
  692. // Qpr: Efficiency factor for the radiation pressure //
  693. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  694. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  695. // S1, S2: Complex scattering amplitudes //
  696. // //
  697. // Return value: //
  698. // Number of multipolar expansion terms used for the calculations //
  699. //**********************************************************************************//
  700. int nMie(int L, std::vector<double> x, std::vector<std::complex<double> > m,
  701. int nTheta, std::vector<double> Theta,
  702. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  703. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  704. return nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  705. }
  706. //**********************************************************************************//
  707. // This function is just a wrapper to call the full 'nMie' function with fewer //
  708. // parameters, it is useful if you want to include a PEC layer but not a limit //
  709. // for the maximum number of terms. //
  710. // //
  711. // Input parameters: //
  712. // L: Number of layers //
  713. // pl: Index of PEC layer. If there is none just send -1 //
  714. // x: Array containing the size parameters of the layers [0..L-1] //
  715. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  716. // nTheta: Number of scattering angles //
  717. // Theta: Array containing all the scattering angles where the scattering //
  718. // amplitudes will be calculated //
  719. // //
  720. // Output parameters: //
  721. // Qext: Efficiency factor for extinction //
  722. // Qsca: Efficiency factor for scattering //
  723. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  724. // Qbk: Efficiency factor for backscattering //
  725. // Qpr: Efficiency factor for the radiation pressure //
  726. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  727. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  728. // S1, S2: Complex scattering amplitudes //
  729. // //
  730. // Return value: //
  731. // Number of multipolar expansion terms used for the calculations //
  732. //**********************************************************************************//
  733. int nMie(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m,
  734. int nTheta, std::vector<double> Theta,
  735. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  736. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  737. return nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  738. }
  739. //**********************************************************************************//
  740. // This function is just a wrapper to call the full 'nMie' function with fewer //
  741. // parameters, it is useful if you want to include a limit for the maximum number //
  742. // of terms but not a PEC layer. //
  743. // //
  744. // Input parameters: //
  745. // L: Number of layers //
  746. // x: Array containing the size parameters of the layers [0..L-1] //
  747. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  748. // nTheta: Number of scattering angles //
  749. // Theta: Array containing all the scattering angles where the scattering //
  750. // amplitudes will be calculated //
  751. // nmax: Maximum number of multipolar expansion terms to be used for the //
  752. // calculations. Only use it if you know what you are doing, otherwise //
  753. // set this parameter to -1 and the function will calculate it //
  754. // //
  755. // Output parameters: //
  756. // Qext: Efficiency factor for extinction //
  757. // Qsca: Efficiency factor for scattering //
  758. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  759. // Qbk: Efficiency factor for backscattering //
  760. // Qpr: Efficiency factor for the radiation pressure //
  761. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  762. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  763. // S1, S2: Complex scattering amplitudes //
  764. // //
  765. // Return value: //
  766. // Number of multipolar expansion terms used for the calculations //
  767. //**********************************************************************************//
  768. int nMie(int L, std::vector<double> x, std::vector<std::complex<double> > m,
  769. int nTheta, std::vector<double> Theta, int nmax,
  770. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  771. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  772. return nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  773. }
  774. //**********************************************************************************//
  775. // This function calculates complex electric and magnetic field in the surroundings //
  776. // and inside (TODO) the particle. //
  777. // //
  778. // Input parameters: //
  779. // L: Number of layers //
  780. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  781. // x: Array containing the size parameters of the layers [0..L-1] //
  782. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  783. // nmax: Maximum number of multipolar expansion terms to be used for the //
  784. // calculations. Only use it if you know what you are doing, otherwise //
  785. // set this parameter to 0 (zero) and the function will calculate it. //
  786. // ncoord: Number of coordinate points //
  787. // Coords: Array containing all coordinates where the complex electric and //
  788. // magnetic fields will be calculated //
  789. // //
  790. // Output parameters: //
  791. // E, H: Complex electric and magnetic field at the provided coordinates //
  792. // //
  793. // Return value: //
  794. // Number of multipolar expansion terms used for the calculations //
  795. //**********************************************************************************//
  796. int nField(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int nmax, int ncoord, std::vector<double> Xp, std::vector<double> Yp, std::vector<double> Zp, std::vector<std::vector<std::complex<double> > >& E, std::vector<std::vector<std::complex<double> > >& H) {
  797. int i, c;
  798. double Rho, Phi, Theta;
  799. std::vector<std::complex<double> > an, bn;
  800. // This array contains the fields in spherical coordinates
  801. std::vector<std::complex<double> > Es, Hs;
  802. Es.resize(3);
  803. Hs.resize(3);
  804. // Calculate scattering coefficients
  805. nmax = ScattCoeffs(L, pl, x, m, nmax, an, bn);
  806. std::vector<double> Pi, Tau;
  807. Pi.resize(nmax);
  808. Tau.resize(nmax);
  809. for (c = 0; c < ncoord; c++) {
  810. // Convert to spherical coordinates
  811. Rho = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  812. // Avoid convergence problems due to Rho too small
  813. if (Rho < 1e-5) {
  814. Rho = 1e-5;
  815. }
  816. //If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  817. if (Rho == 0.0) {
  818. Theta = 0.0;
  819. } else {
  820. Theta = acos(Zp[c]/Rho);
  821. }
  822. //If Xp=Yp=0 then Phi is undefined. Just set it to zero to zero to avoid problems
  823. if ((Xp[c] == 0.0) and (Yp[c] == 0.0)) {
  824. Phi = 0.0;
  825. } else {
  826. Phi = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  827. }
  828. calcPiTau(nmax, Theta, Pi, Tau);
  829. //*******************************************************//
  830. // external scattering field = incident + scattered //
  831. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  832. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  833. //*******************************************************//
  834. // Firstly the easiest case: the field outside the particle
  835. if (Rho >= x[L - 1]) {
  836. fieldExt(nmax, Rho, Phi, Theta, Pi, Tau, an, bn, Es, Hs);
  837. } else {
  838. // TODO, for now just set all the fields to zero
  839. for (i = 0; i < 3; i++) {
  840. Es[i] = std::complex<double>(0.0, 0.0);
  841. Hs[i] = std::complex<double>(0.0, 0.0);
  842. }
  843. }
  844. //Now, convert the fields back to cartesian coordinates
  845. E[c][0] = std::sin(Theta)*std::cos(Phi)*Es[0] + std::cos(Theta)*std::cos(Phi)*Es[1] - std::sin(Phi)*Es[2];
  846. E[c][1] = std::sin(Theta)*std::sin(Phi)*Es[0] + std::cos(Theta)*std::sin(Phi)*Es[1] + std::cos(Phi)*Es[2];
  847. E[c][2] = std::cos(Theta)*Es[0] - std::sin(Theta)*Es[1];
  848. H[c][0] = std::sin(Theta)*std::cos(Phi)*Hs[0] + std::cos(Theta)*std::cos(Phi)*Hs[1] - std::sin(Phi)*Hs[2];
  849. H[c][1] = std::sin(Theta)*std::sin(Phi)*Hs[0] + std::cos(Theta)*std::sin(Phi)*Hs[1] + std::cos(Phi)*Hs[2];
  850. H[c][2] = std::cos(Theta)*Hs[0] - std::sin(Theta)*Hs[1];
  851. }
  852. return nmax;
  853. }