nmie.cc 76 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536
  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: Complex 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, std::vector<std::complex<double> >& jnp,
  482. std::vector<std::complex<double> >& h1n, std::vector<std::complex<double> >& h1np) {
  483. const int limit = 20000;
  484. const double accur = 1.0e-12;
  485. const double tm30 = 1e-30;
  486. const std::complex<double> c_i(0.0, 1.0);
  487. double absc;
  488. std::complex<double> zi, w;
  489. std::complex<double> pl, f, b, d, c, del, jn0;
  490. absc = std::abs(std::real(z)) + std::abs(std::imag(z));
  491. if ((absc < accur) || (std::imag(z) < -3.0)) {
  492. throw std::invalid_argument("TODO add error description for condition if ((absc < accur) || (std::imag(z) < -3.0))");
  493. }
  494. zi = 1.0/z;
  495. w = zi + zi;
  496. pl = double(nmax_)*zi;
  497. f = pl + zi;
  498. b = f + f + zi;
  499. d = 0.0;
  500. c = f;
  501. for (int l = 0; l < limit; l++) {
  502. d = b - d;
  503. c = b - 1.0/c;
  504. absc = std::abs(std::real(d)) + std::abs(std::imag(d));
  505. if (absc < tm30) {
  506. d = tm30;
  507. }
  508. absc = std::abs(std::real(c)) + std::abs(std::imag(c));
  509. if (absc < tm30) {
  510. c = tm30;
  511. }
  512. d = 1.0/d;
  513. del = d*c;
  514. f = f*del;
  515. b += w;
  516. absc = std::abs(std::real(del - 1.0)) + std::abs(std::imag(del - 1.0));
  517. // We have obtained the desired accuracy
  518. if (absc < accur) break;
  519. }
  520. if (absc > accur) {
  521. throw std::invalid_argument("We were not able to obtain the desired accuracy (MultiLayerMie::sbesjh)");
  522. }
  523. jn[nmax_ - 1] = tm30;
  524. jnp[nmax_ - 1] = f*jn[nmax_ - 1];
  525. // Downward recursion to n=0 (N.B. Coulomb Functions)
  526. for (int n = nmax_ - 2; n >= 0; n--) {
  527. jn[n] = pl*jn[n + 1] + jnp[n + 1];
  528. jnp[n] = pl*jn[n] - jn[n + 1];
  529. pl = pl - zi;
  530. }
  531. // Calculate the n=0 Bessel Functions
  532. jn0 = zi*std::sin(z);
  533. h1n[0] = -c_i*zi*std::exp(c_i*z);
  534. h1np[0] = h1n[0]*(c_i - zi);
  535. // Rescale j[n], j'[n], converting to spherical Bessel functions.
  536. // Recur h1[n], h1'[n] as spherical Bessel functions.
  537. w = 1.0/jn[0];
  538. pl = zi;
  539. for (int n = 0; n < nmax_; n++) {
  540. jn[n] = jn0*w*jn[n];
  541. jnp[n] = jn0*w*jnp[n] - zi*jn[n];
  542. if (n > 0) {
  543. h1n[n] = (pl - zi)*h1n[n - 1] - h1np[n - 1];
  544. // check if hankel is increasing (upward stable)
  545. if (std::abs(h1n[n]) < std::abs(h1n[n - 1])) {
  546. throw std::invalid_argument("Error: Hankel not increasing! (MultiLayerMie::sbesjh)");
  547. }
  548. pl += zi;
  549. h1np[n] = -pl*h1n[n] + h1n[n - 1];
  550. }
  551. }
  552. }
  553. // ********************************************************************** //
  554. // Calculate an - equation (5) //
  555. // ********************************************************************** //
  556. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  557. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  558. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  559. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  560. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  561. return Num/Denom;
  562. }
  563. // ********************************************************************** //
  564. // Calculate bn - equation (6) //
  565. // ********************************************************************** //
  566. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  567. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  568. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  569. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  570. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  571. return Num/Denom;
  572. }
  573. // ********************************************************************** //
  574. // Calculates S1 - equation (25a) //
  575. // ********************************************************************** //
  576. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  577. double Pi, double Tau) {
  578. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  579. }
  580. // ********************************************************************** //
  581. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  582. // Pi and Tau) //
  583. // ********************************************************************** //
  584. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  585. double Pi, double Tau) {
  586. return calc_S1(n, an, bn, Tau, Pi);
  587. }
  588. //**********************************************************************************//
  589. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  590. // real argument (x). //
  591. // Equations (20a) - (21b) //
  592. // //
  593. // Input parameters: //
  594. // x: Real argument to evaluate Psi and Zeta //
  595. // nmax: Maximum number of terms to calculate Psi and Zeta //
  596. // //
  597. // Output parameters: //
  598. // Psi, Zeta: Riccati-Bessel functions //
  599. //**********************************************************************************//
  600. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  601. std::vector<std::complex<double> > D1,
  602. std::vector<std::complex<double> > D3,
  603. std::vector<std::complex<double> >& Psi,
  604. std::vector<std::complex<double> >& Zeta) {
  605. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  606. std::complex<double> c_i(0.0, 1.0);
  607. Psi[0] = std::sin(z);
  608. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  609. for (int n = 1; n <= nmax_; n++) {
  610. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  611. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  612. }
  613. }
  614. //**********************************************************************************//
  615. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  616. // functions (D1 and D3) for a complex argument (z). //
  617. // Equations (16a), (16b) and (18a) - (18d) //
  618. // //
  619. // Input parameters: //
  620. // z: Complex argument to evaluate D1 and D3 //
  621. // nmax_: Maximum number of terms to calculate D1 and D3 //
  622. // //
  623. // Output parameters: //
  624. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  625. //**********************************************************************************//
  626. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  627. std::vector<std::complex<double> >& D1,
  628. std::vector<std::complex<double> >& D3) {
  629. // Downward recurrence for D1 - equations (16a) and (16b)
  630. D1[nmax_] = std::complex<double>(0.0, 0.0);
  631. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  632. for (int n = nmax_; n > 0; n--) {
  633. D1[n - 1] = double(n)*zinv - 1.0/(D1[n] + double(n)*zinv);
  634. }
  635. if (std::abs(D1[0]) > 100000.0)
  636. throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  637. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  638. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  639. *std::exp(-2.0*z.imag()));
  640. D3[0] = std::complex<double>(0.0, 1.0);
  641. for (int n = 1; n <= nmax_; n++) {
  642. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  643. *(static_cast<double>(n)*zinv- D3[n - 1]);
  644. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  645. }
  646. }
  647. //**********************************************************************************//
  648. // This function calculates Pi and Tau for a given value of cos(Theta). //
  649. // Equations (26a) - (26c) //
  650. // //
  651. // Input parameters: //
  652. // nmax_: Maximum number of terms to calculate Pi and Tau //
  653. // nTheta: Number of scattering angles //
  654. // Theta: Array containing all the scattering angles where the scattering //
  655. // amplitudes will be calculated //
  656. // //
  657. // Output parameters: //
  658. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  659. //**********************************************************************************//
  660. void MultiLayerMie::calcPiTau(const double& costheta,
  661. std::vector<double>& Pi, std::vector<double>& Tau) {
  662. int n;
  663. //****************************************************//
  664. // Equations (26a) - (26c) //
  665. //****************************************************//
  666. // Initialize Pi and Tau
  667. Pi[0] = 1.0;
  668. Tau[0] = costheta;
  669. // Calculate the actual values
  670. if (nmax_ > 1) {
  671. Pi[1] = 3*costheta*Pi[0];
  672. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  673. for (n = 2; n < nmax_; n++) {
  674. Pi[n] = ((n + n + 1)*costheta*Pi[n - 1] - (n + 1)*Pi[n - 2])/n;
  675. Tau[n] = (n + 1)*costheta*Pi[n] - (n + 2)*Pi[n - 1];
  676. }
  677. }
  678. } // end of MultiLayerMie::calcPiTau(...)
  679. void MultiLayerMie::calcSpherHarm(const double Rho, const double Phi, const double Theta,
  680. const std::complex<double>& zn, const std::complex<double>& dzn,
  681. const double& Pi, const double& Tau, const double& n,
  682. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  683. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  684. // using eq 4.50 in BH
  685. std::complex<double> c_zero(0.0, 0.0);
  686. std::complex<double> deriv = Rho*dzn + zn;
  687. using std::sin;
  688. using std::cos;
  689. Mo1n[0] = c_zero;
  690. Mo1n[1] = cos(Phi)*Pi*zn;
  691. Mo1n[2] = -sin(Phi)*Tau*zn;
  692. Me1n[0] = c_zero;
  693. Me1n[1] = -sin(Phi)*Pi*zn;
  694. Me1n[2] = -cos(Phi)*Tau*zn;
  695. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  696. No1n[1] = sin(Phi)*Tau*deriv/Rho;
  697. No1n[2] = cos(Phi)*Pi*deriv/Rho;
  698. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  699. Ne1n[1] = cos(Phi)*Tau*deriv/Rho;
  700. Ne1n[2] = -sin(Phi)*Pi*deriv/Rho;
  701. } // end of MultiLayerMie::calcSpherHarm(...)
  702. //**********************************************************************************//
  703. // This function calculates the scattering coefficients required to calculate //
  704. // both the near- and far-field parameters. //
  705. // //
  706. // Input parameters: //
  707. // L: Number of layers //
  708. // pl: Index of PEC layer. If there is none just send -1 //
  709. // x: Array containing the size parameters of the layers [0..L-1] //
  710. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  711. // nmax: Maximum number of multipolar expansion terms to be used for the //
  712. // calculations. Only use it if you know what you are doing, otherwise //
  713. // set this parameter to -1 and the function will calculate it. //
  714. // //
  715. // Output parameters: //
  716. // an, bn: Complex scattering amplitudes //
  717. // //
  718. // Return value: //
  719. // Number of multipolar expansion terms used for the calculations //
  720. //**********************************************************************************//
  721. void MultiLayerMie::ExtScattCoeffs() {
  722. areExtCoeffsCalc_ = false;
  723. const std::vector<double>& x = size_param_;
  724. const std::vector<std::complex<double> >& m = refr_index_;
  725. const int& pl = PEC_layer_position_;
  726. const int L = refr_index_.size();
  727. //************************************************************************//
  728. // Calculate the index of the first layer. It can be either 0 (default) //
  729. // or the index of the outermost PEC layer. In the latter case all layers //
  730. // below the PEC are discarded. //
  731. // ***********************************************************************//
  732. // TODO, is it possible for PEC to have a zero index? If yes than
  733. // is should be:
  734. // int fl = (pl > - 1) ? pl : 0;
  735. // This will give the same result, however, it corresponds the
  736. // logic - if there is PEC, than first layer is PEC.
  737. // Well, I followed the logic: First layer is always zero unless it has
  738. // an upper PEC layer.
  739. int fl = (pl > 0) ? pl : 0;
  740. if (nmax_preset_ <= 0) calcNmax(fl);
  741. else nmax_ = nmax_preset_;
  742. std::complex<double> z1, z2;
  743. //**************************************************************************//
  744. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  745. // means that index = layer number - 1 or index = n - 1. The only exception //
  746. // are the arrays for representing D1, D3 and Q because they need a value //
  747. // for the index 0 (zero), hence it is important to consider this shift //
  748. // between different arrays. The change was done to optimize memory usage. //
  749. //**************************************************************************//
  750. // Allocate memory to the arrays
  751. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  752. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  753. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  754. for (int l = 0; l < L; l++) {
  755. Q[l].resize(nmax_ + 1);
  756. Ha[l].resize(nmax_);
  757. Hb[l].resize(nmax_);
  758. }
  759. an_.resize(nmax_);
  760. bn_.resize(nmax_);
  761. PsiZeta_.resize(nmax_ + 1);
  762. std::vector<std::complex<double> > D1XL(nmax_ + 1), D3XL(nmax_ + 1),
  763. PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  764. //*************************************************//
  765. // Calculate D1 and D3 for z1 in the first layer //
  766. //*************************************************//
  767. if (fl == pl) { // PEC layer
  768. for (int n = 0; n <= nmax_; n++) {
  769. D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
  770. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  771. }
  772. } else { // Regular layer
  773. z1 = x[fl]* m[fl];
  774. // Calculate D1 and D3
  775. calcD1D3(z1, D1_mlxl, D3_mlxl);
  776. }
  777. //******************************************************************//
  778. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  779. //******************************************************************//
  780. for (int n = 0; n < nmax_; n++) {
  781. Ha[fl][n] = D1_mlxl[n + 1];
  782. Hb[fl][n] = D1_mlxl[n + 1];
  783. }
  784. //*****************************************************//
  785. // Iteration from the second layer to the last one (L) //
  786. //*****************************************************//
  787. std::complex<double> Temp, Num, Denom;
  788. std::complex<double> G1, G2;
  789. for (int l = fl + 1; l < L; l++) {
  790. //************************************************************//
  791. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  792. //************************************************************//
  793. z1 = x[l]*m[l];
  794. z2 = x[l - 1]*m[l];
  795. //Calculate D1 and D3 for z1
  796. calcD1D3(z1, D1_mlxl, D3_mlxl);
  797. //Calculate D1 and D3 for z2
  798. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  799. //*********************************************//
  800. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  801. //*********************************************//
  802. // Upward recurrence for Q - equations (19a) and (19b)
  803. Num = std::exp(-2.0*(z1.imag() - z2.imag()))
  804. *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
  805. Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
  806. Q[l][0] = Num/Denom;
  807. for (int n = 1; n <= nmax_; n++) {
  808. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  809. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  810. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  811. }
  812. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  813. for (int n = 1; n <= nmax_; n++) {
  814. //Ha
  815. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  816. G1 = -D1_mlxlM1[n];
  817. G2 = -D3_mlxlM1[n];
  818. } else {
  819. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  820. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  821. } // end of if PEC
  822. Temp = Q[l][n]*G1;
  823. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  824. Denom = G2 - Temp;
  825. Ha[l][n - 1] = Num/Denom;
  826. //Hb
  827. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  828. G1 = Hb[l - 1][n - 1];
  829. G2 = Hb[l - 1][n - 1];
  830. } else {
  831. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  832. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  833. } // end of if PEC
  834. Temp = Q[l][n]*G1;
  835. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  836. Denom = (G2- Temp);
  837. Hb[l][n - 1] = (Num/ Denom);
  838. } // end of for Ha and Hb terms
  839. } // end of for layers iteration
  840. //**************************************//
  841. //Calculate D1, D3, Psi and Zeta for XL //
  842. //**************************************//
  843. // Calculate D1XL and D3XL
  844. calcD1D3(x[L - 1], D1XL, D3XL);
  845. // Calculate PsiXL and ZetaXL
  846. calcPsiZeta(x[L - 1], D1XL, D3XL, PsiXL, ZetaXL);
  847. //*********************************************************************//
  848. // Finally, we calculate the scattering coefficients (an and bn) and //
  849. // the angular functions (Pi and Tau). Note that for these arrays the //
  850. // first layer is 0 (zero), in future versions all arrays will follow //
  851. // this convention to save memory. (13 Nov, 2014) //
  852. //*********************************************************************//
  853. for (int n = 0; n < nmax_; n++) {
  854. //********************************************************************//
  855. //Expressions for calculating an and bn coefficients are not valid if //
  856. //there is only one PEC layer (ie, for a simple PEC sphere). //
  857. //********************************************************************//
  858. if (pl < (L - 1)) {
  859. 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]);
  860. 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]);
  861. } else {
  862. 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]);
  863. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  864. }
  865. } // end of for an and bn terms
  866. areExtCoeffsCalc_ = true;
  867. } // end of void MultiLayerMie::ExtScattCoeffs(...)
  868. //**********************************************************************************//
  869. // This function calculates the actual scattering parameters and amplitudes //
  870. // //
  871. // Input parameters: //
  872. // L: Number of layers //
  873. // pl: Index of PEC layer. If there is none just send -1 //
  874. // x: Array containing the size parameters of the layers [0..L-1] //
  875. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  876. // nTheta: Number of scattering angles //
  877. // Theta: Array containing all the scattering angles where the scattering //
  878. // amplitudes will be calculated //
  879. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  880. // calculations. Only use it if you know what you are doing, otherwise //
  881. // set this parameter to -1 and the function will calculate it //
  882. // //
  883. // Output parameters: //
  884. // Qext: Efficiency factor for extinction //
  885. // Qsca: Efficiency factor for scattering //
  886. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  887. // Qbk: Efficiency factor for backscattering //
  888. // Qpr: Efficiency factor for the radiation pressure //
  889. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  890. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  891. // S1, S2: Complex scattering amplitudes //
  892. // //
  893. // Return value: //
  894. // Number of multipolar expansion terms used for the calculations //
  895. //**********************************************************************************//
  896. void MultiLayerMie::RunMieCalculation() {
  897. if (size_param_.size() != refr_index_.size())
  898. throw std::invalid_argument("Each size parameter should have only one index!");
  899. if (size_param_.size() == 0)
  900. throw std::invalid_argument("Initialize model first!");
  901. const std::vector<double>& x = size_param_;
  902. areIntCoeffsCalc_ = false;
  903. areExtCoeffsCalc_ = false;
  904. isMieCalculated_ = false;
  905. // Calculate scattering coefficients
  906. ExtScattCoeffs();
  907. // for (int i = 0; i < nmax_; i++) {
  908. // printf("a[%i] = %g, %g; b[%i] = %g, %g\n", i, an_[i].real(), an_[i].imag(), i, bn_[i].real(), bn_[i].imag());
  909. // }
  910. if (!areExtCoeffsCalc_)
  911. throw std::invalid_argument("Calculation of scattering coefficients failed!");
  912. // Initialize the scattering parameters
  913. Qext_ = 0;
  914. Qsca_ = 0;
  915. Qabs_ = 0;
  916. Qbk_ = 0;
  917. Qpr_ = 0;
  918. asymmetry_factor_ = 0;
  919. albedo_ = 0;
  920. Qsca_ch_.clear();
  921. Qext_ch_.clear();
  922. Qabs_ch_.clear();
  923. Qbk_ch_.clear();
  924. Qpr_ch_.clear();
  925. Qsca_ch_.resize(nmax_ - 1);
  926. Qext_ch_.resize(nmax_ - 1);
  927. Qabs_ch_.resize(nmax_ - 1);
  928. Qbk_ch_.resize(nmax_ - 1);
  929. Qpr_ch_.resize(nmax_ - 1);
  930. Qsca_ch_norm_.resize(nmax_ - 1);
  931. Qext_ch_norm_.resize(nmax_ - 1);
  932. Qabs_ch_norm_.resize(nmax_ - 1);
  933. Qbk_ch_norm_.resize(nmax_ - 1);
  934. Qpr_ch_norm_.resize(nmax_ - 1);
  935. // Initialize the scattering amplitudes
  936. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  937. S1_.swap(tmp1);
  938. S2_ = S1_;
  939. std::vector<double> Pi(nmax_), Tau(nmax_);
  940. std::complex<double> Qbktmp(0.0, 0.0);
  941. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  942. // By using downward recurrence we avoid loss of precision due to float rounding errors
  943. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  944. // http://en.wikipedia.org/wiki/Loss_of_significance
  945. for (int i = nmax_ - 2; i >= 0; i--) {
  946. const int n = i + 1;
  947. // Equation (27)
  948. Qext_ch_norm_[i] = (an_[i].real() + bn_[i].real());
  949. Qext_ch_[i] = (n + n + 1.0)*Qext_ch_norm_[i];
  950. //Qext_ch_[i] = (n + n + 1)*(an_[i].real() + bn_[i].real());
  951. Qext_ += Qext_ch_[i];
  952. // Equation (28)
  953. Qsca_ch_norm_[i] = (an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  954. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  955. Qsca_ch_[i] = (n + n + 1.0)*Qsca_ch_norm_[i];
  956. Qsca_ += Qsca_ch_[i];
  957. // Qsca_ch_[i] += (n + n + 1)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  958. // + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  959. // Equation (29) TODO We must check carefully this equation. If we
  960. // remove the typecast to double then the result changes. Which is
  961. // the correct one??? Ovidio (2014/12/10) With cast ratio will
  962. // give double, without cast (n + n + 1)/(n*(n + 1)) will be
  963. // rounded to integer. Tig (2015/02/24)
  964. Qpr_ch_[i]=((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  965. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  966. Qpr_ += Qpr_ch_[i];
  967. // Equation (33)
  968. Qbktmp_ch[i] = (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  969. Qbktmp += Qbktmp_ch[i];
  970. // Calculate the scattering amplitudes (S1 and S2) //
  971. // Equations (25a) - (25b) //
  972. for (int t = 0; t < theta_.size(); t++) {
  973. calcPiTau(std::cos(theta_[t]), Pi, Tau);
  974. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
  975. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
  976. }
  977. }
  978. double x2 = pow2(x.back());
  979. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  980. for (double& Q : Qext_ch_) Q = 2.0*Q/x2;
  981. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  982. for (double& Q : Qsca_ch_) Q = 2.0*Q/x2;
  983. //for (double& Q : Qsca_ch_norm_) Q = 2.0*Q/x2;
  984. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  985. for (int i = 0; i < nmax_ - 1; ++i) Qpr_ch_[i] = Qext_ch_[i] - 4.0*Qpr_ch_[i]/x2;
  986. Qabs_ = Qext_ - Qsca_; // Equation (30)
  987. for (int i = 0; i < nmax_ - 1; ++i) {
  988. Qabs_ch_[i] = Qext_ch_[i] - Qsca_ch_[i];
  989. Qabs_ch_norm_[i] = Qext_ch_norm_[i] - Qsca_ch_norm_[i];
  990. }
  991. albedo_ = Qsca_/Qext_; // Equation (31)
  992. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  993. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  994. isMieCalculated_ = true;
  995. }
  996. // ********************************************************************** //
  997. // ********************************************************************** //
  998. // ********************************************************************** //
  999. void MultiLayerMie::IntScattCoeffs() {
  1000. if (!areExtCoeffsCalc_)
  1001. throw std::invalid_argument("(IntScattCoeffs) You should calculate external coefficients first!");
  1002. areIntCoeffsCalc_ = false;
  1003. std::complex<double> c_one(1.0, 0.0);
  1004. std::complex<double> c_zero(0.0, 0.0);
  1005. const int L = refr_index_.size();
  1006. // we need to fill
  1007. // std::vector< std::vector<std::complex<double> > > anl_, bnl_, cnl_, dnl_;
  1008. // for n = [0..nmax_) and for l=[L..0)
  1009. // TODO: to decrease cache miss outer loop is with n and inner with reversed l
  1010. // at the moment outer is forward l and inner in n
  1011. anl_.resize(L + 1);
  1012. bnl_.resize(L + 1);
  1013. cnl_.resize(L + 1);
  1014. dnl_.resize(L + 1);
  1015. for (auto& element:anl_) element.resize(nmax_);
  1016. for (auto& element:bnl_) element.resize(nmax_);
  1017. for (auto& element:cnl_) element.resize(nmax_);
  1018. for (auto& element:dnl_) element.resize(nmax_);
  1019. // Yang, paragraph under eq. A3
  1020. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  1021. for (int i = 0; i < nmax_; ++i) {
  1022. anl_[L][i] = an_[i];
  1023. bnl_[L][i] = bn_[i];
  1024. cnl_[L][i] = c_one;
  1025. dnl_[L][i] = c_one;
  1026. }
  1027. std::vector<std::complex<double> > z(L), z1(L);
  1028. for (int i = 0; i < L - 1; ++i) {
  1029. z[i] = size_param_[i]*refr_index_[i];
  1030. z1[i] = size_param_[i]*refr_index_[i + 1];
  1031. }
  1032. z[L - 1] = size_param_[L - 1]*refr_index_[L - 1];
  1033. z1[L - 1] = size_param_[L - 1];
  1034. std::vector< std::vector<std::complex<double> > > D1z(L), D1z1(L), D3z(L), D3z1(L);
  1035. std::vector< std::vector<std::complex<double> > > Psiz(L), Psiz1(L), Zetaz(L), Zetaz1(L);
  1036. for (int l = 0; l < L; ++l) {
  1037. D1z[l].resize(nmax_ + 1);
  1038. D1z1[l].resize(nmax_ + 1);
  1039. D3z[l].resize(nmax_ + 1);
  1040. D3z1[l].resize(nmax_ + 1);
  1041. Psiz[l].resize(nmax_ + 1);
  1042. Psiz1[l].resize(nmax_ + 1);
  1043. Zetaz[l].resize(nmax_ + 1);
  1044. Zetaz1[l].resize(nmax_ + 1);
  1045. }
  1046. for (int l = 0; l < L; ++l) {
  1047. calcD1D3(z[l], D1z[l], D3z[l]);
  1048. calcD1D3(z1[l], D1z1[l], D3z1[l]);
  1049. calcPsiZeta(z[l], D1z[l], D3z[l], Psiz[l], Zetaz[l]);
  1050. calcPsiZeta(z1[l], D1z1[l], D3z1[l], Psiz1[l], Zetaz1[l]);
  1051. }
  1052. auto& m = refr_index_;
  1053. std::vector< std::complex<double> > m1(L);
  1054. for (int l = 0; l < L - 1; ++l) m1[l] = m[l + 1];
  1055. m1[L - 1] = std::complex<double> (1.0, 0.0);
  1056. // for (auto zz : m) printf ("m[i]=%g \n\n ", zz.real());
  1057. for (int l = L - 1; l >= 0; l--) {
  1058. for (int n = nmax_ - 2; n >= 0; n--) {
  1059. auto denomZeta = m1[l]*Zetaz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1060. auto denomPsi = m1[l]*Psiz[l][n + 1]*(D1z[l][n + 1] - D3z[l][n + 1]);
  1061. auto T1 = anl_[l + 1][n]*Zetaz1[l][n + 1] - dnl_[l + 1][n]*Psiz1[l][n + 1];
  1062. auto T2 = bnl_[l + 1][n]*Zetaz1[l][n + 1] - cnl_[l + 1][n]*Psiz1[l][n + 1];
  1063. auto T3 = -D1z1[l][n + 1]*dnl_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*anl_[l + 1][n]*Zetaz1[l][n + 1];
  1064. auto T4 = -D1z1[l][n + 1]*cnl_[l + 1][n]*Psiz1[l][n + 1] + D3z1[l][n + 1]*bnl_[l + 1][n]*Zetaz1[l][n + 1];
  1065. // anl
  1066. anl_[l][n] = (D1z[l][n + 1]*m1[l]*T1 - m[l]*T3)/denomZeta;
  1067. // bnl
  1068. bnl_[l][n] = (D1z[l][n + 1]*m[l]*T2 - m1[l]*T4)/denomZeta;
  1069. // cnl
  1070. cnl_[l][n] = (D3z[l][n + 1]*m[l]*T2 - m1[l]*T4)/denomPsi;
  1071. // dnl
  1072. dnl_[l][n] = (D3z[l][n + 1]*m1[l]*T1 - m[l]*T3)/denomPsi;
  1073. } // end of all n
  1074. } // end of all l
  1075. // Check the result and change an__0 and bn__0 for exact zero
  1076. for (int n = 0; n < nmax_; ++n) {
  1077. if (std::abs(anl_[0][n]) < 1e-10) anl_[0][n] = 0.0;
  1078. else throw std::invalid_argument("Unstable calculation of a__0_n!");
  1079. if (std::abs(bnl_[0][n]) < 1e-10) bnl_[0][n] = 0.0;
  1080. else throw std::invalid_argument("Unstable calculation of b__0_n!");
  1081. }
  1082. // for (int l = 0; l < L; ++l) {
  1083. // printf("l=%d --> ", l);
  1084. // for (int n = 0; n < nmax_ + 1; ++n) {
  1085. // if (n < 20) continue;
  1086. // printf("n=%d --> D1zn=%g, D3zn=%g, D1zn=%g, D3zn=%g || ",
  1087. // n,
  1088. // D1z[l][n].real(), D3z[l][n].real(),
  1089. // D1z1[l][n].real(), D3z1[l][n].real());
  1090. // }
  1091. // printf("\n\n");
  1092. // }
  1093. // for (int l = 0; l < L; ++l) {
  1094. // printf("l=%d --> ", l);
  1095. // for (int n = 0; n < nmax_ + 1; ++n) {
  1096. // printf("n=%d --> D1zn=%g, D3zn=%g, D1zn=%g, D3zn=%g || ",
  1097. // n,
  1098. // D1z[l][n].real(), D3z[l][n].real(),
  1099. // D1z1[l][n].real(), D3z1[l][n].real());
  1100. // }
  1101. // printf("\n\n");
  1102. // }
  1103. //for (int i = 0; i < L + 1; ++i) {
  1104. // printf("Layer =%d ---> \n", i);
  1105. // for (int n = 0; n < nmax_; ++n) {
  1106. // if (n < 20) continue;
  1107. // printf(" || n=%d --> a=%g,%g b=%g,%g c=%g,%g d=%g,%g\n",
  1108. // n,
  1109. // anl_[i][n].real(), anl_[i][n].imag(),
  1110. // bnl_[i][n].real(), bnl_[i][n].imag(),
  1111. // cnl_[i][n].real(), cnl_[i][n].imag(),
  1112. // dnl_[i][n].real(), dnl_[i][n].imag());
  1113. // }
  1114. // printf("\n\n");
  1115. //}
  1116. areIntCoeffsCalc_ = true;
  1117. }
  1118. // ********************************************************************** //
  1119. // ********************************************************************** //
  1120. // ********************************************************************** //
  1121. // external scattering field = incident + scattered //
  1122. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1123. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1124. // ********************************************************************** //
  1125. void MultiLayerMie::fieldExt(const double Rho, const double Phi, const double Theta,
  1126. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1127. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1128. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1129. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1130. std::vector<std::complex<double> > Ei(3, c_zero), Hi(3, c_zero), Es(3, c_zero), Hs(3, c_zero);
  1131. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1132. std::vector<double> Pi(nmax_), Tau(nmax_);
  1133. // Calculate spherical Bessel and Hankel functions
  1134. sbesjh(Rho, jn, jnp, h1n, h1np);
  1135. // Calculate angular functions Pi and Tau
  1136. calcPiTau(std::cos(Theta), Pi, Tau);
  1137. for (int n = 0; n < nmax_; n++) {
  1138. int n1 = n + 1;
  1139. double rn = static_cast<double>(n1);
  1140. // using BH 4.12 and 4.50
  1141. calcSpherHarm(Rho, Phi, Theta, h1n[n1], h1np[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1142. // scattered field: BH p.94 (4.45)
  1143. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1144. for (int i = 0; i < 3; i++) {
  1145. Es[i] = Es[i] + En*(c_i*an_[n]*N3e1n[i] - bn_[n]*M3o1n[i]);
  1146. Hs[i] = Hs[i] + En*(c_i*bn_[n]*N3o1n[i] + an_[n]*M3e1n[i]);
  1147. }
  1148. }
  1149. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1150. // basis unit vectors = er, etheta, ephi
  1151. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  1152. {
  1153. using std::sin;
  1154. using std::cos;
  1155. Ei[0] = eifac*sin(Theta)*cos(Phi);
  1156. Ei[1] = eifac*cos(Theta)*cos(Phi);
  1157. Ei[2] = -eifac*sin(Phi);
  1158. }
  1159. // magnetic field
  1160. double hffact = 1.0/(cc_*mu_);
  1161. for (int i = 0; i < 3; i++) {
  1162. Hs[i] = hffact*Hs[i];
  1163. }
  1164. // incident H field: BH p.26 (2.43), p.89 (4.21)
  1165. std::complex<double> hffacta = hffact;
  1166. std::complex<double> hifac = eifac*hffacta;
  1167. {
  1168. using std::sin;
  1169. using std::cos;
  1170. Hi[0] = hifac*sin(Theta)*sin(Phi);
  1171. Hi[1] = hifac*cos(Theta)*sin(Phi);
  1172. Hi[2] = hifac*cos(Phi);
  1173. }
  1174. for (int i = 0; i < 3; i++) {
  1175. // electric field E [V m - 1] = EF*E0
  1176. E[i] = Ei[i] + Es[i];
  1177. H[i] = Hi[i] + Hs[i];
  1178. // printf("ext E[%d]=%g",i,std::abs(E[i]));
  1179. }
  1180. } // end of MultiLayerMie::fieldExt(...)
  1181. // ********************************************************************** //
  1182. // ********************************************************************** //
  1183. // ********************************************************************** //
  1184. void MultiLayerMie::fieldInt(const double Rho, const double Phi, const double Theta,
  1185. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1186. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1187. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1188. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1189. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  1190. std::vector<std::complex<double> > El(3, c_zero), Ei(3, c_zero), Eic(3, c_zero), Hl(3, c_zero);
  1191. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1192. std::vector<double> Pi(nmax_), Tau(nmax_);
  1193. int l = 0; // Layer number
  1194. std::complex<double> ml;
  1195. if (Rho > size_param_.back()) {
  1196. l = size_param_.size();
  1197. ml = c_one;
  1198. } else {
  1199. for (int i = 0; i < size_param_.size() - 1; ++i) {
  1200. if (size_param_[i] < Rho && Rho <= size_param_[i + 1]) {
  1201. l = i;
  1202. }
  1203. }
  1204. ml = refr_index_[l];
  1205. }
  1206. // for (int i = 0; i < size_param_.size(); i++) {
  1207. // printf("x[%i] = %g; m[%i] = %g, %g; ", i, size_param_[i], i, refr_index_[i].real(), refr_index_[i].imag());
  1208. // }
  1209. // printf("\nRho = %g; Phi = %g; Theta = %g; x[%i] = %g; m[%i] = %g, %g\n", Rho, Phi, Theta, l, size_param_[l], l, ml.real(), ml.imag());
  1210. // Calculate spherical Bessel and Hankel functions
  1211. sbesjh(Rho*ml, jn, jnp, h1n, h1np);
  1212. // Calculate angular functions Pi and Tau
  1213. calcPiTau(std::cos(Theta), Pi, Tau);
  1214. for (int n = nmax_ - 1; n >= 0; n--) {
  1215. int n1 = n + 1;
  1216. double rn = static_cast<double>(n1);
  1217. // using BH 4.12 and 4.50
  1218. calcSpherHarm(Rho, Phi, Theta, jn[n1], jnp[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1219. calcSpherHarm(Rho, Phi, Theta, h1n[n1], h1np[n], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1220. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1221. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1222. for (int i = 0; i < 3; i++) {
  1223. Ei[i] = Ei[i] + En*(M1o1n[i] - c_i*N1e1n[i]);
  1224. El[i] = El[i] + En*(cnl_[l][n]*M1o1n[i] - c_i*dnl_[l][n]*N1e1n[i]
  1225. + c_i*anl_[l][n]*N3e1n[i] - bnl_[l][n]*M3o1n[i]);
  1226. Hl[i] = Hl[i] + En*(-dnl_[l][n]*M1e1n[i] - c_i*cnl_[l][n]*N1o1n[i]
  1227. + c_i*bnl_[l][n]*N3o1n[i] + anl_[l][n]*M3e1n[i]);
  1228. }
  1229. } // end of for all n
  1230. // Debug field calculation outside the particle
  1231. if (l == size_param_.size()) {
  1232. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1233. // basis unit vectors = er, etheta, ephi
  1234. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  1235. {
  1236. using std::sin;
  1237. using std::cos;
  1238. Eic[0] = eifac*sin(Theta)*cos(Phi);
  1239. Eic[1] = eifac*cos(Theta)*cos(Phi);
  1240. Eic[2] = -eifac*sin(Phi);
  1241. }
  1242. printf("Rho = %g; Phi = %g; Theta = %g\n", Rho, Phi, Theta);
  1243. for (int i = 0; i < 3; i++) {
  1244. printf("Ei[%i] = %g, %g; Eic[%i] = %g, %g\n", i, Ei[i].real(), Ei[i].imag(), i, Eic[i].real(), Eic[i].imag());
  1245. }
  1246. }
  1247. // magnetic field
  1248. double hffact = 1.0/(cc_*mu_);
  1249. for (int i = 0; i < 3; i++) {
  1250. Hl[i] = hffact*Hl[i];
  1251. }
  1252. for (int i = 0; i < 3; i++) {
  1253. // electric field E [V m - 1] = EF*E0
  1254. E[i] = El[i];
  1255. H[i] = Hl[i];
  1256. }
  1257. } // end of MultiLayerMie::fieldInt(...)
  1258. //**********************************************************************************//
  1259. // This function calculates complex electric and magnetic field in the surroundings //
  1260. // and inside (TODO) the particle. //
  1261. // //
  1262. // Input parameters: //
  1263. // L: Number of layers //
  1264. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1265. // x: Array containing the size parameters of the layers [0..L-1] //
  1266. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1267. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1268. // calculations. Only use it if you know what you are doing, otherwise //
  1269. // set this parameter to 0 (zero) and the function will calculate it. //
  1270. // ncoord: Number of coordinate points //
  1271. // Coords: Array containing all coordinates where the complex electric and //
  1272. // magnetic fields will be calculated //
  1273. // //
  1274. // Output parameters: //
  1275. // E, H: Complex electric and magnetic field at the provided coordinates //
  1276. // //
  1277. // Return value: //
  1278. // Number of multipolar expansion terms used for the calculations //
  1279. //**********************************************************************************//
  1280. void MultiLayerMie::RunFieldCalculation() {
  1281. // Calculate external scattering coefficients an_ and bn_
  1282. ExtScattCoeffs();
  1283. // Calculate internal scattering coefficients anl_ and bnl_
  1284. IntScattCoeffs();
  1285. // for (int i = 0; i < nmax_; i++) {
  1286. // printf("a[%i] = %g, %g; b[%i] = %g, %g\n", i, an_[i].real(), an_[i].imag(), i, bn_[i].real(), bn_[i].imag());
  1287. // }
  1288. long total_points = coords_[0].size();
  1289. E_.resize(total_points);
  1290. H_.resize(total_points);
  1291. for (auto& f : E_) f.resize(3);
  1292. for (auto& f : H_) f.resize(3);
  1293. for (int point = 0; point < total_points; point++) {
  1294. const double& Xp = coords_[0][point];
  1295. const double& Yp = coords_[1][point];
  1296. const double& Zp = coords_[2][point];
  1297. // Convert to spherical coordinates
  1298. double Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1299. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1300. double Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1301. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1302. double Phi = (Xp != 0.0 || Yp != 0.0) ? std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1303. // Avoid convergence problems due to Rho too small
  1304. if (Rho < 1e-5) Rho = 1e-5;
  1305. //*******************************************************//
  1306. // external scattering field = incident + scattered //
  1307. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1308. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1309. //*******************************************************//
  1310. // This array contains the fields in spherical coordinates
  1311. std::vector<std::complex<double> > Es(3), Hs(3);
  1312. // Firstly the easiest case: the field outside the particle
  1313. if (Rho >= GetSizeParameter()) {
  1314. fieldInt(Rho, Phi, Theta, Es, Hs);
  1315. // fieldExt(Rho, Phi, Theta, Es, Hs);
  1316. //printf("\nFin E ext: %g,%g,%g Rho=%g\n", std::abs(Es[0]), std::abs(Es[1]),std::abs(Es[2]), Rho);
  1317. } else {
  1318. fieldInt(Rho, Phi, Theta, Es, Hs);
  1319. // printf("\nFin E int: %g,%g,%g Rho=%g\n", std::abs(Es[0]), std::abs(Es[1]),std::abs(Es[2]), Rho);
  1320. }
  1321. //Now, convert the fields back to cartesian coordinates
  1322. {
  1323. using std::sin;
  1324. using std::cos;
  1325. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1326. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1327. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1328. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1329. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1330. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1331. }
  1332. //printf("Cart E: %g,%g,%g Rho=%g\n", std::abs(Ex), std::abs(Ey),std::abs(Ez), Rho);
  1333. } // end of for all field coordinates
  1334. } // end of MultiLayerMie::RunFieldCalculation()
  1335. } // end of namespace nmie