nmie.cc 56 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278
  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. // TODO should be replaced with std::round() or std::lround()
  44. #define round(x) ((x) >= 0 ? (int)((x) + 0.5):(int)((x) - 0.5))
  45. complex C_ZERO = {0.0, 0.0};
  46. complex C_ONE = {1.0, 0.0};
  47. complex C_I = {0.0, 1.0};
  48. int compare(std::string operation, complex *a, std::vector< std::complex<double> > b) {
  49. for (int i = 0; i < b.size(); ++i) {
  50. //if (i > 50) continue;
  51. double diff_r = std::abs((a[i].r - b[i].real())/a[i].r);
  52. double diff_i = std::abs((a[i].i - b[i].imag())/a[i].i);
  53. double epsilon= 1e-16;
  54. if (diff_r > epsilon ||diff_i > epsilon) {
  55. printf("\n*** WARNING!! Non-zero diff!!! ***\n");
  56. printf("Op: %s at i=%i, diff_r=%g, diff_i=%g, a_r=%g, a_i = %g\n",
  57. operation.c_str(), i, diff_r, diff_i, a[i].r, a[i].i);
  58. }
  59. double factor = 1.0;
  60. if ((diff_r > epsilon * factor && a[i].r/a[0].r > epsilon)
  61. || (diff_i > epsilon*factor && a[i].i/a[0].i > epsilon)) {
  62. printf("\n******** ERROR!! Non-zero diff!!! ********\n");
  63. printf("Op: %s at i=%i, diff_r=%g, diff_i=%g, a_r=%g, a_i = %g\n",
  64. operation.c_str(), i, diff_r, diff_i, a[i].r, a[i].i);
  65. }
  66. }
  67. return 0;
  68. }
  69. int firstLayer(int L, int pl) {
  70. if (pl >= 0) {
  71. return pl;
  72. } else {
  73. return 0;
  74. }
  75. }
  76. // Calculate Nstop - equation (17)
  77. int Nstop(double xL) {
  78. int result;
  79. if (xL <= 8) {
  80. result = std::lround(xL + 4*pow(xL, 1/3) + 1);
  81. } else if (xL <= 4200) {
  82. result = std::lround(xL + 4.05*pow(xL, 1/3) + 2);
  83. } else {
  84. result = std::lround(xL + 4*pow(xL, 1/3) + 2);
  85. }
  86. return result;
  87. }
  88. int Nmax(int L, int fl, int pl, double x[], complex m[]) {
  89. int i, result, ri, riM1;
  90. result = Nstop(x[L - 1]);
  91. for (i = fl; i < L; i++) {
  92. if (i > pl) {
  93. ri = round(Cabs(RCmul(x[i], m[i])));
  94. } else {
  95. ri = 0;
  96. }
  97. if (result < ri) {
  98. result = ri;
  99. }
  100. if ((i > fl) && ((i - 1) > pl)) {
  101. riM1 = round(Cabs(RCmul(x[i - 1], m[i])));
  102. } else {
  103. riM1 = 0;
  104. }
  105. if (result < riM1) {
  106. result = riM1;
  107. }
  108. }
  109. return result + 15;
  110. }
  111. //**********************************************************************************//
  112. int Nmax_std(int L, int fl, int pl, std::vector<double> x,
  113. std::vector<std::complex<double> > m) {
  114. int i, result, ri, riM1;
  115. result = Nstop(x[L - 1]);
  116. for (i = fl; i < L; i++) {
  117. if (i > pl) {
  118. ri = std::lround(std::abs(x[i]*m[i]));
  119. } else {
  120. ri = 0;
  121. }
  122. if (result < ri) {
  123. result = ri;
  124. }
  125. if ((i > fl) && ((i - 1) > pl)) {
  126. riM1 = std::lround(std::abs(x[i - 1]* m[i]));
  127. } else {
  128. riM1 = 0;
  129. }
  130. if (result < riM1) {
  131. result = riM1;
  132. }
  133. }
  134. return result + 15;
  135. }
  136. //**********************************************************************************//
  137. // This function calculates the spherical Bessel functions (jn and hn) for a given //
  138. // value of r. //
  139. // //
  140. // Input parameters: //
  141. // r: Real argument to evaluate jn and hn //
  142. // n_max: Maximum number of terms to calculate jn and hn //
  143. // //
  144. // Output parameters: //
  145. // jn, hn: Spherical Bessel functions (double) //
  146. //**********************************************************************************//
  147. void sphericalBessel(double r, int n_max, double **j, double **h) {
  148. int n;
  149. (*j)[0] = sin(r)/r;
  150. (*j)[1] = sin(r)/r/r - cos(r)/r;
  151. (*h)[0] = -cos(r)/r;
  152. (*h)[1] = -cos(r)/r/r - sin(r)/r;
  153. for (n = 2; n < n_max; n++) {
  154. (*j)[n] = (n + n + 1)*(*j)[n - 1]/r - (*j)[n - 2];
  155. (*h)[n] = (n + n + 1)*(*h)[n - 1]/r - (*h)[n - 2];
  156. }
  157. }
  158. //**********************************************************************************//
  159. // This function calculates the spherical Bessel functions (jn and hn) for a given //
  160. // value of z. //
  161. // //
  162. // Input parameters: //
  163. // z: Real argument to evaluate jn and hn //
  164. // n_max: Maximum number of terms to calculate jn and hn //
  165. // //
  166. // Output parameters: //
  167. // jn, hn: Spherical Bessel functions (complex) //
  168. //**********************************************************************************//
  169. void CsphericalBessel(complex z, int n_max, complex **j, complex **h) {
  170. int n;
  171. (*j)[0] = Cdiv(Csin(z), z);
  172. (*j)[1] = Csub(Cdiv(Cdiv(Csin(z), z), z), Cdiv(Ccos(z), z));
  173. (*h)[0] = Csub(C_ZERO, Cdiv(Ccos(z), z));
  174. (*h)[1] = Csub(C_ZERO, Cadd(Cdiv(Cdiv(Ccos(z), z), z), Cdiv(Csin(z), z)));
  175. for (n = 2; n < n_max; n++) {
  176. (*j)[n] = Csub(RCmul(n + n + 1, Cdiv((*j)[n - 1], z)), (*j)[n - 2]);
  177. (*h)[n] = Csub(RCmul(n + n + 1, Cdiv((*h)[n - 1], z)), (*h)[n - 2]);
  178. }
  179. }
  180. // Calculate an - equation (5)
  181. complex calc_an(int n, double XL, complex Ha, complex mL, complex PsiXL, complex ZetaXL, complex PsiXLM1, complex ZetaXLM1) {
  182. complex Num = Csub(Cmul(Cadd(Cdiv(Ha, mL), Complex(n/XL, 0)), PsiXL), PsiXLM1);
  183. complex Denom = Csub(Cmul(Cadd(Cdiv(Ha, mL), Complex(n/XL, 0)), ZetaXL), ZetaXLM1);
  184. return Cdiv(Num, Denom);
  185. }
  186. /////////////////////////////////////////////
  187. std::complex<double>
  188. calc_an_std(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  189. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  190. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  191. std::complex<double> Num = (( (Ha/mL) + std::complex<double>(n/XL, 0) ) * PsiXL) - PsiXLM1;
  192. std::complex<double> Denom = (( (Ha/mL) + std::complex<double>(n/XL, 0) ) * ZetaXL) - ZetaXLM1;
  193. return Num/Denom;
  194. }
  195. // Calculate bn - equation (6)
  196. complex calc_bn(int n, double XL, complex Hb, complex mL, complex PsiXL, complex ZetaXL, complex PsiXLM1, complex ZetaXLM1) {
  197. complex Num = Csub(Cmul(Cadd(Cmul(Hb, mL), Complex(n/XL, 0)), PsiXL), PsiXLM1);
  198. complex Denom = Csub(Cmul(Cadd(Cmul(Hb, mL), Complex(n/XL, 0)), ZetaXL), ZetaXLM1);
  199. return Cdiv(Num, Denom);
  200. }
  201. // Calculate bn - equation (6)
  202. std::complex<double>
  203. calc_bn_std(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  204. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  205. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  206. std::complex<double> Num = (( (Hb*mL) + std::complex<double>(n/XL,0) ) * PsiXL) - PsiXLM1;
  207. std::complex<double> Denom = (( (Hb*mL) + std::complex<double>(n/XL,0) ) * ZetaXL) - ZetaXLM1;
  208. return Num/Denom;
  209. }
  210. // Calculates S1_n - equation (25a)
  211. complex calc_S1_n(int n, complex an, complex bn, double Pin, double Taun) {
  212. return RCmul((double)(n + n + 1)/(double)(n*n + n), Cadd(RCmul(Pin, an), RCmul(Taun, bn)));
  213. }
  214. // Calculates S2_n - equation (25b) (it's the same as (25a), just switches Pin and Taun)
  215. complex calc_S2_n(int n, complex an, complex bn, double Pin, double Taun) {
  216. return calc_S1_n(n, an, bn, Taun, Pin);
  217. }
  218. //**********************************************************************************//
  219. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  220. // given value of z. //
  221. // //
  222. // Input parameters: //
  223. // z: Complex argument to evaluate Psi and Zeta //
  224. // n_max: Maximum number of terms to calculate Psi and Zeta //
  225. // //
  226. // Output parameters: //
  227. // Psi, Zeta: Riccati-Bessel functions //
  228. //**********************************************************************************//
  229. void calcPsiZeta(complex z, int n_max, complex *D1, complex *D3, complex **Psi, complex **Zeta) {
  230. int n;
  231. complex cn;
  232. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  233. (*Psi)[0] = Complex(sin(z.r), 0);
  234. (*Zeta)[0] = Complex(sin(z.r), -cos(z.r));
  235. for (n = 1; n <= n_max; n++) {
  236. cn = Complex(n, 0);
  237. (*Psi)[n] = Cmul((*Psi)[n - 1], Csub(Cdiv(cn, z), D1[n - 1]));
  238. (*Zeta)[n] = Cmul((*Zeta)[n - 1], Csub(Cdiv(cn, z), D3[n - 1]));
  239. }
  240. }
  241. //**********************************************************************************//
  242. void calcPsiZeta_std(std::complex<double> z, int n_max,
  243. std::vector< std::complex<double> > D1,
  244. std::vector< std::complex<double> > D3,
  245. std::vector< std::complex<double> > &Psi,
  246. std::vector< std::complex<double> > &Zeta) {
  247. int n;
  248. std::complex<double> cn;
  249. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  250. Psi[0] = std::complex<double>(sin(z.real()), 0);
  251. Zeta[0] = std::complex<double>(sin(z.real()), -cos(z.real()));
  252. for (n = 1; n <= n_max; n++) {
  253. cn = std::complex<double>(n, 0);
  254. Psi[n] = Psi[n-1] * ( (cn/z) - D1[n-1] );
  255. Zeta[n] = Zeta[n-1] * ( (cn/z) - D3[n-1] );
  256. }
  257. }
  258. //**********************************************************************************//
  259. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  260. // functions (D1 and D3) for a given value of z. //
  261. // //
  262. // Input parameters: //
  263. // z: Complex argument to evaluate D1 and D3 //
  264. // n_max: Maximum number of terms to calculate D1 and D3 //
  265. // //
  266. // Output parameters: //
  267. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  268. //**********************************************************************************//
  269. void calcD1D3(complex z, int n_max, complex **D1, complex **D3) {
  270. int n;
  271. complex cn;
  272. complex *PsiZeta = (complex *) malloc((n_max + 1)*sizeof(complex));
  273. // Downward recurrence for D1 - equations (16a) and (16b)
  274. (*D1)[n_max] = C_ZERO;
  275. for (n = n_max; n > 0; n--) {
  276. cn = Complex(n, 0);
  277. (*D1)[n - 1] = Csub(Cdiv(cn, z), Cdiv(C_ONE, Cadd((*D1)[n], Cdiv(cn, z))));
  278. }
  279. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  280. PsiZeta[0] = RCmul(0.5, Csub(C_ONE, Cmul(Complex(cos(2*z.r), sin(2*z.r)), Complex(exp(-2*z.i), 0))));
  281. (*D3)[0] = C_I;
  282. for (n = 1; n <= n_max; n++) {
  283. cn = Complex(n, 0);
  284. PsiZeta[n] = Cmul(PsiZeta[n - 1], Cmul(Csub(Cdiv(cn, z), (*D1)[n - 1]), Csub(Cdiv(cn, z), (*D3)[n - 1])));
  285. (*D3)[n] = Cadd((*D1)[n], Cdiv(C_I, PsiZeta[n]));
  286. }
  287. free(PsiZeta);
  288. }
  289. //**********************************************************************************//
  290. void calcD1D3_std(std::complex<double> z, int n_max,
  291. std::vector< std::complex<double> > &D1,
  292. std::vector< std::complex<double> > &D3) {
  293. int n;
  294. std::complex<double> cn;
  295. //complex *PsiZeta = (complex *) malloc((n_max + 1)*sizeof(complex));
  296. std::vector< std::complex<double> > PsiZeta;
  297. PsiZeta.resize(n_max + 1);
  298. // Downward recurrence for D1 - equations (16a) and (16b)
  299. D1[n_max] = std::complex<double>(0.0, 0.0);
  300. for (n = n_max; n > 0; n--) {
  301. cn = std::complex<double>(n, 0.0);
  302. D1[n - 1] = cn/z - 1.0/(D1[n] + (cn/z));
  303. }
  304. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  305. PsiZeta[0] = 0.5 * (1.0 - (std::complex<double>(cos(2*z.real()), sin(2*z.real()))
  306. * std::complex<double>(exp(-2*z.imag()), 0)));
  307. D3[0] = std::complex<double>(0.0, 1.0);
  308. for (n = 1; n <= n_max; n++) {
  309. cn = std::complex<double>(n, 0.0);
  310. PsiZeta[n] = PsiZeta[n-1] * (((cn/ z) - D1[n - 1]) * ((cn/ z)- D3[n - 1]));
  311. D3[n] = D1[n]+ (std::complex<double>(0.0, 1.0)/ PsiZeta[n]);
  312. }
  313. // free(PsiZeta);
  314. }
  315. //**********************************************************************************//
  316. // This function calculates Pi and Tau for all values of Theta. //
  317. // //
  318. // Input parameters: //
  319. // n_max: Maximum number of terms to calculate Pi and Tau //
  320. // nTheta: Number of scattering angles //
  321. // Theta: Array containing all the scattering angles where the scattering //
  322. // amplitudes will be calculated //
  323. // //
  324. // Output parameters: //
  325. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  326. //**********************************************************************************//
  327. void calcPiTau(int n_max, int nTheta, double Theta[], double ***Pi, double ***Tau) {
  328. int n, t;
  329. for (n = 0; n < n_max; n++) {
  330. //****************************************************//
  331. // Equations (26a) - (26c) //
  332. //****************************************************//
  333. for (t = 0; t < nTheta; t++) {
  334. if (n == 0) {
  335. // Initialize Pi and Tau
  336. (*Pi)[n][t] = 1.0;
  337. (*Tau)[n][t] = (n + 1)*cos(Theta[t]);
  338. } else {
  339. // Calculate the actual values
  340. (*Pi)[n][t] = ((n == 1) ? ((n + n + 1)*cos(Theta[t])*(*Pi)[n - 1][t]/n)
  341. : (((n + n + 1)*cos(Theta[t])*(*Pi)[n - 1][t] - (n + 1)*(*Pi)[n - 2][t])/n));
  342. (*Tau)[n][t] = (n + 1)*cos(Theta[t])*(*Pi)[n][t] - (n + 2)*(*Pi)[n - 1][t];
  343. }
  344. }
  345. }
  346. }
  347. //**********************************************************************************//
  348. // This function calculates the scattering coefficients required to calculate //
  349. // both the near- and far-field parameters. //
  350. // //
  351. // Input parameters: //
  352. // L: Number of layers //
  353. // pl: Index of PEC layer. If there is none just send -1 //
  354. // x: Array containing the size parameters of the layers [0..L-1] //
  355. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  356. // n_max: Maximum number of multipolar expansion terms to be used for the //
  357. // calculations. Only used if you know what you are doing, otherwise set //
  358. // this parameter to -1 and the function will calculate it. //
  359. // //
  360. // Output parameters: //
  361. // an, bn: Complex scattering amplitudes //
  362. // //
  363. // Return value: //
  364. // Number of multipolar expansion terms used for the calculations //
  365. //**********************************************************************************//
  366. int ScattCoeff(int L, int pl, double x[], complex m[], int n_max, complex **an, complex **bn){
  367. //************************************************************************//
  368. // Calculate the index of the first layer. It can be either 0 (default) //
  369. // or the index of the outermost PEC layer. In the latter case all layers //
  370. // below the PEC are discarded. //
  371. //************************************************************************//
  372. int fl = firstLayer(L, pl);
  373. if (n_max <= 0) {
  374. n_max = Nmax(L, fl, pl, x, m);
  375. }
  376. complex z1, z2, cn;
  377. complex Num, Denom;
  378. complex G1, G2;
  379. complex Temp;
  380. double Tmp;
  381. int n, l, t;
  382. //**************************************************************************//
  383. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  384. // means that index = layer number - 1 or index = n - 1. The only exception //
  385. // are the arrays for representing D1, D3 and Q because they need a value //
  386. // for the index 0 (zero), hence it is important to consider this shift //
  387. // between different arrays. The change was done to optimize memory usage. //
  388. //**************************************************************************//
  389. // Allocate memory to the arrays
  390. complex **D1_mlxl = (complex **) malloc(L*sizeof(complex *));
  391. complex **D1_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  392. complex **D3_mlxl = (complex **) malloc(L*sizeof(complex *));
  393. complex **D3_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  394. complex **Q = (complex **) malloc(L*sizeof(complex *));
  395. complex **Ha = (complex **) malloc(L*sizeof(complex *));
  396. complex **Hb = (complex **) malloc(L*sizeof(complex *));
  397. for (l = 0; l < L; l++) {
  398. D1_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  399. D1_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  400. D3_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  401. D3_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  402. Q[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  403. Ha[l] = (complex *) malloc(n_max*sizeof(complex));
  404. Hb[l] = (complex *) malloc(n_max*sizeof(complex));
  405. }
  406. (*an) = (complex *) malloc(n_max*sizeof(complex));
  407. (*bn) = (complex *) malloc(n_max*sizeof(complex));
  408. complex *D1XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  409. complex *D3XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  410. complex *PsiXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  411. complex *ZetaXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  412. //*************************************************//
  413. // Calculate D1 and D3 for z1 in the first layer //
  414. //*************************************************//
  415. if (fl == pl) { // PEC layer
  416. for (n = 0; n <= n_max; n++) {
  417. D1_mlxl[fl][n] = Complex(0, -1);
  418. D3_mlxl[fl][n] = C_I;
  419. }
  420. } else { // Regular layer
  421. z1 = RCmul(x[fl], m[fl]);
  422. // Calculate D1 and D3
  423. calcD1D3(z1, n_max, &(D1_mlxl[fl]), &(D3_mlxl[fl]));
  424. }
  425. //******************************************************************//
  426. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  427. //******************************************************************//
  428. for (n = 0; n < n_max; n++) {
  429. Ha[fl][n] = D1_mlxl[fl][n + 1];
  430. Hb[fl][n] = D1_mlxl[fl][n + 1];
  431. }
  432. //*****************************************************//
  433. // Iteration from the second layer to the last one (L) //
  434. //*****************************************************//
  435. for (l = fl + 1; l < L; l++) {
  436. //************************************************************//
  437. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  438. //************************************************************//
  439. z1 = RCmul(x[l], m[l]);
  440. z2 = RCmul(x[l - 1], m[l]);
  441. //Calculate D1 and D3 for z1
  442. calcD1D3(z1, n_max, &(D1_mlxl[l]), &(D3_mlxl[l]));
  443. //Calculate D1 and D3 for z2
  444. calcD1D3(z2, n_max, &(D1_mlxlM1[l]), &(D3_mlxlM1[l]));
  445. //*********************************************//
  446. //Calculate Q, Ha and Hb in the layers fl+1..L //
  447. //*********************************************//
  448. // Upward recurrence for Q - equations (19a) and (19b)
  449. Num = RCmul(exp(-2*(z1.i - z2.i)), Complex(cos(-2*z2.r) - exp(-2*z2.i), sin(-2*z2.r)));
  450. Denom = Complex(cos(-2*z1.r) - exp(-2*z1.i), sin(-2*z1.r));
  451. Q[l][0] = Cdiv(Num, Denom);
  452. for (n = 1; n <= n_max; n++) {
  453. cn = Complex(n, 0);
  454. Num = Cmul(Cadd(Cmul(z1, D1_mlxl[l][n]), cn), Csub(cn, Cmul(z1, D3_mlxl[l][n - 1])));
  455. Denom = Cmul(Cadd(Cmul(z2, D1_mlxlM1[l][n]), cn), Csub(cn, Cmul(z2, D3_mlxlM1[l][n - 1])));
  456. Q[l][n] = Cdiv(Cmul(RCmul((x[l - 1]*x[l - 1])/(x[l]*x[l]), Q[l][n - 1]), Num), Denom);
  457. }
  458. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  459. for (n = 1; n <= n_max; n++) {
  460. //Ha
  461. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  462. G1 = RCmul(-1.0, D1_mlxlM1[l][n]);
  463. G2 = RCmul(-1.0, D3_mlxlM1[l][n]);
  464. } else {
  465. G1 = Csub(Cmul(m[l], Ha[l - 1][n - 1]), Cmul(m[l - 1], D1_mlxlM1[l][n]));
  466. G2 = Csub(Cmul(m[l], Ha[l - 1][n - 1]), Cmul(m[l - 1], D3_mlxlM1[l][n]));
  467. }
  468. Temp = Cmul(Q[l][n], G1);
  469. Num = Csub(Cmul(G2, D1_mlxl[l][n]), Cmul(Temp, D3_mlxl[l][n]));
  470. Denom = Csub(G2, Temp);
  471. Ha[l][n - 1] = Cdiv(Num, Denom);
  472. //Hb
  473. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  474. G1 = Hb[l - 1][n - 1];
  475. G2 = Hb[l - 1][n - 1];
  476. } else {
  477. G1 = Csub(Cmul(m[l - 1], Hb[l - 1][n - 1]), Cmul(m[l], D1_mlxlM1[l][n]));
  478. G2 = Csub(Cmul(m[l - 1], Hb[l - 1][n - 1]), Cmul(m[l], D3_mlxlM1[l][n]));
  479. }
  480. Temp = Cmul(Q[l][n], G1);
  481. Num = Csub(Cmul(G2, D1_mlxl[l][n]), Cmul(Temp, D3_mlxl[l][n]));
  482. Denom = Csub(G2, Temp);
  483. Hb[l][n - 1] = Cdiv(Num, Denom);
  484. }
  485. }
  486. //**************************************//
  487. //Calculate D1, D3, Psi and Zeta for XL //
  488. //**************************************//
  489. z1 = Complex(x[L - 1], 0);
  490. // Calculate D1XL and D3XL
  491. calcD1D3(z1, n_max, &D1XL, &D3XL);
  492. // Calculate PsiXL and ZetaXL
  493. calcPsiZeta(z1, n_max, D1XL, D3XL, &PsiXL, &ZetaXL);
  494. //*********************************************************************//
  495. // Finally, we calculate the scattering coefficients (an and bn) and //
  496. // the angular functions (Pi and Tau). Note that for these arrays the //
  497. // first layer is 0 (zero), in future versions all arrays will follow //
  498. // this convention to save memory. (13 Nov, 2014) //
  499. //*********************************************************************//
  500. for (n = 0; n < n_max; n++) {
  501. //********************************************************************//
  502. //Expressions for calculating an and bn coefficients are not valid if //
  503. //there is only one PEC layer (ie, for a simple PEC sphere). //
  504. //********************************************************************//
  505. if (pl < (L - 1)) {
  506. (*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]);
  507. (*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]);
  508. } else {
  509. (*an)[n] = calc_an(n + 1, x[L - 1], C_ZERO, C_ONE, PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  510. (*bn)[n] = Cdiv(PsiXL[n + 1], ZetaXL[n + 1]);
  511. }
  512. }
  513. // Free the memory used for the arrays
  514. for (l = 0; l < L; l++) {
  515. free(D1_mlxl[l]);
  516. free(D1_mlxlM1[l]);
  517. free(D3_mlxl[l]);
  518. free(D3_mlxlM1[l]);
  519. free(Q[l]);
  520. free(Ha[l]);
  521. free(Hb[l]);
  522. }
  523. free(D1_mlxl);
  524. free(D1_mlxlM1);
  525. free(D3_mlxl);
  526. free(D3_mlxlM1);
  527. free(Q);
  528. free(Ha);
  529. free(Hb);
  530. free(D1XL);
  531. free(D3XL);
  532. free(PsiXL);
  533. free(ZetaXL);
  534. return n_max;
  535. }
  536. //**********************************************************************************//
  537. //**********************************************************************************//
  538. //**********************************************************************************//
  539. int ScattCoeff_std(double x[], complex m[],
  540. int L, int pl, std::vector<double> x_std,
  541. std::vector<std::complex<double> > m_std, int n_max,
  542. std::vector< std::complex<double> > &an_std,
  543. std::vector< std::complex<double> > &bn_std){
  544. //************************************************************************//
  545. // Calculate the index of the first layer. It can be either 0 (default) //
  546. // or the index of the outermost PEC layer. In the latter case all layers //
  547. // below the PEC are discarded. //
  548. //************************************************************************//
  549. //TODO: Why?
  550. // int fl = firstLayer(L, pl);
  551. // instead of
  552. int fl = (pl > 0) ? pl : 0;
  553. // fl - first layer, pl - pec layer.
  554. if (n_max <= 0) {
  555. int tmp_n_max = Nmax(L, fl, pl, x, m);
  556. n_max = Nmax_std(L, fl, pl, x_std, m_std);
  557. if (n_max != tmp_n_max) printf("n_max mismatch 2\n");
  558. }
  559. // complex z1, z2, cn;
  560. // complex Num, Denom;
  561. // complex G1, G2;
  562. // complex Temp;
  563. std::complex<double> z1, z2, cn;
  564. std::complex<double> Num, Denom;
  565. std::complex<double> G1, G2;
  566. std::complex<double> Temp;
  567. double Tmp;
  568. int n, l, t;
  569. //**************************************************************************//
  570. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  571. // means that index = layer number - 1 or index = n - 1. The only exception //
  572. // are the arrays for representing D1, D3 and Q because they need a value //
  573. // for the index 0 (zero), hence it is important to consider this shift //
  574. // between different arrays. The change was done to optimize memory usage. //
  575. //**************************************************************************//
  576. // Allocate memory to the arrays
  577. //complex **D1_mlxl = (complex **) malloc(L*sizeof(complex *));
  578. //complex **D1_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  579. std::vector< std::vector< std::complex<double> > > D1_mlxl;
  580. D1_mlxl.resize(L);
  581. std::vector< std::vector< std::complex<double> > >D1_mlxlM1;
  582. D1_mlxlM1.resize(L);
  583. // complex **D3_mlxl = (complex **) malloc(L*sizeof(complex *));
  584. // complex **D3_mlxlM1 = (complex **) malloc(L*sizeof(complex *));
  585. std::vector< std::vector< std::complex<double> > > D3_mlxl;
  586. D3_mlxl.resize(L);
  587. std::vector< std::vector< std::complex<double> > > D3_mlxlM1;
  588. D3_mlxlM1.resize(L);
  589. //complex **Q = (complex **) malloc(L*sizeof(complex *));
  590. std::vector< std::vector< std::complex<double> > > Q;
  591. Q.resize(L);
  592. //complex **Ha = (complex **) malloc(L*sizeof(complex *));
  593. std::vector< std::vector< std::complex<double> > > Ha;
  594. Ha.resize(L);
  595. //complex **Hb = (complex **) malloc(L*sizeof(complex *));
  596. std::vector< std::vector< std::complex<double> > > Hb;
  597. Hb.resize(L);
  598. for (l = 0; l < L; l++) {
  599. // D1_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  600. // D1_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  601. D1_mlxl[l].resize(n_max +1);
  602. D1_mlxlM1[l].resize(n_max +1);
  603. // D3_mlxl[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  604. // D3_mlxlM1[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  605. D3_mlxl[l].resize(n_max +1);
  606. D3_mlxlM1[l].resize(n_max +1);
  607. //Q[l] = (complex *) malloc((n_max + 1)*sizeof(complex));
  608. Q[l].resize(n_max + 1);
  609. // Ha[l] = (complex *) malloc(n_max*sizeof(complex));
  610. // Hb[l] = (complex *) malloc(n_max*sizeof(complex));
  611. Ha[l].resize(n_max);
  612. Hb[l].resize(n_max);
  613. }
  614. // (*an) = (complex *) malloc(n_max*sizeof(complex));
  615. // (*bn) = (complex *) malloc(n_max*sizeof(complex));
  616. an_std.resize(n_max);
  617. bn_std.resize(n_max);
  618. // complex *D1XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  619. // complex *D3XL = (complex *) malloc((n_max + 1)*sizeof(complex));
  620. std::vector< std::complex<double> > D1XL;
  621. D1XL.resize(n_max+1);
  622. std::vector< std::complex<double> > D3XL;
  623. D3XL.resize(n_max+1);
  624. // complex *PsiXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  625. // complex *ZetaXL = (complex *) malloc((n_max + 1)*sizeof(complex));
  626. std::vector< std::complex<double> > PsiXL;
  627. PsiXL.resize(n_max+1);
  628. std::vector< std::complex<double> > ZetaXL;
  629. ZetaXL.resize(n_max+1);
  630. //*************************************************//
  631. // Calculate D1 and D3 for z1 in the first layer //
  632. //*************************************************//
  633. if (fl == pl) { // PEC layer
  634. for (n = 0; n <= n_max; n++) {
  635. D1_mlxl[fl][n] = std::complex<double>(0, -1);
  636. D3_mlxl[fl][n] = std::complex<double>(0, 1);
  637. }
  638. } else { // Regular layer
  639. z1 = x_std[fl]* m_std[fl];
  640. // Calculate D1 and D3
  641. calcD1D3_std(z1, n_max, D1_mlxl[fl], D3_mlxl[fl]);
  642. }
  643. //******************************************************************//
  644. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  645. //******************************************************************//
  646. for (n = 0; n < n_max; n++) {
  647. Ha[fl][n] = D1_mlxl[fl][n + 1];
  648. Hb[fl][n] = D1_mlxl[fl][n + 1];
  649. }
  650. //*****************************************************//
  651. // Iteration from the second layer to the last one (L) //
  652. //*****************************************************//
  653. for (l = fl + 1; l < L; l++) {
  654. //************************************************************//
  655. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  656. //************************************************************//
  657. z1 = x_std[l] * m_std[l];
  658. z2 = x_std[l - 1] * m_std[l];
  659. //Calculate D1 and D3 for z1
  660. calcD1D3_std(z1, n_max, D1_mlxl[l], D3_mlxl[l]);
  661. //Calculate D1 and D3 for z2
  662. calcD1D3_std(z2, n_max, D1_mlxlM1[l], D3_mlxlM1[l]);
  663. //*********************************************//
  664. //Calculate Q, Ha and Hb in the layers fl+1..L //
  665. //*********************************************//
  666. // Upward recurrence for Q - equations (19a) and (19b)
  667. //Num = RCmul(exp(-2*(z1.i - z2.i)), Complex(cos(-2*z2.r) - exp(-2*z2.i), sin(-2*z2.r)));
  668. Num = exp( -2.0 * ( z1.imag() - z2.imag() ) ) *
  669. std::complex<double>(cos(-2.0*z2.real()) - exp(-2.0*z2.imag()), sin(-2.0*z2.real()));
  670. Denom = std::complex<double>(cos(-2.0*z1.real()) - exp(-2.0*z1.imag()),
  671. sin(-2.0*z1.real()));
  672. Q[l][0] = Num/Denom;
  673. for (n = 1; n <= n_max; n++) {
  674. cn = std::complex<double>(n, 0);
  675. Num = ( (z1*D1_mlxl[l][n]) + cn) * (cn - (z1*D3_mlxl[l][n - 1]) );
  676. Denom = ( (z2*D1_mlxlM1[l][n]) + cn) * (cn- (z2* D3_mlxlM1[l][n - 1]) );
  677. Q[l][n] = (((x[l - 1]*x[l - 1])/(x[l]*x[l])* Q[l][n - 1]) * Num)/Denom;
  678. }
  679. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  680. for (n = 1; n <= n_max; n++) {
  681. //Ha
  682. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  683. G1 = -1.0 * D1_mlxlM1[l][n];
  684. G2 = -1.0 * D3_mlxlM1[l][n];
  685. } else {
  686. G1 = (m_std[l] * Ha[l - 1][n - 1]) - (m_std[l - 1] * D1_mlxlM1[l][n]);
  687. G2 = (m_std[l] * Ha[l - 1][n - 1]) - (m_std[l - 1] * D3_mlxlM1[l][n]);
  688. }
  689. Temp = Q[l][n] * G1;
  690. Num = (G2*D1_mlxl[l][n]) - (Temp*D3_mlxl[l][n]);
  691. Denom = G2 - Temp;
  692. Ha[l][n - 1] = Num / Denom;
  693. //Hb
  694. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  695. G1 = Hb[l - 1][n - 1];
  696. G2 = Hb[l - 1][n - 1];
  697. } else {
  698. G1 = (m_std[l - 1] * Hb[l - 1][n - 1]) - (m_std[l] * D1_mlxlM1[l][n]);
  699. G2 = (m_std[l - 1] * Hb[l - 1][n - 1]) - (m_std[l] * D3_mlxlM1[l][n]);
  700. }
  701. Temp = Q[l][n] * G1;
  702. Num = (G2*D1_mlxl[l][n]) - (Temp* D3_mlxl[l][n]);
  703. Denom = (G2- Temp);
  704. Hb[l][n - 1] = (Num/ Denom);
  705. }
  706. }
  707. //**************************************//
  708. //Calculate D1, D3, Psi and Zeta for XL //
  709. //**************************************//
  710. z1 = std::complex<double>(x_std[L - 1], 0);
  711. // Calculate D1XL and D3XL
  712. calcD1D3_std(z1, n_max, D1XL, D3XL);
  713. // Calculate PsiXL and ZetaXL
  714. calcPsiZeta_std(z1, n_max, D1XL, D3XL, PsiXL, ZetaXL);
  715. //*********************************************************************//
  716. // Finally, we calculate the scattering coefficients (an and bn) and //
  717. // the angular functions (Pi and Tau). Note that for these arrays the //
  718. // first layer is 0 (zero), in future versions all arrays will follow //
  719. // this convention to save memory. (13 Nov, 2014) //
  720. //*********************************************************************//
  721. for (n = 0; n < n_max; n++) {
  722. //********************************************************************//
  723. //Expressions for calculating an and bn coefficients are not valid if //
  724. //there is only one PEC layer (ie, for a simple PEC sphere). //
  725. //********************************************************************//
  726. if (pl < (L - 1)) {
  727. an_std[n] = calc_an_std(n + 1, x_std[L-1], Ha[L-1][n], m_std[L-1],
  728. PsiXL[n + 1], ZetaXL[n+1], PsiXL[n], ZetaXL[n]);
  729. bn_std[n] = calc_bn_std(n + 1, x_std[L-1], Hb[L-1][n], m_std[L-1],
  730. PsiXL[n+1], ZetaXL[n+1], PsiXL[n], ZetaXL[n]);
  731. } else {
  732. an_std[n] = calc_an_std(n+1, x_std[L-1], std::complex<double>(0.0,0.0),
  733. std::complex<double>(1.0,0.0),
  734. PsiXL[n+1], ZetaXL[n+1], PsiXL[n], ZetaXL[n]);
  735. bn_std[n] = PsiXL[n+1] / ZetaXL[n+1];
  736. }
  737. }
  738. return n_max;
  739. }
  740. //**********************************************************************************//
  741. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  742. // parameters, it is here mainly for compatibility with older versions of the //
  743. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  744. // any limit for the maximum number of terms. //
  745. //**********************************************************************************//
  746. 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[]) {
  747. return nMieScatt(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  748. }
  749. int nMie_std(double x[], complex m[], double Theta[], complex S1[], complex S2[],
  750. 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) {
  751. 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);
  752. }
  753. //**********************************************************************************//
  754. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  755. // parameters, it is useful if you want to include a PEC layer but not a limit //
  756. // for the maximum number of terms. //
  757. // //
  758. // Input parameters: //
  759. // L: Number of layers //
  760. // pl: Index of PEC layer. If there is none just send -1 //
  761. // x: Array containing the size parameters of the layers [0..L-1] //
  762. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  763. // nTheta: Number of scattering angles //
  764. // Theta: Array containing all the scattering angles where the scattering //
  765. // amplitudes will be calculated //
  766. // //
  767. // Output parameters: //
  768. // Qext: Efficiency factor for extinction //
  769. // Qsca: Efficiency factor for scattering //
  770. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  771. // Qbk: Efficiency factor for backscattering //
  772. // Qpr: Efficiency factor for the radiation pressure //
  773. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  774. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  775. // S1, S2: Complex scattering amplitudes //
  776. // //
  777. // Return value: //
  778. // Number of multipolar expansion terms used for the calculations //
  779. //**********************************************************************************//
  780. 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[]) {
  781. return nMieScatt(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  782. }
  783. //**********************************************************************************//
  784. // This function is just a wrapper to call the function 'nMieScatt' with fewer //
  785. // parameters, it is useful if you want to include a limit for the maximum number //
  786. // of terms but not a PEC layer. //
  787. // //
  788. // Input parameters: //
  789. // L: Number of layers //
  790. // x: Array containing the size parameters of the layers [0..L-1] //
  791. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  792. // nTheta: Number of scattering angles //
  793. // Theta: Array containing all the scattering angles where the scattering //
  794. // amplitudes will be calculated //
  795. // n_max: Maximum number of multipolar expansion terms to be used for the //
  796. // calculations. Only used if you know what you are doing, otherwise set //
  797. // this parameter to -1 and the function will calculate it //
  798. // //
  799. // Output parameters: //
  800. // Qext: Efficiency factor for extinction //
  801. // Qsca: Efficiency factor for scattering //
  802. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  803. // Qbk: Efficiency factor for backscattering //
  804. // Qpr: Efficiency factor for the radiation pressure //
  805. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  806. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  807. // S1, S2: Complex scattering amplitudes //
  808. // //
  809. // Return value: //
  810. // Number of multipolar expansion terms used for the calculations //
  811. //**********************************************************************************//
  812. 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[]) {
  813. return nMieScatt(L, -1, x, m, nTheta, Theta, n_max, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  814. }
  815. //**********************************************************************************//
  816. // This function calculates the actual scattering parameters and amplitudes //
  817. // //
  818. // Input parameters: //
  819. // L: Number of layers //
  820. // pl: Index of PEC layer. If there is none just send -1 //
  821. // x: Array containing the size parameters of the layers [0..L-1] //
  822. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  823. // nTheta: Number of scattering angles //
  824. // Theta: Array containing all the scattering angles where the scattering //
  825. // amplitudes will be calculated //
  826. // n_max: Maximum number of multipolar expansion terms to be used for the //
  827. // calculations. Only used if you know what you are doing, otherwise set //
  828. // this parameter to -1 and the function will calculate it //
  829. // //
  830. // Output parameters: //
  831. // Qext: Efficiency factor for extinction //
  832. // Qsca: Efficiency factor for scattering //
  833. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  834. // Qbk: Efficiency factor for backscattering //
  835. // Qpr: Efficiency factor for the radiation pressure //
  836. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  837. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  838. // S1, S2: Complex scattering amplitudes //
  839. // //
  840. // Return value: //
  841. // Number of multipolar expansion terms used for the calculations //
  842. //**********************************************************************************//
  843. 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[]) {
  844. int i, n, t;
  845. double **Pi, **Tau;
  846. complex *an, *bn;
  847. complex Qbktmp;
  848. n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  849. Pi = (double **) malloc(n_max*sizeof(double *));
  850. Tau = (double **) malloc(n_max*sizeof(double *));
  851. for (n = 0; n < n_max; n++) {
  852. Pi[n] = (double *) malloc(nTheta*sizeof(double));
  853. Tau[n] = (double *) malloc(nTheta*sizeof(double));
  854. }
  855. calcPiTau(n_max, nTheta, Theta, &Pi, &Tau);
  856. double x2 = x[L - 1]*x[L - 1];
  857. // Initialize the scattering parameters
  858. *Qext = 0;
  859. *Qsca = 0;
  860. *Qabs = 0;
  861. *Qbk = 0;
  862. Qbktmp = C_ZERO;
  863. *Qpr = 0;
  864. *g = 0;
  865. *Albedo = 0;
  866. // Initialize the scattering amplitudes
  867. for (t = 0; t < nTheta; t++) {
  868. S1[t] = C_ZERO;
  869. S2[t] = C_ZERO;
  870. }
  871. // By using downward recurrence we avoid loss of precision due to float rounding errors
  872. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  873. // http://en.wikipedia.org/wiki/Loss_of_significance
  874. for (i = n_max - 2; i >= 0; i--) {
  875. n = i + 1;
  876. // Equation (27)
  877. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  878. // Equation (28)
  879. *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);
  880. // Equation (29)
  881. *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));
  882. // Equation (33)
  883. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  884. //****************************************************//
  885. // Calculate the scattering amplitudes (S1 and S2) //
  886. // Equations (25a) - (25b) //
  887. //****************************************************//
  888. for (t = 0; t < nTheta; t++) {
  889. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  890. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  891. }
  892. }
  893. *Qext = 2*(*Qext)/x2; // Equation (27)
  894. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  895. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  896. *Qabs = *Qext - *Qsca; // Equation (30)
  897. *Albedo = *Qsca / *Qext; // Equation (31)
  898. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  899. *Qbk = (Qbktmp.r*Qbktmp.r + Qbktmp.i*Qbktmp.i)/x2; // Equation (33)
  900. // Free the memory used for the arrays
  901. for (n = 0; n < n_max; n++) {
  902. free(Pi[n]);
  903. free(Tau[n]);
  904. }
  905. free(Pi);
  906. free(Tau);
  907. free(an);
  908. free(bn);
  909. return n_max;
  910. }
  911. int nMieScatt_std(double x[], complex m[], double Theta[], complex S1[], complex S2[],
  912. int L, int pl,
  913. std::vector<double> &x_std, std::vector<std::complex<double> > &m_std,
  914. int nTheta, std::vector<double> &Theta_std,
  915. int n_max, double *Qext, double *Qsca, double *Qabs, double *Qbk,
  916. double *Qpr, double *g, double *Albedo,
  917. std::vector< std::complex<double> > &S1_std,
  918. std::vector< std::complex<double> > &S2_std) {
  919. int i, n, t;
  920. double **Pi, **Tau;
  921. std::vector< std::vector<double> > Pi_std, Tau_std;
  922. complex *an, *bn;
  923. std::vector< std::complex<double> > an_std, bn_std;
  924. complex Qbktmp;
  925. std::complex<double> Qbktmp_std;
  926. {
  927. int tmp_n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  928. n_max = ScattCoeff_std(x, m, L, pl, x_std, m_std, n_max, an_std, bn_std);
  929. if (n_max != tmp_n_max) printf("n_max mismatch\n");
  930. // TODO: Check an, bn and an_std, bn_std are equal
  931. compare("an vs an_std: ", an, an_std);
  932. }
  933. Pi = (double **) malloc(n_max*sizeof(double *));
  934. Tau = (double **) malloc(n_max*sizeof(double *));
  935. for (n = 0; n < n_max; n++) {
  936. Pi[n] = (double *) malloc(nTheta*sizeof(double));
  937. Tau[n] = (double *) malloc(nTheta*sizeof(double));
  938. }
  939. calcPiTau(n_max, nTheta, Theta, &Pi, &Tau);
  940. double x2 = x[L - 1]*x[L - 1];
  941. // Initialize the scattering parameters
  942. *Qext = 0;
  943. *Qsca = 0;
  944. *Qabs = 0;
  945. *Qbk = 0;
  946. Qbktmp = C_ZERO;
  947. *Qpr = 0;
  948. *g = 0;
  949. *Albedo = 0;
  950. // Initialize the scattering amplitudes
  951. for (t = 0; t < nTheta; t++) {
  952. S1[t] = C_ZERO;
  953. S2[t] = C_ZERO;
  954. }
  955. // By using downward recurrence we avoid loss of precision due to float rounding errors
  956. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  957. // http://en.wikipedia.org/wiki/Loss_of_significance
  958. for (i = n_max - 2; i >= 0; i--) {
  959. n = i + 1;
  960. // Equation (27)
  961. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  962. // Equation (28)
  963. *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);
  964. // Equation (29)
  965. *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));
  966. // Equation (33)
  967. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  968. //****************************************************//
  969. // Calculate the scattering amplitudes (S1 and S2) //
  970. // Equations (25a) - (25b) //
  971. //****************************************************//
  972. for (t = 0; t < nTheta; t++) {
  973. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  974. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  975. }
  976. }
  977. *Qext = 2*(*Qext)/x2; // Equation (27)
  978. *Qsca = 2*(*Qsca)/x2; // Equation (28)
  979. *Qpr = *Qext - 4*(*Qpr)/x2; // Equation (29)
  980. *Qabs = *Qext - *Qsca; // Equation (30)
  981. *Albedo = *Qsca / *Qext; // Equation (31)
  982. *g = (*Qext - *Qpr) / *Qsca; // Equation (32)
  983. *Qbk = (Qbktmp.r*Qbktmp.r + Qbktmp.i*Qbktmp.i)/x2; // Equation (33)
  984. // Free the memory used for the arrays
  985. for (n = 0; n < n_max; n++) {
  986. free(Pi[n]);
  987. free(Tau[n]);
  988. }
  989. free(Pi);
  990. free(Tau);
  991. free(an);
  992. free(bn);
  993. return n_max;
  994. }
  995. //**********************************************************************************//
  996. // This function calculates complex electric and magnetic field in the surroundings //
  997. // and inside (TODO) the particle. //
  998. // //
  999. // Input parameters: //
  1000. // L: Number of layers //
  1001. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1002. // x: Array containing the size parameters of the layers [0..L-1] //
  1003. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1004. // n_max: Maximum number of multipolar expansion terms to be used for the //
  1005. // calculations. Only used if you know what you are doing, otherwise set //
  1006. // this parameter to 0 (zero) and the function will calculate it. //
  1007. // nCoords: Number of coordinate points //
  1008. // Coords: Array containing all coordinates where the complex electric and //
  1009. // magnetic fields will be calculated //
  1010. // //
  1011. // Output parameters: //
  1012. // E, H: Complex electric and magnetic field at the provided coordinates //
  1013. // //
  1014. // Return value: //
  1015. // Number of multipolar expansion terms used for the calculations //
  1016. //**********************************************************************************//
  1017. 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[]){
  1018. int i, n, c;
  1019. double **Pi, **Tau;
  1020. complex *an, *bn;
  1021. double *Rho = (double *) malloc(nCoords*sizeof(double));
  1022. double *Phi = (double *) malloc(nCoords*sizeof(double));
  1023. double *Theta = (double *) malloc(nCoords*sizeof(double));
  1024. for (c = 0; c < nCoords; c++) {
  1025. Rho[c] = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  1026. if (Rho[c] < 1e-3) {
  1027. Rho[c] = 1e-3;
  1028. }
  1029. Phi[c] = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  1030. Theta[c] = acos(Xp[c]/Rho[c]);
  1031. }
  1032. n_max = ScattCoeff(L, pl, x, m, n_max, &an, &bn);
  1033. Pi = (double **) malloc(n_max*sizeof(double *));
  1034. Tau = (double **) malloc(n_max*sizeof(double *));
  1035. for (n = 0; n < n_max; n++) {
  1036. Pi[n] = (double *) malloc(nCoords*sizeof(double));
  1037. Tau[n] = (double *) malloc(nCoords*sizeof(double));
  1038. }
  1039. calcPiTau(n_max, nCoords, Theta, &Pi, &Tau);
  1040. double x2 = x[L - 1]*x[L - 1];
  1041. // Initialize the fields
  1042. for (c = 0; c < nCoords; c++) {
  1043. E[c] = C_ZERO;
  1044. H[c] = C_ZERO;
  1045. }
  1046. //*******************************************************//
  1047. // external scattering field = incident + scattered //
  1048. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1049. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1050. //*******************************************************//
  1051. // Firstly the easiest case, we want the field outside the particle
  1052. if (Rho[c] >= x[L - 1]) {
  1053. }
  1054. // for (i = 1; i < (n_max - 1); i++) {
  1055. // n = i - 1;
  1056. /* // Equation (27)
  1057. *Qext = *Qext + (double)(n + n + 1)*(an[i].r + bn[i].r);
  1058. // Equation (28)
  1059. *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);
  1060. // Equation (29)
  1061. *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));
  1062. // Equation (33)
  1063. Qbktmp = Cadd(Qbktmp, RCmul((double)((n + n + 1)*(1 - 2*(n % 2))), Csub(an[i], bn[i])));
  1064. */
  1065. //****************************************************//
  1066. // Calculate the scattering amplitudes (S1 and S2) //
  1067. // Equations (25a) - (25b) //
  1068. //****************************************************//
  1069. /* for (t = 0; t < nTheta; t++) {
  1070. S1[t] = Cadd(S1[t], calc_S1_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  1071. S2[t] = Cadd(S2[t], calc_S2_n(n, an[i], bn[i], Pi[i][t], Tau[i][t]));
  1072. }*/
  1073. // }
  1074. // Free the memory used for the arrays
  1075. for (n = 0; n < n_max; n++) {
  1076. free(Pi[n]);
  1077. free(Tau[n]);
  1078. }
  1079. free(Pi);
  1080. free(Tau);
  1081. free(an);
  1082. free(bn);
  1083. free(Rho);
  1084. free(Phi);
  1085. free(Theta);
  1086. return n_max;
  1087. }