nmie.cc 74 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404
  1. //**********************************************************************************//
  2. // Copyright (C) 2009-2015 Ovidio Pena <ovidio@bytesfall.com> //
  3. // Copyright (C) 2013-2015 Konstantin Ladutenko <kostyfisik@gmail.com> //
  4. // //
  5. // This file is part of scattnlay //
  6. // //
  7. // This program is free software: you can redistribute it and/or modify //
  8. // it under the terms of the GNU General Public License as published by //
  9. // the Free Software Foundation, either version 3 of the License, or //
  10. // (at your option) any later version. //
  11. // //
  12. // This program is distributed in the hope that it will be useful, //
  13. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  14. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  15. // GNU General Public License for more details. //
  16. // //
  17. // The only additional remark is that we expect that all publications //
  18. // describing work using this software, or all commercial products //
  19. // using it, cite the following reference: //
  20. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  21. // a multilayered sphere," Computer Physics Communications, //
  22. // vol. 180, Nov. 2009, pp. 2348-2354. //
  23. // //
  24. // You should have received a copy of the GNU General Public License //
  25. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  26. //**********************************************************************************//
  27. //**********************************************************************************//
  28. // This class implements the algorithm for a multilayered sphere described by: //
  29. // [1] W. Yang, "Improved recursive algorithm for light scattering by a //
  30. // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. //
  31. // //
  32. // You can find the description of all the used equations in: //
  33. // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  34. // a multilayered sphere," Computer Physics Communications, //
  35. // vol. 180, Nov. 2009, pp. 2348-2354. //
  36. // //
  37. // Hereinafter all equations numbers refer to [2] //
  38. //**********************************************************************************//
  39. #include "nmie.h"
  40. #include <array>
  41. #include <algorithm>
  42. #include <cstdio>
  43. #include <cstdlib>
  44. #include <stdexcept>
  45. #include <vector>
  46. namespace nmie {
  47. //helpers
  48. template<class T> inline T pow2(const T value) {return value*value;}
  49. int round(double x) {
  50. return x >= 0 ? (int)(x + 0.5):(int)(x - 0.5);
  51. }
  52. //**********************************************************************************//
  53. // This function emulates a C call to calculate the actual scattering parameters //
  54. // and amplitudes. //
  55. // //
  56. // Input parameters: //
  57. // L: Number of layers //
  58. // pl: Index of PEC layer. If there is none just send -1 //
  59. // x: Array containing the size parameters of the layers [0..L-1] //
  60. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  61. // nTheta: Number of scattering angles //
  62. // Theta: Array containing all the scattering angles where the scattering //
  63. // amplitudes will be calculated //
  64. // nmax: Maximum number of multipolar expansion terms to be used for the //
  65. // calculations. Only use it if you know what you are doing, otherwise //
  66. // set this parameter to -1 and the function will calculate it //
  67. // //
  68. // Output parameters: //
  69. // Qext: Efficiency factor for extinction //
  70. // Qsca: Efficiency factor for scattering //
  71. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  72. // Qbk: Efficiency factor for backscattering //
  73. // Qpr: Efficiency factor for the radiation pressure //
  74. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  75. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  76. // S1, S2: Complex scattering amplitudes //
  77. // //
  78. // Return value: //
  79. // Number of multipolar expansion terms used for the calculations //
  80. //**********************************************************************************//
  81. int nMie(const 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) {
  82. if (x.size() != L || m.size() != L)
  83. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  84. if (Theta.size() != nTheta)
  85. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  86. try {
  87. MultiLayerMie multi_layer_mie;
  88. multi_layer_mie.SetLayersSize(x);
  89. multi_layer_mie.SetLayersIndex(m);
  90. multi_layer_mie.SetAngles(Theta);
  91. multi_layer_mie.SetPECLayer(pl);
  92. multi_layer_mie.SetMaxTerms(nmax);
  93. multi_layer_mie.RunMieCalculation();
  94. *Qext = multi_layer_mie.GetQext();
  95. *Qsca = multi_layer_mie.GetQsca();
  96. *Qabs = multi_layer_mie.GetQabs();
  97. *Qbk = multi_layer_mie.GetQbk();
  98. *Qpr = multi_layer_mie.GetQpr();
  99. *g = multi_layer_mie.GetAsymmetryFactor();
  100. *Albedo = multi_layer_mie.GetAlbedo();
  101. S1 = multi_layer_mie.GetS1();
  102. S2 = multi_layer_mie.GetS2();
  103. } catch(const std::invalid_argument& ia) {
  104. // Will catch if multi_layer_mie fails or other errors.
  105. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  106. throw std::invalid_argument(ia);
  107. return -1;
  108. }
  109. return 0;
  110. }
  111. //**********************************************************************************//
  112. // This function is just a wrapper to call the full 'nMie' function with fewer //
  113. // parameters, it is here mainly for compatibility with older versions of the //
  114. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  115. // any limit for the maximum number of terms. //
  116. // //
  117. // Input parameters: //
  118. // L: Number of layers //
  119. // x: Array containing the size parameters of the layers [0..L-1] //
  120. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  121. // nTheta: Number of scattering angles //
  122. // Theta: Array containing all the scattering angles where the scattering //
  123. // amplitudes will be calculated //
  124. // //
  125. // Output parameters: //
  126. // Qext: Efficiency factor for extinction //
  127. // Qsca: Efficiency factor for scattering //
  128. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  129. // Qbk: Efficiency factor for backscattering //
  130. // Qpr: Efficiency factor for the radiation pressure //
  131. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  132. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  133. // S1, S2: Complex scattering amplitudes //
  134. // //
  135. // Return value: //
  136. // Number of multipolar expansion terms used for the calculations //
  137. //**********************************************************************************//
  138. int nMie(const 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) {
  139. return nmie::nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  140. }
  141. //**********************************************************************************//
  142. // This function is just a wrapper to call the full 'nMie' function with fewer //
  143. // parameters, it is useful if you want to include a PEC layer but not a limit //
  144. // for the maximum number of terms. //
  145. // //
  146. // Input parameters: //
  147. // L: Number of layers //
  148. // pl: Index of PEC layer. If there is none just send -1 //
  149. // x: Array containing the size parameters of the layers [0..L-1] //
  150. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  151. // nTheta: Number of scattering angles //
  152. // Theta: Array containing all the scattering angles where the scattering //
  153. // amplitudes will be calculated //
  154. // //
  155. // Output parameters: //
  156. // Qext: Efficiency factor for extinction //
  157. // Qsca: Efficiency factor for scattering //
  158. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  159. // Qbk: Efficiency factor for backscattering //
  160. // Qpr: Efficiency factor for the radiation pressure //
  161. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  162. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  163. // S1, S2: Complex scattering amplitudes //
  164. // //
  165. // Return value: //
  166. // Number of multipolar expansion terms used for the calculations //
  167. //**********************************************************************************//
  168. int nMie(const 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) {
  169. return nmie::nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  170. }
  171. //**********************************************************************************//
  172. // This function is just a wrapper to call the full 'nMie' function with fewer //
  173. // parameters, it is useful if you want to include a limit for the maximum number //
  174. // of terms but not a PEC layer. //
  175. // //
  176. // Input parameters: //
  177. // L: Number of layers //
  178. // x: Array containing the size parameters of the layers [0..L-1] //
  179. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  180. // nTheta: Number of scattering angles //
  181. // Theta: Array containing all the scattering angles where the scattering //
  182. // amplitudes will be calculated //
  183. // nmax: Maximum number of multipolar expansion terms to be used for the //
  184. // calculations. Only use it if you know what you are doing, otherwise //
  185. // set this parameter to -1 and the function will calculate it //
  186. // //
  187. // Output parameters: //
  188. // Qext: Efficiency factor for extinction //
  189. // Qsca: Efficiency factor for scattering //
  190. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  191. // Qbk: Efficiency factor for backscattering //
  192. // Qpr: Efficiency factor for the radiation pressure //
  193. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  194. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  195. // S1, S2: Complex scattering amplitudes //
  196. // //
  197. // Return value: //
  198. // Number of multipolar expansion terms used for the calculations //
  199. //**********************************************************************************//
  200. int nMie(const 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) {
  201. return nmie::nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  202. }
  203. //**********************************************************************************//
  204. // This function emulates a C call to calculate complex electric and magnetic field //
  205. // in the surroundings and inside (TODO) the particle. //
  206. // //
  207. // Input parameters: //
  208. // L: Number of layers //
  209. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  210. // x: Array containing the size parameters of the layers [0..L-1] //
  211. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  212. // nmax: Maximum number of multipolar expansion terms to be used for the //
  213. // calculations. Only use it if you know what you are doing, otherwise //
  214. // set this parameter to 0 (zero) and the function will calculate it. //
  215. // ncoord: Number of coordinate points //
  216. // Coords: Array containing all coordinates where the complex electric and //
  217. // magnetic fields will be calculated //
  218. // //
  219. // Output parameters: //
  220. // E, H: Complex electric and magnetic field at the provided coordinates //
  221. // //
  222. // Return value: //
  223. // Number of multipolar expansion terms used for the calculations //
  224. //**********************************************************************************//
  225. int nField(const 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) {
  226. if (x.size() != L || m.size() != L)
  227. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  228. if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
  229. || E.size() != ncoord || H.size() != ncoord)
  230. throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
  231. for (auto f:E)
  232. if (f.size() != 3)
  233. throw std::invalid_argument("Field E is not 3D!");
  234. for (auto f:H)
  235. if (f.size() != 3)
  236. throw std::invalid_argument("Field H is not 3D!");
  237. try {
  238. MultiLayerMie multi_layer_mie;
  239. //multi_layer_mie.SetPECLayer(pl); // TODO add PEC layer to field plotting
  240. multi_layer_mie.SetLayersSize(x);
  241. multi_layer_mie.SetLayersIndex(m);
  242. multi_layer_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
  243. multi_layer_mie.RunFieldCalculation();
  244. E = multi_layer_mie.GetFieldE();
  245. H = multi_layer_mie.GetFieldH();
  246. } catch(const std::invalid_argument& ia) {
  247. // Will catch if multi_layer_mie fails or other errors.
  248. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  249. throw std::invalid_argument(ia);
  250. return - 1;
  251. }
  252. return 0;
  253. }
  254. // ********************************************************************** //
  255. // Returns previously calculated Qext //
  256. // ********************************************************************** //
  257. double MultiLayerMie::GetQext() {
  258. if (!isMieCalculated_)
  259. throw std::invalid_argument("You should run calculations before result request!");
  260. return Qext_;
  261. }
  262. // ********************************************************************** //
  263. // Returns previously calculated Qabs //
  264. // ********************************************************************** //
  265. double MultiLayerMie::GetQabs() {
  266. if (!isMieCalculated_)
  267. throw std::invalid_argument("You should run calculations before result request!");
  268. return Qabs_;
  269. }
  270. // ********************************************************************** //
  271. // Returns previously calculated Qsca //
  272. // ********************************************************************** //
  273. double MultiLayerMie::GetQsca() {
  274. if (!isMieCalculated_)
  275. throw std::invalid_argument("You should run calculations before result request!");
  276. return Qsca_;
  277. }
  278. // ********************************************************************** //
  279. // Returns previously calculated Qbk //
  280. // ********************************************************************** //
  281. double MultiLayerMie::GetQbk() {
  282. if (!isMieCalculated_)
  283. throw std::invalid_argument("You should run calculations before result request!");
  284. return Qbk_;
  285. }
  286. // ********************************************************************** //
  287. // Returns previously calculated Qpr //
  288. // ********************************************************************** //
  289. double MultiLayerMie::GetQpr() {
  290. if (!isMieCalculated_)
  291. throw std::invalid_argument("You should run calculations before result request!");
  292. return Qpr_;
  293. }
  294. // ********************************************************************** //
  295. // Returns previously calculated assymetry factor //
  296. // ********************************************************************** //
  297. double MultiLayerMie::GetAsymmetryFactor() {
  298. if (!isMieCalculated_)
  299. throw std::invalid_argument("You should run calculations before result request!");
  300. return asymmetry_factor_;
  301. }
  302. // ********************************************************************** //
  303. // Returns previously calculated Albedo //
  304. // ********************************************************************** //
  305. double MultiLayerMie::GetAlbedo() {
  306. if (!isMieCalculated_)
  307. throw std::invalid_argument("You should run calculations before result request!");
  308. return albedo_;
  309. }
  310. // ********************************************************************** //
  311. // Returns previously calculated S1 //
  312. // ********************************************************************** //
  313. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  314. if (!isMieCalculated_)
  315. throw std::invalid_argument("You should run calculations before result request!");
  316. return S1_;
  317. }
  318. // ********************************************************************** //
  319. // Returns previously calculated S2 //
  320. // ********************************************************************** //
  321. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  322. if (!isMieCalculated_)
  323. throw std::invalid_argument("You should run calculations before result request!");
  324. return S2_;
  325. }
  326. // ********************************************************************** //
  327. // Modify scattering (theta) angles //
  328. // ********************************************************************** //
  329. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  330. isIntCoeffsCalc_ = false;
  331. isExtCoeffsCalc_ = false;
  332. isMieCalculated_ = false;
  333. theta_ = angles;
  334. }
  335. // ********************************************************************** //
  336. // Modify size of all layers //
  337. // ********************************************************************** //
  338. void MultiLayerMie::SetLayersSize(const std::vector<double>& layer_size) {
  339. isIntCoeffsCalc_ = false;
  340. isExtCoeffsCalc_ = false;
  341. isMieCalculated_ = false;
  342. size_param_.clear();
  343. double prev_layer_size = 0.0;
  344. for (auto curr_layer_size : layer_size) {
  345. if (curr_layer_size <= 0.0)
  346. throw std::invalid_argument("Size parameter should be positive!");
  347. if (prev_layer_size > curr_layer_size)
  348. throw std::invalid_argument
  349. ("Size parameter for next layer should be larger than the previous one!");
  350. prev_layer_size = curr_layer_size;
  351. size_param_.push_back(curr_layer_size);
  352. }
  353. }
  354. // ********************************************************************** //
  355. // Modify refractive index of all layers //
  356. // ********************************************************************** //
  357. void MultiLayerMie::SetLayersIndex(const std::vector< std::complex<double> >& index) {
  358. isIntCoeffsCalc_ = false;
  359. isExtCoeffsCalc_ = false;
  360. isMieCalculated_ = false;
  361. refractive_index_ = index;
  362. }
  363. // ********************************************************************** //
  364. // Modify coordinates for field calculation //
  365. // ********************************************************************** //
  366. void MultiLayerMie::SetFieldCoords(const std::vector< std::vector<double> >& coords) {
  367. if (coords.size() != 3)
  368. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  369. if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
  370. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  371. coords_ = coords;
  372. }
  373. // ********************************************************************** //
  374. // ********************************************************************** //
  375. // ********************************************************************** //
  376. void MultiLayerMie::SetPECLayer(int layer_position) {
  377. isIntCoeffsCalc_ = false;
  378. isExtCoeffsCalc_ = false;
  379. isMieCalculated_ = false;
  380. if (layer_position < 0 && layer_position != -1)
  381. throw std::invalid_argument("Error! Layers are numbered from 0!");
  382. PEC_layer_position_ = layer_position;
  383. }
  384. // ********************************************************************** //
  385. // Set maximun number of terms to be used //
  386. // ********************************************************************** //
  387. void MultiLayerMie::SetMaxTerms(int nmax) {
  388. isIntCoeffsCalc_ = false;
  389. isExtCoeffsCalc_ = false;
  390. isMieCalculated_ = false;
  391. nmax_preset_ = nmax;
  392. }
  393. // ********************************************************************** //
  394. // ********************************************************************** //
  395. // ********************************************************************** //
  396. double MultiLayerMie::GetSizeParameter() {
  397. if (size_param_.size() > 0)
  398. return size_param_.back();
  399. else
  400. return 0;
  401. }
  402. // ********************************************************************** //
  403. // Clear layer information //
  404. // ********************************************************************** //
  405. void MultiLayerMie::ClearLayers() {
  406. isIntCoeffsCalc_ = false;
  407. isExtCoeffsCalc_ = false;
  408. isMieCalculated_ = false;
  409. size_param_.clear();
  410. refractive_index_.clear();
  411. }
  412. // ********************************************************************** //
  413. // ********************************************************************** //
  414. // ********************************************************************** //
  415. // Computational core
  416. // ********************************************************************** //
  417. // ********************************************************************** //
  418. // ********************************************************************** //
  419. // ********************************************************************** //
  420. // Calculate calcNstop - equation (17) //
  421. // ********************************************************************** //
  422. void MultiLayerMie::calcNstop() {
  423. const double& xL = size_param_.back();
  424. if (xL <= 8) {
  425. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  426. } else if (xL <= 4200) {
  427. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  428. } else {
  429. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  430. }
  431. }
  432. // ********************************************************************** //
  433. // Maximum number of terms required for the calculation //
  434. // ********************************************************************** //
  435. void MultiLayerMie::calcNmax(unsigned int first_layer) {
  436. int ri, riM1;
  437. const std::vector<double>& x = size_param_;
  438. const std::vector<std::complex<double> >& m = refractive_index_;
  439. calcNstop(); // Set initial nmax_ value
  440. for (unsigned int i = first_layer; i < x.size(); i++) {
  441. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  442. ri = round(std::abs(x[i]*m[i]));
  443. else
  444. ri = 0;
  445. nmax_ = std::max(nmax_, ri);
  446. // first layer is pec, if pec is present
  447. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  448. riM1 = round(std::abs(x[i - 1]* m[i]));
  449. else
  450. riM1 = 0;
  451. nmax_ = std::max(nmax_, riM1);
  452. }
  453. nmax_ += 15; // Final nmax_ value
  454. }
  455. // ********************************************************************** //
  456. // Calculate an - equation (5) //
  457. // ********************************************************************** //
  458. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  459. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  460. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  461. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  462. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  463. return Num/Denom;
  464. }
  465. // ********************************************************************** //
  466. // Calculate bn - equation (6) //
  467. // ********************************************************************** //
  468. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  469. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  470. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  471. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  472. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  473. return Num/Denom;
  474. }
  475. // ********************************************************************** //
  476. // Calculates S1 - equation (25a) //
  477. // ********************************************************************** //
  478. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  479. double Pi, double Tau) {
  480. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  481. }
  482. // ********************************************************************** //
  483. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  484. // Pi and Tau) //
  485. // ********************************************************************** //
  486. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  487. double Pi, double Tau) {
  488. return calc_S1(n, an, bn, Tau, Pi);
  489. }
  490. //**********************************************************************************//
  491. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  492. // functions (D1 and D3) for a complex argument (z). //
  493. // Equations (16a), (16b) and (18a) - (18d) //
  494. // //
  495. // Input parameters: //
  496. // z: Complex argument to evaluate D1 and D3 //
  497. // nmax_: Maximum number of terms to calculate D1 and D3 //
  498. // //
  499. // Output parameters: //
  500. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  501. //**********************************************************************************//
  502. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  503. std::vector<std::complex<double> >& D1,
  504. std::vector<std::complex<double> >& D3) {
  505. // Downward recurrence for D1 - equations (16a) and (16b)
  506. D1[nmax_] = std::complex<double>(0.0, 0.0);
  507. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  508. for (int n = nmax_; n > 0; n--) {
  509. D1[n - 1] = double(n)*zinv - 1.0/(D1[n] + double(n)*zinv);
  510. }
  511. if (std::abs(D1[0]) > 100000.0)
  512. throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  513. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  514. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  515. *std::exp(-2.0*z.imag()));
  516. D3[0] = std::complex<double>(0.0, 1.0);
  517. for (int n = 1; n <= nmax_; n++) {
  518. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  519. *(static_cast<double>(n)*zinv- D3[n - 1]);
  520. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  521. }
  522. }
  523. //**********************************************************************************//
  524. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  525. // complex argument (z). //
  526. // Equations (20a) - (21b) //
  527. // //
  528. // Input parameters: //
  529. // z: Complex argument to evaluate Psi and Zeta //
  530. // nmax: Maximum number of terms to calculate Psi and Zeta //
  531. // //
  532. // Output parameters: //
  533. // Psi, Zeta: Riccati-Bessel functions //
  534. //**********************************************************************************//
  535. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  536. std::vector<std::complex<double> >& Psi,
  537. std::vector<std::complex<double> >& Zeta) {
  538. std::complex<double> c_i(0.0, 1.0);
  539. std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
  540. // First, calculate the logarithmic derivatives
  541. calcD1D3(z, D1, D3);
  542. // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  543. Psi[0] = std::sin(z);
  544. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  545. for (int n = 1; n <= nmax_; n++) {
  546. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  547. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  548. }
  549. }
  550. //**********************************************************************************//
  551. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  552. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  553. // //
  554. // Input parameters: //
  555. // z: Complex argument to evaluate jn and h1n //
  556. // nmax_: Maximum number of terms to calculate jn and h1n //
  557. // //
  558. // Output parameters: //
  559. // jn, h1n: Spherical Bessel and Hankel functions //
  560. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  561. // //
  562. // What we actually calculate are the Ricatti-Bessel fucntions and then simply //
  563. // evaluate the spherical Bessel and Hankel functions and their derivatives //
  564. // using the relations: //
  565. // //
  566. // j[n] = Psi[n]/z //
  567. // j'[n] = 0.5*(Psi[n-1]-Psi[n+1]-jn[n])/z //
  568. // h1[n] = Zeta[n]/z //
  569. // h1'[n] = 0.5*(Zeta[n-1]-Zeta[n+1]-h1n[n])/z //
  570. // //
  571. //**********************************************************************************//
  572. void MultiLayerMie::sbesjh(std::complex<double> z,
  573. std::vector<std::complex<double> >& jn, std::vector<std::complex<double> >& jnp,
  574. std::vector<std::complex<double> >& h1n, std::vector<std::complex<double> >& h1np) {
  575. std::vector<std::complex<double> > Psi(nmax_ + 1), Zeta(nmax_ + 1);
  576. // First, calculate the Riccati-Bessel functions
  577. calcPsiZeta(z, Psi, Zeta);
  578. // Now, calculate Spherical Bessel and Hankel functions and their derivatives
  579. for (int n = 0; n < nmax_; n++) {
  580. jn[n] = Psi[n]/z;
  581. h1n[n] = Zeta[n]/z;
  582. if (n == 0) {
  583. jnp[n] = -Psi[1]/z - 0.5*jn[n]/z;
  584. h1np[n] = -Zeta[1]/z - 0.5*h1n[n]/z;
  585. } else {
  586. jnp[n] = 0.5*(Psi[n - 1] - Psi[n + 1] - jn[n])/z;
  587. h1np[n] = 0.5*(Zeta[n - 1] - Zeta[n + 1] - h1n[n])/z;
  588. }
  589. }
  590. }
  591. //**********************************************************************************//
  592. // This function calculates Pi and Tau for a given value of cos(Theta). //
  593. // Equations (26a) - (26c) //
  594. // //
  595. // Input parameters: //
  596. // nmax_: Maximum number of terms to calculate Pi and Tau //
  597. // nTheta: Number of scattering angles //
  598. // Theta: Array containing all the scattering angles where the scattering //
  599. // amplitudes will be calculated //
  600. // //
  601. // Output parameters: //
  602. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  603. //**********************************************************************************//
  604. void MultiLayerMie::calcPiTau(const double& costheta,
  605. std::vector<double>& Pi, std::vector<double>& Tau) {
  606. int n;
  607. //****************************************************//
  608. // Equations (26a) - (26c) //
  609. //****************************************************//
  610. // Initialize Pi and Tau
  611. Pi[0] = 1.0;
  612. Tau[0] = costheta;
  613. // Calculate the actual values
  614. if (nmax_ > 1) {
  615. Pi[1] = 3*costheta*Pi[0];
  616. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  617. for (n = 2; n < nmax_; n++) {
  618. Pi[n] = ((n + n + 1)*costheta*Pi[n - 1] - (n + 1)*Pi[n - 2])/n;
  619. Tau[n] = (n + 1)*costheta*Pi[n] - (n + 2)*Pi[n - 1];
  620. }
  621. }
  622. } // end of MultiLayerMie::calcPiTau(...)
  623. //**********************************************************************************//
  624. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  625. // required to calculate the near-field parameters. //
  626. // //
  627. // Input parameters: //
  628. // Rho: Radial distance //
  629. // Phi: Azimuthal angle //
  630. // Theta: Polar angle //
  631. // zn: Either the spherical Bessel or Hankel function of first kind //
  632. // dzn: Derivative of zn //
  633. // Pi, Tau: Angular functions Pi and Tau //
  634. // n: Order of vector spherical harmonics //
  635. // //
  636. // Output parameters: //
  637. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  638. //**********************************************************************************//
  639. void MultiLayerMie::calcSpherHarm(const double Rho, const double Phi, const double Theta,
  640. const std::complex<double>& zn, const std::complex<double>& dzn,
  641. const double& Pi, const double& Tau, const double& n,
  642. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  643. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  644. // using eq 4.50 in BH
  645. std::complex<double> c_zero(0.0, 0.0);
  646. std::complex<double> deriv = Rho*dzn + zn;
  647. using std::sin;
  648. using std::cos;
  649. Mo1n[0] = c_zero;
  650. Mo1n[1] = cos(Phi)*Pi*zn;
  651. Mo1n[2] = -sin(Phi)*Tau*zn;
  652. Me1n[0] = c_zero;
  653. Me1n[1] = -sin(Phi)*Pi*zn;
  654. Me1n[2] = -cos(Phi)*Tau*zn;
  655. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  656. No1n[1] = sin(Phi)*Tau*deriv/Rho;
  657. No1n[2] = cos(Phi)*Pi*deriv/Rho;
  658. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  659. Ne1n[1] = cos(Phi)*Tau*deriv/Rho;
  660. Ne1n[2] = -sin(Phi)*Pi*deriv/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::ExternalScattCoeffs() {
  682. isExtCoeffsCalc_ = 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. isExtCoeffsCalc_ = true;
  817. } // end of MultiLayerMie::ExternalScattCoeffs(...)
  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. isIntCoeffsCalc_ = false;
  853. isExtCoeffsCalc_ = false;
  854. isMieCalculated_ = false;
  855. // Calculate scattering coefficients
  856. ExternalScattCoeffs();
  857. if (!isExtCoeffsCalc_) // 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 scattering 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::InternalScattCoeffs() {
  927. if (!isExtCoeffsCalc_)
  928. throw std::invalid_argument("(InternalScattCoeffs) You should calculate external coefficients first!");
  929. isIntCoeffsCalc_ = false;
  930. std::complex<double> c_one(1.0, 0.0);
  931. std::complex<double> c_zero(0.0, 0.0);
  932. const int L = refractive_index_.size();
  933. aln_.resize(L + 1);
  934. bln_.resize(L + 1);
  935. cln_.resize(L + 1);
  936. dln_.resize(L + 1);
  937. for (auto& element:aln_) element.resize(nmax_);
  938. for (auto& element:bln_) element.resize(nmax_);
  939. for (auto& element:cln_) element.resize(nmax_);
  940. for (auto& element:dln_) element.resize(nmax_);
  941. // Yang, paragraph under eq. A3
  942. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  943. for (int i = 0; i < nmax_; ++i) {
  944. aln_[L][i] = an_[i];
  945. bln_[L][i] = bn_[i];
  946. cln_[L][i] = c_one;
  947. dln_[L][i] = c_one;
  948. }
  949. std::vector<std::complex<double> > z(L), z1(L);
  950. for (int i = 0; i < L - 1; ++i) {
  951. z[i] = size_param_[i]*refractive_index_[i];
  952. z1[i] = size_param_[i]*refractive_index_[i + 1];
  953. }
  954. z[L - 1] = size_param_[L - 1]*refractive_index_[L - 1];
  955. z1[L - 1] = size_param_[L - 1];
  956. std::vector< std::vector<std::complex<double> > > D1z(L), D1z1(L), D3z(L), D3z1(L);
  957. std::vector< std::vector<std::complex<double> > > Psiz(L), Psiz1(L), Zetaz(L), Zetaz1(L);
  958. for (int l = 0; l < L; ++l) {
  959. D1z[l].resize(nmax_ + 1);
  960. D1z1[l].resize(nmax_ + 1);
  961. D3z[l].resize(nmax_ + 1);
  962. D3z1[l].resize(nmax_ + 1);
  963. Psiz[l].resize(nmax_ + 1);
  964. Psiz1[l].resize(nmax_ + 1);
  965. Zetaz[l].resize(nmax_ + 1);
  966. Zetaz1[l].resize(nmax_ + 1);
  967. }
  968. for (int l = 0; l < L; ++l) {
  969. calcD1D3(z[l], D1z[l], D3z[l]);
  970. calcD1D3(z1[l], D1z1[l], D3z1[l]);
  971. calcPsiZeta(z[l], Psiz[l], Zetaz[l]);
  972. calcPsiZeta(z1[l], Psiz1[l], Zetaz1[l]);
  973. }
  974. auto& m = refractive_index_;
  975. std::vector< std::complex<double> > m1(L);
  976. for (int l = 0; l < L - 1; ++l) m1[l] = m[l + 1];
  977. m1[L - 1] = std::complex<double> (1.0, 0.0);
  978. for (int l = L - 1; l >= 0; l--) {
  979. for (int n = nmax_ - 2; n >= 0; n--) {
  980. auto denomZeta = m1[l]*Zetaz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  981. auto denomPsi = m1[l]*Psiz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  982. auto T1 = aln_[l + 1][n]*Zetaz1[l][n + 1] - dln_[l + 1][n]*Psiz1[l][n + 1];
  983. auto T2 = bln_[l + 1][n]*Zetaz1[l][n + 1] - cln_[l + 1][n]*Psiz1[l][n + 1];
  984. auto T3 = -D1z1[l][n + 1]*dln_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*aln_[l + 1][n]*Zetaz1[l][n + 1];
  985. auto T4 = -D1z1[l][n + 1]*cln_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*bln_[l + 1][n]*Zetaz1[l][n + 1];
  986. // anl
  987. aln_[l][n] = (D1z[l][n + 1]*m1[l]*T1 - m[l]*T3)/denomZeta;
  988. // bnl
  989. bln_[l][n] = (D1z[l][n + 1]*m[l]*T2 - m1[l]*T4)/denomZeta;
  990. // cnl
  991. cln_[l][n] = (D3z[l][n + 1]*m[l]*T2 - m1[l]*T4)/denomPsi;
  992. // dnl
  993. dln_[l][n] = (D3z[l][n + 1]*m1[l]*T1 - m[l]*T3)/denomPsi;
  994. } // end of all n
  995. } // end of all l
  996. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  997. for (int n = 0; n < nmax_; ++n) {
  998. if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  999. else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  1000. if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  1001. else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  1002. }
  1003. isIntCoeffsCalc_ = true;
  1004. } // end of void MultiLayerMie::InternalScattCoeffs()
  1005. // ********************************************************************** //
  1006. // external scattering field = incident + scattered //
  1007. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1008. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1009. // ********************************************************************** //
  1010. void MultiLayerMie::fieldExt(const double Rho, const double Phi, const double Theta,
  1011. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1012. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1013. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1014. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1015. std::vector<std::complex<double> > Ei(3, c_zero), Hi(3, c_zero), Es(3, c_zero), Hs(3, c_zero);
  1016. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1017. std::vector<double> Pi(nmax_), Tau(nmax_);
  1018. // Calculate spherical Bessel and Hankel functions
  1019. sbesjh(Rho, jn, jnp, h1n, h1np);
  1020. // Calculate angular functions Pi and Tau
  1021. calcPiTau(std::cos(Theta), Pi, Tau);
  1022. for (int n = 0; n < nmax_; n++) {
  1023. int n1 = n + 1;
  1024. double rn = static_cast<double>(n1);
  1025. // using BH 4.12 and 4.50
  1026. calcSpherHarm(Rho, Phi, Theta, h1n[n1], h1np[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1027. // scattered field: BH p.94 (4.45)
  1028. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1029. for (int i = 0; i < 3; i++) {
  1030. Es[i] = Es[i] + En*(c_i*an_[n]*N3e1n[i] - bn_[n]*M3o1n[i]);
  1031. Hs[i] = Hs[i] + En*(c_i*bn_[n]*N3o1n[i] + an_[n]*M3e1n[i]);
  1032. }
  1033. }
  1034. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1035. // basis unit vectors = er, etheta, ephi
  1036. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  1037. {
  1038. using std::sin;
  1039. using std::cos;
  1040. Ei[0] = eifac*sin(Theta)*cos(Phi);
  1041. Ei[1] = eifac*cos(Theta)*cos(Phi);
  1042. Ei[2] = -eifac*sin(Phi);
  1043. }
  1044. // magnetic field
  1045. double hffact = 1.0/(cc_*mu_);
  1046. for (int i = 0; i < 3; i++) {
  1047. Hs[i] = hffact*Hs[i];
  1048. }
  1049. // incident H field: BH p.26 (2.43), p.89 (4.21)
  1050. std::complex<double> hffacta = hffact;
  1051. std::complex<double> hifac = eifac*hffacta;
  1052. {
  1053. using std::sin;
  1054. using std::cos;
  1055. Hi[0] = hifac*sin(Theta)*sin(Phi);
  1056. Hi[1] = hifac*cos(Theta)*sin(Phi);
  1057. Hi[2] = hifac*cos(Phi);
  1058. }
  1059. for (int i = 0; i < 3; i++) {
  1060. // electric field E [V m - 1] = EF*E0
  1061. E[i] = Ei[i] + Es[i];
  1062. H[i] = Hi[i] + Hs[i];
  1063. }
  1064. } // end of MultiLayerMie::fieldExt(...)
  1065. //**********************************************************************************//
  1066. // This function calculates the electric (E) and magnetic (H) fields inside and //
  1067. // around the particle. //
  1068. // //
  1069. // Input parameters (coordinates of the point): //
  1070. // Rho: Radial distance //
  1071. // Phi: Azimuthal angle //
  1072. // Theta: Polar angle //
  1073. // //
  1074. // Output parameters: //
  1075. // E, H: Complex electric and magnetic fields //
  1076. //**********************************************************************************//
  1077. void MultiLayerMie::calcField(const double Rho, const double Phi, const double Theta,
  1078. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1079. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1080. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1081. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1082. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  1083. std::vector<std::complex<double> > El(3, c_zero), Hl(3, c_zero);
  1084. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1085. std::vector<double> Pi(nmax_), Tau(nmax_);
  1086. int l = 0; // Layer number
  1087. std::complex<double> ml;
  1088. if (Rho > size_param_.back()) {
  1089. l = size_param_.size();
  1090. ml = c_one;
  1091. } else {
  1092. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  1093. if (Rho <= size_param_[i]) {
  1094. l = i;
  1095. }
  1096. }
  1097. ml = refractive_index_[l];
  1098. }
  1099. // std::complex<double> z1(1.0, 0.0);
  1100. // sbesjh(z1, jn, jnp, h1n, h1np);
  1101. // printf("nmax_ = %d\n", nmax_);
  1102. // printf("$$$$ REAL $$$$$$ Ricatti-Bessel: Calculate and compare against Wolfram Alpha\n");
  1103. // printf("j(0,1) = %16.18f\n", real(jn[0]));
  1104. // printf("WA j() = 0.841470984807896506652502321630\n");
  1105. // printf("j(1,1) = %16.18f\n", real(jn[1]));
  1106. // printf("WA j() = 0.301168678939756789251565714187322395890252640\n");
  1107. // printf("j(10,1) = %.14g\n", real(jn[10]));
  1108. // printf("WA j() = 7.116552640047313023966191737248811458533809572 × 10^-11\n");
  1109. // printf("j(20,1) = %.14g\n", real(jn[16]));
  1110. // printf("WA j() = 7.537795722236872993957557741584960855614358030 × 10^-26\n");
  1111. // std::complex<double> z(1.0, 2.0);
  1112. // sbesjh(z, jn, jnp, h1n, h1np);
  1113. // auto result = jn;
  1114. // printf("===========Calculate and compare against Wolfram Alpha\n");
  1115. // printf("j(0,1+i2) = re(%16.18f)\n im(%16.18f)\n",
  1116. // real(result[0]),
  1117. // imag(result[0]));
  1118. // printf("WA j() = re(1.4169961192118759)\n im(-0.874391197002)\n");
  1119. // printf("j(1,1+i2) = re(%16.18f)\n im(%16.18f)\n",
  1120. // real(result[1]),
  1121. // imag(result[1]));
  1122. // printf("WA j() = re(0.74785726329830368)\n im(0.68179207555304)\n");
  1123. // printf("j(4,1+i2) = re(%16.18f)\n im(%16.18f)\n",
  1124. // real(result[4]),
  1125. // imag(result[4]));
  1126. // printf("WA j() = re(-0.01352410550046)\n im(-0.027169663050653)\n");
  1127. // Calculate spherical Bessel and Hankel functions and their derivatives
  1128. sbesjh(Rho*ml, jn, jnp, h1n, h1np);
  1129. // Calculate angular functions Pi and Tau
  1130. calcPiTau(std::cos(Theta), Pi, Tau);
  1131. for (int n = nmax_ - 1; n >= 0; n--) {
  1132. int n1 = n + 1;
  1133. double rn = static_cast<double>(n1);
  1134. // using BH 4.12 and 4.50
  1135. calcSpherHarm(Rho, Phi, Theta, jn[n1], jnp[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1136. calcSpherHarm(Rho, Phi, Theta, h1n[n1], h1np[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1137. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1138. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1139. for (int i = 0; i < 3; i++) {
  1140. El[i] = El[i] + En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  1141. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  1142. Hl[i] = Hl[i] + En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  1143. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  1144. }
  1145. } // end of for all n
  1146. // magnetic field
  1147. double hffact = 1.0/(cc_*mu_);
  1148. for (int i = 0; i < 3; i++) {
  1149. Hl[i] = hffact*Hl[i];
  1150. }
  1151. for (int i = 0; i < 3; i++) {
  1152. // electric field E [V m - 1] = EF*E0
  1153. E[i] = El[i];
  1154. H[i] = Hl[i];
  1155. }
  1156. } // end of MultiLayerMie::calcField(...)
  1157. //**********************************************************************************//
  1158. // This function calculates complex electric and magnetic field in the surroundings //
  1159. // and inside the particle. //
  1160. // //
  1161. // Input parameters: //
  1162. // L: Number of layers //
  1163. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1164. // x: Array containing the size parameters of the layers [0..L-1] //
  1165. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1166. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1167. // calculations. Only use it if you know what you are doing, otherwise //
  1168. // set this parameter to 0 (zero) and the function will calculate it. //
  1169. // ncoord: Number of coordinate points //
  1170. // Coords: Array containing all coordinates where the complex electric and //
  1171. // magnetic fields will be calculated //
  1172. // //
  1173. // Output parameters: //
  1174. // E, H: Complex electric and magnetic field at the provided coordinates //
  1175. // //
  1176. // Return value: //
  1177. // Number of multipolar expansion terms used for the calculations //
  1178. //**********************************************************************************//
  1179. void MultiLayerMie::RunFieldCalculation() {
  1180. // Calculate external scattering coefficients an_ and bn_
  1181. ExternalScattCoeffs();
  1182. // Calculate internal scattering coefficients aln_, bln_, cln_, and dln_
  1183. InternalScattCoeffs();
  1184. long total_points = coords_[0].size();
  1185. E_.resize(total_points);
  1186. H_.resize(total_points);
  1187. for (auto& f : E_) f.resize(3);
  1188. for (auto& f : H_) f.resize(3);
  1189. for (int point = 0; point < total_points; point++) {
  1190. const double& Xp = coords_[0][point];
  1191. const double& Yp = coords_[1][point];
  1192. const double& Zp = coords_[2][point];
  1193. // Convert to spherical coordinates
  1194. double Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1195. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1196. double Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1197. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1198. double Phi = (Xp != 0.0 || Yp != 0.0) ? std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1199. // Avoid convergence problems due to Rho too small
  1200. if (Rho < 1e-5) Rho = 1e-5;
  1201. //*******************************************************//
  1202. // external scattering field = incident + scattered //
  1203. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1204. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1205. //*******************************************************//
  1206. // This array contains the fields in spherical coordinates
  1207. std::vector<std::complex<double> > Es(3), Hs(3);
  1208. // Firstly the easiest case: the field outside the particle
  1209. // if (Rho >= GetSizeParameter()) {
  1210. // fieldExt(Rho, Phi, Theta, Es, Hs);
  1211. // } else {
  1212. calcField(Rho, Phi, Theta, Es, Hs); //Should work fine both: inside and outside the particle
  1213. //}
  1214. { //Now, convert the fields back to cartesian coordinates
  1215. using std::sin;
  1216. using std::cos;
  1217. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1218. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1219. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1220. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1221. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1222. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1223. }
  1224. } // end of for all field coordinates
  1225. } // end of MultiLayerMie::RunFieldCalculation()
  1226. } // end of namespace nmie