nmie.cc 40 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803
  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. // Calculate Nstop - equation (17)
  44. int Nstop(double xL) {
  45. int result;
  46. if (xL <= 8) {
  47. result = round(xL + 4*pow(xL, 1/3) + 1);
  48. } else if (xL <= 4200) {
  49. result = round(xL + 4.05*pow(xL, 1/3) + 2);
  50. } else {
  51. result = round(xL + 4*pow(xL, 1/3) + 2);
  52. }
  53. return result;
  54. }
  55. //**********************************************************************************//
  56. int Nmax(int L, int fl, int pl,
  57. std::vector<double> x,
  58. std::vector<std::complex<double> > m) {
  59. int i, result, ri, riM1;
  60. result = Nstop(x[L - 1]);
  61. for (i = fl; i < L; i++) {
  62. if (i > pl) {
  63. ri = round(std::abs(x[i]*m[i]));
  64. } else {
  65. ri = 0;
  66. }
  67. if (result < ri) {
  68. result = ri;
  69. }
  70. if ((i > fl) && ((i - 1) > pl)) {
  71. riM1 = round(std::abs(x[i - 1]* m[i]));
  72. } else {
  73. riM1 = 0;
  74. }
  75. if (result < riM1) {
  76. result = riM1;
  77. }
  78. }
  79. return result + 15;
  80. }
  81. //**********************************************************************************//
  82. // This function calculates the spherical Bessel functions (jn and hn) for a given //
  83. // value of z. //
  84. // //
  85. // Input parameters: //
  86. // z: Real argument to evaluate jn and hn //
  87. // n_max: Maximum number of terms to calculate jn and hn //
  88. // //
  89. // Output parameters: //
  90. // jn, hn: Spherical Bessel functions (complex) //
  91. //**********************************************************************************//
  92. void sphericalBessel(std::complex<double> r, int n_max, std::vector<std::complex<double> > &j, std::vector<std::complex<double> > &h) {
  93. int n;
  94. j[0] = sin(r)/r;
  95. j[1] = sin(r)/r/r - cos(r)/r;
  96. h[0] = -cos(r)/r;
  97. h[1] = -cos(r)/r/r - sin(r)/r;
  98. for (n = 2; n < n_max; n++) {
  99. j[n] = double(n + n + 1)*j[n - 1]/r - j[n - 2];
  100. h[n] = double(n + n + 1)*h[n - 1]/r - h[n - 2];
  101. }
  102. }
  103. //**********************************************************************************//
  104. // This function calculates the spherical Bessel functions (jn and hn) for a given //
  105. // value of r. //
  106. // //
  107. // Input parameters: //
  108. // r: Real argument to evaluate jn and hn //
  109. // n_max: Maximum number of terms to calculate jn and hn //
  110. // //
  111. // Output parameters: //
  112. // jn, hn: Spherical Bessel functions (double) //
  113. //**********************************************************************************//
  114. void sphericalBessel(double r, int n_max, std::vector<double> &j, std::vector<double> &h) {
  115. int n;
  116. j[0] = sin(r)/r;
  117. j[1] = sin(r)/r/r - cos(r)/r;
  118. h[0] = -cos(r)/r;
  119. h[1] = -cos(r)/r/r - sin(r)/r;
  120. for (n = 2; n < n_max; n++) {
  121. j[n] = double(n + n + 1)*j[n - 1]/r - j[n - 2];
  122. h[n] = double(n + n + 1)*h[n - 1]/r - h[n - 2];
  123. }
  124. }
  125. // Calculate an - equation (5)
  126. std::complex<double> calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  127. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  128. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  129. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  130. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  131. return Num/Denom;
  132. }
  133. // Calculate bn - equation (6)
  134. std::complex<double> calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  135. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  136. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  137. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  138. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  139. return Num/Denom;
  140. }
  141. // Calculates S1 - equation (25a)
  142. std::complex<double> calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  143. double Pi, double Tau) {
  144. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  145. }
  146. // Calculates S2 - equation (25b) (it's the same as (25a), just switches Pi and Tau)
  147. std::complex<double> calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  148. double Pi, double Tau) {
  149. return calc_S1(n, an, bn, Tau, Pi);
  150. }
  151. //**********************************************************************************//
  152. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  153. // real argument (x). //
  154. // Equations (20a) - (21b) //
  155. // //
  156. // Input parameters: //
  157. // x: Real argument to evaluate Psi and Zeta //
  158. // n_max: Maximum number of terms to calculate Psi and Zeta //
  159. // //
  160. // Output parameters: //
  161. // Psi, Zeta: Riccati-Bessel functions //
  162. //**********************************************************************************//
  163. void calcPsiZeta(double x, int n_max,
  164. std::vector<std::complex<double> > D1,
  165. std::vector<std::complex<double> > D3,
  166. std::vector<std::complex<double> > &Psi,
  167. std::vector<std::complex<double> > &Zeta) {
  168. int n;
  169. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  170. Psi[0] = std::complex<double>(sin(x), 0);
  171. Zeta[0] = std::complex<double>(sin(x), -cos(x));
  172. for (n = 1; n <= n_max; n++) {
  173. Psi[n] = Psi[n - 1]*(n/x - D1[n - 1]);
  174. Zeta[n] = Zeta[n - 1]*(n/x - D3[n - 1]);
  175. }
  176. }
  177. //**********************************************************************************//
  178. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  179. // functions (D1 and D3) for a complex argument (z). //
  180. // Equations (16a), (16b) and (18a) - (18d) //
  181. // //
  182. // Input parameters: //
  183. // z: Complex argument to evaluate D1 and D3 //
  184. // n_max: Maximum number of terms to calculate D1 and D3 //
  185. // //
  186. // Output parameters: //
  187. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  188. //**********************************************************************************//
  189. void calcD1D3(std::complex<double> z, int n_max,
  190. std::vector<std::complex<double> > &D1,
  191. std::vector<std::complex<double> > &D3) {
  192. int n;
  193. std::vector<std::complex<double> > PsiZeta;
  194. PsiZeta.resize(n_max + 1);
  195. // Downward recurrence for D1 - equations (16a) and (16b)
  196. D1[n_max] = std::complex<double>(0.0, 0.0);
  197. for (n = n_max; n > 0; n--) {
  198. D1[n - 1] = double(n)/z - 1.0/(D1[n] + double(n)/z);
  199. }
  200. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  201. PsiZeta[0] = 0.5*(1.0 - std::complex<double>(cos(2.0*z.real()), sin(2.0*z.real()))*exp(-2.0*z.imag()));
  202. D3[0] = std::complex<double>(0.0, 1.0);
  203. for (n = 1; n <= n_max; n++) {
  204. PsiZeta[n] = PsiZeta[n - 1]*(double(n)/z - D1[n - 1])*(double(n)/z- D3[n - 1]);
  205. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta[n];
  206. }
  207. }
  208. //**********************************************************************************//
  209. // This function calculates Pi and Tau for all values of Theta. //
  210. // Equations (26a) - (26c) //
  211. // //
  212. // Input parameters: //
  213. // n_max: Maximum number of terms to calculate Pi and Tau //
  214. // nTheta: Number of scattering angles //
  215. // Theta: Array containing all the scattering angles where the scattering //
  216. // amplitudes will be calculated //
  217. // //
  218. // Output parameters: //
  219. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  220. //**********************************************************************************//
  221. void calcPiTau(int n_max, int nTheta, std::vector<double> Theta,
  222. std::vector< std::vector<double> > &Pi,
  223. std::vector< std::vector<double> > &Tau) {
  224. int n, t;
  225. for (n = 0; n < n_max; n++) {
  226. //****************************************************//
  227. // Equations (26a) - (26c) //
  228. //****************************************************//
  229. for (t = 0; t < nTheta; t++) {
  230. if (n == 0) {
  231. // Initialize Pi and Tau
  232. Pi[n][t] = 1.0;
  233. Tau[n][t] = (n + 1)*cos(Theta[t]);
  234. } else {
  235. // Calculate the actual values
  236. Pi[n][t] = ((n == 1) ? ((n + n + 1)*cos(Theta[t])*Pi[n - 1][t]/n)
  237. : (((n + n + 1)*cos(Theta[t])*Pi[n - 1][t] - (n + 1)*Pi[n - 2][t])/n));
  238. Tau[n][t] = (n + 1)*cos(Theta[t])*Pi[n][t] - (n + 2)*Pi[n - 1][t];
  239. }
  240. }
  241. }
  242. }
  243. //**********************************************************************************//
  244. // This function calculates the scattering coefficients required to calculate //
  245. // both the near- and far-field parameters. //
  246. // //
  247. // Input parameters: //
  248. // L: Number of layers //
  249. // pl: Index of PEC layer. If there is none just send -1 //
  250. // x: Array containing the size parameters of the layers [0..L-1] //
  251. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  252. // n_max: Maximum number of multipolar expansion terms to be used for the //
  253. // calculations. Only used if you know what you are doing, otherwise set //
  254. // this parameter to -1 and the function will calculate it. //
  255. // //
  256. // Output parameters: //
  257. // an, bn: Complex scattering amplitudes //
  258. // //
  259. // Return value: //
  260. // Number of multipolar expansion terms used for the calculations //
  261. //**********************************************************************************//
  262. int ScattCoeffs(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int n_max,
  263. std::vector<std::complex<double> > &an, std::vector<std::complex<double> > &bn) {
  264. //************************************************************************//
  265. // Calculate the index of the first layer. It can be either 0 (default) //
  266. // or the index of the outermost PEC layer. In the latter case all layers //
  267. // below the PEC are discarded. //
  268. //************************************************************************//
  269. int fl = (pl > 0) ? pl : 0;
  270. if (n_max <= 0) {
  271. n_max = Nmax(L, fl, pl, x, m);
  272. }
  273. std::complex<double> z1, z2;
  274. std::complex<double> Num, Denom;
  275. std::complex<double> G1, G2;
  276. std::complex<double> Temp;
  277. int n, l;
  278. //**************************************************************************//
  279. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  280. // means that index = layer number - 1 or index = n - 1. The only exception //
  281. // are the arrays for representing D1, D3 and Q because they need a value //
  282. // for the index 0 (zero), hence it is important to consider this shift //
  283. // between different arrays. The change was done to optimize memory usage. //
  284. //**************************************************************************//
  285. // Allocate memory to the arrays
  286. std::vector<std::vector<std::complex<double> > > D1_mlxl, D1_mlxlM1;
  287. D1_mlxl.resize(L);
  288. D1_mlxlM1.resize(L);
  289. std::vector<std::vector<std::complex<double> > > D3_mlxl, D3_mlxlM1;
  290. D3_mlxl.resize(L);
  291. D3_mlxlM1.resize(L);
  292. std::vector<std::vector<std::complex<double> > > Q;
  293. Q.resize(L);
  294. std::vector<std::vector<std::complex<double> > > Ha, Hb;
  295. Ha.resize(L);
  296. Hb.resize(L);
  297. for (l = 0; l < L; l++) {
  298. D1_mlxl[l].resize(n_max + 1);
  299. D1_mlxlM1[l].resize(n_max + 1);
  300. D3_mlxl[l].resize(n_max + 1);
  301. D3_mlxlM1[l].resize(n_max + 1);
  302. Q[l].resize(n_max + 1);
  303. Ha[l].resize(n_max);
  304. Hb[l].resize(n_max);
  305. }
  306. an.resize(n_max);
  307. bn.resize(n_max);
  308. std::vector<std::complex<double> > D1XL, D3XL;
  309. D1XL.resize(n_max + 1);
  310. D3XL.resize(n_max + 1);
  311. std::vector<std::complex<double> > PsiXL, ZetaXL;
  312. PsiXL.resize(n_max + 1);
  313. ZetaXL.resize(n_max + 1);
  314. //*************************************************//
  315. // Calculate D1 and D3 for z1 in the first layer //
  316. //*************************************************//
  317. if (fl == pl) { // PEC layer
  318. for (n = 0; n <= n_max; n++) {
  319. D1_mlxl[fl][n] = std::complex<double>(0.0, -1.0);
  320. D3_mlxl[fl][n] = std::complex<double>(0.0, 1.0);
  321. }
  322. } else { // Regular layer
  323. z1 = x[fl]* m[fl];
  324. // Calculate D1 and D3
  325. calcD1D3(z1, n_max, D1_mlxl[fl], D3_mlxl[fl]);
  326. }
  327. //******************************************************************//
  328. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  329. //******************************************************************//
  330. for (n = 0; n < n_max; n++) {
  331. Ha[fl][n] = D1_mlxl[fl][n + 1];
  332. Hb[fl][n] = D1_mlxl[fl][n + 1];
  333. }
  334. //*****************************************************//
  335. // Iteration from the second layer to the last one (L) //
  336. //*****************************************************//
  337. for (l = fl + 1; l < L; l++) {
  338. //************************************************************//
  339. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  340. //************************************************************//
  341. z1 = x[l]*m[l];
  342. z2 = x[l - 1]*m[l];
  343. //Calculate D1 and D3 for z1
  344. calcD1D3(z1, n_max, D1_mlxl[l], D3_mlxl[l]);
  345. //Calculate D1 and D3 for z2
  346. calcD1D3(z2, n_max, D1_mlxlM1[l], D3_mlxlM1[l]);
  347. //*********************************************//
  348. //Calculate Q, Ha and Hb in the layers fl+1..L //
  349. //*********************************************//
  350. // Upward recurrence for Q - equations (19a) and (19b)
  351. 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()));
  352. Denom = std::complex<double>(cos(-2.0*z1.real()) - exp(-2.0*z1.imag()), sin(-2.0*z1.real()));
  353. Q[l][0] = Num/Denom;
  354. for (n = 1; n <= n_max; n++) {
  355. Num = (z1*D1_mlxl[l][n] + double(n))*(double(n) - z1*D3_mlxl[l][n - 1]);
  356. Denom = (z2*D1_mlxlM1[l][n] + double(n))*(double(n) - z2*D3_mlxlM1[l][n - 1]);
  357. Q[l][n] = (((x[l - 1]*x[l - 1])/(x[l]*x[l])* Q[l][n - 1])*Num)/Denom;
  358. }
  359. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  360. for (n = 1; n <= n_max; n++) {
  361. //Ha
  362. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  363. G1 = -D1_mlxlM1[l][n];
  364. G2 = -D3_mlxlM1[l][n];
  365. } else {
  366. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[l][n]);
  367. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[l][n]);
  368. }
  369. Temp = Q[l][n]*G1;
  370. Num = (G2*D1_mlxl[l][n]) - (Temp*D3_mlxl[l][n]);
  371. Denom = G2 - Temp;
  372. Ha[l][n - 1] = Num/Denom;
  373. //Hb
  374. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  375. G1 = Hb[l - 1][n - 1];
  376. G2 = Hb[l - 1][n - 1];
  377. } else {
  378. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[l][n]);
  379. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[l][n]);
  380. }
  381. Temp = Q[l][n]*G1;
  382. Num = (G2*D1_mlxl[l][n]) - (Temp* D3_mlxl[l][n]);
  383. Denom = (G2- Temp);
  384. Hb[l][n - 1] = (Num/ Denom);
  385. }
  386. }
  387. //**************************************//
  388. //Calculate D1, D3, Psi and Zeta for XL //
  389. //**************************************//
  390. // Calculate D1XL and D3XL
  391. calcD1D3(x[L - 1], n_max, D1XL, D3XL);
  392. // Calculate PsiXL and ZetaXL
  393. calcPsiZeta(x[L - 1], n_max, D1XL, D3XL, PsiXL, ZetaXL);
  394. //*********************************************************************//
  395. // Finally, we calculate the scattering coefficients (an and bn) and //
  396. // the angular functions (Pi and Tau). Note that for these arrays the //
  397. // first layer is 0 (zero), in future versions all arrays will follow //
  398. // this convention to save memory. (13 Nov, 2014) //
  399. //*********************************************************************//
  400. for (n = 0; n < n_max; n++) {
  401. //********************************************************************//
  402. //Expressions for calculating an and bn coefficients are not valid if //
  403. //there is only one PEC layer (ie, for a simple PEC sphere). //
  404. //********************************************************************//
  405. if (pl < (L - 1)) {
  406. 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]);
  407. 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]);
  408. } else {
  409. 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]);
  410. bn[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  411. }
  412. }
  413. return n_max;
  414. }
  415. //**********************************************************************************//
  416. // This function calculates the actual scattering parameters and amplitudes //
  417. // //
  418. // Input parameters: //
  419. // L: Number of layers //
  420. // pl: Index of PEC layer. If there is none just send -1 //
  421. // x: Array containing the size parameters of the layers [0..L-1] //
  422. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  423. // nTheta: Number of scattering angles //
  424. // Theta: Array containing all the scattering angles where the scattering //
  425. // amplitudes will be calculated //
  426. // n_max: Maximum number of multipolar expansion terms to be used for the //
  427. // calculations. Only used if you know what you are doing, otherwise set //
  428. // this parameter to -1 and the function will calculate it //
  429. // //
  430. // Output parameters: //
  431. // Qext: Efficiency factor for extinction //
  432. // Qsca: Efficiency factor for scattering //
  433. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  434. // Qbk: Efficiency factor for backscattering //
  435. // Qpr: Efficiency factor for the radiation pressure //
  436. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  437. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  438. // S1, S2: Complex scattering amplitudes //
  439. // //
  440. // Return value: //
  441. // Number of multipolar expansion terms used for the calculations //
  442. //**********************************************************************************//
  443. int nMie(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m,
  444. int nTheta, std::vector<double> Theta, int n_max,
  445. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  446. std::vector<std::complex<double> > &S1, std::vector<std::complex<double> > &S2) {
  447. int i, n, t;
  448. std::vector<std::complex<double> > an, bn;
  449. std::complex<double> Qbktmp;
  450. // Calculate scattering coefficients
  451. n_max = ScattCoeffs(L, pl, x, m, n_max, an, bn);
  452. std::vector< std::vector<double> > Pi;
  453. Pi.resize(n_max);
  454. std::vector< std::vector<double> > Tau;
  455. Tau.resize(n_max);
  456. for (n = 0; n < n_max; n++) {
  457. Pi[n].resize(nTheta);
  458. Tau[n].resize(nTheta);
  459. }
  460. calcPiTau(n_max, nTheta, Theta, Pi, Tau);
  461. double x2 = x[L - 1]*x[L - 1];
  462. // Initialize the scattering parameters
  463. *Qext = 0;
  464. *Qsca = 0;
  465. *Qabs = 0;
  466. *Qbk = 0;
  467. Qbktmp = std::complex<double>(0.0, 0.0);
  468. *Qpr = 0;
  469. *g = 0;
  470. *Albedo = 0;
  471. // Initialize the scattering amplitudes
  472. for (t = 0; t < nTheta; t++) {
  473. S1[t] = std::complex<double>(0.0, 0.0);
  474. S2[t] = std::complex<double>(0.0, 0.0);
  475. }
  476. // By using downward recurrence we avoid loss of precision due to float rounding errors
  477. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  478. // http://en.wikipedia.org/wiki/Loss_of_significance
  479. for (i = n_max - 2; i >= 0; i--) {
  480. n = i + 1;
  481. // Equation (27)
  482. *Qext += (n + n + 1)*(an[i].real() + bn[i].real());
  483. // Equation (28)
  484. *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());
  485. // Equation (29)
  486. // We must check carefully this equation. If we remove the typecast to double then the result changes. Which is the correct one??? Ovidio (2014/12/10)
  487. *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());
  488. // Equation (33)
  489. Qbktmp = Qbktmp + (double)(n + n + 1)*(1 - 2*(n % 2))*(an[i]- bn[i]);
  490. //****************************************************//
  491. // Calculate the scattering amplitudes (S1 and S2) //
  492. // Equations (25a) - (25b) //
  493. //****************************************************//
  494. for (t = 0; t < nTheta; t++) {
  495. S1[t] += calc_S1(n, an[i], bn[i], Pi[i][t], Tau[i][t]);
  496. S2[t] += calc_S2(n, an[i], bn[i], Pi[i][t], Tau[i][t]);
  497. }
  498. }
  499. *Qext = 2*(*Qext)/x2; // Equation (27)
  500. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  501. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  502. *Qabs = *Qext - *Qsca; // Equation (30)
  503. *Albedo = *Qsca / *Qext; // Equation (31)
  504. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  505. *Qbk = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  506. return n_max;
  507. }
  508. //**********************************************************************************//
  509. // This function is just a wrapper to call the full 'nMie' function with fewer //
  510. // parameters, it is here mainly for compatibility with older versions of the //
  511. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  512. // any limit for the maximum number of terms. //
  513. // //
  514. // Input parameters: //
  515. // L: Number of layers //
  516. // x: Array containing the size parameters of the layers [0..L-1] //
  517. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  518. // nTheta: Number of scattering angles //
  519. // Theta: Array containing all the scattering angles where the scattering //
  520. // amplitudes will be calculated //
  521. // //
  522. // Output parameters: //
  523. // Qext: Efficiency factor for extinction //
  524. // Qsca: Efficiency factor for scattering //
  525. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  526. // Qbk: Efficiency factor for backscattering //
  527. // Qpr: Efficiency factor for the radiation pressure //
  528. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  529. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  530. // S1, S2: Complex scattering amplitudes //
  531. // //
  532. // Return value: //
  533. // Number of multipolar expansion terms used for the calculations //
  534. //**********************************************************************************//
  535. int nMie(int L, std::vector<double> x, std::vector<std::complex<double> > m,
  536. int nTheta, std::vector<double> Theta,
  537. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  538. std::vector<std::complex<double> > &S1, std::vector<std::complex<double> > &S2) {
  539. return nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  540. }
  541. //**********************************************************************************//
  542. // This function is just a wrapper to call the full 'nMie' function with fewer //
  543. // parameters, it is useful if you want to include a PEC layer but not a limit //
  544. // for the maximum number of terms. //
  545. // //
  546. // Input parameters: //
  547. // L: Number of layers //
  548. // pl: Index of PEC layer. If there is none just send -1 //
  549. // x: Array containing the size parameters of the layers [0..L-1] //
  550. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  551. // nTheta: Number of scattering angles //
  552. // Theta: Array containing all the scattering angles where the scattering //
  553. // amplitudes will be calculated //
  554. // //
  555. // Output parameters: //
  556. // Qext: Efficiency factor for extinction //
  557. // Qsca: Efficiency factor for scattering //
  558. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  559. // Qbk: Efficiency factor for backscattering //
  560. // Qpr: Efficiency factor for the radiation pressure //
  561. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  562. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  563. // S1, S2: Complex scattering amplitudes //
  564. // //
  565. // Return value: //
  566. // Number of multipolar expansion terms used for the calculations //
  567. //**********************************************************************************//
  568. int nMie(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m,
  569. int nTheta, std::vector<double> Theta,
  570. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  571. std::vector<std::complex<double> > &S1, std::vector<std::complex<double> > &S2) {
  572. return nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  573. }
  574. //**********************************************************************************//
  575. // This function is just a wrapper to call the full 'nMie' function with fewer //
  576. // parameters, it is useful if you want to include a limit for the maximum number //
  577. // of terms but not a PEC layer. //
  578. // //
  579. // Input parameters: //
  580. // L: Number of layers //
  581. // x: Array containing the size parameters of the layers [0..L-1] //
  582. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  583. // nTheta: Number of scattering angles //
  584. // Theta: Array containing all the scattering angles where the scattering //
  585. // amplitudes will be calculated //
  586. // n_max: Maximum number of multipolar expansion terms to be used for the //
  587. // calculations. Only used if you know what you are doing, otherwise set //
  588. // this parameter to -1 and the function will calculate it //
  589. // //
  590. // Output parameters: //
  591. // Qext: Efficiency factor for extinction //
  592. // Qsca: Efficiency factor for scattering //
  593. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  594. // Qbk: Efficiency factor for backscattering //
  595. // Qpr: Efficiency factor for the radiation pressure //
  596. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  597. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  598. // S1, S2: Complex scattering amplitudes //
  599. // //
  600. // Return value: //
  601. // Number of multipolar expansion terms used for the calculations //
  602. //**********************************************************************************//
  603. int nMie(int L, std::vector<double> x, std::vector<std::complex<double> > m,
  604. int nTheta, std::vector<double> Theta, int n_max,
  605. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  606. std::vector<std::complex<double> > &S1, std::vector<std::complex<double> > &S2) {
  607. return nMie(L, -1, x, m, nTheta, Theta, n_max, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  608. }
  609. //**********************************************************************************//
  610. // This function calculates complex electric and magnetic field in the surroundings //
  611. // and inside (TODO) the particle. //
  612. // //
  613. // Input parameters: //
  614. // L: Number of layers //
  615. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  616. // x: Array containing the size parameters of the layers [0..L-1] //
  617. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  618. // n_max: Maximum number of multipolar expansion terms to be used for the //
  619. // calculations. Only used if you know what you are doing, otherwise set //
  620. // this parameter to 0 (zero) and the function will calculate it. //
  621. // nCoords: Number of coordinate points //
  622. // Coords: Array containing all coordinates where the complex electric and //
  623. // magnetic fields will be calculated //
  624. // //
  625. // Output parameters: //
  626. // E, H: Complex electric and magnetic field at the provided coordinates //
  627. // //
  628. // Return value: //
  629. // Number of multipolar expansion terms used for the calculations //
  630. //**********************************************************************************//
  631. int nField(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int n_max,
  632. int nCoords, std::vector<double> Xp, std::vector<double> Yp, std::vector<double> Zp,
  633. std::vector<std::complex<double> > &E, std::vector<std::complex<double> > &H) {
  634. int i, n, c;
  635. /* double **Pi, **Tau;
  636. complex *an, *bn;
  637. double *Rho = (double *) malloc(nCoords*sizeof(double));
  638. double *Phi = (double *) malloc(nCoords*sizeof(double));
  639. double *Theta = (double *) malloc(nCoords*sizeof(double));
  640. for (c = 0; c < nCoords; c++) {
  641. Rho[c] = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  642. if (Rho[c] < 1e-3) {
  643. Rho[c] = 1e-3;
  644. }
  645. Phi[c] = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  646. Theta[c] = acos(Xp[c]/Rho[c]);
  647. }
  648. n_max = ScattCoeffs(L, pl, x, m, n_max, an, bn);
  649. Pi = (double **) malloc(n_max*sizeof(double *));
  650. Tau = (double **) malloc(n_max*sizeof(double *));
  651. for (n = 0; n < n_max; n++) {
  652. Pi[n] = (double *) malloc(nCoords*sizeof(double));
  653. Tau[n] = (double *) malloc(nCoords*sizeof(double));
  654. }
  655. calcPiTau(n_max, nCoords, Theta, &Pi, &Tau);
  656. double x2 = x[L - 1]*x[L - 1];
  657. // Initialize the fields
  658. for (c = 0; c < nCoords; c++) {
  659. E[c] = C_ZERO;
  660. H[c] = C_ZERO;
  661. }*/
  662. //*******************************************************//
  663. // external scattering field = incident + scattered //
  664. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  665. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  666. //*******************************************************//
  667. // Firstly the easiest case, we want the field outside the particle
  668. // if (Rho[c] >= x[L - 1]) {
  669. // }
  670. for (i = 1; i < (n_max - 1); i++) {
  671. // n = i - 1;
  672. /* // Equation (27)
  673. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  674. // Equation (28)
  675. *Qsca = *Qsca + (double)(n + n + 1)*(an[i].r*an[i].r + an[i].i*an[i].i + bn[i].r*bn[i].r + bn[i].i*bn[i].i);
  676. // Equation (29)
  677. *Qpr = *Qpr + ((n*(n + 2)/(n + 1))*((Cadd(Cmul(an[i], std::conj(an[n])), Cmul(bn[i], std::conj(bn[n])))).r) + ((double)(n + n + 1)/(n*(n + 1)))*(Cmul(an[i], std::conj(bn[i])).r));
  678. // Equation (33)
  679. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  680. */
  681. //****************************************************//
  682. // Calculate the scattering amplitudes (S1 and S2) //
  683. // Equations (25a) - (25b) //
  684. //****************************************************//
  685. /* for (t = 0; t < nTheta; t++) {
  686. S1[t] = Cadd(S1[t], calc_S1(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  687. S2[t] = Cadd(S2[t], calc_S2(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  688. }*/
  689. }
  690. return n_max;
  691. }