nmie.cc 69 KB

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