nmie.cc 66 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246
  1. //**********************************************************************************//
  2. // Copyright (C) 2009-2015 Ovidio Pena <ovidio@bytesfall.com> //
  3. // Copyright (C) 2013-2015 Konstantin Ladutenko <kostyfisik@gmail.com> //
  4. // //
  5. // This file is part of scattnlay //
  6. // //
  7. // This program 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. // This program 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. // The only additional remark is that we expect that all publications //
  18. // describing work using this software, or all commercial products //
  19. // using it, cite the following reference: //
  20. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  21. // a multilayered sphere," Computer Physics Communications, //
  22. // vol. 180, Nov. 2009, pp. 2348-2354. //
  23. // //
  24. // You should have received a copy of the GNU General Public License //
  25. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  26. //**********************************************************************************//
  27. //**********************************************************************************//
  28. // This class implements the algorithm for a multilayered sphere described by: //
  29. // [1] W. Yang, "Improved recursive algorithm for light scattering by a //
  30. // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. //
  31. // //
  32. // You can find the description of all the used equations in: //
  33. // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  34. // a multilayered sphere," Computer Physics Communications, //
  35. // vol. 180, Nov. 2009, pp. 2348-2354. //
  36. // //
  37. // Hereinafter all equations numbers refer to [2] //
  38. //**********************************************************************************//
  39. #include "bessel.h"
  40. #include "nmie.h"
  41. #include <array>
  42. #include <algorithm>
  43. #include <cstdio>
  44. #include <cstdlib>
  45. #include <stdexcept>
  46. #include <vector>
  47. namespace nmie {
  48. //helpers
  49. template<class T> inline T pow2(const T value) {return value*value;}
  50. int round(double x) {
  51. return x >= 0 ? (int)(x + 0.5):(int)(x - 0.5);
  52. }
  53. //**********************************************************************************//
  54. // This function emulates a C call to calculate the actual scattering parameters //
  55. // and amplitudes. //
  56. // //
  57. // Input parameters: //
  58. // L: Number of layers //
  59. // pl: Index of PEC layer. If there is none just send -1 //
  60. // x: Array containing the size parameters of the layers [0..L-1] //
  61. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  62. // nTheta: Number of scattering angles //
  63. // Theta: Array containing all the scattering angles where the scattering //
  64. // amplitudes will be calculated //
  65. // nmax: Maximum number of multipolar expansion terms to be used for the //
  66. // calculations. Only use it if you know what you are doing, otherwise //
  67. // set this parameter to -1 and the function will calculate it //
  68. // //
  69. // Output parameters: //
  70. // Qext: Efficiency factor for extinction //
  71. // Qsca: Efficiency factor for scattering //
  72. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  73. // Qbk: Efficiency factor for backscattering //
  74. // Qpr: Efficiency factor for the radiation pressure //
  75. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  76. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  77. // S1, S2: Complex scattering amplitudes //
  78. // //
  79. // Return value: //
  80. // Number of multipolar expansion terms used for the calculations //
  81. //**********************************************************************************//
  82. int nMie(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, const int nmax, 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) {
  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 (Theta.size() != nTheta)
  86. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  87. try {
  88. MultiLayerMie multi_layer_mie;
  89. multi_layer_mie.SetLayersSize(x);
  90. multi_layer_mie.SetLayersIndex(m);
  91. multi_layer_mie.SetAngles(Theta);
  92. multi_layer_mie.SetPECLayer(pl);
  93. multi_layer_mie.SetMaxTerms(nmax);
  94. multi_layer_mie.RunMieCalculation();
  95. *Qext = multi_layer_mie.GetQext();
  96. *Qsca = multi_layer_mie.GetQsca();
  97. *Qabs = multi_layer_mie.GetQabs();
  98. *Qbk = multi_layer_mie.GetQbk();
  99. *Qpr = multi_layer_mie.GetQpr();
  100. *g = multi_layer_mie.GetAsymmetryFactor();
  101. *Albedo = multi_layer_mie.GetAlbedo();
  102. S1 = multi_layer_mie.GetS1();
  103. S2 = multi_layer_mie.GetS2();
  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. // This function is just a wrapper to call the full 'nMie' function with fewer //
  114. // parameters, it is here mainly for compatibility with older versions of the //
  115. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  116. // any limit for the maximum number of terms. //
  117. // //
  118. // Input parameters: //
  119. // L: Number of layers //
  120. // x: Array containing the size parameters of the layers [0..L-1] //
  121. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  122. // nTheta: Number of scattering angles //
  123. // Theta: Array containing all the scattering angles where the scattering //
  124. // amplitudes will be calculated //
  125. // //
  126. // Output parameters: //
  127. // Qext: Efficiency factor for extinction //
  128. // Qsca: Efficiency factor for scattering //
  129. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  130. // Qbk: Efficiency factor for backscattering //
  131. // Qpr: Efficiency factor for the radiation pressure //
  132. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  133. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  134. // S1, S2: Complex scattering amplitudes //
  135. // //
  136. // Return value: //
  137. // Number of multipolar expansion terms used for the calculations //
  138. //**********************************************************************************//
  139. int nMie(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  140. return nmie::nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  141. }
  142. //**********************************************************************************//
  143. // This function is just a wrapper to call the full 'nMie' function with fewer //
  144. // parameters, it is useful if you want to include a PEC layer but not a limit //
  145. // for the maximum number of terms. //
  146. // //
  147. // Input parameters: //
  148. // L: Number of layers //
  149. // pl: Index of PEC layer. If there is none just send -1 //
  150. // x: Array containing the size parameters of the layers [0..L-1] //
  151. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  152. // nTheta: Number of scattering angles //
  153. // Theta: Array containing all the scattering angles where the scattering //
  154. // amplitudes will be calculated //
  155. // //
  156. // Output parameters: //
  157. // Qext: Efficiency factor for extinction //
  158. // Qsca: Efficiency factor for scattering //
  159. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  160. // Qbk: Efficiency factor for backscattering //
  161. // Qpr: Efficiency factor for the radiation pressure //
  162. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  163. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  164. // S1, S2: Complex scattering amplitudes //
  165. // //
  166. // Return value: //
  167. // Number of multipolar expansion terms used for the calculations //
  168. //**********************************************************************************//
  169. int nMie(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  170. return nmie::nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  171. }
  172. //**********************************************************************************//
  173. // This function is just a wrapper to call the full 'nMie' function with fewer //
  174. // parameters, it is useful if you want to include a limit for the maximum number //
  175. // of terms but not a PEC layer. //
  176. // //
  177. // Input parameters: //
  178. // L: Number of layers //
  179. // x: Array containing the size parameters of the layers [0..L-1] //
  180. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  181. // nTheta: Number of scattering angles //
  182. // Theta: Array containing all the scattering angles where the scattering //
  183. // amplitudes will be calculated //
  184. // nmax: Maximum number of multipolar expansion terms to be used for the //
  185. // calculations. Only use it if you know what you are doing, otherwise //
  186. // set this parameter to -1 and the function will calculate it //
  187. // //
  188. // Output parameters: //
  189. // Qext: Efficiency factor for extinction //
  190. // Qsca: Efficiency factor for scattering //
  191. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  192. // Qbk: Efficiency factor for backscattering //
  193. // Qpr: Efficiency factor for the radiation pressure //
  194. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  195. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  196. // S1, S2: Complex scattering amplitudes //
  197. // //
  198. // Return value: //
  199. // Number of multipolar expansion terms used for the calculations //
  200. //**********************************************************************************//
  201. int nMie(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned int nTheta, std::vector<double>& Theta, const int nmax, 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) {
  202. return nmie::nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  203. }
  204. //**********************************************************************************//
  205. // This function emulates a C call to calculate complex electric and magnetic field //
  206. // in the surroundings and inside (TODO) the particle. //
  207. // //
  208. // Input parameters: //
  209. // L: Number of layers //
  210. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  211. // x: Array containing the size parameters of the layers [0..L-1] //
  212. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  213. // nmax: Maximum number of multipolar expansion terms to be used for the //
  214. // calculations. Only use it if you know what you are doing, otherwise //
  215. // set this parameter to 0 (zero) and the function will calculate it. //
  216. // ncoord: Number of coordinate points //
  217. // Coords: Array containing all coordinates where the complex electric and //
  218. // magnetic fields will be calculated //
  219. // //
  220. // Output parameters: //
  221. // E, H: Complex electric and magnetic field at the provided coordinates //
  222. // //
  223. // Return value: //
  224. // Number of multipolar expansion terms used for the calculations //
  225. //**********************************************************************************//
  226. int nField(const unsigned int L, const int pl, const std::vector<double>& x, const std::vector<std::complex<double> >& m, const int nmax, const unsigned 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) {
  227. if (x.size() != L || m.size() != L)
  228. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  229. if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
  230. || E.size() != ncoord || H.size() != ncoord)
  231. throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
  232. for (auto f:E)
  233. if (f.size() != 3)
  234. throw std::invalid_argument("Field E is not 3D!");
  235. for (auto f:H)
  236. if (f.size() != 3)
  237. throw std::invalid_argument("Field H is not 3D!");
  238. try {
  239. MultiLayerMie multi_layer_mie;
  240. //multi_layer_mie.SetPECLayer(pl); // TODO add PEC layer to field plotting
  241. multi_layer_mie.SetLayersSize(x);
  242. multi_layer_mie.SetLayersIndex(m);
  243. multi_layer_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
  244. multi_layer_mie.RunFieldCalculation();
  245. E = multi_layer_mie.GetFieldE();
  246. H = multi_layer_mie.GetFieldH();
  247. } catch(const std::invalid_argument& ia) {
  248. // Will catch if multi_layer_mie fails or other errors.
  249. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  250. throw std::invalid_argument(ia);
  251. return - 1;
  252. }
  253. return 0;
  254. }
  255. // ********************************************************************** //
  256. // Returns previously calculated Qext //
  257. // ********************************************************************** //
  258. double MultiLayerMie::GetQext() {
  259. if (!isMieCalculated_)
  260. throw std::invalid_argument("You should run calculations before result request!");
  261. return Qext_;
  262. }
  263. // ********************************************************************** //
  264. // Returns previously calculated Qabs //
  265. // ********************************************************************** //
  266. double MultiLayerMie::GetQabs() {
  267. if (!isMieCalculated_)
  268. throw std::invalid_argument("You should run calculations before result request!");
  269. return Qabs_;
  270. }
  271. // ********************************************************************** //
  272. // Returns previously calculated Qsca //
  273. // ********************************************************************** //
  274. double MultiLayerMie::GetQsca() {
  275. if (!isMieCalculated_)
  276. throw std::invalid_argument("You should run calculations before result request!");
  277. return Qsca_;
  278. }
  279. // ********************************************************************** //
  280. // Returns previously calculated Qbk //
  281. // ********************************************************************** //
  282. double MultiLayerMie::GetQbk() {
  283. if (!isMieCalculated_)
  284. throw std::invalid_argument("You should run calculations before result request!");
  285. return Qbk_;
  286. }
  287. // ********************************************************************** //
  288. // Returns previously calculated Qpr //
  289. // ********************************************************************** //
  290. double MultiLayerMie::GetQpr() {
  291. if (!isMieCalculated_)
  292. throw std::invalid_argument("You should run calculations before result request!");
  293. return Qpr_;
  294. }
  295. // ********************************************************************** //
  296. // Returns previously calculated assymetry factor //
  297. // ********************************************************************** //
  298. double MultiLayerMie::GetAsymmetryFactor() {
  299. if (!isMieCalculated_)
  300. throw std::invalid_argument("You should run calculations before result request!");
  301. return asymmetry_factor_;
  302. }
  303. // ********************************************************************** //
  304. // Returns previously calculated Albedo //
  305. // ********************************************************************** //
  306. double MultiLayerMie::GetAlbedo() {
  307. if (!isMieCalculated_)
  308. throw std::invalid_argument("You should run calculations before result request!");
  309. return albedo_;
  310. }
  311. // ********************************************************************** //
  312. // Returns previously calculated S1 //
  313. // ********************************************************************** //
  314. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  315. if (!isMieCalculated_)
  316. throw std::invalid_argument("You should run calculations before result request!");
  317. return S1_;
  318. }
  319. // ********************************************************************** //
  320. // Returns previously calculated S2 //
  321. // ********************************************************************** //
  322. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  323. if (!isMieCalculated_)
  324. throw std::invalid_argument("You should run calculations before result request!");
  325. return S2_;
  326. }
  327. // ********************************************************************** //
  328. // Modify scattering (theta) angles //
  329. // ********************************************************************** //
  330. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  331. isExpCoeffsCalc_ = false;
  332. isScaCoeffsCalc_ = false;
  333. isMieCalculated_ = false;
  334. theta_ = angles;
  335. }
  336. // ********************************************************************** //
  337. // Modify size of all layers //
  338. // ********************************************************************** //
  339. void MultiLayerMie::SetLayersSize(const std::vector<double>& layer_size) {
  340. isExpCoeffsCalc_ = false;
  341. isScaCoeffsCalc_ = false;
  342. isMieCalculated_ = false;
  343. size_param_.clear();
  344. double prev_layer_size = 0.0;
  345. for (auto curr_layer_size : layer_size) {
  346. if (curr_layer_size <= 0.0)
  347. throw std::invalid_argument("Size parameter should be positive!");
  348. if (prev_layer_size > curr_layer_size)
  349. throw std::invalid_argument
  350. ("Size parameter for next layer should be larger than the previous one!");
  351. prev_layer_size = curr_layer_size;
  352. size_param_.push_back(curr_layer_size);
  353. }
  354. }
  355. // ********************************************************************** //
  356. // Modify refractive index of all layers //
  357. // ********************************************************************** //
  358. void MultiLayerMie::SetLayersIndex(const std::vector< std::complex<double> >& index) {
  359. isExpCoeffsCalc_ = false;
  360. isScaCoeffsCalc_ = false;
  361. isMieCalculated_ = false;
  362. refractive_index_ = index;
  363. }
  364. // ********************************************************************** //
  365. // Modify coordinates for field calculation //
  366. // ********************************************************************** //
  367. void MultiLayerMie::SetFieldCoords(const std::vector< std::vector<double> >& coords) {
  368. if (coords.size() != 3)
  369. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  370. if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
  371. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  372. coords_ = coords;
  373. }
  374. // ********************************************************************** //
  375. // ********************************************************************** //
  376. // ********************************************************************** //
  377. void MultiLayerMie::SetPECLayer(int layer_position) {
  378. isExpCoeffsCalc_ = false;
  379. isScaCoeffsCalc_ = false;
  380. isMieCalculated_ = false;
  381. if (layer_position < 0 && layer_position != -1)
  382. throw std::invalid_argument("Error! Layers are numbered from 0!");
  383. PEC_layer_position_ = layer_position;
  384. }
  385. // ********************************************************************** //
  386. // Set maximun number of terms to be used //
  387. // ********************************************************************** //
  388. void MultiLayerMie::SetMaxTerms(int nmax) {
  389. isExpCoeffsCalc_ = false;
  390. isScaCoeffsCalc_ = false;
  391. isMieCalculated_ = false;
  392. nmax_preset_ = nmax;
  393. }
  394. // ********************************************************************** //
  395. // ********************************************************************** //
  396. // ********************************************************************** //
  397. double MultiLayerMie::GetSizeParameter() {
  398. if (size_param_.size() > 0)
  399. return size_param_.back();
  400. else
  401. return 0;
  402. }
  403. // ********************************************************************** //
  404. // Clear layer information //
  405. // ********************************************************************** //
  406. void MultiLayerMie::ClearLayers() {
  407. isExpCoeffsCalc_ = false;
  408. isScaCoeffsCalc_ = false;
  409. isMieCalculated_ = false;
  410. size_param_.clear();
  411. refractive_index_.clear();
  412. }
  413. // ********************************************************************** //
  414. // ********************************************************************** //
  415. // ********************************************************************** //
  416. // Computational core
  417. // ********************************************************************** //
  418. // ********************************************************************** //
  419. // ********************************************************************** //
  420. // ********************************************************************** //
  421. // Calculate calcNstop - equation (17) //
  422. // ********************************************************************** //
  423. void MultiLayerMie::calcNstop() {
  424. const double& xL = size_param_.back();
  425. if (xL <= 8) {
  426. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  427. } else if (xL <= 4200) {
  428. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  429. } else {
  430. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  431. }
  432. }
  433. // ********************************************************************** //
  434. // Maximum number of terms required for the calculation //
  435. // ********************************************************************** //
  436. void MultiLayerMie::calcNmax(unsigned int first_layer) {
  437. int ri, riM1;
  438. const std::vector<double>& x = size_param_;
  439. const std::vector<std::complex<double> >& m = refractive_index_;
  440. calcNstop(); // Set initial nmax_ value
  441. for (unsigned int i = first_layer; i < x.size(); i++) {
  442. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  443. ri = round(std::abs(x[i]*m[i]));
  444. else
  445. ri = 0;
  446. nmax_ = std::max(nmax_, ri);
  447. // first layer is pec, if pec is present
  448. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  449. riM1 = round(std::abs(x[i - 1]* m[i]));
  450. else
  451. riM1 = 0;
  452. nmax_ = std::max(nmax_, riM1);
  453. }
  454. nmax_ += 15; // Final nmax_ value
  455. }
  456. // ********************************************************************** //
  457. // Calculate an - equation (5) //
  458. // ********************************************************************** //
  459. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  460. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  461. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  462. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  463. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  464. return Num/Denom;
  465. }
  466. // ********************************************************************** //
  467. // Calculate bn - equation (6) //
  468. // ********************************************************************** //
  469. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  470. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  471. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  472. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  473. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  474. return Num/Denom;
  475. }
  476. // ********************************************************************** //
  477. // Calculates S1 - equation (25a) //
  478. // ********************************************************************** //
  479. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  480. double Pi, double Tau) {
  481. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  482. }
  483. // ********************************************************************** //
  484. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  485. // Pi and Tau) //
  486. // ********************************************************************** //
  487. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  488. double Pi, double Tau) {
  489. return calc_S1(n, an, bn, Tau, Pi);
  490. }
  491. //**********************************************************************************//
  492. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  493. // functions (D1 and D3) for a complex argument (z). //
  494. // Equations (16a), (16b) and (18a) - (18d) //
  495. // //
  496. // Input parameters: //
  497. // z: Complex argument to evaluate D1 and D3 //
  498. // nmax_: Maximum number of terms to calculate D1 and D3 //
  499. // //
  500. // Output parameters: //
  501. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  502. //**********************************************************************************//
  503. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  504. std::vector<std::complex<double> >& D1,
  505. std::vector<std::complex<double> >& D3) {
  506. // Downward recurrence for D1 - equations (16a) and (16b)
  507. D1[nmax_] = std::complex<double>(0.0, 0.0);
  508. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  509. for (int n = nmax_; n > 0; n--) {
  510. D1[n - 1] = static_cast<double>(n)*zinv - 1.0/(D1[n] + static_cast<double>(n)*zinv);
  511. }
  512. if (std::abs(D1[0]) > 100000.0)
  513. throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  514. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  515. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  516. *std::exp(-2.0*z.imag()));
  517. D3[0] = std::complex<double>(0.0, 1.0);
  518. for (int n = 1; n <= nmax_; n++) {
  519. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  520. *(static_cast<double>(n)*zinv - D3[n - 1]);
  521. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  522. }
  523. }
  524. //**********************************************************************************//
  525. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  526. // complex argument (z). //
  527. // Equations (20a) - (21b) //
  528. // //
  529. // Input parameters: //
  530. // z: Complex argument to evaluate Psi and Zeta //
  531. // nmax: Maximum number of terms to calculate Psi and Zeta //
  532. // //
  533. // Output parameters: //
  534. // Psi, Zeta: Riccati-Bessel functions //
  535. //**********************************************************************************//
  536. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  537. std::vector<std::complex<double> >& Psi,
  538. std::vector<std::complex<double> >& Zeta) {
  539. std::complex<double> c_i(0.0, 1.0);
  540. std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
  541. // First, calculate the logarithmic derivatives
  542. calcD1D3(z, D1, D3);
  543. // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  544. Psi[0] = std::sin(z);
  545. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  546. for (int n = 1; n <= nmax_; n++) {
  547. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  548. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  549. }
  550. }
  551. //**********************************************************************************//
  552. // This function calculates Pi and Tau for a given value of cos(Theta). //
  553. // Equations (26a) - (26c) //
  554. // //
  555. // Input parameters: //
  556. // nmax_: Maximum number of terms to calculate Pi and Tau //
  557. // nTheta: Number of scattering angles //
  558. // Theta: Array containing all the scattering angles where the scattering //
  559. // amplitudes will be calculated //
  560. // //
  561. // Output parameters: //
  562. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  563. //**********************************************************************************//
  564. void MultiLayerMie::calcPiTau(const double& costheta,
  565. std::vector<double>& Pi, std::vector<double>& Tau) {
  566. int i;
  567. //****************************************************//
  568. // Equations (26a) - (26c) //
  569. //****************************************************//
  570. // Initialize Pi and Tau
  571. Pi[0] = 1.0; // n=1
  572. Tau[0] = costheta;
  573. // Calculate the actual values
  574. if (nmax_ > 1) {
  575. Pi[1] = 3*costheta*Pi[0]; //n=2
  576. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  577. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  578. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  579. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  580. }
  581. }
  582. } // end of MultiLayerMie::calcPiTau(...)
  583. //**********************************************************************************//
  584. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  585. // required to calculate the near-field parameters. //
  586. // //
  587. // Input parameters: //
  588. // Rho: Radial distance //
  589. // Phi: Azimuthal angle //
  590. // Theta: Polar angle //
  591. // rn: Either the spherical Ricatti-Bessel function of first or third kind //
  592. // Dn: Logarithmic derivative of rn //
  593. // Pi, Tau: Angular functions Pi and Tau //
  594. // n: Order of vector spherical harmonics //
  595. // //
  596. // Output parameters: //
  597. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  598. //**********************************************************************************//
  599. void MultiLayerMie::calcSpherHarm(const std::complex<double> Rho, const double Theta, const double Phi,
  600. const std::complex<double>& rn, const std::complex<double>& Dn,
  601. const double& Pi, const double& Tau, const double& n,
  602. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  603. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  604. // using eq 4.50 in BH
  605. std::complex<double> c_zero(0.0, 0.0);
  606. using std::sin;
  607. using std::cos;
  608. Mo1n[0] = c_zero;
  609. Mo1n[1] = cos(Phi)*Pi*rn/Rho;
  610. Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
  611. Me1n[0] = c_zero;
  612. Me1n[1] = -sin(Phi)*Pi*rn/Rho;
  613. Me1n[2] = -cos(Phi)*Tau*rn/Rho;
  614. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  615. No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
  616. No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
  617. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  618. Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
  619. Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
  620. } // end of MultiLayerMie::calcSpherHarm(...)
  621. //**********************************************************************************//
  622. // This function calculates the scattering coefficients required to calculate //
  623. // both the near- and far-field parameters. //
  624. // //
  625. // Input parameters: //
  626. // L: Number of layers //
  627. // pl: Index of PEC layer. If there is none just send -1 //
  628. // x: Array containing the size parameters of the layers [0..L-1] //
  629. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  630. // nmax: Maximum number of multipolar expansion terms to be used for the //
  631. // calculations. Only use it if you know what you are doing, otherwise //
  632. // set this parameter to -1 and the function will calculate it. //
  633. // //
  634. // Output parameters: //
  635. // an, bn: Complex scattering amplitudes //
  636. // //
  637. // Return value: //
  638. // Number of multipolar expansion terms used for the calculations //
  639. //**********************************************************************************//
  640. void MultiLayerMie::ScattCoeffs() {
  641. isScaCoeffsCalc_ = false;
  642. const std::vector<double>& x = size_param_;
  643. const std::vector<std::complex<double> >& m = refractive_index_;
  644. const int& pl = PEC_layer_position_;
  645. const int L = refractive_index_.size();
  646. //************************************************************************//
  647. // Calculate the index of the first layer. It can be either 0 (default) //
  648. // or the index of the outermost PEC layer. In the latter case all layers //
  649. // below the PEC are discarded. //
  650. // ***********************************************************************//
  651. int fl = (pl > 0) ? pl : 0;
  652. if (nmax_preset_ <= 0) calcNmax(fl);
  653. else nmax_ = nmax_preset_;
  654. std::complex<double> z1, z2;
  655. //**************************************************************************//
  656. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  657. // means that index = layer number - 1 or index = n - 1. The only exception //
  658. // are the arrays for representing D1, D3 and Q because they need a value //
  659. // for the index 0 (zero), hence it is important to consider this shift //
  660. // between different arrays. The change was done to optimize memory usage. //
  661. //**************************************************************************//
  662. // Allocate memory to the arrays
  663. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  664. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  665. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  666. for (int l = 0; l < L; l++) {
  667. Q[l].resize(nmax_ + 1);
  668. Ha[l].resize(nmax_);
  669. Hb[l].resize(nmax_);
  670. }
  671. an_.resize(nmax_);
  672. bn_.resize(nmax_);
  673. PsiZeta_.resize(nmax_ + 1);
  674. std::vector<std::complex<double> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  675. //*************************************************//
  676. // Calculate D1 and D3 for z1 in the first layer //
  677. //*************************************************//
  678. if (fl == pl) { // PEC layer
  679. for (int n = 0; n <= nmax_; n++) {
  680. D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
  681. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  682. }
  683. } else { // Regular layer
  684. z1 = x[fl]* m[fl];
  685. // Calculate D1 and D3
  686. calcD1D3(z1, D1_mlxl, D3_mlxl);
  687. }
  688. //******************************************************************//
  689. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  690. //******************************************************************//
  691. for (int n = 0; n < nmax_; n++) {
  692. Ha[fl][n] = D1_mlxl[n + 1];
  693. Hb[fl][n] = D1_mlxl[n + 1];
  694. }
  695. //*****************************************************//
  696. // Iteration from the second layer to the last one (L) //
  697. //*****************************************************//
  698. std::complex<double> Temp, Num, Denom;
  699. std::complex<double> G1, G2;
  700. for (int l = fl + 1; l < L; l++) {
  701. //************************************************************//
  702. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  703. //************************************************************//
  704. z1 = x[l]*m[l];
  705. z2 = x[l - 1]*m[l];
  706. //Calculate D1 and D3 for z1
  707. calcD1D3(z1, D1_mlxl, D3_mlxl);
  708. //Calculate D1 and D3 for z2
  709. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  710. //*************************************************//
  711. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  712. //*************************************************//
  713. // Upward recurrence for Q - equations (19a) and (19b)
  714. Num = std::exp(-2.0*(z1.imag() - z2.imag()))
  715. *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
  716. Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
  717. Q[l][0] = Num/Denom;
  718. for (int n = 1; n <= nmax_; n++) {
  719. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  720. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  721. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  722. }
  723. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  724. for (int n = 1; n <= nmax_; n++) {
  725. //Ha
  726. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  727. G1 = -D1_mlxlM1[n];
  728. G2 = -D3_mlxlM1[n];
  729. } else {
  730. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  731. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  732. } // end of if PEC
  733. Temp = Q[l][n]*G1;
  734. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  735. Denom = G2 - Temp;
  736. Ha[l][n - 1] = Num/Denom;
  737. //Hb
  738. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  739. G1 = Hb[l - 1][n - 1];
  740. G2 = Hb[l - 1][n - 1];
  741. } else {
  742. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  743. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  744. } // end of if PEC
  745. Temp = Q[l][n]*G1;
  746. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  747. Denom = (G2- Temp);
  748. Hb[l][n - 1] = (Num/ Denom);
  749. } // end of for Ha and Hb terms
  750. } // end of for layers iteration
  751. //**************************************//
  752. //Calculate Psi and Zeta for XL //
  753. //**************************************//
  754. // Calculate PsiXL and ZetaXL
  755. calcPsiZeta(x[L - 1], PsiXL, ZetaXL);
  756. //*********************************************************************//
  757. // Finally, we calculate the scattering coefficients (an and bn) and //
  758. // the angular functions (Pi and Tau). Note that for these arrays the //
  759. // first layer is 0 (zero), in future versions all arrays will follow //
  760. // this convention to save memory. (13 Nov, 2014) //
  761. //*********************************************************************//
  762. for (int n = 0; n < nmax_; n++) {
  763. //********************************************************************//
  764. //Expressions for calculating an and bn coefficients are not valid if //
  765. //there is only one PEC layer (ie, for a simple PEC sphere). //
  766. //********************************************************************//
  767. if (pl < (L - 1)) {
  768. 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]);
  769. 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]);
  770. } else {
  771. 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]);
  772. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  773. }
  774. } // end of for an and bn terms
  775. isScaCoeffsCalc_ = true;
  776. } // end of MultiLayerMie::ScattCoeffs(...)
  777. //**********************************************************************************//
  778. // This function calculates the actual scattering parameters and amplitudes //
  779. // //
  780. // Input parameters: //
  781. // L: Number of layers //
  782. // pl: Index of PEC layer. If there is none just send -1 //
  783. // x: Array containing the size parameters of the layers [0..L-1] //
  784. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  785. // nTheta: Number of scattering angles //
  786. // Theta: Array containing all the scattering angles where the scattering //
  787. // amplitudes will be calculated //
  788. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  789. // calculations. Only use it if you know what you are doing, otherwise //
  790. // set this parameter to -1 and the function will calculate it //
  791. // //
  792. // Output parameters: //
  793. // Qext: Efficiency factor for extinction //
  794. // Qsca: Efficiency factor for scattering //
  795. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  796. // Qbk: Efficiency factor for backscattering //
  797. // Qpr: Efficiency factor for the radiation pressure //
  798. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  799. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  800. // S1, S2: Complex scattering amplitudes //
  801. // //
  802. // Return value: //
  803. // Number of multipolar expansion terms used for the calculations //
  804. //**********************************************************************************//
  805. void MultiLayerMie::RunMieCalculation() {
  806. if (size_param_.size() != refractive_index_.size())
  807. throw std::invalid_argument("Each size parameter should have only one index!");
  808. if (size_param_.size() == 0)
  809. throw std::invalid_argument("Initialize model first!");
  810. const std::vector<double>& x = size_param_;
  811. isExpCoeffsCalc_ = false;
  812. isScaCoeffsCalc_ = false;
  813. isMieCalculated_ = false;
  814. // Calculate scattering coefficients
  815. ScattCoeffs();
  816. if (!isScaCoeffsCalc_) // TODO seems to be unreachable
  817. throw std::invalid_argument("Calculation of scattering coefficients failed!");
  818. // Initialize the scattering parameters
  819. Qext_ = 0.0;
  820. Qsca_ = 0.0;
  821. Qabs_ = 0.0;
  822. Qbk_ = 0.0;
  823. Qpr_ = 0.0;
  824. asymmetry_factor_ = 0.0;
  825. albedo_ = 0.0;
  826. // Initialize the scattering amplitudes
  827. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  828. S1_.swap(tmp1);
  829. S2_ = S1_;
  830. std::vector<double> Pi(nmax_), Tau(nmax_);
  831. std::complex<double> Qbktmp(0.0, 0.0);
  832. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  833. // By using downward recurrence we avoid loss of precision due to float rounding errors
  834. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  835. // http://en.wikipedia.org/wiki/Loss_of_significance
  836. for (int i = nmax_ - 2; i >= 0; i--) {
  837. const int n = i + 1;
  838. // Equation (27)
  839. Qext_ += (n + n + 1.0)*(an_[i].real() + bn_[i].real());
  840. // Equation (28)
  841. Qsca_ += (n + n + 1.0)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  842. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  843. // Equation (29)
  844. Qpr_ += ((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  845. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  846. // Equation (33)
  847. Qbktmp += (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  848. // Calculate the scattering amplitudes (S1 and S2) //
  849. // Equations (25a) - (25b) //
  850. for (unsigned int t = 0; t < theta_.size(); t++) {
  851. calcPiTau(std::cos(theta_[t]), Pi, Tau);
  852. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
  853. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
  854. }
  855. }
  856. double x2 = pow2(x.back());
  857. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  858. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  859. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  860. Qabs_ = Qext_ - Qsca_; // Equation (30)
  861. albedo_ = Qsca_/Qext_; // Equation (31)
  862. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  863. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  864. isMieCalculated_ = true;
  865. }
  866. //**********************************************************************************//
  867. // This function calculates the expansion coefficients inside the particle, //
  868. // required to calculate the near-field parameters. //
  869. // //
  870. // Input parameters: //
  871. // L: Number of layers //
  872. // pl: Index of PEC layer. If there is none just send -1 //
  873. // x: Array containing the size parameters of the layers [0..L-1] //
  874. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  875. // nmax: Maximum number of multipolar expansion terms to be used for the //
  876. // calculations. Only use it if you know what you are doing, otherwise //
  877. // set this parameter to -1 and the function will calculate it. //
  878. // //
  879. // Output parameters: //
  880. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  881. // //
  882. // Return value: //
  883. // Number of multipolar expansion terms used for the calculations //
  884. //**********************************************************************************//
  885. void MultiLayerMie::ExpanCoeffs() {
  886. if (!isScaCoeffsCalc_)
  887. throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
  888. isExpCoeffsCalc_ = false;
  889. std::complex<double> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  890. const int L = refractive_index_.size();
  891. aln_.resize(L + 1);
  892. bln_.resize(L + 1);
  893. cln_.resize(L + 1);
  894. dln_.resize(L + 1);
  895. for (int l = 0; l <= L; l++) {
  896. aln_[l].resize(nmax_);
  897. bln_[l].resize(nmax_);
  898. cln_[l].resize(nmax_);
  899. dln_[l].resize(nmax_);
  900. }
  901. // Yang, paragraph under eq. A3
  902. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  903. for (int n = 0; n < nmax_; n++) {
  904. aln_[L][n] = an_[n];
  905. bln_[L][n] = bn_[n];
  906. cln_[L][n] = c_one;
  907. dln_[L][n] = c_one;
  908. }
  909. std::vector<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  910. std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  911. std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
  912. auto& m = refractive_index_;
  913. std::vector< std::complex<double> > m1(L);
  914. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  915. m1[L - 1] = std::complex<double> (1.0, 0.0);
  916. std::complex<double> z, z1;
  917. for (int l = L - 1; l >= 0; l--) {
  918. z = size_param_[l]*m[l];
  919. z1 = size_param_[l]*m1[l];
  920. calcD1D3(z, D1z, D3z);
  921. calcD1D3(z1, D1z1, D3z1);
  922. calcPsiZeta(z, Psiz, Zetaz);
  923. calcPsiZeta(z1, Psiz1, Zetaz1);
  924. for (int n = 0; n < nmax_; n++) {
  925. int n1 = n + 1;
  926. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  927. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  928. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  929. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  930. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  931. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  932. // aln
  933. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  934. // bln
  935. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  936. // cln
  937. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  938. // dln
  939. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  940. } // end of all n
  941. } // end of all l
  942. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  943. for (int n = 0; n < nmax_; ++n) {
  944. if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  945. else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  946. if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  947. else throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  948. }
  949. isExpCoeffsCalc_ = true;
  950. } // end of void MultiLayerMie::ExpanCoeffs()
  951. //**********************************************************************************//
  952. // This function calculates the electric (E) and magnetic (H) fields inside and //
  953. // around the particle. //
  954. // //
  955. // Input parameters (coordinates of the point): //
  956. // Rho: Radial distance //
  957. // Phi: Azimuthal angle //
  958. // Theta: Polar angle //
  959. // //
  960. // Output parameters: //
  961. // E, H: Complex electric and magnetic fields //
  962. //**********************************************************************************//
  963. void MultiLayerMie::calcField(const double Rho, const double Theta, const double Phi,
  964. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  965. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  966. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  967. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  968. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  969. std::vector<std::complex<double> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  970. std::vector<double> Pi(nmax_), Tau(nmax_);
  971. int l = 0; // Layer number
  972. std::complex<double> ml;
  973. // Initialize E and H
  974. for (int i = 0; i < 3; i++) {
  975. E[i] = c_zero;
  976. H[i] = c_zero;
  977. }
  978. if (Rho > size_param_.back()) {
  979. l = size_param_.size();
  980. ml = c_one;
  981. } else {
  982. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  983. if (Rho <= size_param_[i]) {
  984. l = i;
  985. }
  986. }
  987. ml = refractive_index_[l];
  988. }
  989. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  990. calcD1D3(Rho*ml, D1n, D3n);
  991. // Calculate Ricatti-Bessel functions
  992. calcPsiZeta(Rho*ml, Psi, Zeta);
  993. // Calculate angular functions Pi and Tau
  994. calcPiTau(std::cos(Theta), Pi, Tau);
  995. for (int n = nmax_ - 2; n >= 0; n--) {
  996. int n1 = n + 1;
  997. double rn = static_cast<double>(n1);
  998. // using BH 4.12 and 4.50
  999. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1000. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1001. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1002. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1003. for (int i = 0; i < 3; i++) {
  1004. // electric field E [V m - 1] = EF*E0
  1005. E[i] += En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  1006. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  1007. H[i] += En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  1008. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  1009. }
  1010. } // end of for all n
  1011. // magnetic field
  1012. std::complex<double> hffact = ml/(cc_*mu_);
  1013. for (int i = 0; i < 3; i++) {
  1014. H[i] = hffact*H[i];
  1015. }
  1016. } // end of MultiLayerMie::calcField(...)
  1017. //**********************************************************************************//
  1018. // This function calculates complex electric and magnetic field in the surroundings //
  1019. // and inside the particle. //
  1020. // //
  1021. // Input parameters: //
  1022. // L: Number of layers //
  1023. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1024. // x: Array containing the size parameters of the layers [0..L-1] //
  1025. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1026. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1027. // calculations. Only use it if you know what you are doing, otherwise //
  1028. // set this parameter to 0 (zero) and the function will calculate it. //
  1029. // ncoord: Number of coordinate points //
  1030. // Coords: Array containing all coordinates where the complex electric and //
  1031. // magnetic fields will be calculated //
  1032. // //
  1033. // Output parameters: //
  1034. // E, H: Complex electric and magnetic field at the provided coordinates //
  1035. // //
  1036. // Return value: //
  1037. // Number of multipolar expansion terms used for the calculations //
  1038. //**********************************************************************************//
  1039. void MultiLayerMie::RunFieldCalculation() {
  1040. double Rho, Theta, Phi;
  1041. // Calculate scattering coefficients an_ and bn_
  1042. ScattCoeffs();
  1043. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  1044. ExpanCoeffs();
  1045. long total_points = coords_[0].size();
  1046. E_.resize(total_points);
  1047. H_.resize(total_points);
  1048. for (auto& f : E_) f.resize(3);
  1049. for (auto& f : H_) f.resize(3);
  1050. for (int point = 0; point < total_points; point++) {
  1051. const double& Xp = coords_[0][point];
  1052. const double& Yp = coords_[1][point];
  1053. const double& Zp = coords_[2][point];
  1054. // Convert to spherical coordinates
  1055. Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1056. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1057. Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1058. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1059. if (Xp == 0.0)
  1060. Phi = (Yp != 0.0) ? std::asin(Yp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1061. else
  1062. Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
  1063. // Avoid convergence problems due to Rho too small
  1064. if (Rho < 1e-5) Rho = 1e-5;
  1065. //*******************************************************//
  1066. // external scattering field = incident + scattered //
  1067. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1068. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1069. //*******************************************************//
  1070. // This array contains the fields in spherical coordinates
  1071. std::vector<std::complex<double> > Es(3), Hs(3);
  1072. // Do the actual calculation of electric and magnetic field
  1073. calcField(Rho, Theta, Phi, Es, Hs);
  1074. { //Now, convert the fields back to cartesian coordinates
  1075. using std::sin;
  1076. using std::cos;
  1077. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1078. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1079. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1080. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1081. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1082. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1083. }
  1084. } // end of for all field coordinates
  1085. } // end of MultiLayerMie::RunFieldCalculation()
  1086. } // end of namespace nmie