nmie.cc 77 KB

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