nmie.cc 69 KB

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