nmie.cc 57 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210
  1. ///
  2. /// @file nmie.cc
  3. /// @author Ladutenko Konstantin <kostyfisik at gmail (.) com>
  4. /// @date Tue Sep 3 00:38:27 2013
  5. /// @copyright 2013 Ladutenko Konstantin
  6. ///
  7. /// nmie is free software: you can redistribute it and/or modify
  8. /// it under the terms of the GNU General Public License as published by
  9. /// the Free Software Foundation, either version 3 of the License, or
  10. /// (at your option) any later version.
  11. ///
  12. /// nmie-wrapper is distributed in the hope that it will be useful,
  13. /// but WITHOUT ANY WARRANTY; without even the implied warranty of
  14. /// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  15. /// GNU General Public License for more details.
  16. ///
  17. /// You should have received a copy of the GNU General Public License
  18. /// along with nmie-wrapper. If not, see <http://www.gnu.org/licenses/>.
  19. ///
  20. /// nmie uses nmie.c from scattnlay by Ovidio Pena
  21. /// <ovidio@bytesfall.com> . He has an additional condition to
  22. /// his library:
  23. // The only additional condition is that we expect that all publications //
  24. // describing work using this software , or all commercial products //
  25. // using it, cite the following reference: //
  26. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  27. // a multilayered sphere," Computer Physics Communications, //
  28. // vol. 180, Nov. 2009, pp. 2348-2354. //
  29. ///
  30. /// @brief Wrapper class around nMie function for ease of use
  31. ///
  32. #include "nmie.h"
  33. #include <array>
  34. #include <algorithm>
  35. #include <cstdio>
  36. #include <cstdlib>
  37. #include <stdexcept>
  38. #include <vector>
  39. namespace nmie {
  40. //helpers
  41. template<class T> inline T pow2(const T value) {return value*value;}
  42. //#define round(x) ((x) >= 0 ? (int)((x) + 0.5):(int)((x) - 0.5))
  43. int round(double x) {
  44. return x >= 0 ? (int)(x + 0.5):(int)(x - 0.5);
  45. }
  46. // ********************************************************************** //
  47. // ********************************************************************** //
  48. // ********************************************************************** //
  49. //emulate C call.
  50. int nMie_wrapper(int L, const std::vector<double>& x, const std::vector<std::complex<double> >& m,
  51. int nTheta, const std::vector<double>& Theta,
  52. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  53. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  54. if (x.size() != L || m.size() != L)
  55. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  56. if (Theta.size() != nTheta)
  57. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  58. try {
  59. MultiLayerMie multi_layer_mie;
  60. multi_layer_mie.SetWidthSP(x);
  61. multi_layer_mie.SetIndexSP(m);
  62. multi_layer_mie.SetAngles(Theta);
  63. multi_layer_mie.RunMieCalculations();
  64. *Qext = multi_layer_mie.GetQext();
  65. *Qsca = multi_layer_mie.GetQsca();
  66. *Qabs = multi_layer_mie.GetQabs();
  67. *Qbk = multi_layer_mie.GetQbk();
  68. *Qpr = multi_layer_mie.GetQpr();
  69. *g = multi_layer_mie.GetAsymmetryFactor();
  70. *Albedo = multi_layer_mie.GetAlbedo();
  71. S1 = multi_layer_mie.GetS1();
  72. S2 = multi_layer_mie.GetS2();
  73. multi_layer_mie.GetFailed();
  74. } catch( const std::invalid_argument& ia ) {
  75. // Will catch if multi_layer_mie fails or other errors.
  76. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  77. throw std::invalid_argument(ia);
  78. return -1;
  79. }
  80. return 0;
  81. }
  82. // ********************************************************************** //
  83. // ********************************************************************** //
  84. // ********************************************************************** //
  85. void MultiLayerMie::GetFailed() {
  86. double faild_x = 9.42477796076938;
  87. //double faild_x = 9.42477796076937;
  88. std::complex<double> z(faild_x, 0.0);
  89. std::vector<int> nmax_local_array = {20, 100, 500, 2500};
  90. for (auto nmax_local : nmax_local_array) {
  91. std::vector<std::complex<double> > D1_failed(nmax_local +1);
  92. // Downward recurrence for D1 - equations (16a) and (16b)
  93. D1_failed[nmax_local] = std::complex<double>(0.0, 0.0);
  94. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  95. for (int n = nmax_local; n > 0; n--) {
  96. D1_failed[n - 1] = double(n)*zinv - 1.0/(D1_failed[n] + double(n)*zinv);
  97. }
  98. printf("Faild D1[0] from reccurence (z = %16.14f, nmax = %d): %g\n",
  99. faild_x, nmax_local, D1_failed[0].real());
  100. }
  101. printf("Faild D1[0] from continued fraction (z = %16.14f): %g\n", faild_x,
  102. calcD1confra(0,z).real());
  103. //D1[nmax_] = calcD1confra(nmax_, z);
  104. }
  105. // ********************************************************************** //
  106. // ********************************************************************** //
  107. // ********************************************************************** //
  108. double MultiLayerMie::GetQext() {
  109. if (!isMieCalculated_)
  110. throw std::invalid_argument("You should run calculations before result reques!");
  111. return Qext_;
  112. }
  113. // ********************************************************************** //
  114. // ********************************************************************** //
  115. // ********************************************************************** //
  116. double MultiLayerMie::GetQabs() {
  117. if (!isMieCalculated_)
  118. throw std::invalid_argument("You should run calculations before result reques!");
  119. return Qabs_;
  120. }
  121. // ********************************************************************** //
  122. // ********************************************************************** //
  123. // ********************************************************************** //
  124. double MultiLayerMie::GetQsca() {
  125. if (!isMieCalculated_)
  126. throw std::invalid_argument("You should run calculations before result reques!");
  127. return Qsca_;
  128. }
  129. // ********************************************************************** //
  130. // ********************************************************************** //
  131. // ********************************************************************** //
  132. double MultiLayerMie::GetQbk() {
  133. if (!isMieCalculated_)
  134. throw std::invalid_argument("You should run calculations before result reques!");
  135. return Qbk_;
  136. }
  137. // ********************************************************************** //
  138. // ********************************************************************** //
  139. // ********************************************************************** //
  140. double MultiLayerMie::GetQpr() {
  141. if (!isMieCalculated_)
  142. throw std::invalid_argument("You should run calculations before result reques!");
  143. return Qpr_;
  144. }
  145. // ********************************************************************** //
  146. // ********************************************************************** //
  147. // ********************************************************************** //
  148. double MultiLayerMie::GetAsymmetryFactor() {
  149. if (!isMieCalculated_)
  150. throw std::invalid_argument("You should run calculations before result reques!");
  151. return asymmetry_factor_;
  152. }
  153. // ********************************************************************** //
  154. // ********************************************************************** //
  155. // ********************************************************************** //
  156. double MultiLayerMie::GetAlbedo() {
  157. if (!isMieCalculated_)
  158. throw std::invalid_argument("You should run calculations before result reques!");
  159. return albedo_;
  160. }
  161. // ********************************************************************** //
  162. // ********************************************************************** //
  163. // ********************************************************************** //
  164. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  165. return S1_;
  166. }
  167. // ********************************************************************** //
  168. // ********************************************************************** //
  169. // ********************************************************************** //
  170. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  171. return S2_;
  172. }
  173. // ********************************************************************** //
  174. // ********************************************************************** //
  175. // ********************************************************************** //
  176. void MultiLayerMie::AddTargetLayer(double width, std::complex<double> layer_index) {
  177. if (width <= 0)
  178. throw std::invalid_argument("Layer width should be positive!");
  179. target_width_.push_back(width);
  180. target_index_.push_back(layer_index);
  181. } // end of void MultiLayerMie::AddTargetLayer(...)
  182. // ********************************************************************** //
  183. // ********************************************************************** //
  184. // ********************************************************************** //
  185. void MultiLayerMie::SetCoatingIndex(std::vector<std::complex<double> > index) {
  186. coating_index_.clear();
  187. for (auto value : index) coating_index_.push_back(value);
  188. } // end of void MultiLayerMie::SetCoatingIndex(std::vector<complex> index);
  189. // ********************************************************************** //
  190. // ********************************************************************** //
  191. // ********************************************************************** //
  192. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  193. isMieCalculated_ = false;
  194. theta_ = angles;
  195. // theta_.clear();
  196. // for (auto value : angles) theta_.push_back(value);
  197. } // end of SetAngles()
  198. // ********************************************************************** //
  199. // ********************************************************************** //
  200. // ********************************************************************** //
  201. void MultiLayerMie::SetCoatingWidth(std::vector<double> width) {
  202. coating_width_.clear();
  203. for (auto w : width)
  204. if (w <= 0)
  205. throw std::invalid_argument("Coating width should be positive!");
  206. else coating_width_.push_back(w);
  207. }
  208. // end of void MultiLayerMie::SetCoatingWidth(...);
  209. // ********************************************************************** //
  210. // ********************************************************************** //
  211. // ********************************************************************** //
  212. void MultiLayerMie::SetWidthSP(const std::vector<double>& size_parameter) {
  213. isMieCalculated_ = false;
  214. size_parameter_.clear();
  215. double prev_size_parameter = 0.0;
  216. for (auto layer_size_parameter : size_parameter) {
  217. if (layer_size_parameter <= 0.0)
  218. throw std::invalid_argument("Size parameter should be positive!");
  219. if (prev_size_parameter > layer_size_parameter)
  220. throw std::invalid_argument
  221. ("Size parameter for next layer should be larger than the previous one!");
  222. prev_size_parameter = layer_size_parameter;
  223. size_parameter_.push_back(layer_size_parameter);
  224. }
  225. }
  226. // end of void MultiLayerMie::SetWidthSP(...);
  227. // ********************************************************************** //
  228. // ********************************************************************** //
  229. // ********************************************************************** //
  230. void MultiLayerMie::SetIndexSP(const std::vector< std::complex<double> >& index) {
  231. isMieCalculated_ = false;
  232. //index_.clear();
  233. index_ = index;
  234. // for (auto value : index) index_.push_back(value);
  235. } // end of void MultiLayerMie::SetIndexSP(...);
  236. // ********************************************************************** //
  237. // ********************************************************************** //
  238. // ********************************************************************** //
  239. void MultiLayerMie::SetPEC(int layer_position) {
  240. if (layer_position < 0)
  241. throw std::invalid_argument("Error! Layers are numbered from 0!");
  242. PEC_layer_position_ = layer_position;
  243. }
  244. // ********************************************************************** //
  245. // ********************************************************************** //
  246. // ********************************************************************** //
  247. void MultiLayerMie::SetMaxTermsNumber(int nmax) {
  248. nmax_ = nmax;
  249. }
  250. // ********************************************************************** //
  251. // ********************************************************************** //
  252. // ********************************************************************** //
  253. void MultiLayerMie::GenerateSizeParameter() {
  254. size_parameter_.clear();
  255. double radius = 0.0;
  256. for (auto width : target_width_) {
  257. radius += width;
  258. size_parameter_.push_back(2*PI*radius / wavelength_);
  259. }
  260. for (auto width : coating_width_) {
  261. radius += width;
  262. size_parameter_.push_back(2*PI*radius / wavelength_);
  263. }
  264. total_radius_ = radius;
  265. } // end of void MultiLayerMie::GenerateSizeParameter();
  266. // ********************************************************************** //
  267. // ********************************************************************** //
  268. // ********************************************************************** //
  269. double MultiLayerMie::GetTotalRadius() {
  270. if (total_radius_ == 0) GenerateSizeParameter();
  271. return total_radius_;
  272. } // end of double MultiLayerMie::GetTotalRadius();
  273. // ********************************************************************** //
  274. // ********************************************************************** //
  275. // ********************************************************************** //
  276. std::vector< std::vector<double> >
  277. MultiLayerMie::GetSpectra(double from_WL, double to_WL, int samples) {
  278. std::vector< std::vector<double> > spectra;
  279. double step_WL = (to_WL - from_WL)/ static_cast<double>(samples);
  280. double wavelength_backup = wavelength_;
  281. long fails = 0;
  282. for (double WL = from_WL; WL < to_WL; WL += step_WL) {
  283. wavelength_ = WL;
  284. try {
  285. RunMieCalculations();
  286. } catch( const std::invalid_argument& ia ) {
  287. fails++;
  288. continue;
  289. }
  290. //printf("%3.1f ",WL);
  291. spectra.push_back(std::vector<double>({wavelength_, Qext_, Qsca_, Qabs_, Qbk_}));
  292. } // end of for each WL in spectra
  293. printf("Spectrum has %li fails\n",fails);
  294. wavelength_ = wavelength_backup;
  295. return spectra;
  296. }
  297. // ********************************************************************** //
  298. // ********************************************************************** //
  299. // ********************************************************************** //
  300. // Calculate Nstop - equation (17)
  301. //
  302. void MultiLayerMie::Nstop() {
  303. const double& xL = size_parameter_.back();
  304. if (xL <= 8) {
  305. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  306. } else if (xL <= 4200) {
  307. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  308. } else {
  309. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  310. }
  311. }
  312. // ********************************************************************** //
  313. // ********************************************************************** //
  314. // ********************************************************************** //
  315. void MultiLayerMie::Nmax(int first_layer) {
  316. int ri, riM1;
  317. const std::vector<double>& x = size_parameter_;
  318. const std::vector<std::complex<double> >& m = index_;
  319. Nstop(); // Set initial nmax_ value
  320. for (int i = first_layer; i < x.size(); i++) {
  321. if (i > PEC_layer_position_)
  322. ri = round(std::abs(x[i]*m[i]));
  323. else
  324. ri = 0;
  325. nmax_ = std::max(nmax_, ri);
  326. // first layer is pec, if pec is present
  327. if ((i > first_layer) && ((i - 1) > PEC_layer_position_))
  328. riM1 = round(std::abs(x[i - 1]* m[i]));
  329. else
  330. riM1 = 0;
  331. nmax_ = std::max(nmax_, riM1);
  332. }
  333. nmax_ += 15; // Final nmax_ value
  334. }
  335. //**********************************************************************************//
  336. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  337. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  338. // //
  339. // Input parameters: //
  340. // z: Real argument to evaluate jn and h1n //
  341. // nmax_: Maximum number of terms to calculate jn and h1n //
  342. // //
  343. // Output parameters: //
  344. // jn, h1n: Spherical Bessel and Hankel functions //
  345. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  346. // //
  347. // The implementation follows the algorithm by I.J. Thompson and A.R. Barnett, //
  348. // Comp. Phys. Comm. 47 (1987) 245-257. //
  349. // //
  350. // Complex spherical Bessel functions from n=0..nmax_-1 for z in the upper half //
  351. // plane (Im(z) > -3). //
  352. // //
  353. // j[n] = j/n(z) Regular solution: j[0]=sin(z)/z //
  354. // j'[n] = d[j/n(z)]/dz //
  355. // h1[n] = h[0]/n(z) Irregular Hankel function: //
  356. // h1'[n] = d[h[0]/n(z)]/dz h1[0] = j0(z) + i*y0(z) //
  357. // = (sin(z)-i*cos(z))/z //
  358. // = -i*exp(i*z)/z //
  359. // Using complex CF1, and trigonometric forms for n=0 solutions. //
  360. //**********************************************************************************//
  361. void MultiLayerMie::sbesjh(std::complex<double> z,
  362. std::vector<std::complex<double> >& jn,
  363. std::vector<std::complex<double> >& jnp,
  364. std::vector<std::complex<double> >& h1n,
  365. std::vector<std::complex<double> >& h1np) {
  366. const int limit = 20000;
  367. const double accur = 1.0e-12;
  368. const double tm30 = 1e-30;
  369. double absc;
  370. std::complex<double> zi, w;
  371. std::complex<double> pl, f, b, d, c, del, jn0, jndb, h1nldb, h1nbdb;
  372. absc = std::abs(std::real(z)) + std::abs(std::imag(z));
  373. if ((absc < accur) || (std::imag(z) < -3.0)) {
  374. throw std::invalid_argument("TODO add error description for condition if ((absc < accur) || (std::imag(z) < -3.0))");
  375. }
  376. zi = 1.0/z;
  377. w = zi + zi;
  378. pl = double(nmax_)*zi;
  379. f = pl + zi;
  380. b = f + f + zi;
  381. d = 0.0;
  382. c = f;
  383. for (int n = 0; n < limit; n++) {
  384. d = b - d;
  385. c = b - 1.0/c;
  386. absc = std::abs(std::real(d)) + std::abs(std::imag(d));
  387. if (absc < tm30) {
  388. d = tm30;
  389. }
  390. absc = std::abs(std::real(c)) + std::abs(std::imag(c));
  391. if (absc < tm30) {
  392. c = tm30;
  393. }
  394. d = 1.0/d;
  395. del = d*c;
  396. f = f*del;
  397. b += w;
  398. absc = std::abs(std::real(del - 1.0)) + std::abs(std::imag(del - 1.0));
  399. if (absc < accur) {
  400. // We have obtained the desired accuracy
  401. break;
  402. }
  403. }
  404. if (absc > accur) {
  405. throw std::invalid_argument("We were not able to obtain the desired accuracy");
  406. }
  407. jn[nmax_ - 1] = tm30;
  408. jnp[nmax_ - 1] = f*jn[nmax_ - 1];
  409. // Downward recursion to n=0 (N.B. Coulomb Functions)
  410. for (int n = nmax_ - 2; n >= 0; n--) {
  411. jn[n] = pl*jn[n + 1] + jnp[n + 1];
  412. jnp[n] = pl*jn[n] - jn[n + 1];
  413. pl = pl - zi;
  414. }
  415. // Calculate the n=0 Bessel Functions
  416. jn0 = zi*std::sin(z);
  417. h1n[0] = std::exp(std::complex<double>(0.0, 1.0)*z)*zi*(-std::complex<double>(0.0, 1.0));
  418. h1np[0] = h1n[0]*(std::complex<double>(0.0, 1.0) - zi);
  419. // Rescale j[n], j'[n], converting to spherical Bessel functions.
  420. // Recur h1[n], h1'[n] as spherical Bessel functions.
  421. w = 1.0/jn[0];
  422. pl = zi;
  423. for (int n = 0; n < nmax_; n++) {
  424. jn[n] = jn0*(w*jn[n]);
  425. jnp[n] = jn0*(w*jnp[n]) - zi*jn[n];
  426. if (n != 0) {
  427. h1n[n] = (pl - zi)*h1n[n - 1] - h1np[n - 1];
  428. // check if hankel is increasing (upward stable)
  429. if (std::abs(h1n[n]) < std::abs(h1n[n - 1])) {
  430. jndb = z;
  431. h1nldb = h1n[n];
  432. h1nbdb = h1n[n - 1];
  433. }
  434. pl += zi;
  435. h1np[n] = -(pl*h1n[n]) + h1n[n - 1];
  436. }
  437. }
  438. }
  439. //**********************************************************************************//
  440. // This function calculates the spherical Bessel functions (bj and by) and the //
  441. // logarithmic derivative (bd) for a given complex value z. See pag. 87 B&H. //
  442. // //
  443. // Input parameters: //
  444. // z: Complex argument to evaluate bj, by and bd //
  445. // nmax_: Maximum number of terms to calculate bj, by and bd //
  446. // //
  447. // Output parameters: //
  448. // bj, by: Spherical Bessel functions //
  449. // bd: Logarithmic derivative //
  450. //**********************************************************************************//
  451. void MultiLayerMie::sphericalBessel(std::complex<double> z,
  452. std::vector<std::complex<double> >& bj,
  453. std::vector<std::complex<double> >& by,
  454. std::vector<std::complex<double> >& bd) {
  455. std::vector<std::complex<double> > jn(nmax_), jnp(nmax_), h1n(nmax_), h1np(nmax_);
  456. sbesjh(z, jn, jnp, h1n, h1np);
  457. for (int n = 0; n < nmax_; n++) {
  458. bj[n] = jn[n];
  459. by[n] = (h1n[n] - jn[n])/std::complex<double>(0.0, 1.0);
  460. bd[n] = jnp[n]/jn[n] + 1.0/z;
  461. }
  462. }
  463. // ********************************************************************** //
  464. // ********************************************************************** //
  465. // ********************************************************************** //
  466. // Calculate an - equation (5)
  467. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  468. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  469. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  470. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  471. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  472. return Num/Denom;
  473. }
  474. // ********************************************************************** //
  475. // ********************************************************************** //
  476. // ********************************************************************** //
  477. // Calculate bn - equation (6)
  478. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  479. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  480. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  481. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  482. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  483. return Num/Denom;
  484. }
  485. // ********************************************************************** //
  486. // ********************************************************************** //
  487. // ********************************************************************** //
  488. // Calculates S1 - equation (25a)
  489. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  490. double Pi, double Tau) {
  491. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  492. }
  493. // ********************************************************************** //
  494. // ********************************************************************** //
  495. // ********************************************************************** //
  496. // Calculates S2 - equation (25b) (it's the same as (25a), just switches Pi and Tau)
  497. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  498. double Pi, double Tau) {
  499. return calc_S1(n, an, bn, Tau, Pi);
  500. }
  501. //**********************************************************************************//
  502. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  503. // real argument (x). //
  504. // Equations (20a) - (21b) //
  505. // //
  506. // Input parameters: //
  507. // x: Real argument to evaluate Psi and Zeta //
  508. // nmax: Maximum number of terms to calculate Psi and Zeta //
  509. // //
  510. // Output parameters: //
  511. // Psi, Zeta: Riccati-Bessel functions //
  512. //**********************************************************************************//
  513. void MultiLayerMie::calcPsiZeta(double x,
  514. std::vector<std::complex<double> > D1,
  515. std::vector<std::complex<double> > D3,
  516. std::vector<std::complex<double> >& Psi,
  517. std::vector<std::complex<double> >& Zeta) {
  518. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  519. Psi[0] = std::complex<double>(sin(x), 0);
  520. Zeta[0] = std::complex<double>(sin(x), -cos(x));
  521. for (int n = 1; n <= nmax_; n++) {
  522. Psi[n] = Psi[n - 1]*(n/x - D1[n - 1]);
  523. Zeta[n] = Zeta[n - 1]*(n/x - D3[n - 1]);
  524. }
  525. }
  526. //**********************************************************************************//
  527. // Function CONFRA ported from MIEV0.f (Wiscombe,1979)
  528. // Ref. to NCAR Technical Notes, Wiscombe, 1979
  529. /*
  530. c Compute Bessel function ratio A-sub-N from its
  531. c continued fraction using Lentz method
  532. c ZINV = Reciprocal of argument of A
  533. c I N T E R N A L V A R I A B L E S
  534. c ------------------------------------
  535. c CAK Term in continued fraction expansion of A (Eq. R25)
  536. c a_k
  537. c CAPT Factor used in Lentz iteration for A (Eq. R27)
  538. c T_k
  539. c CNUMER Numerator in capT ( Eq. R28A )
  540. c N_k
  541. c CDENOM Denominator in capT ( Eq. R28B )
  542. c D_k
  543. c CDTD Product of two successive denominators of capT factors
  544. c ( Eq. R34C )
  545. c xi_1
  546. c CNTN Product of two successive numerators of capT factors
  547. c ( Eq. R34B )
  548. c xi_2
  549. c EPS1 Ill-conditioning criterion
  550. c EPS2 Convergence criterion
  551. c KK Subscript k of cAk ( Eq. R25B )
  552. c k
  553. c KOUNT Iteration counter ( used to prevent infinite looping )
  554. c MAXIT Max. allowed no. of iterations
  555. c MM + 1 and - 1, alternately
  556. */
  557. std::complex<double> MultiLayerMie::calcD1confra(const int N, const std::complex<double> z) {
  558. // NTMR -> nmax_ -1 \\TODO nmax_ ?
  559. //int N = nmax_ - 1;
  560. int KK, KOUNT, MAXIT = 10000, MM;
  561. // double EPS1=1.0e-2;
  562. double EPS2=1.0e-8;
  563. std::complex<double> CAK, CAPT, CDENOM, CDTD, CNTN, CNUMER;
  564. std::complex<double> one = std::complex<double>(1.0,0.0);
  565. std::complex<double> ZINV = one/z;
  566. // c ** Eq. R25a
  567. std::complex<double> CONFRA = static_cast<std::complex<double> >(N+1)*ZINV; //debug ZINV
  568. MM = -1;
  569. KK = 2*N +3; //debug 3
  570. // c ** Eq. R25b, k=2
  571. CAK = static_cast<std::complex<double> >(MM*KK) * ZINV; //debug -3 ZINV
  572. CDENOM = CAK;
  573. CNUMER = CDENOM + one / CONFRA; //-3zinv+z
  574. KOUNT = 1;
  575. //10 CONTINUE
  576. do { ++KOUNT;
  577. if (KOUNT > MAXIT) {
  578. printf("re(%g):im(%g)\t\n", CONFRA.real(), CONFRA.imag());
  579. throw std::invalid_argument("ConFra--Iteration failed to converge!\n");
  580. }
  581. MM *= -1; KK += 2; //debug mm=1 kk=5
  582. CAK = static_cast<std::complex<double> >(MM*KK) * ZINV; // ** Eq. R25b //debug 5zinv
  583. // //c ** Eq. R32 Ill-conditioned case -- stride two terms instead of one
  584. // if (std::abs( CNUMER / CAK ) >= EPS1 || std::abs( CDENOM / CAK ) >= EPS1) {
  585. // //c ** Eq. R34
  586. // CNTN = CAK * CNUMER + 1.0;
  587. // CDTD = CAK * CDENOM + 1.0;
  588. // CONFRA = ( CNTN / CDTD ) * CONFRA; // ** Eq. R33
  589. // MM *= -1; KK += 2;
  590. // CAK = static_cast<std::complex<double> >(MM*KK) * ZINV; // ** Eq. R25b
  591. // //c ** Eq. R35
  592. // CNUMER = CAK + CNUMER / CNTN;
  593. // CDENOM = CAK + CDENOM / CDTD;
  594. // ++KOUNT;
  595. // //GO TO 10
  596. // continue;
  597. // } else { //c *** Well-conditioned case
  598. {
  599. CAPT = CNUMER / CDENOM; // ** Eq. R27 //debug (-3zinv + z)/(-3zinv)
  600. // printf("re(%g):im(%g)**\t", CAPT.real(), CAPT.imag());
  601. CONFRA = CAPT * CONFRA; // ** Eq. R26
  602. //if (N == 0) {output=true;printf(" re:");prn(CONFRA.real());printf(" im:"); prn(CONFRA.imag());output=false;};
  603. //c ** Check for convergence; Eq. R31
  604. if ( std::abs(CAPT.real() - 1.0) >= EPS2 || std::abs(CAPT.imag()) >= EPS2 ) {
  605. //c ** Eq. R30
  606. CNUMER = CAK + one/CNUMER;
  607. CDENOM = CAK + one/CDENOM;
  608. continue;
  609. //GO TO 10
  610. } // end of if < eps2
  611. }
  612. break;
  613. } while(1);
  614. //if (N == 0) printf(" return confra for z=(%g,%g)\n", ZINV.real(), ZINV.imag());
  615. return CONFRA;
  616. }
  617. //**********************************************************************************//
  618. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  619. // functions (D1 and D3) for a complex argument (z). //
  620. // Equations (16a), (16b) and (18a) - (18d) //
  621. // //
  622. // Input parameters: //
  623. // z: Complex argument to evaluate D1 and D3 //
  624. // nmax_: Maximum number of terms to calculate D1 and D3 //
  625. // //
  626. // Output parameters: //
  627. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  628. //**********************************************************************************//
  629. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  630. std::vector<std::complex<double> >& D1,
  631. std::vector<std::complex<double> >& D3) {
  632. // Downward recurrence for D1 - equations (16a) and (16b)
  633. D1[nmax_] = std::complex<double>(0.0, 0.0);
  634. //D1[nmax_] = calcD1confra(nmax_, z);
  635. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  636. // printf(" D:");prn((D1[nmax_]).real()); printf("\t diff:");
  637. // prn((D1[nmax_] + double(nmax_)*zinv).real());
  638. for (int n = nmax_; n > 0; n--) {
  639. D1[n - 1] = double(n)*zinv - 1.0/(D1[n] + double(n)*zinv);
  640. //D1[n-1] = calcD1confra(n-1, z);
  641. // printf(" D:");prn((D1[n-1]).real()); printf("\t diff:");
  642. // prn((D1[n] + double(n)*zinv).real());
  643. }
  644. // printf("\n\n"); iformat=0;
  645. if (std::abs(D1[0]) > 1000.0 )
  646. throw std::invalid_argument
  647. ("Unstable D1! Please, try to change input parameters!\n");
  648. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  649. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(cos(2.0*z.real()), sin(2.0*z.real()))
  650. *exp(-2.0*z.imag()));
  651. D3[0] = std::complex<double>(0.0, 1.0);
  652. for (int n = 1; n <= nmax_; n++) {
  653. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  654. *(static_cast<double>(n)*zinv- D3[n - 1]);
  655. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  656. }
  657. }
  658. //**********************************************************************************//
  659. // This function calculates Pi and Tau for all values of Theta. //
  660. // Equations (26a) - (26c) //
  661. // //
  662. // Input parameters: //
  663. // nmax_: Maximum number of terms to calculate Pi and Tau //
  664. // nTheta: Number of scattering angles //
  665. // Theta: Array containing all the scattering angles where the scattering //
  666. // amplitudes will be calculated //
  667. // //
  668. // Output parameters: //
  669. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  670. //**********************************************************************************//
  671. void MultiLayerMie::calcPiTau(std::vector< std::vector<double> >& Pi,
  672. std::vector< std::vector<double> >& Tau) {
  673. //****************************************************//
  674. // Equations (26a) - (26c) //
  675. //****************************************************//
  676. std::vector<double> costheta(theta_.size(), 0.0);
  677. for (int t = 0; t < theta_.size(); t++) {
  678. costheta[t] = cos(theta_[t]);
  679. }
  680. for (int n = 0; n < nmax_; n++) {
  681. for (int t = 0; t < theta_.size(); t++) {
  682. if (n == 0) {
  683. // Initialize Pi and Tau
  684. Pi[n][t] = 1.0;
  685. Tau[n][t] = (n + 1)*costheta[t];
  686. } else {
  687. // Calculate the actual values
  688. Pi[n][t] = ((n == 1) ? ((n + n + 1)*costheta[t]*Pi[n - 1][t]/n)
  689. : (((n + n + 1)*costheta[t]*Pi[n - 1][t]
  690. - (n + 1)*Pi[n - 2][t])/n));
  691. Tau[n][t] = (n + 1)*costheta[t]*Pi[n][t] - (n + 2)*Pi[n - 1][t];
  692. }
  693. }
  694. }
  695. }
  696. //**********************************************************************************//
  697. // This function calculates the scattering coefficients required to calculate //
  698. // both the near- and far-field parameters. //
  699. // //
  700. // Input parameters: //
  701. // L: Number of layers //
  702. // pl: Index of PEC layer. If there is none just send -1 //
  703. // x: Array containing the size parameters of the layers [0..L-1] //
  704. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  705. // nmax: Maximum number of multipolar expansion terms to be used for the //
  706. // calculations. Only use it if you know what you are doing, otherwise //
  707. // set this parameter to -1 and the function will calculate it. //
  708. // //
  709. // Output parameters: //
  710. // an, bn: Complex scattering amplitudes //
  711. // //
  712. // Return value: //
  713. // Number of multipolar expansion terms used for the calculations //
  714. //**********************************************************************************//
  715. void MultiLayerMie::ScattCoeffs(std::vector<std::complex<double> >& an,
  716. std::vector<std::complex<double> >& bn) {
  717. //************************************************************************//
  718. // Calculate the index of the first layer. It can be either 0 (default) //
  719. // or the index of the outermost PEC layer. In the latter case all layers //
  720. // below the PEC are discarded. //
  721. //************************************************************************//
  722. const std::vector<double>& x = size_parameter_;
  723. const std::vector<std::complex<double> >& m = index_;
  724. const int& pl = PEC_layer_position_;
  725. const int L = index_.size();
  726. // TODO, is it possible for PEC to have a zero index? If yes than is should be:
  727. // int fl = (pl > -1) ? pl : 0;
  728. int fl = (pl > 0) ? pl : 0;
  729. if (nmax_ <= 0) Nmax(fl);
  730. std::complex<double> z1, z2;
  731. //**************************************************************************//
  732. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  733. // means that index = layer number - 1 or index = n - 1. The only exception //
  734. // are the arrays for representing D1, D3 and Q because they need a value //
  735. // for the index 0 (zero), hence it is important to consider this shift //
  736. // between different arrays. The change was done to optimize memory usage. //
  737. //**************************************************************************//
  738. // Allocate memory to the arrays
  739. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  740. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  741. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  742. for (int l = 0; l < L; l++) {
  743. // D1_mlxl[l].resize(nmax_ + 1);
  744. // D1_mlxlM1[l].resize(nmax_ + 1);
  745. // D3_mlxl[l].resize(nmax_ + 1);
  746. // D3_mlxlM1[l].resize(nmax_ + 1);
  747. Q[l].resize(nmax_ + 1);
  748. Ha[l].resize(nmax_);
  749. Hb[l].resize(nmax_);
  750. }
  751. an.resize(nmax_);
  752. bn.resize(nmax_);
  753. PsiZeta_.resize(nmax_ + 1);
  754. std::vector<std::complex<double> > D1XL(nmax_ + 1), D3XL(nmax_ + 1),
  755. PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  756. //*************************************************//
  757. // Calculate D1 and D3 for z1 in the first layer //
  758. //*************************************************//
  759. if (fl == pl) { // PEC layer
  760. for (int n = 0; n <= nmax_; n++) {
  761. D1_mlxl[n] = std::complex<double>(0.0, -1.0);
  762. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  763. }
  764. } else { // Regular layer
  765. z1 = x[fl]* m[fl];
  766. // Calculate D1 and D3
  767. calcD1D3(z1, D1_mlxl, D3_mlxl);
  768. }
  769. // do { \
  770. // ++iformat;\
  771. // if (iformat%5 == 0) printf("%24.16e",z1.real()); \
  772. // } while (false);
  773. //******************************************************************//
  774. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  775. //******************************************************************//
  776. for (int n = 0; n < nmax_; n++) {
  777. Ha[fl][n] = D1_mlxl[n + 1];
  778. Hb[fl][n] = D1_mlxl[n + 1];
  779. }
  780. //*****************************************************//
  781. // Iteration from the second layer to the last one (L) //
  782. //*****************************************************//
  783. std::complex<double> Temp, Num, Denom;
  784. std::complex<double> G1, G2;
  785. for (int l = fl + 1; l < L; l++) {
  786. //************************************************************//
  787. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  788. //************************************************************//
  789. z1 = x[l]*m[l];
  790. z2 = x[l - 1]*m[l];
  791. //Calculate D1 and D3 for z1
  792. calcD1D3(z1, D1_mlxl, D3_mlxl);
  793. //Calculate D1 and D3 for z2
  794. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  795. // prn(z1.real());
  796. // for ( auto i : D1_mlxl) { prn(i.real());
  797. // // prn(i.imag());
  798. // } printf("\n");
  799. //*********************************************//
  800. //Calculate Q, Ha and Hb in the layers fl+1..L //
  801. //*********************************************//
  802. // Upward recurrence for Q - equations (19a) and (19b)
  803. Num = exp(-2.0*(z1.imag() - z2.imag()))
  804. * std::complex<double>(cos(-2.0*z2.real()) - exp(-2.0*z2.imag()), sin(-2.0*z2.real()));
  805. Denom = std::complex<double>(cos(-2.0*z1.real()) - exp(-2.0*z1.imag()), sin(-2.0*z1.real()));
  806. Q[l][0] = Num/Denom;
  807. for (int n = 1; n <= nmax_; n++) {
  808. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  809. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  810. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  811. }
  812. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  813. for (int n = 1; n <= nmax_; n++) {
  814. //Ha
  815. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  816. G1 = -D1_mlxlM1[n];
  817. G2 = -D3_mlxlM1[n];
  818. } else {
  819. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  820. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  821. } // end of if PEC
  822. Temp = Q[l][n]*G1;
  823. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  824. Denom = G2 - Temp;
  825. Ha[l][n - 1] = Num/Denom;
  826. //Hb
  827. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  828. G1 = Hb[l - 1][n - 1];
  829. G2 = Hb[l - 1][n - 1];
  830. } else {
  831. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  832. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  833. } // end of if PEC
  834. Temp = Q[l][n]*G1;
  835. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  836. Denom = (G2- Temp);
  837. Hb[l][n - 1] = (Num/ Denom);
  838. } // end of for Ha and Hb terms
  839. } // end of for layers iteration
  840. //**************************************//
  841. //Calculate D1, D3, Psi and Zeta for XL //
  842. //**************************************//
  843. // Calculate D1XL and D3XL
  844. calcD1D3(x[L - 1], D1XL, D3XL);
  845. //printf("%5.20f\n",Ha[L-1][0].real());
  846. // Calculate PsiXL and ZetaXL
  847. calcPsiZeta(x[L - 1], D1XL, D3XL, PsiXL, ZetaXL);
  848. //*********************************************************************//
  849. // Finally, we calculate the scattering coefficients (an and bn) and //
  850. // the angular functions (Pi and Tau). Note that for these arrays the //
  851. // first layer is 0 (zero), in future versions all arrays will follow //
  852. // this convention to save memory. (13 Nov, 2014) //
  853. //*********************************************************************//
  854. for (int n = 0; n < nmax_; n++) {
  855. //********************************************************************//
  856. //Expressions for calculating an and bn coefficients are not valid if //
  857. //there is only one PEC layer (ie, for a simple PEC sphere). //
  858. //********************************************************************//
  859. if (pl < (L - 1)) {
  860. 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]);
  861. 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]);
  862. } else {
  863. an[n] = calc_an(n + 1, x[L - 1], std::complex<double>(0.0, 0.0), std::complex<double>(1.0, 0.0), PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  864. bn[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  865. }
  866. } // end of for an and bn terms
  867. } // end of void MultiLayerMie::ScattCoeffs(...)
  868. // ********************************************************************** //
  869. // ********************************************************************** //
  870. // ********************************************************************** //
  871. void MultiLayerMie::InitMieCalculations() {
  872. // Initialize the scattering parameters
  873. Qext_ = 0;
  874. Qsca_ = 0;
  875. Qabs_ = 0;
  876. Qbk_ = 0;
  877. Qpr_ = 0;
  878. asymmetry_factor_ = 0;
  879. albedo_ = 0;
  880. // Initialize the scattering amplitudes
  881. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  882. S1_.swap(tmp1);
  883. S2_ = S1_;
  884. }
  885. // ********************************************************************** //
  886. // ********************************************************************** //
  887. // ********************************************************************** //
  888. //**********************************************************************************//
  889. // This function calculates the actual scattering parameters and amplitudes //
  890. // //
  891. // Input parameters: //
  892. // L: Number of layers //
  893. // pl: Index of PEC layer. If there is none just send -1 //
  894. // x: Array containing the size parameters of the layers [0..L-1] //
  895. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  896. // nTheta: Number of scattering angles //
  897. // Theta: Array containing all the scattering angles where the scattering //
  898. // amplitudes will be calculated //
  899. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  900. // calculations. Only use it if you know what you are doing, otherwise //
  901. // set this parameter to -1 and the function will calculate it //
  902. // //
  903. // Output parameters: //
  904. // Qext: Efficiency factor for extinction //
  905. // Qsca: Efficiency factor for scattering //
  906. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  907. // Qbk: Efficiency factor for backscattering //
  908. // Qpr: Efficiency factor for the radiation pressure //
  909. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  910. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  911. // S1, S2: Complex scattering amplitudes //
  912. // //
  913. // Return value: //
  914. // Number of multipolar expansion terms used for the calculations //
  915. //**********************************************************************************//
  916. void MultiLayerMie::RunMieCalculations() {
  917. if (size_parameter_.size() != index_.size())
  918. throw std::invalid_argument("Each size parameter should have only one index!");
  919. if (size_parameter_.size() == 0)
  920. throw std::invalid_argument("Initialize model first!");
  921. std::vector<std::complex<double> > an, bn;
  922. std::complex<double> Qbktmp(0.0, 0.0);
  923. const std::vector<double>& x = size_parameter_;
  924. // Calculate scattering coefficients
  925. ScattCoeffs(an, bn);
  926. // std::vector< std::vector<double> > Pi(nmax_), Tau(nmax_);
  927. std::vector< std::vector<double> > Pi, Tau;
  928. Pi.resize(nmax_);
  929. Tau.resize(nmax_);
  930. for (int i =0; i< nmax_; ++i) {
  931. Pi[i].resize(theta_.size());
  932. Tau[i].resize(theta_.size());
  933. }
  934. calcPiTau(Pi, Tau);
  935. InitMieCalculations();
  936. // By using downward recurrence we avoid loss of precision due to float rounding errors
  937. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  938. // http://en.wikipedia.org/wiki/Loss_of_significance
  939. for (int i = nmax_ - 2; i >= 0; i--) {
  940. int n = i + 1;
  941. // Equation (27)
  942. Qext_ += (n + n + 1)*(an[i].real() + bn[i].real());
  943. // Equation (28)
  944. Qsca_ += (n + n + 1)*(an[i].real()*an[i].real() + an[i].imag()*an[i].imag()
  945. + bn[i].real()*bn[i].real() + bn[i].imag()*bn[i].imag());
  946. // Equation (29) TODO We must check carefully this equation. If we
  947. // remove the typecast to double then the result changes. Which is
  948. // the correct one??? Ovidio (2014/12/10) With cast ratio will
  949. // give double, without cast (n + n + 1)/(n*(n + 1)) will be
  950. // rounded to integer. Tig (2015/02/24)
  951. Qpr_ += ((n*(n + 2)/(n + 1))*((an[i]*std::conj(an[n]) + bn[i]*std::conj(bn[n])).real())
  952. + ((double)(n + n + 1)/(n*(n + 1)))*(an[i]*std::conj(bn[i])).real());
  953. // Equation (33)
  954. Qbktmp = Qbktmp + (double)(n + n + 1)*(1 - 2*(n % 2))*(an[i]- bn[i]);
  955. // Calculate the scattering amplitudes (S1 and S2) //
  956. // Equations (25a) - (25b) //
  957. for (int t = 0; t < theta_.size(); t++) {
  958. S1_[t] += calc_S1(n, an[i], bn[i], Pi[i][t], Tau[i][t]);
  959. S2_[t] += calc_S2(n, an[i], bn[i], Pi[i][t], Tau[i][t]);
  960. }
  961. }
  962. double x2 = pow2(x.back());
  963. Qext_ = 2*(Qext_)/x2; // Equation (27)
  964. Qsca_ = 2*(Qsca_)/x2; // Equation (28)
  965. Qpr_ = Qext_ - 4*(Qpr_)/x2; // Equation (29)
  966. Qabs_ = Qext_ - Qsca_; // Equation (30)
  967. albedo_ = Qsca_ / Qext_; // Equation (31)
  968. asymmetry_factor_ = (Qext_ - Qpr_) / Qsca_; // Equation (32)
  969. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  970. isMieCalculated_ = true;
  971. //return nmax;
  972. }
  973. // ********************************************************************** //
  974. // ********************************************************************** //
  975. // ********************************************************************** //
  976. // external scattering field = incident + scattered
  977. // BH p.92 (4.37), 94 (4.45), 95 (4.50)
  978. // assume: medium is non-absorbing; refim = 0; Uabs = 0
  979. void MultiLayerMie::fieldExt(double Rho, double Phi, double Theta, std::vector<double> Pi, std::vector<double> Tau,
  980. std::vector<std::complex<double> > an, std::vector<std::complex<double> > bn,
  981. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  982. double rn = 0.0;
  983. std::complex<double> zn, xxip, encap;
  984. std::vector<std::complex<double> > vm3o1n, vm3e1n, vn3o1n, vn3e1n;
  985. vm3o1n.resize(3);
  986. vm3e1n.resize(3);
  987. vn3o1n.resize(3);
  988. vn3e1n.resize(3);
  989. std::vector<std::complex<double> > Ei, Hi, Es, Hs;
  990. Ei.resize(3);
  991. Hi.resize(3);
  992. Es.resize(3);
  993. Hs.resize(3);
  994. for (int i = 0; i < 3; i++) {
  995. Ei[i] = std::complex<double>(0.0, 0.0);
  996. Hi[i] = std::complex<double>(0.0, 0.0);
  997. Es[i] = std::complex<double>(0.0, 0.0);
  998. Hs[i] = std::complex<double>(0.0, 0.0);
  999. }
  1000. std::vector<std::complex<double> > bj, by, bd;
  1001. bj.resize(nmax_);
  1002. by.resize(nmax_);
  1003. bd.resize(nmax_);
  1004. // Calculate spherical Bessel and Hankel functions
  1005. sphericalBessel(Rho, bj, by, bd);
  1006. for (int n = 0; n < nmax_; n++) {
  1007. rn = double(n + 1);
  1008. zn = bj[n] + std::complex<double>(0.0, 1.0)*by[n];
  1009. xxip = Rho*(bj[n] + std::complex<double>(0.0, 1.0)*by[n]) - rn*zn;
  1010. vm3o1n[0] = std::complex<double>(0.0, 0.0);
  1011. vm3o1n[1] = std::cos(Phi)*Pi[n]*zn;
  1012. vm3o1n[2] = -(std::sin(Phi)*Tau[n]*zn);
  1013. vm3e1n[0] = std::complex<double>(0.0, 0.0);
  1014. vm3e1n[1] = -(std::sin(Phi)*Pi[n]*zn);
  1015. vm3e1n[2] = -(std::cos(Phi)*Tau[n]*zn);
  1016. vn3o1n[0] = std::sin(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  1017. vn3o1n[1] = std::sin(Phi)*Tau[n]*xxip/Rho;
  1018. vn3o1n[2] = std::cos(Phi)*Pi[n]*xxip/Rho;
  1019. vn3e1n[0] = std::cos(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  1020. vn3e1n[1] = std::cos(Phi)*Tau[n]*xxip/Rho;
  1021. vn3e1n[2] = -(std::sin(Phi)*Pi[n]*xxip/Rho);
  1022. // scattered field: BH p.94 (4.45)
  1023. encap = std::pow(std::complex<double>(0.0, 1.0), rn)*(2.0*rn + 1.0)/(rn*(rn + 1.0));
  1024. for (int i = 0; i < 3; i++) {
  1025. Es[i] = Es[i] + encap*(std::complex<double>(0.0, 1.0)*an[n]*vn3e1n[i] - bn[n]*vm3o1n[i]);
  1026. Hs[i] = Hs[i] + encap*(std::complex<double>(0.0, 1.0)*bn[n]*vn3o1n[i] + an[n]*vm3e1n[i]);
  1027. }
  1028. }
  1029. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1030. // basis unit vectors = er, etheta, ephi
  1031. std::complex<double> eifac = std::exp(std::complex<double>(0.0, 1.0)*Rho*std::cos(Theta));
  1032. Ei[0] = eifac*std::sin(Theta)*std::cos(Phi);
  1033. Ei[1] = eifac*std::cos(Theta)*std::cos(Phi);
  1034. Ei[2] = -(eifac*std::sin(Phi));
  1035. // magnetic field
  1036. double hffact = 1.0/(cc*mu);
  1037. for (int i = 0; i < 3; i++) {
  1038. Hs[i] = hffact*Hs[i];
  1039. }
  1040. // incident H field: BH p.26 (2.43), p.89 (4.21)
  1041. std::complex<double> hffacta = hffact;
  1042. std::complex<double> hifac = eifac*hffacta;
  1043. Hi[0] = hifac*std::sin(Theta)*std::sin(Phi);
  1044. Hi[1] = hifac*std::cos(Theta)*std::sin(Phi);
  1045. Hi[2] = hifac*std::cos(Phi);
  1046. for (int i = 0; i < 3; i++) {
  1047. // electric field E [V m-1] = EF*E0
  1048. E[i] = Ei[i] + Es[i];
  1049. H[i] = Hi[i] + Hs[i];
  1050. }
  1051. }
  1052. // ********************************************************************** //
  1053. // ********************************************************************** //
  1054. // ********************************************************************** //
  1055. //**********************************************************************************//
  1056. // This function calculates complex electric and magnetic field in the surroundings //
  1057. // and inside (TODO) the particle. //
  1058. // //
  1059. // Input parameters: //
  1060. // L: Number of layers //
  1061. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1062. // x: Array containing the size parameters of the layers [0..L-1] //
  1063. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1064. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1065. // calculations. Only use it if you know what you are doing, otherwise //
  1066. // set this parameter to 0 (zero) and the function will calculate it. //
  1067. // ncoord: Number of coordinate points //
  1068. // Coords: Array containing all coordinates where the complex electric and //
  1069. // magnetic fields will be calculated //
  1070. // //
  1071. // Output parameters: //
  1072. // E, H: Complex electric and magnetic field at the provided coordinates //
  1073. // //
  1074. // Return value: //
  1075. // Number of multipolar expansion terms used for the calculations //
  1076. //**********************************************************************************//
  1077. // int MultiLayerMie::nField(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int nmax,
  1078. // int ncoord, std::vector<double> Xp, std::vector<double> Yp, std::vector<double> Zp,
  1079. // std::vector<std::vector<std::complex<double> > >& E, std::vector<std::vector<std::complex<double> > >& H) {
  1080. // double Rho, Phi, Theta;
  1081. // std::vector<std::complex<double> > an, bn;
  1082. // // This array contains the fields in spherical coordinates
  1083. // std::vector<std::complex<double> > Es, Hs;
  1084. // Es.resize(3);
  1085. // Hs.resize(3);
  1086. // // Calculate scattering coefficients
  1087. // ScattCoeffs(L, pl, an, bn);
  1088. // std::vector<double> Pi, Tau;
  1089. // Pi.resize(nmax_);
  1090. // Tau.resize(nmax_);
  1091. // for (int c = 0; c < ncoord; c++) {
  1092. // // Convert to spherical coordinates
  1093. // Rho = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  1094. // if (Rho < 1e-3) {
  1095. // // Avoid convergence problems
  1096. // Rho = 1e-3;
  1097. // }
  1098. // Phi = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  1099. // Theta = acos(Xp[c]/Rho);
  1100. // calcPiTau(Theta, Pi, Tau);
  1101. // //*******************************************************//
  1102. // // external scattering field = incident + scattered //
  1103. // // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1104. // // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1105. // //*******************************************************//
  1106. // // Firstly the easiest case: the field outside the particle
  1107. // if (Rho >= x[L - 1]) {
  1108. // fieldExt(Rho, Phi, Theta, Pi, Tau, an, bn, Es, Hs);
  1109. // } else {
  1110. // // TODO, for now just set all the fields to zero
  1111. // for (int i = 0; i < 3; i++) {
  1112. // Es[i] = std::complex<double>(0.0, 0.0);
  1113. // Hs[i] = std::complex<double>(0.0, 0.0);
  1114. // }
  1115. // }
  1116. // //Now, convert the fields back to cartesian coordinates
  1117. // E[c][0] = std::sin(Theta)*std::cos(Phi)*Es[0] + std::cos(Theta)*std::cos(Phi)*Es[1] - std::sin(Phi)*Es[2];
  1118. // E[c][1] = std::sin(Theta)*std::sin(Phi)*Es[0] + std::cos(Theta)*std::sin(Phi)*Es[1] + std::cos(Phi)*Es[2];
  1119. // E[c][2] = std::cos(Theta)*Es[0] - std::sin(Theta)*Es[1];
  1120. // H[c][0] = std::sin(Theta)*std::cos(Phi)*Hs[0] + std::cos(Theta)*std::cos(Phi)*Hs[1] - std::sin(Phi)*Hs[2];
  1121. // H[c][1] = std::sin(Theta)*std::sin(Phi)*Hs[0] + std::cos(Theta)*std::sin(Phi)*Hs[1] + std::cos(Phi)*Hs[2];
  1122. // H[c][2] = std::cos(Theta)*Hs[0] - std::sin(Theta)*Hs[1];
  1123. // }
  1124. // return nmax;
  1125. // } // end of int nField()
  1126. } // end of namespace nmie