nmie.cc 81 KB

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