nmie.cc 38 KB

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