nmie.cc 63 KB

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