nmie.cc 50 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141
  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. //**********************************************************************************//
  430. //**********************************************************************************//
  431. int ScattCoeff_std(double x[], complex m[], complex **an, complex **bn,
  432. int L, int pl, std::vector<double> x_std,
  433. std::vector<std::complex<double> > m_std, int n_max,
  434. std::vector< std::complex<double> > an_std,
  435. std::vector< std::complex<double> > bn_std){
  436. //************************************************************************//
  437. // Calculate the index of the first layer. It can be either 0 (default) //
  438. // or the index of the outermost PEC layer. In the latter case all layers //
  439. // below the PEC are discarded. //
  440. //************************************************************************//
  441. int fl = firstLayer(L, pl);
  442. if (n_max <= 0) {
  443. n_max = Nmax(L, fl, pl, x, m);
  444. }
  445. complex z1, z2, cn;
  446. complex Num, Denom;
  447. complex G1, G2;
  448. complex Temp;
  449. double Tmp;
  450. int n, l, t;
  451. //**************************************************************************//
  452. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  453. // means that index = layer number - 1 or index = n - 1. The only exception //
  454. // are the arrays for representing D1, D3 and Q because they need a value //
  455. // for the index 0 (zero), hence it is important to consider this shift //
  456. // between different arrays. The change was done to optimize memory usage. //
  457. //**************************************************************************//
  458. // Allocate memory to the arrays
  459. complex **D1_mlxl = (complex **) malloc(L*sizeof(complex *));
  460. complex **D1_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  461. complex **D3_mlxl = (complex **) malloc(L*sizeof(complex *));
  462. complex **D3_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  463. complex **Q = (complex **) malloc(L*sizeof(complex *));
  464. complex **Ha = (complex **) malloc(L*sizeof(complex *));
  465. complex **Hb = (complex **) malloc(L*sizeof(complex *));
  466. for (l = 0; l < L; l++) {
  467. D1_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  468. D1_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  469. D3_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  470. D3_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  471. Q[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  472. Ha[l] = (complex *) malloc(n_max*sizeof(complex));
  473. Hb[l] = (complex *) malloc(n_max*sizeof(complex));
  474. }
  475. (*an) = (complex *) malloc(n_max*sizeof(complex));
  476. (*bn) = (complex *) malloc(n_max*sizeof(complex));
  477. complex *D1XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  478. complex *D3XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  479. complex *PsiXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  480. complex *ZetaXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  481. //*************************************************//
  482. // Calculate D1 and D3 for z1 in the first layer //
  483. //*************************************************//
  484. if (fl == pl) { // PEC layer
  485. for (n = 0; n <= n_max; n++) {
  486. D1_mlxl[fl][n] = Complex(0, -1);
  487. D3_mlxl[fl][n] = C_I;
  488. }
  489. } else { // Regular layer
  490. z1 = RCmul(x[fl], m[fl]);
  491. // Calculate D1 and D3
  492. calcD1D3(z1, n_max, &(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 < n_max; 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 = RCmul(x[l], m[l]);
  509. z2 = RCmul(x[l - 1], m[l]);
  510. //Calculate D1 and D3 for z1
  511. calcD1D3(z1, n_max, &(D1_mlxl[l]), &(D3_mlxl[l]));
  512. //Calculate D1 and D3 for z2
  513. calcD1D3(z2, n_max, &(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 = RCmul(exp(-2*(z1.i - z2.i)), Complex(cos(-2*z2.r) - exp(-2*z2.i), sin(-2*z2.r)));
  519. Denom = Complex(cos(-2*z1.r) - exp(-2*z1.i), sin(-2*z1.r));
  520. Q[l][0] = Cdiv(Num, Denom);
  521. for (n = 1; n <= n_max; n++) {
  522. cn = Complex(n, 0);
  523. Num = Cmul(Cadd(Cmul(z1, D1_mlxl[l][n]), cn), Csub(cn, Cmul(z1, D3_mlxl[l][n - 1])));
  524. Denom = Cmul(Cadd(Cmul(z2, D1_mlxlM1[l][n]), cn), Csub(cn, Cmul(z2, D3_mlxlM1[l][n - 1])));
  525. Q[l][n] = Cdiv(Cmul(RCmul((x[l - 1]*x[l - 1])/(x[l]*x[l]), Q[l][n - 1]), Num), Denom);
  526. }
  527. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  528. for (n = 1; n <= n_max; n++) {
  529. //Ha
  530. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  531. G1 = RCmul(-1.0, D1_mlxlM1[l][n]);
  532. G2 = RCmul(-1.0, D3_mlxlM1[l][n]);
  533. } else {
  534. G1 = Csub(Cmul(m[l], Ha[l - 1][n - 1]), Cmul(m[l - 1], D1_mlxlM1[l][n]));
  535. G2 = Csub(Cmul(m[l], Ha[l - 1][n - 1]), Cmul(m[l - 1], D3_mlxlM1[l][n]));
  536. }
  537. Temp = Cmul(Q[l][n], G1);
  538. Num = Csub(Cmul(G2, D1_mlxl[l][n]), Cmul(Temp, D3_mlxl[l][n]));
  539. Denom = Csub(G2, Temp);
  540. Ha[l][n - 1] = Cdiv(Num, Denom);
  541. //Hb
  542. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  543. G1 = Hb[l - 1][n - 1];
  544. G2 = Hb[l - 1][n - 1];
  545. } else {
  546. G1 = Csub(Cmul(m[l - 1], Hb[l - 1][n - 1]), Cmul(m[l], D1_mlxlM1[l][n]));
  547. G2 = Csub(Cmul(m[l - 1], Hb[l - 1][n - 1]), Cmul(m[l], D3_mlxlM1[l][n]));
  548. }
  549. Temp = Cmul(Q[l][n], G1);
  550. Num = Csub(Cmul(G2, D1_mlxl[l][n]), Cmul(Temp, D3_mlxl[l][n]));
  551. Denom = Csub(G2, Temp);
  552. Hb[l][n - 1] = Cdiv(Num, Denom);
  553. }
  554. }
  555. //**************************************//
  556. //Calculate D1, D3, Psi and Zeta for XL //
  557. //**************************************//
  558. z1 = Complex(x[L - 1], 0);
  559. // Calculate D1XL and D3XL
  560. calcD1D3(z1, n_max, &D1XL, &D3XL);
  561. // Calculate PsiXL and ZetaXL
  562. calcPsiZeta(z1, n_max, D1XL, D3XL, &PsiXL, &ZetaXL);
  563. //*********************************************************************//
  564. // Finally, we calculate the scattering coefficients (an and bn) and //
  565. // the angular functions (Pi and Tau). Note that for these arrays the //
  566. // first layer is 0 (zero), in future versions all arrays will follow //
  567. // this convention to save memory. (13 Nov, 2014) //
  568. //*********************************************************************//
  569. for (n = 0; n < n_max; n++) {
  570. //********************************************************************//
  571. //Expressions for calculating an and bn coefficients are not valid if //
  572. //there is only one PEC layer (ie, for a simple PEC sphere). //
  573. //********************************************************************//
  574. if (pl < (L - 1)) {
  575. (*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]);
  576. (*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]);
  577. } else {
  578. (*an)[n] = calc_an(n + 1, x[L - 1], C_ZERO, C_ONE, PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  579. (*bn)[n] = Cdiv(PsiXL[n + 1], ZetaXL[n + 1]);
  580. }
  581. }
  582. // Free the memory used for the arrays
  583. for (l = 0; l < L; l++) {
  584. free(D1_mlxl[l]);
  585. free(D1_mlxlM1[l]);
  586. free(D3_mlxl[l]);
  587. free(D3_mlxlM1[l]);
  588. free(Q[l]);
  589. free(Ha[l]);
  590. free(Hb[l]);
  591. }
  592. free(D1_mlxl);
  593. free(D1_mlxlM1);
  594. free(D3_mlxl);
  595. free(D3_mlxlM1);
  596. free(Q);
  597. free(Ha);
  598. free(Hb);
  599. free(D1XL);
  600. free(D3XL);
  601. free(PsiXL);
  602. free(ZetaXL);
  603. return n_max;
  604. }
  605. //**********************************************************************************//
  606. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  607. // parameters, it is here mainly for compatibility with older versions of the //
  608. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  609. // any limit for the maximum number of terms. //
  610. //**********************************************************************************//
  611. 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[]) {
  612. return nMieScatt(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  613. }
  614. int nMie_std(double x[], complex m[], double Theta[], complex S1[], complex S2[],
  615. int L, std::vector<double> &x_std, std::vector<std::complex<double> > &m_std, int nTheta, std::vector<double> &Theta_std, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector< std::complex<double> > &S1_std, std::vector< std::complex<double> > &S2_std) {
  616. return nMieScatt_std(x, m, Theta, S1, S2, L, -1, x_std, m_std, nTheta, Theta_std, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1_std, S2_std);
  617. }
  618. //**********************************************************************************//
  619. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  620. // parameters, it is useful if you want to include a PEC layer but not a limit //
  621. // for the maximum number of terms. //
  622. // //
  623. // Input parameters: //
  624. // L: Number of layers //
  625. // pl: Index of PEC layer. If there is none just send -1 //
  626. // x: Array containing the size parameters of the layers [0..L-1] //
  627. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  628. // nTheta: Number of scattering angles //
  629. // Theta: Array containing all the scattering angles where the scattering //
  630. // amplitudes will be calculated //
  631. // //
  632. // Output parameters: //
  633. // Qext: Efficiency factor for extinction //
  634. // Qsca: Efficiency factor for scattering //
  635. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  636. // Qbk: Efficiency factor for backscattering //
  637. // Qpr: Efficiency factor for the radiation pressure //
  638. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  639. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  640. // S1, S2: Complex scattering amplitudes //
  641. // //
  642. // Return value: //
  643. // Number of multipolar expansion terms used for the calculations //
  644. //**********************************************************************************//
  645. 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[]) {
  646. return nMieScatt(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  647. }
  648. //**********************************************************************************//
  649. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  650. // parameters, it is useful if you want to include a limit for the maximum number //
  651. // of terms but not a PEC layer. //
  652. // //
  653. // Input parameters: //
  654. // L: Number of layers //
  655. // x: Array containing the size parameters of the layers [0..L-1] //
  656. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  657. // nTheta: Number of scattering angles //
  658. // Theta: Array containing all the scattering angles where the scattering //
  659. // amplitudes will be calculated //
  660. // n_max: Maximum number of multipolar expansion terms to be used for the //
  661. // calculations. Only used if you know what you are doing, otherwise set //
  662. // this parameter to -1 and the function will calculate it //
  663. // //
  664. // Output parameters: //
  665. // Qext: Efficiency factor for extinction //
  666. // Qsca: Efficiency factor for scattering //
  667. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  668. // Qbk: Efficiency factor for backscattering //
  669. // Qpr: Efficiency factor for the radiation pressure //
  670. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  671. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  672. // S1, S2: Complex scattering amplitudes //
  673. // //
  674. // Return value: //
  675. // Number of multipolar expansion terms used for the calculations //
  676. //**********************************************************************************//
  677. 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[]) {
  678. return nMieScatt(L, -1, x, m, nTheta, Theta, n_max, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  679. }
  680. //**********************************************************************************//
  681. // This function calculates the actual scattering parameters and amplitudes //
  682. // //
  683. // Input parameters: //
  684. // L: Number of layers //
  685. // pl: Index of PEC layer. If there is none just send -1 //
  686. // x: Array containing the size parameters of the layers [0..L-1] //
  687. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  688. // nTheta: Number of scattering angles //
  689. // Theta: Array containing all the scattering angles where the scattering //
  690. // amplitudes will be calculated //
  691. // n_max: Maximum number of multipolar expansion terms to be used for the //
  692. // calculations. Only used if you know what you are doing, otherwise set //
  693. // this parameter to -1 and the function will calculate it //
  694. // //
  695. // Output parameters: //
  696. // Qext: Efficiency factor for extinction //
  697. // Qsca: Efficiency factor for scattering //
  698. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  699. // Qbk: Efficiency factor for backscattering //
  700. // Qpr: Efficiency factor for the radiation pressure //
  701. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  702. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  703. // S1, S2: Complex scattering amplitudes //
  704. // //
  705. // Return value: //
  706. // Number of multipolar expansion terms used for the calculations //
  707. //**********************************************************************************//
  708. 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[]) {
  709. int i, n, t;
  710. double **Pi, **Tau;
  711. complex *an, *bn;
  712. complex Qbktmp;
  713. n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  714. Pi = (double **) malloc(n_max*sizeof(double *));
  715. Tau = (double **) malloc(n_max*sizeof(double *));
  716. for (n = 0; n < n_max; n++) {
  717. Pi[n] = (double *) malloc(nTheta*sizeof(double));
  718. Tau[n] = (double *) malloc(nTheta*sizeof(double));
  719. }
  720. calcPiTau(n_max, nTheta, Theta, &Pi, &Tau);
  721. double x2 = x[L - 1]*x[L - 1];
  722. // Initialize the scattering parameters
  723. *Qext = 0;
  724. *Qsca = 0;
  725. *Qabs = 0;
  726. *Qbk = 0;
  727. Qbktmp = C_ZERO;
  728. *Qpr = 0;
  729. *g = 0;
  730. *Albedo = 0;
  731. // Initialize the scattering amplitudes
  732. for (t = 0; t < nTheta; t++) {
  733. S1[t] = C_ZERO;
  734. S2[t] = C_ZERO;
  735. }
  736. // By using downward recurrence we avoid loss of precision due to float rounding errors
  737. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  738. // http://en.wikipedia.org/wiki/Loss_of_significance
  739. for (i = n_max - 2; i >= 0; i--) {
  740. n = i + 1;
  741. // Equation (27)
  742. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  743. // Equation (28)
  744. *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);
  745. // Equation (29)
  746. *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));
  747. // Equation (33)
  748. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  749. //****************************************************//
  750. // Calculate the scattering amplitudes (S1 and S2) //
  751. // Equations (25a) - (25b) //
  752. //****************************************************//
  753. for (t = 0; t < nTheta; t++) {
  754. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  755. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  756. }
  757. }
  758. *Qext = 2*(*Qext)/x2; // Equation (27)
  759. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  760. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  761. *Qabs = *Qext - *Qsca; // Equation (30)
  762. *Albedo = *Qsca / *Qext; // Equation (31)
  763. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  764. *Qbk = (Qbktmp.r*Qbktmp.r + Qbktmp.i*Qbktmp.i)/x2; // Equation (33)
  765. // Free the memory used for the arrays
  766. for (n = 0; n < n_max; n++) {
  767. free(Pi[n]);
  768. free(Tau[n]);
  769. }
  770. free(Pi);
  771. free(Tau);
  772. free(an);
  773. free(bn);
  774. return n_max;
  775. }
  776. int nMieScatt_std(double x[], complex m[], double Theta[], complex S1[], complex S2[],
  777. int L, int pl,
  778. std::vector<double> &x_std, std::vector<std::complex<double> > &m_std,
  779. int nTheta, std::vector<double> &Theta_std,
  780. int n_max, double *Qext, double *Qsca, double *Qabs, double *Qbk,
  781. double *Qpr, double *g, double *Albedo,
  782. std::vector< std::complex<double> > &S1_std,
  783. std::vector< std::complex<double> > &S2_std) {
  784. int i, n, t;
  785. double **Pi, **Tau;
  786. std::vector< std::vector<double> > Pi_std, Tau_std;
  787. complex *an, *bn;
  788. std::vector< std::complex<double> > an_std, bn_std;
  789. complex Qbktmp;
  790. std::complex<double> Qbktmp_std;
  791. {
  792. int tmp_n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  793. n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  794. }
  795. Pi = (double **) malloc(n_max*sizeof(double *));
  796. Tau = (double **) malloc(n_max*sizeof(double *));
  797. for (n = 0; n < n_max; n++) {
  798. Pi[n] = (double *) malloc(nTheta*sizeof(double));
  799. Tau[n] = (double *) malloc(nTheta*sizeof(double));
  800. }
  801. calcPiTau(n_max, nTheta, Theta, &Pi, &Tau);
  802. double x2 = x[L - 1]*x[L - 1];
  803. // Initialize the scattering parameters
  804. *Qext = 0;
  805. *Qsca = 0;
  806. *Qabs = 0;
  807. *Qbk = 0;
  808. Qbktmp = C_ZERO;
  809. *Qpr = 0;
  810. *g = 0;
  811. *Albedo = 0;
  812. // Initialize the scattering amplitudes
  813. for (t = 0; t < nTheta; t++) {
  814. S1[t] = C_ZERO;
  815. S2[t] = C_ZERO;
  816. }
  817. // By using downward recurrence we avoid loss of precision due to float rounding errors
  818. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  819. // http://en.wikipedia.org/wiki/Loss_of_significance
  820. for (i = n_max - 2; i >= 0; i--) {
  821. n = i + 1;
  822. // Equation (27)
  823. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  824. // Equation (28)
  825. *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);
  826. // Equation (29)
  827. *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));
  828. // Equation (33)
  829. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  830. //****************************************************//
  831. // Calculate the scattering amplitudes (S1 and S2) //
  832. // Equations (25a) - (25b) //
  833. //****************************************************//
  834. for (t = 0; t < nTheta; t++) {
  835. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  836. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  837. }
  838. }
  839. *Qext = 2*(*Qext)/x2; // Equation (27)
  840. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  841. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  842. *Qabs = *Qext - *Qsca; // Equation (30)
  843. *Albedo = *Qsca / *Qext; // Equation (31)
  844. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  845. *Qbk = (Qbktmp.r*Qbktmp.r + Qbktmp.i*Qbktmp.i)/x2; // Equation (33)
  846. // Free the memory used for the arrays
  847. for (n = 0; n < n_max; n++) {
  848. free(Pi[n]);
  849. free(Tau[n]);
  850. }
  851. free(Pi);
  852. free(Tau);
  853. free(an);
  854. free(bn);
  855. return n_max;
  856. }
  857. //**********************************************************************************//
  858. // This function calculates complex electric and magnetic field in the surroundings //
  859. // and inside (TODO) the particle. //
  860. // //
  861. // Input parameters: //
  862. // L: Number of layers //
  863. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  864. // x: Array containing the size parameters of the layers [0..L-1] //
  865. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  866. // n_max: Maximum number of multipolar expansion terms to be used for the //
  867. // calculations. Only used if you know what you are doing, otherwise set //
  868. // this parameter to 0 (zero) and the function will calculate it. //
  869. // nCoords: Number of coordinate points //
  870. // Coords: Array containing all coordinates where the complex electric and //
  871. // magnetic fields will be calculated //
  872. // //
  873. // Output parameters: //
  874. // E, H: Complex electric and magnetic field at the provided coordinates //
  875. // //
  876. // Return value: //
  877. // Number of multipolar expansion terms used for the calculations //
  878. //**********************************************************************************//
  879. 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[]){
  880. int i, n, c;
  881. double **Pi, **Tau;
  882. complex *an, *bn;
  883. double *Rho = (double *) malloc(nCoords*sizeof(double));
  884. double *Phi = (double *) malloc(nCoords*sizeof(double));
  885. double *Theta = (double *) malloc(nCoords*sizeof(double));
  886. for (c = 0; c < nCoords; c++) {
  887. Rho[c] = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  888. if (Rho[c] < 1e-3) {
  889. Rho[c] = 1e-3;
  890. }
  891. Phi[c] = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  892. Theta[c] = acos(Xp[c]/Rho[c]);
  893. }
  894. n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  895. Pi = (double **) malloc(n_max*sizeof(double *));
  896. Tau = (double **) malloc(n_max*sizeof(double *));
  897. for (n = 0; n < n_max; n++) {
  898. Pi[n] = (double *) malloc(nCoords*sizeof(double));
  899. Tau[n] = (double *) malloc(nCoords*sizeof(double));
  900. }
  901. calcPiTau(n_max, nCoords, Theta, &Pi, &Tau);
  902. double x2 = x[L - 1]*x[L - 1];
  903. // Initialize the fields
  904. for (c = 0; c < nCoords; c++) {
  905. E[c] = C_ZERO;
  906. H[c] = C_ZERO;
  907. }
  908. //*******************************************************//
  909. // external scattering field = incident + scattered //
  910. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  911. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  912. //*******************************************************//
  913. // Firstly the easiest case, we want the field outside the particle
  914. if (Rho[c] >= x[L - 1]) {
  915. }
  916. // for (i = 1; i < (n_max - 1); i++) {
  917. // n = i - 1;
  918. /* // Equation (27)
  919. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  920. // Equation (28)
  921. *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);
  922. // Equation (29)
  923. *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));
  924. // Equation (33)
  925. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  926. */
  927. //****************************************************//
  928. // Calculate the scattering amplitudes (S1 and S2) //
  929. // Equations (25a) - (25b) //
  930. //****************************************************//
  931. /* for (t = 0; t < nTheta; t++) {
  932. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  933. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  934. }*/
  935. // }
  936. // Free the memory used for the arrays
  937. for (n = 0; n < n_max; n++) {
  938. free(Pi[n]);
  939. free(Tau[n]);
  940. }
  941. free(Pi);
  942. free(Tau);
  943. free(an);
  944. free(bn);
  945. free(Rho);
  946. free(Phi);
  947. free(Theta);
  948. return n_max;
  949. }