nmie.cc 69 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302
  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.SetLayersSize(x);
  77. ml_mie.SetLayersIndex(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.SetLayersSize(x);
  129. ml_mie.SetLayersIndex(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.SetLayersSize(x);
  282. ml_mie.SetLayersIndex(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::SetLayersSize(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::SetLayersIndex(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]) > 100000.0)
  554. throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  555. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  556. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  557. *std::exp(-2.0*z.imag()));
  558. D3[0] = std::complex<double>(0.0, 1.0);
  559. for (int n = 1; n <= nmax_; n++) {
  560. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  561. *(static_cast<double>(n)*zinv - D3[n - 1]);
  562. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  563. }
  564. }
  565. //**********************************************************************************//
  566. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  567. // complex argument (z). //
  568. // Equations (20a) - (21b) //
  569. // //
  570. // Input parameters: //
  571. // z: Complex argument to evaluate Psi and Zeta //
  572. // nmax: Maximum number of terms to calculate Psi and Zeta //
  573. // //
  574. // Output parameters: //
  575. // Psi, Zeta: Riccati-Bessel functions //
  576. //**********************************************************************************//
  577. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  578. std::vector<std::complex<double> >& Psi,
  579. std::vector<std::complex<double> >& Zeta) {
  580. std::complex<double> c_i(0.0, 1.0);
  581. std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
  582. // First, calculate the logarithmic derivatives
  583. calcD1D3(z, D1, D3);
  584. // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  585. Psi[0] = std::sin(z);
  586. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  587. for (int n = 1; n <= nmax_; n++) {
  588. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  589. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  590. }
  591. }
  592. //**********************************************************************************//
  593. // This function calculates Pi and Tau for a given value of cos(Theta). //
  594. // Equations (26a) - (26c) //
  595. // //
  596. // Input parameters: //
  597. // nmax_: Maximum number of terms to calculate Pi and Tau //
  598. // nTheta: Number of scattering angles //
  599. // Theta: Array containing all the scattering angles where the scattering //
  600. // amplitudes will be calculated //
  601. // //
  602. // Output parameters: //
  603. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  604. //**********************************************************************************//
  605. void MultiLayerMie::calcPiTau(const double& costheta,
  606. std::vector<double>& Pi, std::vector<double>& Tau) {
  607. int i;
  608. //****************************************************//
  609. // Equations (26a) - (26c) //
  610. //****************************************************//
  611. // Initialize Pi and Tau
  612. Pi[0] = 1.0; // n=1
  613. Tau[0] = costheta;
  614. // Calculate the actual values
  615. if (nmax_ > 1) {
  616. Pi[1] = 3*costheta*Pi[0]; //n=2
  617. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  618. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  619. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  620. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  621. }
  622. }
  623. } // end of MultiLayerMie::calcPiTau(...)
  624. //**********************************************************************************//
  625. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  626. // required to calculate the near-field parameters. //
  627. // //
  628. // Input parameters: //
  629. // Rho: Radial distance //
  630. // Phi: Azimuthal angle //
  631. // Theta: Polar angle //
  632. // rn: Either the spherical Ricatti-Bessel function of first or third kind //
  633. // Dn: Logarithmic derivative of rn //
  634. // Pi, Tau: Angular functions Pi and Tau //
  635. // n: Order of vector spherical harmonics //
  636. // //
  637. // Output parameters: //
  638. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  639. //**********************************************************************************//
  640. void MultiLayerMie::calcSpherHarm(const std::complex<double> Rho, const double Theta, const double Phi,
  641. const std::complex<double>& rn, const std::complex<double>& Dn,
  642. const double& Pi, const double& Tau, const double& n,
  643. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  644. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  645. // using eq 4.50 in BH
  646. std::complex<double> c_zero(0.0, 0.0);
  647. using std::sin;
  648. using std::cos;
  649. Mo1n[0] = c_zero;
  650. Mo1n[1] = cos(Phi)*Pi*rn/Rho;
  651. Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
  652. Me1n[0] = c_zero;
  653. Me1n[1] = -sin(Phi)*Pi*rn/Rho;
  654. Me1n[2] = -cos(Phi)*Tau*rn/Rho;
  655. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  656. No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
  657. No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
  658. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  659. Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
  660. Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
  661. } // end of MultiLayerMie::calcSpherHarm(...)
  662. //**********************************************************************************//
  663. // This function calculates the scattering coefficients required to calculate //
  664. // both the near- and far-field parameters. //
  665. // //
  666. // Input parameters: //
  667. // L: Number of layers //
  668. // pl: Index of PEC layer. If there is none just send -1 //
  669. // x: Array containing the size parameters of the layers [0..L-1] //
  670. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  671. // nmax: Maximum number of multipolar expansion terms to be used for the //
  672. // calculations. Only use it if you know what you are doing, otherwise //
  673. // set this parameter to -1 and the function will calculate it. //
  674. // //
  675. // Output parameters: //
  676. // an, bn: Complex scattering amplitudes //
  677. // //
  678. // Return value: //
  679. // Number of multipolar expansion terms used for the calculations //
  680. //**********************************************************************************//
  681. void MultiLayerMie::calcScattCoeffs() {
  682. isScaCoeffsCalc_ = false;
  683. const std::vector<double>& x = size_param_;
  684. const std::vector<std::complex<double> >& m = refractive_index_;
  685. const int& pl = PEC_layer_position_;
  686. const int L = refractive_index_.size();
  687. //************************************************************************//
  688. // Calculate the index of the first layer. It can be either 0 (default) //
  689. // or the index of the outermost PEC layer. In the latter case all layers //
  690. // below the PEC are discarded. //
  691. // ***********************************************************************//
  692. int fl = (pl > 0) ? pl : 0;
  693. if (nmax_preset_ <= 0) calcNmax(fl);
  694. else nmax_ = nmax_preset_;
  695. std::complex<double> z1, z2;
  696. //**************************************************************************//
  697. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  698. // means that index = layer number - 1 or index = n - 1. The only exception //
  699. // are the arrays for representing D1, D3 and Q because they need a value //
  700. // for the index 0 (zero), hence it is important to consider this shift //
  701. // between different arrays. The change was done to optimize memory usage. //
  702. //**************************************************************************//
  703. // Allocate memory to the arrays
  704. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  705. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  706. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  707. for (int l = 0; l < L; l++) {
  708. Q[l].resize(nmax_ + 1);
  709. Ha[l].resize(nmax_);
  710. Hb[l].resize(nmax_);
  711. }
  712. an_.resize(nmax_);
  713. bn_.resize(nmax_);
  714. PsiZeta_.resize(nmax_ + 1);
  715. std::vector<std::complex<double> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  716. //*************************************************//
  717. // Calculate D1 and D3 for z1 in the first layer //
  718. //*************************************************//
  719. if (fl == pl) { // PEC layer
  720. for (int n = 0; n <= nmax_; n++) {
  721. D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
  722. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  723. }
  724. } else { // Regular layer
  725. z1 = x[fl]* m[fl];
  726. // Calculate D1 and D3
  727. calcD1D3(z1, D1_mlxl, D3_mlxl);
  728. }
  729. //******************************************************************//
  730. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  731. //******************************************************************//
  732. for (int n = 0; n < nmax_; n++) {
  733. Ha[fl][n] = D1_mlxl[n + 1];
  734. Hb[fl][n] = D1_mlxl[n + 1];
  735. }
  736. //*****************************************************//
  737. // Iteration from the second layer to the last one (L) //
  738. //*****************************************************//
  739. std::complex<double> Temp, Num, Denom;
  740. std::complex<double> G1, G2;
  741. for (int l = fl + 1; l < L; l++) {
  742. //************************************************************//
  743. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  744. //************************************************************//
  745. z1 = x[l]*m[l];
  746. z2 = x[l - 1]*m[l];
  747. //Calculate D1 and D3 for z1
  748. calcD1D3(z1, D1_mlxl, D3_mlxl);
  749. //Calculate D1 and D3 for z2
  750. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  751. //*************************************************//
  752. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  753. //*************************************************//
  754. // Upward recurrence for Q - equations (19a) and (19b)
  755. Num = std::exp(-2.0*(z1.imag() - z2.imag()))
  756. *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
  757. Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
  758. Q[l][0] = Num/Denom;
  759. for (int n = 1; n <= nmax_; n++) {
  760. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  761. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  762. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  763. }
  764. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  765. for (int n = 1; n <= nmax_; n++) {
  766. //Ha
  767. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  768. G1 = -D1_mlxlM1[n];
  769. G2 = -D3_mlxlM1[n];
  770. } else {
  771. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  772. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  773. } // end of if PEC
  774. Temp = Q[l][n]*G1;
  775. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  776. Denom = G2 - Temp;
  777. Ha[l][n - 1] = Num/Denom;
  778. //Hb
  779. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  780. G1 = Hb[l - 1][n - 1];
  781. G2 = Hb[l - 1][n - 1];
  782. } else {
  783. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  784. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  785. } // end of if PEC
  786. Temp = Q[l][n]*G1;
  787. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  788. Denom = (G2- Temp);
  789. Hb[l][n - 1] = (Num/ Denom);
  790. } // end of for Ha and Hb terms
  791. } // end of for layers iteration
  792. //**************************************//
  793. //Calculate Psi and Zeta for XL //
  794. //**************************************//
  795. // Calculate PsiXL and ZetaXL
  796. calcPsiZeta(x[L - 1], PsiXL, ZetaXL);
  797. //*********************************************************************//
  798. // Finally, we calculate the scattering coefficients (an and bn) and //
  799. // the angular functions (Pi and Tau). Note that for these arrays the //
  800. // first layer is 0 (zero), in future versions all arrays will follow //
  801. // this convention to save memory. (13 Nov, 2014) //
  802. //*********************************************************************//
  803. for (int n = 0; n < nmax_; n++) {
  804. //********************************************************************//
  805. //Expressions for calculating an and bn coefficients are not valid if //
  806. //there is only one PEC layer (ie, for a simple PEC sphere). //
  807. //********************************************************************//
  808. if (pl < (L - 1)) {
  809. 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]);
  810. 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]);
  811. } else {
  812. 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]);
  813. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  814. }
  815. } // end of for an and bn terms
  816. isScaCoeffsCalc_ = true;
  817. } // end of MultiLayerMie::calcScattCoeffs()
  818. //**********************************************************************************//
  819. // This function calculates the actual scattering parameters and amplitudes //
  820. // //
  821. // Input parameters: //
  822. // L: Number of layers //
  823. // pl: Index of PEC layer. If there is none just send -1 //
  824. // x: Array containing the size parameters of the layers [0..L-1] //
  825. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  826. // nTheta: Number of scattering angles //
  827. // Theta: Array containing all the scattering angles where the scattering //
  828. // amplitudes will be calculated //
  829. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  830. // calculations. Only use it if you know what you are doing, otherwise //
  831. // set this parameter to -1 and the function will calculate it //
  832. // //
  833. // Output parameters: //
  834. // Qext: Efficiency factor for extinction //
  835. // Qsca: Efficiency factor for scattering //
  836. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  837. // Qbk: Efficiency factor for backscattering //
  838. // Qpr: Efficiency factor for the radiation pressure //
  839. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  840. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  841. // S1, S2: Complex scattering amplitudes //
  842. // //
  843. // Return value: //
  844. // Number of multipolar expansion terms used for the calculations //
  845. //**********************************************************************************//
  846. void MultiLayerMie::RunMieCalculation() {
  847. if (size_param_.size() != refractive_index_.size())
  848. throw std::invalid_argument("Each size parameter should have only one index!");
  849. if (size_param_.size() == 0)
  850. throw std::invalid_argument("Initialize model first!");
  851. const std::vector<double>& x = size_param_;
  852. isExpCoeffsCalc_ = false;
  853. isScaCoeffsCalc_ = false;
  854. isMieCalculated_ = false;
  855. // Calculate scattering coefficients
  856. calcScattCoeffs();
  857. if (!isScaCoeffsCalc_) // TODO seems to be unreachable
  858. throw std::invalid_argument("Calculation of scattering coefficients failed!");
  859. // Initialize the scattering parameters
  860. Qext_ = 0.0;
  861. Qsca_ = 0.0;
  862. Qabs_ = 0.0;
  863. Qbk_ = 0.0;
  864. Qpr_ = 0.0;
  865. asymmetry_factor_ = 0.0;
  866. albedo_ = 0.0;
  867. // Initialize the scattering amplitudes
  868. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  869. S1_.swap(tmp1);
  870. S2_ = S1_;
  871. std::vector<double> Pi(nmax_), Tau(nmax_);
  872. std::complex<double> Qbktmp(0.0, 0.0);
  873. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  874. // By using downward recurrence we avoid loss of precision due to float rounding errors
  875. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  876. // http://en.wikipedia.org/wiki/Loss_of_significance
  877. for (int i = nmax_ - 2; i >= 0; i--) {
  878. const int n = i + 1;
  879. // Equation (27)
  880. Qext_ += (n + n + 1.0)*(an_[i].real() + bn_[i].real());
  881. // Equation (28)
  882. Qsca_ += (n + n + 1.0)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  883. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  884. // Equation (29)
  885. Qpr_ += ((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  886. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  887. // Equation (33)
  888. Qbktmp += (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  889. // Calculate the scattering amplitudes (S1 and S2) //
  890. // Equations (25a) - (25b) //
  891. for (unsigned int t = 0; t < theta_.size(); t++) {
  892. calcPiTau(std::cos(theta_[t]), Pi, Tau);
  893. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
  894. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
  895. }
  896. }
  897. double x2 = pow2(x.back());
  898. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  899. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  900. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  901. Qabs_ = Qext_ - Qsca_; // Equation (30)
  902. albedo_ = Qsca_/Qext_; // Equation (31)
  903. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  904. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  905. isMieCalculated_ = true;
  906. }
  907. //**********************************************************************************//
  908. // This function calculates the expansion coefficients inside the particle, //
  909. // required to calculate the near-field parameters. //
  910. // //
  911. // Input parameters: //
  912. // L: Number of layers //
  913. // pl: Index of PEC layer. If there is none just send -1 //
  914. // x: Array containing the size parameters of the layers [0..L-1] //
  915. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  916. // nmax: Maximum number of multipolar expansion terms to be used for the //
  917. // calculations. Only use it if you know what you are doing, otherwise //
  918. // set this parameter to -1 and the function will calculate it. //
  919. // //
  920. // Output parameters: //
  921. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  922. // //
  923. // Return value: //
  924. // Number of multipolar expansion terms used for the calculations //
  925. //**********************************************************************************//
  926. void MultiLayerMie::calcExpanCoeffs() {
  927. if (!isScaCoeffsCalc_)
  928. throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
  929. isExpCoeffsCalc_ = false;
  930. std::complex<double> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  931. const int L = refractive_index_.size();
  932. aln_.resize(L + 1);
  933. bln_.resize(L + 1);
  934. cln_.resize(L + 1);
  935. dln_.resize(L + 1);
  936. for (int l = 0; l <= L; l++) {
  937. aln_[l].resize(nmax_);
  938. bln_[l].resize(nmax_);
  939. cln_[l].resize(nmax_);
  940. dln_[l].resize(nmax_);
  941. }
  942. // Yang, paragraph under eq. A3
  943. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  944. for (int n = 0; n < nmax_; n++) {
  945. aln_[L][n] = an_[n];
  946. bln_[L][n] = bn_[n];
  947. cln_[L][n] = c_one;
  948. dln_[L][n] = c_one;
  949. }
  950. std::vector<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  951. std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  952. std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
  953. auto& m = refractive_index_;
  954. std::vector< std::complex<double> > m1(L);
  955. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  956. m1[L - 1] = std::complex<double> (1.0, 0.0);
  957. std::complex<double> z, z1;
  958. for (int l = L - 1; l >= 0; l--) {
  959. z = size_param_[l]*m[l];
  960. z1 = size_param_[l]*m1[l];
  961. calcD1D3(z, D1z, D3z);
  962. calcD1D3(z1, D1z1, D3z1);
  963. calcPsiZeta(z, Psiz, Zetaz);
  964. calcPsiZeta(z1, Psiz1, Zetaz1);
  965. for (int n = 0; n < nmax_; n++) {
  966. int n1 = n + 1;
  967. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  968. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  969. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  970. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  971. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  972. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  973. // aln
  974. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  975. // bln
  976. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  977. // cln
  978. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  979. // dln
  980. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  981. } // end of all n
  982. } // end of all l
  983. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  984. for (int n = 0; n < nmax_; ++n) {
  985. if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  986. else {
  987. //throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  988. printf("Warning: Potentially unstable calculation of aln (aln[0][%i] = %g, %gi)\n", n, aln_[0][n].real(), aln_[0][n].imag());
  989. aln_[0][n] = 0.0;
  990. }
  991. if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  992. else {
  993. //throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  994. printf("Warning: Potentially unstable calculation of bln (bln[0][%i] = %g, %gi)\n", n, bln_[0][n].real(), bln_[0][n].imag());
  995. bln_[0][n] = 0.0;
  996. }
  997. }
  998. isExpCoeffsCalc_ = true;
  999. } // end of void MultiLayerMie::calcExpanCoeffs()
  1000. //**********************************************************************************//
  1001. // This function calculates the electric (E) and magnetic (H) fields inside and //
  1002. // around the particle. //
  1003. // //
  1004. // Input parameters (coordinates of the point): //
  1005. // Rho: Radial distance //
  1006. // Phi: Azimuthal angle //
  1007. // Theta: Polar angle //
  1008. // //
  1009. // Output parameters: //
  1010. // E, H: Complex electric and magnetic fields //
  1011. //**********************************************************************************//
  1012. void MultiLayerMie::calcField(const double Rho, const double Theta, const double Phi,
  1013. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1014. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1015. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1016. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1017. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  1018. std::vector<std::complex<double> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  1019. std::vector<double> Pi(nmax_), Tau(nmax_);
  1020. int l = 0; // Layer number
  1021. std::complex<double> ml;
  1022. // Initialize E and H
  1023. for (int i = 0; i < 3; i++) {
  1024. E[i] = c_zero;
  1025. H[i] = c_zero;
  1026. }
  1027. if (Rho > size_param_.back()) {
  1028. l = size_param_.size();
  1029. ml = c_one;
  1030. } else {
  1031. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  1032. if (Rho <= size_param_[i]) {
  1033. l = i;
  1034. }
  1035. }
  1036. ml = refractive_index_[l];
  1037. }
  1038. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  1039. calcD1D3(Rho*ml, D1n, D3n);
  1040. // Calculate Ricatti-Bessel functions
  1041. calcPsiZeta(Rho*ml, Psi, Zeta);
  1042. // Calculate angular functions Pi and Tau
  1043. calcPiTau(std::cos(Theta), Pi, Tau);
  1044. for (int n = nmax_ - 2; n >= 0; n--) {
  1045. int n1 = n + 1;
  1046. double rn = static_cast<double>(n1);
  1047. // using BH 4.12 and 4.50
  1048. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1049. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1050. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1051. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1052. for (int i = 0; i < 3; i++) {
  1053. // electric field E [V m - 1] = EF*E0
  1054. E[i] += En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  1055. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  1056. H[i] += En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  1057. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  1058. }
  1059. } // end of for all n
  1060. // magnetic field
  1061. std::complex<double> hffact = ml/(cc_*mu_);
  1062. for (int i = 0; i < 3; i++) {
  1063. H[i] = hffact*H[i];
  1064. }
  1065. } // end of MultiLayerMie::calcField(...)
  1066. //**********************************************************************************//
  1067. // This function calculates complex electric and magnetic field in the surroundings //
  1068. // and inside the particle. //
  1069. // //
  1070. // Input parameters: //
  1071. // L: Number of layers //
  1072. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1073. // x: Array containing the size parameters of the layers [0..L-1] //
  1074. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1075. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1076. // calculations. Only use it if you know what you are doing, otherwise //
  1077. // set this parameter to 0 (zero) and the function will calculate it. //
  1078. // ncoord: Number of coordinate points //
  1079. // Coords: Array containing all coordinates where the complex electric and //
  1080. // magnetic fields will be calculated //
  1081. // //
  1082. // Output parameters: //
  1083. // E, H: Complex electric and magnetic field at the provided coordinates //
  1084. // //
  1085. // Return value: //
  1086. // Number of multipolar expansion terms used for the calculations //
  1087. //**********************************************************************************//
  1088. void MultiLayerMie::RunFieldCalculation() {
  1089. double Rho, Theta, Phi;
  1090. // Calculate scattering coefficients an_ and bn_
  1091. calcScattCoeffs();
  1092. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  1093. calcExpanCoeffs();
  1094. long total_points = coords_[0].size();
  1095. E_.resize(total_points);
  1096. H_.resize(total_points);
  1097. for (auto& f : E_) f.resize(3);
  1098. for (auto& f : H_) f.resize(3);
  1099. for (int point = 0; point < total_points; point++) {
  1100. const double& Xp = coords_[0][point];
  1101. const double& Yp = coords_[1][point];
  1102. const double& Zp = coords_[2][point];
  1103. // Convert to spherical coordinates
  1104. Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1105. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1106. Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1107. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1108. if (Xp == 0.0)
  1109. Phi = (Yp != 0.0) ? std::asin(Yp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1110. else
  1111. Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
  1112. // Avoid convergence problems due to Rho too small
  1113. if (Rho < 1e-5) Rho = 1e-5;
  1114. //*******************************************************//
  1115. // external scattering field = incident + scattered //
  1116. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1117. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1118. //*******************************************************//
  1119. // This array contains the fields in spherical coordinates
  1120. std::vector<std::complex<double> > Es(3), Hs(3);
  1121. // Do the actual calculation of electric and magnetic field
  1122. calcField(Rho, Theta, Phi, Es, Hs);
  1123. { //Now, convert the fields back to cartesian coordinates
  1124. using std::sin;
  1125. using std::cos;
  1126. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1127. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1128. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1129. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1130. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1131. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1132. }
  1133. } // end of for all field coordinates
  1134. } // end of MultiLayerMie::RunFieldCalculation()
  1135. } // end of namespace nmie