nmie.cc 82 KB

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