nmie-wrapper.cc 66 KB

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