nmie-wrapper.cc 85 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735
  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,2014,2015 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, std::vector<double>& x, std::vector<std::complex<double> >& m, int nTheta, std::vector<double>& Theta, double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo, std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  51. if (x.size() != L || m.size() != L)
  52. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  53. if (Theta.size() != nTheta)
  54. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  55. try {
  56. MultiLayerMie multi_layer_mie;
  57. multi_layer_mie.SetWidthSP(x);
  58. multi_layer_mie.SetIndexSP(m);
  59. multi_layer_mie.SetAngles(Theta);
  60. multi_layer_mie.RunMieCalculations();
  61. *Qext = multi_layer_mie.GetQext();
  62. *Qsca = multi_layer_mie.GetQsca();
  63. *Qabs = multi_layer_mie.GetQabs();
  64. *Qbk = multi_layer_mie.GetQbk();
  65. *Qpr = multi_layer_mie.GetQpr();
  66. *g = multi_layer_mie.GetAsymmetryFactor();
  67. *Albedo = multi_layer_mie.GetAlbedo();
  68. S1 = multi_layer_mie.GetS1();
  69. S2 = multi_layer_mie.GetS2();
  70. //multi_layer_mie.GetFailed();
  71. } catch(const std::invalid_argument& ia) {
  72. // Will catch if multi_layer_mie fails or other errors.
  73. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  74. throw std::invalid_argument(ia);
  75. return -1;
  76. }
  77. return 0;
  78. }
  79. // ********************************************************************** //
  80. // ********************************************************************** //
  81. // ********************************************************************** //
  82. int nField(const int L, const int pl, const std::vector<double>& x, const std::vector<std::complex<double> >& m, const int nmax, const int ncoord, const std::vector<double>& Xp_vec, const std::vector<double>& Yp_vec, const std::vector<double>& Zp_vec, std::vector<std::vector<std::complex<double> > >& E, std::vector<std::vector<std::complex<double> > >& H) {
  83. if (x.size() != L || m.size() != L)
  84. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  85. if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
  86. || E.size() != ncoord || H.size() != ncoord)
  87. throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
  88. for (auto f:E)
  89. if (f.size() != 3)
  90. throw std::invalid_argument("Field E is not 3D!");
  91. for (auto f:H)
  92. if (f.size() != 3)
  93. throw std::invalid_argument("Field H is not 3D!");
  94. try {
  95. MultiLayerMie multi_layer_mie;
  96. //multi_layer_mie.SetPEC(pl);
  97. multi_layer_mie.SetWidthSP(x);
  98. multi_layer_mie.SetIndexSP(m);
  99. multi_layer_mie.SetFieldPointsSP({Xp_vec, Yp_vec, Zp_vec});
  100. multi_layer_mie.RunFieldCalculations();
  101. E = multi_layer_mie.GetFieldE();
  102. H = multi_layer_mie.GetFieldH();
  103. //multi_layer_mie.GetFailed();
  104. } catch(const std::invalid_argument& ia) {
  105. // Will catch if multi_layer_mie fails or other errors.
  106. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  107. throw std::invalid_argument(ia);
  108. return - 1;
  109. }
  110. return 0;
  111. }
  112. // ********************************************************************** //
  113. // ********************************************************************** //
  114. // ********************************************************************** //
  115. void MultiLayerMie::GetFailed() {
  116. double faild_x = 9.42477796076938;
  117. //double faild_x = 9.42477796076937;
  118. std::complex<double> z(faild_x, 0.0);
  119. std::vector<int> nmax_local_array = {20, 100, 500, 2500};
  120. for (auto nmax_local : nmax_local_array) {
  121. std::vector<std::complex<double> > D1_failed(nmax_local + 1);
  122. // Downward recurrence for D1 - equations (16a) and (16b)
  123. D1_failed[nmax_local] = std::complex<double>(0.0, 0.0);
  124. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  125. for (int n = nmax_local; n > 0; n--) {
  126. D1_failed[n - 1] = double(n)*zinv - 1.0/(D1_failed[n] + double(n)*zinv);
  127. }
  128. printf("Faild D1[0] from reccurence (z = %16.14f, nmax = %d): %g\n",
  129. faild_x, nmax_local, D1_failed[0].real());
  130. }
  131. printf("Faild D1[0] from continued fraction (z = %16.14f): %g\n", faild_x,
  132. calcD1confra(0,z).real());
  133. //D1[nmax_] = calcD1confra(nmax_, z);
  134. }
  135. // ********************************************************************** //
  136. // ********************************************************************** //
  137. // ********************************************************************** //
  138. double MultiLayerMie::GetQext() {
  139. if (!isMieCalculated_)
  140. throw std::invalid_argument("You should run calculations before result request!");
  141. return Qext_;
  142. }
  143. // ********************************************************************** //
  144. // ********************************************************************** //
  145. // ********************************************************************** //
  146. double MultiLayerMie::GetQabs() {
  147. if (!isMieCalculated_)
  148. throw std::invalid_argument("You should run calculations before result request!");
  149. return Qabs_;
  150. }
  151. // ********************************************************************** //
  152. // ********************************************************************** //
  153. // ********************************************************************** //
  154. std::vector<double> MultiLayerMie::GetQabs_channel() {
  155. if (!isMieCalculated_)
  156. throw std::invalid_argument("You should run calculations before result request!");
  157. return Qabs_ch_;
  158. }
  159. // ********************************************************************** //
  160. // ********************************************************************** //
  161. // ********************************************************************** //
  162. std::vector<double> MultiLayerMie::GetQabs_channel_normalized() {
  163. if (!isMieCalculated_)
  164. throw std::invalid_argument("You should run calculations before result request!");
  165. // std::vector<double> NACS(nmax_ - 1, 0.0);
  166. // double x2 = pow2(size_parameter_.back());
  167. // for (int i = 0; i < nmax_ - 1; ++i) {
  168. // const int n = i + 1;
  169. // NACS[i] = Qabs_ch_[i]*x2/(2.0*(2.0*static_cast<double>(n) + 1));
  170. // // if (NACS[i] > 0.250000001)
  171. // // throw std::invalid_argument("Unexpected normalized absorption cross-section value!");
  172. // }
  173. //return NACS;
  174. return Qabs_ch_norm_;
  175. }
  176. // ********************************************************************** //
  177. // ********************************************************************** //
  178. // ********************************************************************** //
  179. double MultiLayerMie::GetQsca() {
  180. if (!isMieCalculated_)
  181. throw std::invalid_argument("You should run calculations before result request!");
  182. return Qsca_;
  183. }
  184. // ********************************************************************** //
  185. // ********************************************************************** //
  186. // ********************************************************************** //
  187. std::vector<double> MultiLayerMie::GetQsca_channel() {
  188. if (!isMieCalculated_)
  189. throw std::invalid_argument("You should run calculations before result request!");
  190. return Qsca_ch_;
  191. }
  192. // ********************************************************************** //
  193. // ********************************************************************** //
  194. // ********************************************************************** //
  195. std::vector<double> MultiLayerMie::GetQsca_channel_normalized() {
  196. if (!isMieCalculated_)
  197. throw std::invalid_argument("You should run calculations before result request!");
  198. // std::vector<double> NACS(nmax_ - 1, 0.0);
  199. // double x2 = pow2(size_parameter_.back());
  200. // for (int i = 0; i < nmax_ - 1; ++i) {
  201. // const int n = i + 1;
  202. // NACS[i] = Qsca_ch_[i]*x2/(2.0*(2.0*static_cast<double>(n) + 1.0));
  203. // }
  204. // return NACS;
  205. return Qsca_ch_norm_;
  206. }
  207. // ********************************************************************** //
  208. // ********************************************************************** //
  209. // ********************************************************************** //
  210. double MultiLayerMie::GetQbk() {
  211. if (!isMieCalculated_)
  212. throw std::invalid_argument("You should run calculations before result request!");
  213. return Qbk_;
  214. }
  215. // ********************************************************************** //
  216. // ********************************************************************** //
  217. // ********************************************************************** //
  218. double MultiLayerMie::GetQpr() {
  219. if (!isMieCalculated_)
  220. throw std::invalid_argument("You should run calculations before result request!");
  221. return Qpr_;
  222. }
  223. // ********************************************************************** //
  224. // ********************************************************************** //
  225. // ********************************************************************** //
  226. double MultiLayerMie::GetAsymmetryFactor() {
  227. if (!isMieCalculated_)
  228. throw std::invalid_argument("You should run calculations before result request!");
  229. return asymmetry_factor_;
  230. }
  231. // ********************************************************************** //
  232. // ********************************************************************** //
  233. // ********************************************************************** //
  234. double MultiLayerMie::GetAlbedo() {
  235. if (!isMieCalculated_)
  236. throw std::invalid_argument("You should run calculations before result request!");
  237. return albedo_;
  238. }
  239. // ********************************************************************** //
  240. // ********************************************************************** //
  241. // ********************************************************************** //
  242. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  243. if (!isMieCalculated_)
  244. throw std::invalid_argument("You should run calculations before result request!");
  245. return S1_;
  246. }
  247. // ********************************************************************** //
  248. // ********************************************************************** //
  249. // ********************************************************************** //
  250. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  251. if (!isMieCalculated_)
  252. throw std::invalid_argument("You should run calculations before result request!");
  253. return S2_;
  254. }
  255. // ********************************************************************** //
  256. // ********************************************************************** //
  257. // ********************************************************************** //
  258. void MultiLayerMie::AddTargetLayer(double width, std::complex<double> layer_index) {
  259. isMieCalculated_ = false;
  260. if (width <= 0)
  261. throw std::invalid_argument("Layer width should be positive!");
  262. target_width_.push_back(width);
  263. target_index_.push_back(layer_index);
  264. } // end of void MultiLayerMie::AddTargetLayer(...)
  265. // ********************************************************************** //
  266. // ********************************************************************** //
  267. // ********************************************************************** //
  268. void MultiLayerMie::SetTargetPEC(double radius) {
  269. isMieCalculated_ = false;
  270. if (target_width_.size() != 0 || target_index_.size() != 0)
  271. throw std::invalid_argument("Error! Define PEC target radius before any other layers!");
  272. // Add layer of any index...
  273. AddTargetLayer(radius, std::complex<double>(0.0, 0.0));
  274. // ... and mark it as PEC
  275. SetPEC(0.0);
  276. }
  277. // ********************************************************************** //
  278. // ********************************************************************** //
  279. // ********************************************************************** //
  280. void MultiLayerMie::SetCoatingIndex(std::vector<std::complex<double> > index) {
  281. isMieCalculated_ = false;
  282. coating_index_.clear();
  283. for (auto value : index) coating_index_.push_back(value);
  284. } // end of void MultiLayerMie::SetCoatingIndex(std::vector<complex> index);
  285. // ********************************************************************** //
  286. // ********************************************************************** //
  287. // ********************************************************************** //
  288. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  289. isMieCalculated_ = false;
  290. theta_ = angles;
  291. // theta_.clear();
  292. // for (auto value : angles) theta_.push_back(value);
  293. } // end of SetAngles()
  294. // ********************************************************************** //
  295. // ********************************************************************** //
  296. // ********************************************************************** //
  297. void MultiLayerMie::SetCoatingWidth(std::vector<double> width) {
  298. isMieCalculated_ = false;
  299. coating_width_.clear();
  300. for (auto w : width)
  301. if (w <= 0)
  302. throw std::invalid_argument("Coating width should be positive!");
  303. else coating_width_.push_back(w);
  304. }
  305. // end of void MultiLayerMie::SetCoatingWidth(...);
  306. // ********************************************************************** //
  307. // ********************************************************************** //
  308. // ********************************************************************** //
  309. void MultiLayerMie::SetWidthSP(const std::vector<double>& size_parameter) {
  310. isMieCalculated_ = false;
  311. size_parameter_.clear();
  312. double prev_size_parameter = 0.0;
  313. for (auto layer_size_parameter : size_parameter) {
  314. if (layer_size_parameter <= 0.0)
  315. throw std::invalid_argument("Size parameter should be positive!");
  316. if (prev_size_parameter > layer_size_parameter)
  317. throw std::invalid_argument
  318. ("Size parameter for next layer should be larger than the previous one!");
  319. prev_size_parameter = layer_size_parameter;
  320. size_parameter_.push_back(layer_size_parameter);
  321. }
  322. }
  323. // end of void MultiLayerMie::SetWidthSP(...);
  324. // ********************************************************************** //
  325. // ********************************************************************** //
  326. // ********************************************************************** //
  327. void MultiLayerMie::SetIndexSP(const std::vector< std::complex<double> >& index) {
  328. isMieCalculated_ = false;
  329. //index_.clear();
  330. index_ = index;
  331. // for (auto value : index) index_.push_back(value);
  332. } // end of void MultiLayerMie::SetIndexSP(...);
  333. // ********************************************************************** //
  334. // ********************************************************************** //
  335. // ********************************************************************** //
  336. void MultiLayerMie::SetFieldPointsSP(const std::vector< std::vector<double> >& coords_sp) {
  337. if (coords_sp.size() != 3)
  338. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  339. if (coords_sp[0].size() != coords_sp[1].size() || coords_sp[0].size() != coords_sp[2].size())
  340. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  341. coords_sp_ = coords_sp;
  342. // for (int i = 0; i < coords_sp_[0].size(); ++i) {
  343. // printf("%g, %g, %g\n", coords_sp_[0][i], coords_sp_[1][i], coords_sp_[2][i]);
  344. // }
  345. } // end of void MultiLayerMie::SetFieldPointsSP(...)
  346. // ********************************************************************** //
  347. // ********************************************************************** //
  348. // ********************************************************************** //
  349. void MultiLayerMie::SetPEC(int layer_position) {
  350. isMieCalculated_ = false;
  351. if (layer_position < 0)
  352. throw std::invalid_argument("Error! Layers are numbered from 0!");
  353. PEC_layer_position_ = layer_position;
  354. }
  355. // ********************************************************************** //
  356. // ********************************************************************** //
  357. // ********************************************************************** //
  358. void MultiLayerMie::SetMaxTermsNumber(int nmax) {
  359. isMieCalculated_ = false;
  360. nmax_preset_ = nmax;
  361. //debug
  362. printf("Setting max terms: %d\n", nmax_preset_);
  363. }
  364. // ********************************************************************** //
  365. // ********************************************************************** //
  366. // ********************************************************************** //
  367. void MultiLayerMie::GenerateSizeParameter() {
  368. isMieCalculated_ = false;
  369. size_parameter_.clear();
  370. double radius = 0.0;
  371. for (auto width : target_width_) {
  372. radius += width;
  373. size_parameter_.push_back(2*PI_*radius/wavelength_);
  374. }
  375. for (auto width : coating_width_) {
  376. radius += width;
  377. size_parameter_.push_back(2*PI_*radius/wavelength_);
  378. }
  379. total_radius_ = radius;
  380. } // end of void MultiLayerMie::GenerateSizeParameter();
  381. // ********************************************************************** //
  382. // ********************************************************************** //
  383. // ********************************************************************** //
  384. void MultiLayerMie::GenerateIndex() {
  385. isMieCalculated_ = false;
  386. index_.clear();
  387. for (auto index : target_index_) index_.push_back(index);
  388. for (auto index : coating_index_) index_.push_back(index);
  389. } // end of void MultiLayerMie::GenerateIndex();
  390. // ********************************************************************** //
  391. // ********************************************************************** //
  392. // ********************************************************************** //
  393. double MultiLayerMie::GetTotalRadius() {
  394. if (!isMieCalculated_)
  395. throw std::invalid_argument("You should run calculations before result request!");
  396. if (total_radius_ == 0) GenerateSizeParameter();
  397. return total_radius_;
  398. } // end of double MultiLayerMie::GetTotalRadius();
  399. // ********************************************************************** //
  400. // ********************************************************************** //
  401. // ********************************************************************** //
  402. std::vector< std::vector<double> >
  403. MultiLayerMie::GetSpectra(double from_WL, double to_WL, int samples) {
  404. if (!isMieCalculated_)
  405. throw std::invalid_argument("You should run calculations before result request!");
  406. std::vector< std::vector<double> > spectra;
  407. double step_WL = (to_WL - from_WL)/static_cast<double>(samples);
  408. double wavelength_backup = wavelength_;
  409. long fails = 0;
  410. for (double WL = from_WL; WL < to_WL; WL += step_WL) {
  411. wavelength_ = WL;
  412. try {
  413. RunMieCalculations();
  414. } catch(const std::invalid_argument& ia) {
  415. fails++;
  416. continue;
  417. }
  418. //printf("%3.1f ",WL);
  419. spectra.push_back(std::vector<double>({wavelength_, Qext_, Qsca_, Qabs_, Qbk_}));
  420. } // end of for each WL in spectra
  421. printf("Spectrum has %li fails\n",fails);
  422. wavelength_ = wavelength_backup;
  423. return spectra;
  424. }
  425. // ********************************************************************** //
  426. // ********************************************************************** //
  427. // ********************************************************************** //
  428. void MultiLayerMie::ClearTarget() {
  429. isMieCalculated_ = false;
  430. target_width_.clear();
  431. target_index_.clear();
  432. }
  433. // ********************************************************************** //
  434. // ********************************************************************** //
  435. // ********************************************************************** //
  436. void MultiLayerMie::ClearCoating() {
  437. isMieCalculated_ = false;
  438. coating_width_.clear();
  439. coating_index_.clear();
  440. }
  441. // ********************************************************************** //
  442. // ********************************************************************** //
  443. // ********************************************************************** //
  444. void MultiLayerMie::ClearLayers() {
  445. isMieCalculated_ = false;
  446. ClearTarget();
  447. ClearCoating();
  448. }
  449. // ********************************************************************** //
  450. // ********************************************************************** //
  451. // ********************************************************************** //
  452. void MultiLayerMie::ClearAllDesign() {
  453. isMieCalculated_ = false;
  454. ClearLayers();
  455. size_parameter_.clear();
  456. index_.clear();
  457. }
  458. // ********************************************************************** //
  459. // ********************************************************************** //
  460. // ********************************************************************** //
  461. // Computational core
  462. // ********************************************************************** //
  463. // ********************************************************************** //
  464. // ********************************************************************** //
  465. // Calculate Nstop - equation (17)
  466. //
  467. void MultiLayerMie::Nstop() {
  468. const double& xL = size_parameter_.back();
  469. if (xL <= 8) {
  470. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  471. } else if (xL <= 4200) {
  472. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  473. } else {
  474. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  475. }
  476. }
  477. // ********************************************************************** //
  478. // ********************************************************************** //
  479. // ********************************************************************** //
  480. void MultiLayerMie::Nmax(int first_layer) {
  481. int ri, riM1;
  482. const std::vector<double>& x = size_parameter_;
  483. const std::vector<std::complex<double> >& m = index_;
  484. Nstop(); // Set initial nmax_ value
  485. for (int i = first_layer; i < x.size(); i++) {
  486. if (i > PEC_layer_position_)
  487. ri = round(std::abs(x[i]*m[i]));
  488. else
  489. ri = 0;
  490. nmax_ = std::max(nmax_, ri);
  491. // first layer is pec, if pec is present
  492. if ((i > first_layer) && ((i - 1) > PEC_layer_position_))
  493. riM1 = round(std::abs(x[i - 1]* m[i]));
  494. else
  495. riM1 = 0;
  496. nmax_ = std::max(nmax_, riM1);
  497. }
  498. nmax_ += 15; // Final nmax_ value
  499. }
  500. //**********************************************************************************//
  501. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  502. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  503. // //
  504. // Input parameters: //
  505. // z: Real argument to evaluate jn and h1n //
  506. // nmax_: Maximum number of terms to calculate jn and h1n //
  507. // //
  508. // Output parameters: //
  509. // jn, h1n: Spherical Bessel and Hankel functions //
  510. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  511. // //
  512. // The implementation follows the algorithm by I.J. Thompson and A.R. Barnett, //
  513. // Comp. Phys. Comm. 47 (1987) 245-257. //
  514. // //
  515. // Complex spherical Bessel functions from n=0..nmax_ - 1 for z in the upper half //
  516. // plane (Im(z) > -3). //
  517. // //
  518. // j[n] = j/n(z) Regular solution: j[0]=sin(z)/z //
  519. // j'[n] = d[j/n(z)]/dz //
  520. // h1[n] = h[0]/n(z) Irregular Hankel function: //
  521. // h1'[n] = d[h[0]/n(z)]/dz h1[0] = j0(z) + i*y0(z) //
  522. // = (sin(z)-i*cos(z))/z //
  523. // = -i*exp(i*z)/z //
  524. // Using complex CF1, and trigonometric forms for n=0 solutions. //
  525. //**********************************************************************************//
  526. void MultiLayerMie::sbesjh(std::complex<double> z,
  527. std::vector<std::complex<double> >& jn,
  528. std::vector<std::complex<double> >& jnp,
  529. std::vector<std::complex<double> >& h1n,
  530. std::vector<std::complex<double> >& h1np) {
  531. const int limit = 20000;
  532. const double accur = 1.0e-12;
  533. const double tm30 = 1e-30;
  534. double absc;
  535. std::complex<double> zi, w;
  536. std::complex<double> pl, f, b, d, c, del, jn0, jndb, h1nldb, h1nbdb;
  537. absc = std::abs(std::real(z)) + std::abs(std::imag(z));
  538. if ((absc < accur) || (std::imag(z) < -3.0)) {
  539. throw std::invalid_argument("TODO add error description for condition if ((absc < accur) || (std::imag(z) < -3.0))");
  540. }
  541. zi = 1.0/z;
  542. w = zi + zi;
  543. pl = double(nmax_)*zi;
  544. f = pl + zi;
  545. b = f + f + zi;
  546. d = 0.0;
  547. c = f;
  548. for (int n = 0; n < limit; n++) {
  549. d = b - d;
  550. c = b - 1.0/c;
  551. absc = std::abs(std::real(d)) + std::abs(std::imag(d));
  552. if (absc < tm30) {
  553. d = tm30;
  554. }
  555. absc = std::abs(std::real(c)) + std::abs(std::imag(c));
  556. if (absc < tm30) {
  557. c = tm30;
  558. }
  559. d = 1.0/d;
  560. del = d*c;
  561. f = f*del;
  562. b += w;
  563. absc = std::abs(std::real(del - 1.0)) + std::abs(std::imag(del - 1.0));
  564. if (absc < accur) {
  565. // We have obtained the desired accuracy
  566. break;
  567. }
  568. }
  569. if (absc > accur) {
  570. throw std::invalid_argument("We were not able to obtain the desired accuracy");
  571. }
  572. jn[nmax_ - 1] = tm30;
  573. jnp[nmax_ - 1] = f*jn[nmax_ - 1];
  574. // Downward recursion to n=0 (N.B. Coulomb Functions)
  575. for (int n = nmax_ - 2; n >= 0; n--) {
  576. jn[n] = pl*jn[n + 1] + jnp[n + 1];
  577. jnp[n] = pl*jn[n] - jn[n + 1];
  578. pl = pl - zi;
  579. }
  580. // Calculate the n=0 Bessel Functions
  581. jn0 = zi*std::sin(z);
  582. h1n[0] = std::exp(std::complex<double>(0.0, 1.0)*z)*zi*(-std::complex<double>(0.0, 1.0));
  583. h1np[0] = h1n[0]*(std::complex<double>(0.0, 1.0) - zi);
  584. // Rescale j[n], j'[n], converting to spherical Bessel functions.
  585. // Recur h1[n], h1'[n] as spherical Bessel functions.
  586. w = 1.0/jn[0];
  587. pl = zi;
  588. for (int n = 0; n < nmax_; n++) {
  589. jn[n] = jn0*(w*jn[n]);
  590. jnp[n] = jn0*(w*jnp[n]) - zi*jn[n];
  591. if (n != 0) {
  592. h1n[n] = (pl - zi)*h1n[n - 1] - h1np[n - 1];
  593. // check if hankel is increasing (upward stable)
  594. if (std::abs(h1n[n]) < std::abs(h1n[n - 1])) {
  595. jndb = z;
  596. h1nldb = h1n[n];
  597. h1nbdb = h1n[n - 1];
  598. }
  599. pl += zi;
  600. h1np[n] = -(pl*h1n[n]) + h1n[n - 1];
  601. }
  602. }
  603. }
  604. //**********************************************************************************//
  605. // This function calculates the spherical Bessel functions (bj and by) and the //
  606. // logarithmic derivative (bd) for a given complex value z. See pag. 87 B&H. //
  607. // //
  608. // Input parameters: //
  609. // z: Complex argument to evaluate bj, by and bd //
  610. // nmax_: Maximum number of terms to calculate bj, by and bd //
  611. // //
  612. // Output parameters: //
  613. // bj, by: Spherical Bessel functions //
  614. // bd: Logarithmic derivative //
  615. //**********************************************************************************//
  616. void MultiLayerMie::sphericalBessel(std::complex<double> z,
  617. std::vector<std::complex<double> >& bj,
  618. std::vector<std::complex<double> >& by,
  619. std::vector<std::complex<double> >& bd) {
  620. std::vector<std::complex<double> > jn(nmax_), jnp(nmax_), h1n(nmax_), h1np(nmax_);
  621. sbesjh(z, jn, jnp, h1n, h1np);
  622. for (int n = 0; n < nmax_; n++) {
  623. bj[n] = jn[n];
  624. by[n] = (h1n[n] - jn[n])/std::complex<double>(0.0, 1.0);
  625. bd[n] = jnp[n]/jn[n] + 1.0/z;
  626. }
  627. // std::complex<double> besselj_0 = std::sin(z)/z;
  628. // std::complex<double> bessely_0 = -std::cos(z)/z;
  629. // if (nmax_>0) {
  630. // bj[0] = std::sin(z)/pow2(z)-std::cos(z)/z; //bj1
  631. // by[0] = std::cos(z)/pow2(z)-std::sin(z)/z; //by1
  632. // }
  633. // if (nmax_>1) {
  634. // bj[1] = bj[0]*3.0/z-besselj_0;//bj2
  635. // by[1] = by[0]*3.0/z-bessely_0;//bj2
  636. // }
  637. // for (int n = 2; n < nmax_; n++) {
  638. // bj[n] = (2.0*n - 1.0)/z*bj[n - 1] - bj[n];
  639. // by[n] = (2.0*n - 1.0)/z*by[n - 1] - by[n];
  640. // }
  641. }
  642. // ********************************************************************** //
  643. // ********************************************************************** //
  644. // ********************************************************************** //
  645. // Calculate an - equation (5)
  646. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  647. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  648. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  649. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  650. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  651. return Num/Denom;
  652. }
  653. // ********************************************************************** //
  654. // ********************************************************************** //
  655. // ********************************************************************** //
  656. // Calculate bn - equation (6)
  657. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  658. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  659. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  660. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  661. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  662. return Num/Denom;
  663. }
  664. // ********************************************************************** //
  665. // ********************************************************************** //
  666. // ********************************************************************** //
  667. // Calculates S1 - equation (25a)
  668. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  669. double Pi, double Tau) {
  670. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  671. }
  672. // ********************************************************************** //
  673. // ********************************************************************** //
  674. // ********************************************************************** //
  675. // Calculates S2 - equation (25b) (it's the same as (25a), just switches Pi and Tau)
  676. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  677. double Pi, double Tau) {
  678. return calc_S1(n, an, bn, Tau, Pi);
  679. }
  680. //**********************************************************************************//
  681. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  682. // real argument (x). //
  683. // Equations (20a) - (21b) //
  684. // //
  685. // Input parameters: //
  686. // x: Real argument to evaluate Psi and Zeta //
  687. // nmax: Maximum number of terms to calculate Psi and Zeta //
  688. // //
  689. // Output parameters: //
  690. // Psi, Zeta: Riccati-Bessel functions //
  691. //**********************************************************************************//
  692. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  693. std::vector<std::complex<double> > D1,
  694. std::vector<std::complex<double> > D3,
  695. std::vector<std::complex<double> >& Psi,
  696. std::vector<std::complex<double> >& Zeta) {
  697. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  698. //Psi[0] = std::complex<double>(std::sin(x), 0);
  699. std::complex<double> c_i(0.0, 1.0);
  700. Psi[0] = std::sin(z);
  701. //Zeta[0] = std::complex<double>(std::sin(x), -std::cos(x));
  702. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  703. for (int n = 1; n <= nmax_; n++) {
  704. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  705. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  706. }
  707. }
  708. //**********************************************************************************//
  709. // Function CONFRA ported from MIEV0.f (Wiscombe,1979)
  710. // Ref. to NCAR Technical Notes, Wiscombe, 1979
  711. /*
  712. c Compute Bessel function ratio A-sub-N from its
  713. c continued fraction using Lentz method
  714. c ZINV = Reciprocal of argument of A
  715. c I N T E R N A L V A R I A B L E S
  716. c ------------------------------------
  717. c CAK Term in continued fraction expansion of A (Eq. R25)
  718. c a_k
  719. c CAPT Factor used in Lentz iteration for A (Eq. R27)
  720. c T_k
  721. c CNUMER Numerator in capT (Eq. R28A)
  722. c N_k
  723. c CDENOM Denominator in capT (Eq. R28B)
  724. c D_k
  725. c CDTD Product of two successive denominators of capT factors
  726. c (Eq. R34C)
  727. c xi_1
  728. c CNTN Product of two successive numerators of capT factors
  729. c (Eq. R34B)
  730. c xi_2
  731. c EPS1 Ill-conditioning criterion
  732. c EPS2 Convergence criterion
  733. c KK Subscript k of cAk (Eq. R25B)
  734. c k
  735. c KOUNT Iteration counter (used to prevent infinite looping)
  736. c MAXIT Max. allowed no. of iterations
  737. c MM + 1 and - 1, alternately
  738. */
  739. std::complex<double> MultiLayerMie::calcD1confra(const int N, const std::complex<double> z) {
  740. // NTMR -> nmax_ - 1 \\TODO nmax_ ?
  741. //int N = nmax_ - 1;
  742. int KK, KOUNT, MAXIT = 10000, MM;
  743. // double EPS1=1.0e-2;
  744. double EPS2=1.0e-8;
  745. std::complex<double> CAK, CAPT, CDENOM, CDTD, CNTN, CNUMER;
  746. std::complex<double> one = std::complex<double>(1.0,0.0);
  747. std::complex<double> ZINV = one/z;
  748. // c ** Eq. R25a
  749. std::complex<double> CONFRA = static_cast<std::complex<double> >(N + 1)*ZINV; //debug ZINV
  750. MM = - 1;
  751. KK = 2*N +3; //debug 3
  752. // c ** Eq. R25b, k=2
  753. CAK = static_cast<std::complex<double> >(MM*KK)*ZINV; //debug -3 ZINV
  754. CDENOM = CAK;
  755. CNUMER = CDENOM + one/CONFRA; //-3zinv+z
  756. KOUNT = 1;
  757. //10 CONTINUE
  758. do { ++KOUNT;
  759. if (KOUNT > MAXIT) {
  760. printf("re(%g):im(%g)\t\n", CONFRA.real(), CONFRA.imag());
  761. throw std::invalid_argument("ConFra--Iteration failed to converge!\n");
  762. }
  763. MM *= - 1; KK += 2; //debug mm=1 kk=5
  764. CAK = static_cast<std::complex<double> >(MM*KK)*ZINV; // ** Eq. R25b //debug 5zinv
  765. // //c ** Eq. R32 Ill-conditioned case -- stride two terms instead of one
  766. // if (std::abs(CNUMER/CAK) >= EPS1 || std::abs(CDENOM/CAK) >= EPS1) {
  767. // //c ** Eq. R34
  768. // CNTN = CAK*CNUMER + 1.0;
  769. // CDTD = CAK*CDENOM + 1.0;
  770. // CONFRA = (CNTN/CDTD)*CONFRA; // ** Eq. R33
  771. // MM *= - 1; KK += 2;
  772. // CAK = static_cast<std::complex<double> >(MM*KK)*ZINV; // ** Eq. R25b
  773. // //c ** Eq. R35
  774. // CNUMER = CAK + CNUMER/CNTN;
  775. // CDENOM = CAK + CDENOM/CDTD;
  776. // ++KOUNT;
  777. // //GO TO 10
  778. // continue;
  779. // } else { //c *** Well-conditioned case
  780. {
  781. CAPT = CNUMER/CDENOM; // ** Eq. R27 //debug (-3zinv + z)/(-3zinv)
  782. // printf("re(%g):im(%g)**\t", CAPT.real(), CAPT.imag());
  783. CONFRA = CAPT*CONFRA; // ** Eq. R26
  784. //if (N == 0) {output=true;printf(" re:");prn(CONFRA.real());printf(" im:"); prn(CONFRA.imag());output=false;};
  785. //c ** Check for convergence; Eq. R31
  786. if (std::abs(CAPT.real() - 1.0) >= EPS2 || std::abs(CAPT.imag()) >= EPS2) {
  787. //c ** Eq. R30
  788. CNUMER = CAK + one/CNUMER;
  789. CDENOM = CAK + one/CDENOM;
  790. continue;
  791. //GO TO 10
  792. } // end of if < eps2
  793. }
  794. break;
  795. } while(1);
  796. //if (N == 0) printf(" return confra for z=(%g,%g)\n", ZINV.real(), ZINV.imag());
  797. return CONFRA;
  798. }
  799. //**********************************************************************************//
  800. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  801. // functions (D1 and D3) for a complex argument (z). //
  802. // Equations (16a), (16b) and (18a) - (18d) //
  803. // //
  804. // Input parameters: //
  805. // z: Complex argument to evaluate D1 and D3 //
  806. // nmax_: Maximum number of terms to calculate D1 and D3 //
  807. // //
  808. // Output parameters: //
  809. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  810. //**********************************************************************************//
  811. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  812. std::vector<std::complex<double> >& D1,
  813. std::vector<std::complex<double> >& D3) {
  814. // Downward recurrence for D1 - equations (16a) and (16b)
  815. D1[nmax_] = std::complex<double>(0.0, 0.0);
  816. //D1[nmax_] = calcD1confra(nmax_, z);
  817. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  818. // printf(" D:");prn((D1[nmax_]).real()); printf("\t diff:");
  819. // prn((D1[nmax_] + double(nmax_)*zinv).real());
  820. for (int n = nmax_; n > 0; n--) {
  821. D1[n - 1] = double(n)*zinv - 1.0/(D1[n] + double(n)*zinv);
  822. //D1[n - 1] = calcD1confra(n - 1, z);
  823. // printf(" D:");prn((D1[n - 1]).real()); printf("\t diff:");
  824. // prn((D1[n] + double(n)*zinv).real());
  825. }
  826. // printf("\n\n"); iformat=0;
  827. if (std::abs(D1[0]) > 100000.0)
  828. throw std::invalid_argument
  829. ("Unstable D1! Please, try to change input parameters!\n");
  830. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  831. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  832. *std::exp(-2.0*z.imag()));
  833. D3[0] = std::complex<double>(0.0, 1.0);
  834. for (int n = 1; n <= nmax_; n++) {
  835. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  836. *(static_cast<double>(n)*zinv- D3[n - 1]);
  837. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  838. }
  839. }
  840. //**********************************************************************************//
  841. // This function calculates Pi and Tau for all values of Theta. //
  842. // Equations (26a) - (26c) //
  843. // //
  844. // Input parameters: //
  845. // nmax_: Maximum number of terms to calculate Pi and Tau //
  846. // nTheta: Number of scattering angles //
  847. // Theta: Array containing all the scattering angles where the scattering //
  848. // amplitudes will be calculated //
  849. // //
  850. // Output parameters: //
  851. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  852. //**********************************************************************************//
  853. void MultiLayerMie::calcSinglePiTau(const double& costheta, std::vector<double>& Pi,
  854. std::vector<double>& Tau) {
  855. //****************************************************//
  856. // Equations (26a) - (26c) //
  857. //****************************************************//
  858. for (int n = 0; n < nmax_; n++) {
  859. if (n == 0) {
  860. // Initialize Pi and Tau
  861. Pi[n] = 1.0;
  862. Tau[n] = (n + 1)*costheta;
  863. } else {
  864. // Calculate the actual values
  865. Pi[n] = ((n == 1) ? ((n + n + 1)*costheta*Pi[n - 1]/n)
  866. : (((n + n + 1)*costheta*Pi[n - 1]
  867. - (n + 1)*Pi[n - 2])/n));
  868. Tau[n] = (n + 1)*costheta*Pi[n] - (n + 2)*Pi[n - 1];
  869. }
  870. }
  871. } // end of void MultiLayerMie::calcPiTau(...)
  872. void MultiLayerMie::calcAllPiTau(std::vector< std::vector<double> >& Pi,
  873. std::vector< std::vector<double> >& Tau) {
  874. std::vector<double> costheta(theta_.size(), 0.0);
  875. for (int t = 0; t < theta_.size(); t++) {
  876. costheta[t] = std::cos(theta_[t]);
  877. }
  878. // Do not join upper and lower 'for' to a single one! It will slow
  879. // down the code!!! (For about 0.5-2.0% of runtime, it is probably
  880. // due to increased cache missing rate originated from the
  881. // recurrence in calcPiTau...)
  882. for (int t = 0; t < theta_.size(); t++) {
  883. calcSinglePiTau(costheta[t], Pi[t], Tau[t]);
  884. //calcSinglePiTau(std::cos(theta_[t]), Pi[t], Tau[t]); // It is slow!!
  885. }
  886. } // end of void MultiLayerMie::calcAllPiTau(...)
  887. //**********************************************************************************//
  888. // This function calculates the scattering coefficients required to calculate //
  889. // both the near- and far-field parameters. //
  890. // //
  891. // Input parameters: //
  892. // L: Number of layers //
  893. // pl: Index of PEC layer. If there is none just send -1 //
  894. // x: Array containing the size parameters of the layers [0..L-1] //
  895. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  896. // nmax: Maximum number of multipolar expansion terms to be used for the //
  897. // calculations. Only use it if you know what you are doing, otherwise //
  898. // set this parameter to -1 and the function will calculate it. //
  899. // //
  900. // Output parameters: //
  901. // an, bn: Complex scattering amplitudes //
  902. // //
  903. // Return value: //
  904. // Number of multipolar expansion terms used for the calculations //
  905. //**********************************************************************************//
  906. void MultiLayerMie::ExtScattCoeffs(std::vector<std::complex<double> >& an,
  907. std::vector<std::complex<double> >& bn) {
  908. const std::vector<double>& x = size_parameter_;
  909. const std::vector<std::complex<double> >& m = index_;
  910. const int& pl = PEC_layer_position_;
  911. const int L = index_.size();
  912. //************************************************************************//
  913. // Calculate the index of the first layer. It can be either 0 (default) //
  914. // or the index of the outermost PEC layer. In the latter case all layers //
  915. // below the PEC are discarded. //
  916. // ***********************************************************************//
  917. // TODO, is it possible for PEC to have a zero index? If yes than
  918. // is should be:
  919. // int fl = (pl > - 1) ? pl : 0;
  920. // This will give the same result, however, it corresponds the
  921. // logic - if there is PEC, than first layer is PEC.
  922. // Well, I followed the logic: First layer is always zero unless it has
  923. // an upper PEC layer.
  924. int fl = (pl > 0) ? pl : 0;
  925. if (nmax_ <= 0) Nmax(fl);
  926. std::complex<double> z1, z2;
  927. //**************************************************************************//
  928. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  929. // means that index = layer number - 1 or index = n - 1. The only exception //
  930. // are the arrays for representing D1, D3 and Q because they need a value //
  931. // for the index 0 (zero), hence it is important to consider this shift //
  932. // between different arrays. The change was done to optimize memory usage. //
  933. //**************************************************************************//
  934. // Allocate memory to the arrays
  935. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  936. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  937. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  938. for (int l = 0; l < L; l++) {
  939. Q[l].resize(nmax_ + 1);
  940. Ha[l].resize(nmax_);
  941. Hb[l].resize(nmax_);
  942. }
  943. an.resize(nmax_);
  944. bn.resize(nmax_);
  945. PsiZeta_.resize(nmax_ + 1);
  946. std::vector<std::complex<double> > D1XL(nmax_ + 1), D3XL(nmax_ + 1),
  947. PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  948. //*************************************************//
  949. // Calculate D1 and D3 for z1 in the first layer //
  950. //*************************************************//
  951. if (fl == pl) { // PEC layer
  952. for (int n = 0; n <= nmax_; n++) {
  953. D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
  954. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  955. }
  956. } else { // Regular layer
  957. z1 = x[fl]* m[fl];
  958. // Calculate D1 and D3
  959. calcD1D3(z1, D1_mlxl, D3_mlxl);
  960. }
  961. // do { \
  962. // ++iformat;\
  963. // if (iformat%5 == 0) printf("%24.16e",z1.real());
  964. // } while (false);
  965. //******************************************************************//
  966. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  967. //******************************************************************//
  968. for (int n = 0; n < nmax_; n++) {
  969. Ha[fl][n] = D1_mlxl[n + 1];
  970. Hb[fl][n] = D1_mlxl[n + 1];
  971. }
  972. //*****************************************************//
  973. // Iteration from the second layer to the last one (L) //
  974. //*****************************************************//
  975. std::complex<double> Temp, Num, Denom;
  976. std::complex<double> G1, G2;
  977. for (int l = fl + 1; l < L; l++) {
  978. //************************************************************//
  979. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  980. //************************************************************//
  981. z1 = x[l]*m[l];
  982. z2 = x[l - 1]*m[l];
  983. //Calculate D1 and D3 for z1
  984. calcD1D3(z1, D1_mlxl, D3_mlxl);
  985. //Calculate D1 and D3 for z2
  986. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  987. // prn(z1.real());
  988. // for (auto i : D1_mlxl) { prn(i.real());
  989. // // prn(i.imag());
  990. // } printf("\n");
  991. //*********************************************//
  992. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  993. //*********************************************//
  994. // Upward recurrence for Q - equations (19a) and (19b)
  995. Num = std::exp(-2.0*(z1.imag() - z2.imag()))
  996. *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
  997. Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
  998. Q[l][0] = Num/Denom;
  999. for (int n = 1; n <= nmax_; n++) {
  1000. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  1001. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  1002. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  1003. }
  1004. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  1005. for (int n = 1; n <= nmax_; n++) {
  1006. //Ha
  1007. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  1008. G1 = -D1_mlxlM1[n];
  1009. G2 = -D3_mlxlM1[n];
  1010. } else {
  1011. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  1012. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  1013. } // end of if PEC
  1014. Temp = Q[l][n]*G1;
  1015. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  1016. Denom = G2 - Temp;
  1017. Ha[l][n - 1] = Num/Denom;
  1018. //Hb
  1019. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  1020. G1 = Hb[l - 1][n - 1];
  1021. G2 = Hb[l - 1][n - 1];
  1022. } else {
  1023. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  1024. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  1025. } // end of if PEC
  1026. Temp = Q[l][n]*G1;
  1027. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  1028. Denom = (G2- Temp);
  1029. Hb[l][n - 1] = (Num/ Denom);
  1030. } // end of for Ha and Hb terms
  1031. } // end of for layers iteration
  1032. //**************************************//
  1033. //Calculate D1, D3, Psi and Zeta for XL //
  1034. //**************************************//
  1035. // Calculate D1XL and D3XL
  1036. calcD1D3(x[L - 1], D1XL, D3XL);
  1037. //printf("%5.20f\n",Ha[L - 1][0].real());
  1038. // Calculate PsiXL and ZetaXL
  1039. calcPsiZeta(x[L - 1], D1XL, D3XL, PsiXL, ZetaXL);
  1040. //*********************************************************************//
  1041. // Finally, we calculate the scattering coefficients (an and bn) and //
  1042. // the angular functions (Pi and Tau). Note that for these arrays the //
  1043. // first layer is 0 (zero), in future versions all arrays will follow //
  1044. // this convention to save memory. (13 Nov, 2014) //
  1045. //*********************************************************************//
  1046. for (int n = 0; n < nmax_; n++) {
  1047. //********************************************************************//
  1048. //Expressions for calculating an and bn coefficients are not valid if //
  1049. //there is only one PEC layer (ie, for a simple PEC sphere). //
  1050. //********************************************************************//
  1051. if (pl < (L - 1)) {
  1052. 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]);
  1053. 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]);
  1054. } else {
  1055. 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]);
  1056. bn[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  1057. }
  1058. } // end of for an and bn terms
  1059. } // end of void MultiLayerMie::ExtScattCoeffs(...)
  1060. // ********************************************************************** //
  1061. // ********************************************************************** //
  1062. // ********************************************************************** //
  1063. void MultiLayerMie::InitMieCalculations() {
  1064. isMieCalculated_ = false;
  1065. // Initialize the scattering parameters
  1066. Qext_ = 0;
  1067. Qsca_ = 0;
  1068. Qabs_ = 0;
  1069. Qbk_ = 0;
  1070. Qpr_ = 0;
  1071. asymmetry_factor_ = 0;
  1072. albedo_ = 0;
  1073. Qsca_ch_.clear();
  1074. Qext_ch_.clear();
  1075. Qabs_ch_.clear();
  1076. Qbk_ch_.clear();
  1077. Qpr_ch_.clear();
  1078. Qsca_ch_.resize(nmax_ - 1);
  1079. Qext_ch_.resize(nmax_ - 1);
  1080. Qabs_ch_.resize(nmax_ - 1);
  1081. Qbk_ch_.resize(nmax_ - 1);
  1082. Qpr_ch_.resize(nmax_ - 1);
  1083. Qsca_ch_norm_.resize(nmax_ - 1);
  1084. Qext_ch_norm_.resize(nmax_ - 1);
  1085. Qabs_ch_norm_.resize(nmax_ - 1);
  1086. Qbk_ch_norm_.resize(nmax_ - 1);
  1087. Qpr_ch_norm_.resize(nmax_ - 1);
  1088. // Initialize the scattering amplitudes
  1089. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  1090. S1_.swap(tmp1);
  1091. S2_ = S1_;
  1092. }
  1093. // ********************************************************************** //
  1094. // ********************************************************************** //
  1095. // ********************************************************************** //
  1096. void MultiLayerMie::ConvertToSP() {
  1097. isMieCalculated_ = false;
  1098. if (target_width_.size() + coating_width_.size() == 0)
  1099. return; // Nothing to convert, we suppose that SP was set directly
  1100. GenerateSizeParameter();
  1101. GenerateIndex();
  1102. if (size_parameter_.size() != index_.size())
  1103. throw std::invalid_argument("Ivalid conversion of width to size parameter units!/n");
  1104. }
  1105. // ********************************************************************** //
  1106. // ********************************************************************** //
  1107. // ********************************************************************** //
  1108. //**********************************************************************************//
  1109. // This function calculates the actual scattering parameters and amplitudes //
  1110. // //
  1111. // Input parameters: //
  1112. // L: Number of layers //
  1113. // pl: Index of PEC layer. If there is none just send -1 //
  1114. // x: Array containing the size parameters of the layers [0..L-1] //
  1115. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1116. // nTheta: Number of scattering angles //
  1117. // Theta: Array containing all the scattering angles where the scattering //
  1118. // amplitudes will be calculated //
  1119. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  1120. // calculations. Only use it if you know what you are doing, otherwise //
  1121. // set this parameter to -1 and the function will calculate it //
  1122. // //
  1123. // Output parameters: //
  1124. // Qext: Efficiency factor for extinction //
  1125. // Qsca: Efficiency factor for scattering //
  1126. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  1127. // Qbk: Efficiency factor for backscattering //
  1128. // Qpr: Efficiency factor for the radiation pressure //
  1129. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  1130. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  1131. // S1, S2: Complex scattering amplitudes //
  1132. // //
  1133. // Return value: //
  1134. // Number of multipolar expansion terms used for the calculations //
  1135. //**********************************************************************************//
  1136. void MultiLayerMie::RunMieCalculations() {
  1137. isMieCalculated_ = false;
  1138. ConvertToSP();
  1139. nmax_ = nmax_preset_;
  1140. if (size_parameter_.size() != index_.size())
  1141. throw std::invalid_argument("Each size parameter should have only one index!");
  1142. if (size_parameter_.size() == 0)
  1143. throw std::invalid_argument("Initialize model first!");
  1144. const std::vector<double>& x = size_parameter_;
  1145. // Calculate scattering coefficients
  1146. ExtScattCoeffs(an_, bn_);
  1147. // std::vector< std::vector<double> > Pi(nmax_), Tau(nmax_);
  1148. std::vector< std::vector<double> > Pi, Tau;
  1149. Pi.resize(theta_.size());
  1150. Tau.resize(theta_.size());
  1151. for (int i =0; i< theta_.size(); ++i) {
  1152. Pi[i].resize(nmax_);
  1153. Tau[i].resize(nmax_);
  1154. }
  1155. calcAllPiTau(Pi, Tau);
  1156. InitMieCalculations(); //
  1157. std::complex<double> Qbktmp(0.0, 0.0);
  1158. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  1159. // By using downward recurrence we avoid loss of precision due to float rounding errors
  1160. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  1161. // http://en.wikipedia.org/wiki/Loss_of_significance
  1162. for (int i = nmax_ - 2; i >= 0; i--) {
  1163. const int n = i + 1;
  1164. // Equation (27)
  1165. Qext_ch_norm_[i] = (an_[i].real() + bn_[i].real());
  1166. Qext_ch_[i] = (n + n + 1.0)*Qext_ch_norm_[i];
  1167. //Qext_ch_[i] = (n + n + 1)*(an_[i].real() + bn_[i].real());
  1168. Qext_ += Qext_ch_[i];
  1169. // Equation (28)
  1170. Qsca_ch_norm_[i] = (an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  1171. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  1172. Qsca_ch_[i] = (n + n + 1.0)*Qsca_ch_norm_[i];
  1173. Qsca_ += Qsca_ch_[i];
  1174. // Qsca_ch_[i] += (n + n + 1)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  1175. // + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  1176. // Equation (29) TODO We must check carefully this equation. If we
  1177. // remove the typecast to double then the result changes. Which is
  1178. // the correct one??? Ovidio (2014/12/10) With cast ratio will
  1179. // give double, without cast (n + n + 1)/(n*(n + 1)) will be
  1180. // rounded to integer. Tig (2015/02/24)
  1181. Qpr_ch_[i]=((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  1182. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  1183. Qpr_ += Qpr_ch_[i];
  1184. // Equation (33)
  1185. Qbktmp_ch[i] = (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  1186. Qbktmp += Qbktmp_ch[i];
  1187. // Calculate the scattering amplitudes (S1 and S2) //
  1188. // Equations (25a) - (25b) //
  1189. for (int t = 0; t < theta_.size(); t++) {
  1190. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[t][i], Tau[t][i]);
  1191. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[t][i], Tau[t][i]);
  1192. }
  1193. }
  1194. double x2 = pow2(x.back());
  1195. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  1196. for (double& Q : Qext_ch_) Q = 2.0*Q/x2;
  1197. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  1198. for (double& Q : Qsca_ch_) Q = 2.0*Q/x2;
  1199. //for (double& Q : Qsca_ch_norm_) Q = 2.0*Q/x2;
  1200. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  1201. for (int i = 0; i < nmax_ - 1; ++i) Qpr_ch_[i] = Qext_ch_[i] - 4.0*Qpr_ch_[i]/x2;
  1202. Qabs_ = Qext_ - Qsca_; // Equation (30)
  1203. for (int i = 0; i < nmax_ - 1; ++i) {
  1204. Qabs_ch_[i] = Qext_ch_[i] - Qsca_ch_[i];
  1205. Qabs_ch_norm_[i] = Qext_ch_norm_[i] - Qsca_ch_norm_[i];
  1206. }
  1207. albedo_ = Qsca_/Qext_; // Equation (31)
  1208. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  1209. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  1210. isMieCalculated_ = true;
  1211. nmax_used_ = nmax_;
  1212. // printf("Run Mie result: Qext = %g, Qsca = %g, Qabs = %g, Qbk = %g \n",
  1213. // GetQext(), GetQsca(), GetQabs(), GetQbk());
  1214. //return nmax;
  1215. }
  1216. // ********************************************************************** //
  1217. // ********************************************************************** //
  1218. // ********************************************************************** //
  1219. void MultiLayerMie::IntScattCoeffsInit() {
  1220. const int L = index_.size();
  1221. // we need to fill
  1222. // std::vector< std::vector<std::complex<double> > > al_n_, bl_n_, cl_n_, dl_n_;
  1223. // for n = [0..nmax_) and for l=[L..0)
  1224. // TODO: to decrease cache miss outer loop is with n and inner with reversed l
  1225. // at the moment outer is forward l and inner in n
  1226. al_n_.resize(L + 1);
  1227. bl_n_.resize(L + 1);
  1228. cl_n_.resize(L + 1);
  1229. dl_n_.resize(L + 1);
  1230. for (auto& element:al_n_) element.resize(nmax_);
  1231. for (auto& element:bl_n_) element.resize(nmax_);
  1232. for (auto& element:cl_n_) element.resize(nmax_);
  1233. for (auto& element:dl_n_) element.resize(nmax_);
  1234. std::complex<double> c_one(1.0, 0.0);
  1235. std::complex<double> c_zero(0.0, 0.0);
  1236. // Yang, paragraph under eq. A3
  1237. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  1238. for (int i = 0; i < nmax_; ++i) {
  1239. al_n_[L][i] = an_[i];
  1240. bl_n_[L][i] = bn_[i];
  1241. cl_n_[L][i] = c_one;
  1242. dl_n_[L][i] = c_one;
  1243. if (i < 3) printf(" (%g) ", std::abs(an_[i]));
  1244. }
  1245. }
  1246. // ********************************************************************** //
  1247. // ********************************************************************** //
  1248. // ********************************************************************** //
  1249. void MultiLayerMie::IntScattCoeffs() {
  1250. if (!isMieCalculated_)
  1251. throw std::invalid_argument("(IntScattCoeffs) You should run calculations first!");
  1252. IntScattCoeffsInit();
  1253. const int L = index_.size();
  1254. std::vector<std::complex<double> > z(L), z1(L);
  1255. for (int i = 0; i < L - 1; ++i) {
  1256. z[i] =size_parameter_[i]*index_[i];
  1257. z1[i]=size_parameter_[i]*index_[i + 1];
  1258. }
  1259. z[L - 1] = size_parameter_[L - 1]*index_[L - 1];
  1260. z1[L - 1] = size_parameter_[L - 1];
  1261. std::vector< std::vector<std::complex<double> > > D1z(L), D1z1(L), D3z(L), D3z1(L);
  1262. std::vector< std::vector<std::complex<double> > > Psiz(L), Psiz1(L), Zetaz(L), Zetaz1(L);
  1263. for (int l = 0; l < L; ++l) {
  1264. D1z[l].resize(nmax_ + 1);
  1265. D1z1[l].resize(nmax_ + 1);
  1266. D3z[l].resize(nmax_ + 1);
  1267. D3z1[l].resize(nmax_ + 1);
  1268. Psiz[l].resize(nmax_ + 1);
  1269. Psiz1[l].resize(nmax_ + 1);
  1270. Zetaz[l].resize(nmax_ + 1);
  1271. Zetaz1[l].resize(nmax_ + 1);
  1272. }
  1273. for (int l = 0; l < L; ++l) {
  1274. calcD1D3(z[l],D1z[l],D3z[l]);
  1275. calcD1D3(z1[l],D1z1[l],D3z1[l]);
  1276. calcPsiZeta(z[l],D1z[l],D3z[l], Psiz[l],Zetaz[l]);
  1277. calcPsiZeta(z1[l],D1z1[l],D3z1[l], Psiz1[l],Zetaz1[l]);
  1278. }
  1279. auto& m = index_;
  1280. std::vector< std::complex<double> > m1(L);
  1281. for (int l = 0; l < L - 1; ++l) m1[l] = m[l + 1];
  1282. m1[L - 1] = std::complex<double> (1.0, 0.0);
  1283. // for (auto zz : m) printf ("m[i]=%g \n\n ", zz.real());
  1284. for (int l = L - 1; l >= 0; --l) {
  1285. for (int n = 0; n < nmax_; ++n) {
  1286. // al_n
  1287. auto denom = m1[l]*Zetaz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1288. al_n_[l][n] = D1z[l][n + 1]*m1[l]*(al_n_[l + 1][n]*Zetaz1[l][n + 1] - dl_n_[l + 1][n]*Psiz1[l][n + 1])
  1289. - m[l]*(-D1z1[l][n + 1]*dl_n_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*al_n_[l + 1][n]*Zetaz1[l][n + 1]);
  1290. al_n_[l][n] /= denom;
  1291. // dl_n
  1292. denom = m1[l]*Psiz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1293. dl_n_[l][n] = D3z[l][n + 1]*m1[l]*(al_n_[l + 1][n]*Zetaz1[l][n + 1] - dl_n_[l + 1][n]*Psiz1[l][n + 1])
  1294. - m[l]*(-D1z1[l][n + 1]*dl_n_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*al_n_[l + 1][n]*Zetaz1[l][n + 1]);
  1295. dl_n_[l][n] /= denom;
  1296. // bl_n
  1297. denom = m1[l]*Zetaz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1298. bl_n_[l][n] = D1z[l][n + 1]*m[l]*(bl_n_[l + 1][n]*Zetaz1[l][n + 1] - cl_n_[l + 1][n]*Psiz1[l][n + 1])
  1299. - m1[l]*(-D1z1[l][n + 1]*cl_n_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*bl_n_[l + 1][n]*Zetaz1[l][n + 1]);
  1300. bl_n_[l][n] /= denom;
  1301. // cl_n
  1302. denom = m1[l]*Psiz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1303. cl_n_[l][n] = D3z[l][n + 1]*m[l]*(bl_n_[l + 1][n]*Zetaz1[l][n + 1] - cl_n_[l + 1][n]*Psiz1[l][n + 1])
  1304. - m1[l]*(-D1z1[l][n + 1]*cl_n_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*bl_n_[l + 1][n]*Zetaz1[l][n + 1]);
  1305. cl_n_[l][n] /= denom;
  1306. } // end of all n
  1307. } // end of for all l
  1308. // Check the result and change an__0 and bn__0 for exact zero
  1309. for (int n = 0; n < nmax_; ++n) {
  1310. if (std::abs(al_n_[0][n]) < 1e-10) al_n_[0][n] = 0.0;
  1311. else throw std::invalid_argument("Unstable calculation of a__0_n!");
  1312. if (std::abs(bl_n_[0][n]) < 1e-10) bl_n_[0][n] = 0.0;
  1313. else throw std::invalid_argument("Unstable calculation of b__0_n!");
  1314. }
  1315. // for (int l = 0; l < L; ++l) {
  1316. // printf("l=%d --> ", l);
  1317. // for (int n = 0; n < nmax_ + 1; ++n) {
  1318. // if (n < 20) continue;
  1319. // printf("n=%d --> D1zn=%g, D3zn=%g, D1zn=%g, D3zn=%g || ",
  1320. // n,
  1321. // D1z[l][n].real(), D3z[l][n].real(),
  1322. // D1z1[l][n].real(), D3z1[l][n].real());
  1323. // }
  1324. // printf("\n\n");
  1325. // }
  1326. // for (int l = 0; l < L; ++l) {
  1327. // printf("l=%d --> ", l);
  1328. // for (int n = 0; n < nmax_ + 1; ++n) {
  1329. // printf("n=%d --> D1zn=%g, D3zn=%g, D1zn=%g, D3zn=%g || ",
  1330. // n,
  1331. // D1z[l][n].real(), D3z[l][n].real(),
  1332. // D1z1[l][n].real(), D3z1[l][n].real());
  1333. // }
  1334. // printf("\n\n");
  1335. // }
  1336. for (int i = 0; i < L + 1; ++i) {
  1337. printf("Layer =%d ---> ", i);
  1338. for (int n = 0; n < nmax_; ++n) {
  1339. // if (n < 20) continue;
  1340. printf(" || n=%d --> a=%g,%g b=%g,%g c=%g,%g d=%g,%g",
  1341. n,
  1342. al_n_[i][n].real(), al_n_[i][n].imag(),
  1343. bl_n_[i][n].real(), bl_n_[i][n].imag(),
  1344. cl_n_[i][n].real(), cl_n_[i][n].imag(),
  1345. dl_n_[i][n].real(), dl_n_[i][n].imag());
  1346. }
  1347. printf("\n\n");
  1348. }
  1349. }
  1350. // ********************************************************************** //
  1351. // ********************************************************************** //
  1352. // ********************************************************************** //
  1353. // external scattering field = incident + scattered
  1354. // BH p.92 (4.37), 94 (4.45), 95 (4.50)
  1355. // assume: medium is non-absorbing; refim = 0; Uabs = 0
  1356. void MultiLayerMie::fieldExt(const double Rho, const double Phi, const double Theta, const std::vector<double>& Pi, const std::vector<double>& Tau, std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1357. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0);
  1358. std::vector<std::complex<double> > vm3o1n(3), vm3e1n(3), vn3o1n(3), vn3e1n(3);
  1359. std::vector<std::complex<double> > Ei(3,c_zero), Hi(3,c_zero), Es(3,c_zero), Hs(3,c_zero);
  1360. std::vector<std::complex<double> > bj(nmax_ + 1), by(nmax_ + 1), bd(nmax_ + 1);
  1361. // Calculate spherical Bessel and Hankel functions
  1362. printf("########## layer OUT ############\n");
  1363. sphericalBessel(Rho,bj, by, bd);
  1364. for (int n = 0; n < nmax_; n++) {
  1365. double rn = static_cast<double>(n + 1);
  1366. std::complex<double> zn = bj[n + 1] + c_i*by[n + 1];
  1367. // using BH 4.12 and 4.50
  1368. std::complex<double> xxip = Rho*(bj[n] + c_i*by[n]) - rn*zn;
  1369. using std::sin;
  1370. using std::cos;
  1371. vm3o1n[0] = c_zero;
  1372. vm3o1n[1] = cos(Phi)*Pi[n]*zn;
  1373. vm3o1n[2] = -sin(Phi)*Tau[n]*zn;
  1374. vm3e1n[0] = c_zero;
  1375. vm3e1n[1] = -sin(Phi)*Pi[n]*zn;
  1376. vm3e1n[2] = -cos(Phi)*Tau[n]*zn;
  1377. vn3o1n[0] = sin(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1378. vn3o1n[1] = sin(Phi)*Tau[n]*xxip/Rho;
  1379. vn3o1n[2] = cos(Phi)*Pi[n]*xxip/Rho;
  1380. vn3e1n[0] = cos(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1381. vn3e1n[1] = cos(Phi)*Tau[n]*xxip/Rho;
  1382. vn3e1n[2] = -sin(Phi)*Pi[n]*xxip/Rho;
  1383. // scattered field: BH p.94 (4.45)
  1384. std::complex<double> encap = std::pow(c_i, rn)*(2.0*rn + 1.0)/(rn*rn + rn);
  1385. for (int i = 0; i < 3; i++) {
  1386. Es[i] = Es[i] + encap*(c_i*an_[n]*vn3e1n[i] - bn_[n]*vm3o1n[i]);
  1387. Hs[i] = Hs[i] + encap*(c_i*bn_[n]*vn3o1n[i] + an_[n]*vm3e1n[i]);
  1388. //if (n<3) printf(" E[%d]=%g ", i,std::abs(Es[i]));
  1389. if (n<3) printf(" !!=%d=== %g ", i,std::abs(Es[i]));
  1390. }
  1391. }
  1392. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1393. // basis unit vectors = er, etheta, ephi
  1394. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  1395. {
  1396. using std::sin;
  1397. using std::cos;
  1398. Ei[0] = eifac*sin(Theta)*cos(Phi);
  1399. Ei[1] = eifac*cos(Theta)*cos(Phi);
  1400. Ei[2] = -eifac*sin(Phi);
  1401. }
  1402. // magnetic field
  1403. double hffact = 1.0/(cc_*mu_);
  1404. for (int i = 0; i < 3; i++) {
  1405. Hs[i] = hffact*Hs[i];
  1406. }
  1407. // incident H field: BH p.26 (2.43), p.89 (4.21)
  1408. std::complex<double> hffacta = hffact;
  1409. std::complex<double> hifac = eifac*hffacta;
  1410. {
  1411. using std::sin;
  1412. using std::cos;
  1413. Hi[0] = hifac*sin(Theta)*sin(Phi);
  1414. Hi[1] = hifac*cos(Theta)*sin(Phi);
  1415. Hi[2] = hifac*cos(Phi);
  1416. }
  1417. for (int i = 0; i < 3; i++) {
  1418. // electric field E [V m - 1] = EF*E0
  1419. E[i] = Ei[i] + Es[i];
  1420. H[i] = Hi[i] + Hs[i];
  1421. // printf("ext E[%d]=%g",i,std::abs(E[i]));
  1422. }
  1423. } // end of void fieldExt(...)
  1424. // ********************************************************************** //
  1425. // ********************************************************************** //
  1426. // ********************************************************************** //
  1427. void MultiLayerMie::fieldInt(const double Rho, const double Phi, const double Theta, const std::vector<double>& Pi, const std::vector<double>& Tau, std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1428. // printf("field int Qext = %g, Qsca = %g, Qabs = %g, Qbk = %g, \n",
  1429. // GetQext(), GetQsca(), GetQabs(), GetQbk());
  1430. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1431. std::vector<std::complex<double> > vm3o1n(3), vm3e1n(3), vn3o1n(3), vn3e1n(3);
  1432. std::vector<std::complex<double> > vm1o1n(3), vm1e1n(3), vn1o1n(3), vn1e1n(3);
  1433. std::vector<std::complex<double> > El(3,c_zero),Ei(3,c_zero), Hl(3,c_zero);
  1434. std::vector<std::complex<double> > bj(nmax_ + 1), by(nmax_ + 1), bd(nmax_ + 1);
  1435. int layer=0; // layer number
  1436. std::complex<double> index_l;
  1437. for (int i = 0; i < size_parameter_.size() - 1; ++i) {
  1438. if (size_parameter_[i] < Rho && Rho <= size_parameter_[i + 1]) {
  1439. layer=i;
  1440. }
  1441. }
  1442. if (Rho > size_parameter_.back()) {
  1443. layer = size_parameter_.size();
  1444. index_l = c_one;
  1445. } else {
  1446. index_l = index_[layer];
  1447. }
  1448. std::complex<double> bessel_arg = Rho*index_l;
  1449. std::complex<double>& rh = bessel_arg;
  1450. std::complex<double> besselj_1 = std::sin(rh)/pow2(rh)-std::cos(rh)/rh;
  1451. printf("bessel arg = %g,%g index=%g,%g besselj[1]=%g,%g\n", bessel_arg.real(), bessel_arg.imag(), index_l.real(), index_l.imag(), besselj_1.real(), besselj_1.imag());
  1452. const int& l = layer;
  1453. printf("########## layer %d ############\n",l);
  1454. // Calculate spherical Bessel and Hankel functions
  1455. sphericalBessel(bessel_arg,bj, by, bd);
  1456. printf("besselj[1]=%g,%g\n", bj[1].real(), bj[1].imag());
  1457. printf("bessely[1]=%g,%g\n", by[1].real(), by[1].imag());
  1458. for (int n = 0; n < nmax_; n++) {
  1459. double rn = static_cast<double>(n + 1);
  1460. std::complex<double> znm1 = bj[n] + c_i*by[n];
  1461. std::complex<double> zn = bj[n + 1] + c_i*by[n + 1];
  1462. //if (n<3) printf("\nbesselh = %g,%g", zn.real(), zn.imag()); //!
  1463. // using BH 4.12 and 4.50
  1464. std::complex<double> xxip = Rho*(bj[n] + c_i*by[n]) - rn*zn;
  1465. //if (n<3) printf("\nxxip = %g,%g", xxip.real(), xxip.imag()); //!
  1466. using std::sin;
  1467. using std::cos;
  1468. vm3o1n[0] = c_zero;
  1469. vm3o1n[1] = cos(Phi)*Pi[n]*zn;
  1470. vm3o1n[2] = -sin(Phi)*Tau[n]*zn;
  1471. // if (n<3) printf("\nRE vm3o1n[0]%g vm3o1n[1]%g vm3o1n[2]%g \nIM vm3o1n[0]%g vm3o1n[1]%g vm3o1n[2]%g",
  1472. // vm3o1n[0].real(), vm3o1n[1].real(), vm3o1n[2].real(),
  1473. // vm3o1n[0].imag(), vm3o1n[1].imag(), vm3o1n[2].imag());
  1474. vm3e1n[0] = c_zero;
  1475. vm3e1n[1] = -sin(Phi)*Pi[n]*zn;
  1476. vm3e1n[2] = -cos(Phi)*Tau[n]*zn;
  1477. vn3o1n[0] = sin(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1478. vn3o1n[1] = sin(Phi)*Tau[n]*xxip/Rho;
  1479. vn3o1n[2] = cos(Phi)*Pi[n]*xxip/Rho;
  1480. vn3e1n[0] = cos(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1481. vn3e1n[1] = cos(Phi)*Tau[n]*xxip/Rho;
  1482. vn3e1n[2] = -sin(Phi)*Pi[n]*xxip/Rho;
  1483. // if (n<3) printf("\nRE vn3e1n[0]%g vn3e1n[1]%g vn3e1n[2]%g \nIM vn3e1n[0]%g vn3e1n[1]%g vn3e1n[2]%g",
  1484. // vn3e1n[0].real(), vn3e1n[1].real(), vn3e1n[2].real(),
  1485. // vn3e1n[0].imag(), vn3e1n[1].imag(), vn3e1n[2].imag());
  1486. znm1 = bj[n];
  1487. zn = bj[n + 1];
  1488. // znm1 = (bj[n] + c_i*by[n]).real();
  1489. // zn = (bj[n + 1] + c_i*by[n + 1]).real();
  1490. xxip = Rho*(bj[n]) - rn*zn;
  1491. if (n<3)printf("\nbesselj = %g,%g", zn.real(), zn.imag()); //!
  1492. vm1o1n[0] = c_zero;
  1493. vm1o1n[1] = cos(Phi)*Pi[n]*zn;
  1494. vm1o1n[2] = -sin(Phi)*Tau[n]*zn;
  1495. vm1e1n[0] = c_zero;
  1496. vm1e1n[1] = -sin(Phi)*Pi[n]*zn;
  1497. vm1e1n[2] = -cos(Phi)*Tau[n]*zn;
  1498. vn1o1n[0] = sin(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1499. vn1o1n[1] = sin(Phi)*Tau[n]*xxip/Rho;
  1500. vn1o1n[2] = cos(Phi)*Pi[n]*xxip/Rho;
  1501. // if (n<3) printf("\nvn1o1n[2](%g) = cos(Phi)(%g)*Pi[n](%g)*xxip(%g)/Rho(%g)",
  1502. // std::abs(vn1o1n[2]), cos(Phi),Pi[n],std::abs(xxip),Rho);
  1503. vn1e1n[0] = cos(Phi)*rn*(rn + 1.0)*sin(Theta)*Pi[n]*zn/Rho;
  1504. vn1e1n[1] = cos(Phi)*Tau[n]*xxip/Rho;
  1505. vn1e1n[2] = -sin(Phi)*Pi[n]*xxip/Rho;
  1506. // if (n<3) printf("\nRE vm3o1n[0]%g vm3o1n[1]%g vm3o1n[2]%g \nIM vm3o1n[0]%g vm3o1n[1]%g vm3o1n[2]%g",
  1507. // vm3o1n[0].real(), vm3o1n[1].real(), vm3o1n[2].real(),
  1508. // vm3o1n[0].imag(), vm3o1n[1].imag(), vm3o1n[2].imag());
  1509. // scattered field: BH p.94 (4.45)
  1510. std::complex<double> encap = std::pow(c_i, rn)*(2.0*rn + 1.0)/(rn*rn + rn);
  1511. // if (n<3) printf("\n===== n=%d ======\n",n);
  1512. for (int i = 0; i < 3; i++) {
  1513. // if (n<3 && i==0) printf("\nn=%d",n);
  1514. // if (n<3) printf("\nbefore !El[%d]=%g,%g! ", i, El[i].real(), El[i].imag());
  1515. Ei[i] = encap*(cl_n_[l][n]*vm1o1n[i] - c_i*dl_n_[l][n]*vn1e1n[i]
  1516. + c_i*al_n_[l][n]*vn3e1n[i] - bl_n_[l][n]*vm3o1n[i]);
  1517. El[i] = El[i] + encap*(cl_n_[l][n]*vm1o1n[i] - c_i*dl_n_[l][n]*vn1e1n[i]
  1518. + c_i*al_n_[l][n]*vn3e1n[i] - bl_n_[l][n]*vm3o1n[i]);
  1519. Hl[i] = Hl[i] + encap*(-dl_n_[l][n]*vm1e1n[i] - c_i*cl_n_[l][n]*vn1o1n[i]
  1520. + c_i*bl_n_[l][n]*vn3o1n[i] + al_n_[l][n]*vm3e1n[i]);
  1521. // printf("\n !Ei[%d]=%g,%g! ", i, Ei[i].real(), Ei[i].imag());
  1522. // if (n<3) printf("\n !El[%d]=%g,%g! ", i, El[i].real(), El[i].imag());
  1523. // //printf(" ===%d=== %g ", i,std::abs(cl_n_[l][n]*vm1o1n[i] - c_i*dl_n_[l][n]*vn1e1n[i]));
  1524. // if (n<3) printf(" ===%d=== %g ", i,std::abs(//-dl_n_[l][n]*vm1e1n[i]
  1525. // //- c_i*cl_n_[l][n]*
  1526. // vn1o1n[i]
  1527. // // + c_i*bl_n_[l][n]*vn3o1n[i]
  1528. // // + al_n_[l][n]*vm3e1n[i]
  1529. // ));
  1530. // if (n<3) printf(" --- Ei[%d]=%g! ", i,std::abs(encap*(vm1o1n[i] - c_i*vn1e1n[i])));
  1531. }
  1532. //if (n<3) printf(" bj=%g \n", std::abs(bj[n]));
  1533. } // end of for all n
  1534. // magnetic field
  1535. double hffact = 1.0/(cc_*mu_);
  1536. for (int i = 0; i < 3; i++) {
  1537. Hl[i] = hffact*Hl[i];
  1538. }
  1539. for (int i = 0; i < 3; i++) {
  1540. // electric field E [V m - 1] = EF*E0
  1541. E[i] = El[i];
  1542. H[i] = Hl[i];
  1543. printf("\n !El[%d]=%g,%g! ", i, El[i].real(), El[i].imag());
  1544. //printf(" E[%d]=%g",i,std::abs(El[i]));
  1545. }
  1546. } // end of void fieldExt(...)
  1547. // ********************************************************************** //
  1548. // ********************************************************************** //
  1549. // ********************************************************************** //
  1550. //**********************************************************************************//
  1551. // This function calculates complex electric and magnetic field in the surroundings //
  1552. // and inside (TODO) the particle. //
  1553. // //
  1554. // Input parameters: //
  1555. // L: Number of layers //
  1556. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1557. // x: Array containing the size parameters of the layers [0..L - 1] //
  1558. // m: Array containing the relative refractive indexes of the layers [0..L - 1] //
  1559. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1560. // calculations. Only use it if you know what you are doing, otherwise //
  1561. // set this parameter to 0 (zero) and the function will calculate it. //
  1562. // ncoord: Number of coordinate points //
  1563. // Coords: Array containing all coordinates where the complex electric and //
  1564. // magnetic fields will be calculated //
  1565. // //
  1566. // Output parameters: //
  1567. // E, H: Complex electric and magnetic field at the provided coordinates //
  1568. // //
  1569. // Return value: //
  1570. // Number of multipolar expansion terms used for the calculations //
  1571. //**********************************************************************************//
  1572. void MultiLayerMie::RunFieldCalculations() {
  1573. // Calculate scattering coefficients an_ and bn_
  1574. RunMieCalculations();
  1575. //nmax_=10;
  1576. IntScattCoeffs();
  1577. std::vector<double> Pi(nmax_), Tau(nmax_);
  1578. long total_points = coords_sp_[0].size();
  1579. E_field_.resize(total_points);
  1580. H_field_.resize(total_points);
  1581. for (auto& f:E_field_) f.resize(3);
  1582. for (auto& f:H_field_) f.resize(3);
  1583. for (int point = 0; point < total_points; ++point) {
  1584. const double& Xp = coords_sp_[0][point];
  1585. const double& Yp = coords_sp_[1][point];
  1586. const double& Zp = coords_sp_[2][point];
  1587. printf("X=%g, Y=%g, Z=%g\n", Xp, Yp, Zp);
  1588. // Convert to spherical coordinates
  1589. double Rho, Phi, Theta;
  1590. Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1591. // printf("Rho=%g\n", Rho);
  1592. // Avoid convergence problems due to Rho too small
  1593. if (Rho < 1e-10) Rho = 1e-10;
  1594. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1595. if (Rho == 0.0) Theta = 0.0;
  1596. else Theta = std::acos(Zp/Rho);
  1597. // printf("Theta=%g\n", Theta);
  1598. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to zero to avoid problems
  1599. if (Xp == 0.0 && Yp == 0.0) Phi = 0.0;
  1600. else Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
  1601. // printf("Phi=%g\n", Phi);
  1602. calcSinglePiTau(std::cos(Theta), Pi, Tau);
  1603. //*******************************************************//
  1604. // external scattering field = incident + scattered //
  1605. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1606. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1607. //*******************************************************//
  1608. // This array contains the fields in spherical coordinates
  1609. std::vector<std::complex<double> > Es(3), Hs(3);
  1610. const double outer_size = size_parameter_.back();
  1611. // Firstly the easiest case: the field outside the particle
  1612. printf("rho=%g, outer=%g ", Rho, outer_size);
  1613. if (Rho >= outer_size) {
  1614. fieldExt(Rho, Phi, Theta, Pi, Tau, Es, Hs);
  1615. printf("\nFin E ext: %g,%g,%g Rho=%g\n", std::abs(Es[0]), std::abs(Es[1]),std::abs(Es[2]), Rho);
  1616. } else {
  1617. fieldInt(Rho, Phi, Theta, Pi, Tau, Es, Hs);
  1618. printf("\nFin E int: %g,%g,%g Rho=%g\n", std::abs(Es[0]), std::abs(Es[1]),std::abs(Es[2]), Rho);
  1619. }
  1620. std::complex<double>& Ex = E_field_[point][0];
  1621. std::complex<double>& Ey = E_field_[point][1];
  1622. std::complex<double>& Ez = E_field_[point][2];
  1623. std::complex<double>& Hx = H_field_[point][0];
  1624. std::complex<double>& Hy = H_field_[point][1];
  1625. std::complex<double>& Hz = H_field_[point][2];
  1626. //Now, convert the fields back to cartesian coordinates
  1627. {
  1628. using std::sin;
  1629. using std::cos;
  1630. Ex = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1631. Ey = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1632. Ez = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1633. Hx = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1634. Hy = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1635. Hz = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1636. }
  1637. printf("Cart E: %g,%g,%g Rho=%g\n", std::abs(Ex), std::abs(Ey),std::abs(Ez),
  1638. Rho);
  1639. } // end of for all field coordinates
  1640. } // end of void MultiLayerMie::RunFieldCalculations()
  1641. } // end of namespace nmie