nmie-impl.hpp 59 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191
  1. #ifndef SRC_NMIE_IMPL_HPP_
  2. #define SRC_NMIE_IMPL_HPP_
  3. //**********************************************************************************//
  4. // Copyright (C) 2009-2018 Ovidio Pena <ovidio@bytesfall.com> //
  5. // Copyright (C) 2013-2018 Konstantin Ladutenko <kostyfisik@gmail.com> //
  6. // //
  7. // This file is part of scattnlay //
  8. // //
  9. // This program is free software: you can redistribute it and/or modify //
  10. // it under the terms of the GNU General Public License as published by //
  11. // the Free Software Foundation, either version 3 of the License, or //
  12. // (at your option) any later version. //
  13. // //
  14. // This program is distributed in the hope that it will be useful, //
  15. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  16. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  17. // GNU General Public License for more details. //
  18. // //
  19. // The only additional remark is that we expect that all publications //
  20. // describing work using this software, or all commercial products //
  21. // using it, cite at least one of the following references: //
  22. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  23. // a multilayered sphere," Computer Physics Communications, //
  24. // vol. 180, Nov. 2009, pp. 2348-2354. //
  25. // [2] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie //
  26. // calculation of electromagnetic near-field for a multilayered //
  27. // sphere," Computer Physics Communications, vol. 214, May 2017, //
  28. // pp. 225-230. //
  29. // //
  30. // You should have received a copy of the GNU General Public License //
  31. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  32. //**********************************************************************************//
  33. //**********************************************************************************//
  34. // This class implements the algorithm for a multilayered sphere described by: //
  35. // [1] W. Yang, "Improved recursive algorithm for light scattering by a //
  36. // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. //
  37. // //
  38. // You can find the description of all the used equations in: //
  39. // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  40. // a multilayered sphere," Computer Physics Communications, //
  41. // vol. 180, Nov. 2009, pp. 2348-2354. //
  42. // [3] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie //
  43. // calculation of electromagnetic near-field for a multilayered //
  44. // sphere," Computer Physics Communications, vol. 214, May 2017, //
  45. // pp. 225-230. //
  46. // //
  47. // Hereinafter all equations numbers refer to [2] //
  48. //**********************************************************************************//
  49. #include <iostream>
  50. #include <iomanip>
  51. #include <stdexcept>
  52. #include <vector>
  53. #include "special-functions-impl.hpp"
  54. #include "nmie.hpp"
  55. #include "nmie-precision.hpp"
  56. namespace nmie {
  57. //class implementation
  58. // ********************************************************************** //
  59. // Returns previously calculated Qext //
  60. // ********************************************************************** //
  61. template <typename FloatType>
  62. FloatType MultiLayerMie<FloatType>::GetQext() {
  63. if (!isMieCalculated_)
  64. throw std::invalid_argument("You should run calculations before result request!");
  65. return Qext_;
  66. }
  67. // ********************************************************************** //
  68. // Returns previously calculated Qabs //
  69. // ********************************************************************** //
  70. template <typename FloatType>
  71. FloatType MultiLayerMie<FloatType>::GetQabs() {
  72. if (!isMieCalculated_)
  73. throw std::invalid_argument("You should run calculations before result request!");
  74. return Qabs_;
  75. }
  76. // ********************************************************************** //
  77. // Returns previously calculated Qsca //
  78. // ********************************************************************** //
  79. template <typename FloatType>
  80. FloatType MultiLayerMie<FloatType>::GetQsca() {
  81. if (!isMieCalculated_)
  82. throw std::invalid_argument("You should run calculations before result request!");
  83. return Qsca_;
  84. }
  85. // ********************************************************************** //
  86. // Returns previously calculated Qbk //
  87. // ********************************************************************** //
  88. template <typename FloatType>
  89. FloatType MultiLayerMie<FloatType>::GetQbk() {
  90. if (!isMieCalculated_)
  91. throw std::invalid_argument("You should run calculations before result request!");
  92. return Qbk_;
  93. }
  94. // ********************************************************************** //
  95. // Returns previously calculated Qpr //
  96. // ********************************************************************** //
  97. template <typename FloatType>
  98. FloatType MultiLayerMie<FloatType>::GetQpr() {
  99. if (!isMieCalculated_)
  100. throw std::invalid_argument("You should run calculations before result request!");
  101. return Qpr_;
  102. }
  103. // ********************************************************************** //
  104. // Returns previously calculated assymetry factor //
  105. // ********************************************************************** //
  106. template <typename FloatType>
  107. FloatType MultiLayerMie<FloatType>::GetAsymmetryFactor() {
  108. if (!isMieCalculated_)
  109. throw std::invalid_argument("You should run calculations before result request!");
  110. return asymmetry_factor_;
  111. }
  112. // ********************************************************************** //
  113. // Returns previously calculated Albedo //
  114. // ********************************************************************** //
  115. template <typename FloatType>
  116. FloatType MultiLayerMie<FloatType>::GetAlbedo() {
  117. if (!isMieCalculated_)
  118. throw std::invalid_argument("You should run calculations before result request!");
  119. return albedo_;
  120. }
  121. // ********************************************************************** //
  122. // Returns previously calculated S1 //
  123. // ********************************************************************** //
  124. template <typename FloatType>
  125. std::vector<std::complex<FloatType> > MultiLayerMie<FloatType>::GetS1() {
  126. if (!isMieCalculated_)
  127. throw std::invalid_argument("You should run calculations before result request!");
  128. return S1_;
  129. }
  130. // ********************************************************************** //
  131. // Returns previously calculated S2 //
  132. // ********************************************************************** //
  133. template <typename FloatType>
  134. std::vector<std::complex<FloatType> > MultiLayerMie<FloatType>::GetS2() {
  135. if (!isMieCalculated_)
  136. throw std::invalid_argument("You should run calculations before result request!");
  137. return S2_;
  138. }
  139. // ********************************************************************** //
  140. // Modify scattering (theta) angles //
  141. // ********************************************************************** //
  142. template <typename FloatType>
  143. void MultiLayerMie<FloatType>::SetAngles(const std::vector<FloatType>& angles) {
  144. MarkUncalculated();
  145. theta_ = angles;
  146. }
  147. // ********************************************************************** //
  148. // Modify size of all layers //
  149. // ********************************************************************** //
  150. template <typename FloatType>
  151. void MultiLayerMie<FloatType>::SetLayersSize(const std::vector<FloatType>& layer_size) {
  152. MarkUncalculated();
  153. size_param_.clear();
  154. FloatType prev_layer_size = 0.0;
  155. for (auto curr_layer_size : layer_size) {
  156. if (curr_layer_size <= 0.0)
  157. throw std::invalid_argument("Size parameter should be positive!");
  158. if (prev_layer_size > curr_layer_size)
  159. throw std::invalid_argument
  160. ("Size parameter for next layer should be larger than the previous one!");
  161. prev_layer_size = curr_layer_size;
  162. size_param_.push_back(curr_layer_size);
  163. }
  164. }
  165. // ********************************************************************** //
  166. // Modify refractive index of all layers //
  167. // ********************************************************************** //
  168. template <typename FloatType>
  169. void MultiLayerMie<FloatType>::SetLayersIndex(const std::vector< std::complex<FloatType> >& index) {
  170. MarkUncalculated();
  171. refractive_index_ = index;
  172. }
  173. // ********************************************************************** //
  174. // Modify coordinates for field calculation //
  175. // ********************************************************************** //
  176. template <typename FloatType>
  177. void MultiLayerMie<FloatType>::SetFieldCoords(const std::vector< std::vector<FloatType> >& coords) {
  178. if (coords.size() != 3)
  179. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  180. if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
  181. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  182. coords_ = coords;
  183. }
  184. // ********************************************************************** //
  185. // Modify index of PEC layer //
  186. // ********************************************************************** //
  187. template <typename FloatType>
  188. void MultiLayerMie<FloatType>::SetPECLayer(int layer_position) {
  189. MarkUncalculated();
  190. if (layer_position < 0 && layer_position != -1)
  191. throw std::invalid_argument("Error! Layers are numbered from 0!");
  192. PEC_layer_position_ = layer_position;
  193. }
  194. // ********************************************************************** //
  195. // Set maximun number of terms to be used //
  196. // ********************************************************************** //
  197. template <typename FloatType>
  198. void MultiLayerMie<FloatType>::SetMaxTerms(int nmax) {
  199. MarkUncalculated();
  200. nmax_preset_ = nmax;
  201. }
  202. // ********************************************************************** //
  203. // Get total size parameter of particle //
  204. // ********************************************************************** //
  205. template <typename FloatType>
  206. FloatType MultiLayerMie<FloatType>::GetSizeParameter() {
  207. if (size_param_.size() > 0)
  208. return size_param_.back();
  209. else
  210. return 0;
  211. }
  212. // ********************************************************************** //
  213. // Mark uncalculated //
  214. // ********************************************************************** //
  215. template <typename FloatType>
  216. void MultiLayerMie<FloatType>::MarkUncalculated() {
  217. isExpCoeffsCalc_ = false;
  218. isScaCoeffsCalc_ = false;
  219. isMieCalculated_ = false;
  220. }
  221. // ********************************************************************** //
  222. // Clear layer information //
  223. // ********************************************************************** //
  224. template <typename FloatType>
  225. void MultiLayerMie<FloatType>::ClearLayers() {
  226. MarkUncalculated();
  227. size_param_.clear();
  228. refractive_index_.clear();
  229. }
  230. // ********************************************************************** //
  231. // ********************************************************************** //
  232. // ********************************************************************** //
  233. // Computational core
  234. // ********************************************************************** //
  235. // ********************************************************************** //
  236. // ********************************************************************** //
  237. int LeRu_cutoff(std::complex<double> z) {
  238. auto x = std::abs(z);
  239. return std::round(x + 11 * std::pow(x, (1.0 / 3.0)) + 1);
  240. // return 10000;
  241. }
  242. // ********************************************************************** //
  243. // Calculate calcNstop - equation (17) //
  244. // ********************************************************************** //
  245. template <typename FloatType>
  246. void MultiLayerMie<FloatType>::calcNstop() {
  247. //Wiscombe
  248. const FloatType& xL = size_param_.back();
  249. if (xL <= 8) {
  250. nmax_ = newround(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  251. } else if (xL <= 4200) {
  252. nmax_ = newround(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  253. } else {
  254. nmax_ = newround(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  255. }
  256. //Le Ru
  257. auto Nstop = nmie::LeRu_cutoff(static_cast<double>(xL))+1;
  258. if (Nstop > nmax_) nmax_ = Nstop;
  259. }
  260. // ********************************************************************** //
  261. // Maximum number of terms required for the calculation //
  262. // ********************************************************************** //
  263. template <typename FloatType>
  264. void MultiLayerMie<FloatType>::calcNmax(unsigned int first_layer) {
  265. int ri, riM1;
  266. const std::vector<FloatType>& x = size_param_;
  267. const std::vector<std::complex<FloatType> >& m = refractive_index_;
  268. calcNstop(); // Set initial nmax_ value
  269. for (unsigned int i = first_layer; i < x.size(); i++) {
  270. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  271. ri = newround(cabs(x[i]*m[i]));
  272. else
  273. ri = 0;
  274. nmax_ = std::max(nmax_, ri);
  275. // first layer is pec, if pec is present
  276. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  277. riM1 = newround(cabs(x[i - 1]* m[i]));
  278. else
  279. riM1 = 0;
  280. nmax_ = std::max(nmax_, riM1);
  281. }
  282. nmax_ += 100; // Final nmax_ value
  283. #ifdef MULTI_PRECISION
  284. nmax_ += MULTI_PRECISION; //TODO we may need to use more terms that this for MP computations.
  285. #endif
  286. // nmax_ *= nmax_;
  287. // printf("using nmax %i\n", nmax_);
  288. }
  289. // ********************************************************************** //
  290. // Calculate an - equation (5) //
  291. // ********************************************************************** //
  292. template <typename FloatType>
  293. std::complex<FloatType> MultiLayerMie<FloatType>::
  294. calc_an(int n, FloatType XL, std::complex<FloatType> Ha, std::complex<FloatType> mL,
  295. std::complex<FloatType> PsiXL, std::complex<FloatType> ZetaXL,
  296. std::complex<FloatType> PsiXLM1, std::complex<FloatType> ZetaXLM1) {
  297. std::complex<FloatType> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  298. std::complex<FloatType> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  299. // std::cout<< std::setprecision(100)
  300. // << "Ql " << PsiXL
  301. // << std::endl;
  302. return Num/Denom;
  303. }
  304. // ********************************************************************** //
  305. // Calculate bn - equation (6) //
  306. // ********************************************************************** //
  307. template <typename FloatType>
  308. std::complex<FloatType> MultiLayerMie<FloatType>::calc_bn(int n, FloatType XL, std::complex<FloatType> Hb, std::complex<FloatType> mL,
  309. std::complex<FloatType> PsiXL, std::complex<FloatType> ZetaXL,
  310. std::complex<FloatType> PsiXLM1, std::complex<FloatType> ZetaXLM1) {
  311. std::complex<FloatType> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  312. std::complex<FloatType> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  313. return Num/Denom;
  314. }
  315. // ********************************************************************** //
  316. // Calculates S1 - equation (25a) //
  317. // ********************************************************************** //
  318. template <typename FloatType>
  319. std::complex<FloatType> MultiLayerMie<FloatType>::calc_S1(int n, std::complex<FloatType> an, std::complex<FloatType> bn,
  320. FloatType Pi, FloatType Tau) {
  321. return FloatType(n + n + 1)*(Pi*an + Tau*bn)/FloatType(n*n + n);
  322. }
  323. // ********************************************************************** //
  324. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  325. // Pi and Tau) //
  326. // ********************************************************************** //
  327. template <typename FloatType>
  328. std::complex<FloatType> MultiLayerMie<FloatType>::calc_S2(int n, std::complex<FloatType> an, std::complex<FloatType> bn,
  329. FloatType Pi, FloatType Tau) {
  330. return calc_S1(n, an, bn, Tau, Pi);
  331. }
  332. //**********************************************************************************//
  333. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  334. // functions (D1 and D3) for a complex argument (z). //
  335. // Equations (16a), (16b) and (18a) - (18d) //
  336. // //
  337. // Input parameters: //
  338. // z: Complex argument to evaluate D1 and D3 //
  339. // nmax_: Maximum number of terms to calculate D1 and D3 //
  340. // //
  341. // Output parameters: //
  342. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  343. //**********************************************************************************//
  344. template <typename FloatType>
  345. void MultiLayerMie<FloatType>::calcD1D3(const std::complex<FloatType> z,
  346. std::vector<std::complex<FloatType> >& D1,
  347. std::vector<std::complex<FloatType> >& D3) {
  348. std::vector<std::complex<FloatType> > PsiZeta(nmax_+1);
  349. evalDownwardD1(z, D1);
  350. evalUpwardD3 (z, D1, D3, PsiZeta);
  351. }
  352. //**********************************************************************************//
  353. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  354. // complex argument (z). //
  355. // Equations (20a) - (21b) //
  356. // //
  357. // Input parameters: //
  358. // z: Complex argument to evaluate Psi and Zeta //
  359. // nmax: Maximum number of terms to calculate Psi and Zeta //
  360. // //
  361. // Output parameters: //
  362. // Psi, Zeta: Riccati-Bessel functions //
  363. //**********************************************************************************//
  364. template <typename FloatType>
  365. void MultiLayerMie<FloatType>::calcPsiZeta(std::complex<FloatType> z,
  366. std::vector<std::complex<FloatType> >& Psi,
  367. std::vector<std::complex<FloatType> >& Zeta) {
  368. std::vector<std::complex<FloatType> > D1(nmax_ + 1), D3(nmax_ + 1),
  369. PsiZeta(nmax_+1);
  370. // First, calculate the logarithmic derivatives
  371. evalDownwardD1(z, D1);
  372. // Now, use the upward recurrence to calculate Psi equations (20ab)
  373. evalUpwardPsi(z, D1, Psi);
  374. // Now, use the upward recurrence to calculate Psi*Zeta equations (18ad)
  375. evalUpwardD3 (z, D1, D3, PsiZeta);
  376. for (unsigned int i = 0; i < Zeta.size(); i++) {
  377. Zeta[i] = PsiZeta[i]/Psi[i];
  378. }
  379. // evalUpwardZeta(z, D3, Zeta);
  380. }
  381. //**********************************************************************************//
  382. // This function calculates Pi and Tau for a given value of cos(Theta). //
  383. // Equations (26a) - (26c) //
  384. // //
  385. // Input parameters: //
  386. // nmax_: Maximum number of terms to calculate Pi and Tau //
  387. // nTheta: Number of scattering angles //
  388. // Theta: Array containing all the scattering angles where the scattering //
  389. // amplitudes will be calculated //
  390. // //
  391. // Output parameters: //
  392. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  393. //**********************************************************************************//
  394. template <typename FloatType>
  395. void MultiLayerMie<FloatType>::calcPiTau(const FloatType& costheta,
  396. std::vector<FloatType>& Pi, std::vector<FloatType>& Tau) {
  397. int i;
  398. //****************************************************//
  399. // Equations (26a) - (26c) //
  400. //****************************************************//
  401. // Initialize Pi and Tau
  402. Pi[0] = 1.0; // n=1
  403. Tau[0] = costheta;
  404. // Calculate the actual values
  405. if (nmax_ > 1) {
  406. Pi[1] = 3*costheta*Pi[0]; //n=2
  407. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  408. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  409. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  410. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  411. }
  412. }
  413. } // end of MultiLayerMie::calcPiTau(...)
  414. //**********************************************************************************//
  415. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  416. // required to calculate the near-field parameters. //
  417. // //
  418. // Input parameters: //
  419. // Rho: Radial distance //
  420. // Phi: Azimuthal angle //
  421. // Theta: Polar angle //
  422. // rn: Either the spherical Ricatti-Bessel function of first or third kind //
  423. // Dn: Logarithmic derivative of rn //
  424. // Pi, Tau: Angular functions Pi and Tau //
  425. // n: Order of vector spherical harmonics //
  426. // //
  427. // Output parameters: //
  428. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  429. //**********************************************************************************//
  430. template <typename FloatType>
  431. void MultiLayerMie<FloatType>::calcSpherHarm(const std::complex<FloatType> Rho, const FloatType Theta, const FloatType Phi,
  432. const std::complex<FloatType>& rn, const std::complex<FloatType>& Dn,
  433. const FloatType& Pi, const FloatType& Tau, const FloatType& n,
  434. std::vector<std::complex<FloatType> >& Mo1n, std::vector<std::complex<FloatType> >& Me1n,
  435. std::vector<std::complex<FloatType> >& No1n, std::vector<std::complex<FloatType> >& Ne1n) {
  436. // using eq 4.50 in BH
  437. std::complex<FloatType> c_zero(0.0, 0.0);
  438. using nmm::sin;
  439. using nmm::cos;
  440. Mo1n[0] = c_zero;
  441. Mo1n[1] = cos(Phi)*Pi*rn/Rho;
  442. Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
  443. Me1n[0] = c_zero;
  444. Me1n[1] = -sin(Phi)*Pi*rn/Rho;
  445. Me1n[2] = -cos(Phi)*Tau*rn/Rho;
  446. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  447. No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
  448. No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
  449. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  450. Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
  451. Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
  452. } // end of MultiLayerMie::calcSpherHarm(...)
  453. //**********************************************************************************//
  454. // This function calculates the scattering coefficients required to calculate //
  455. // both the near- and far-field parameters. //
  456. // //
  457. // Input parameters: //
  458. // L: Number of layers //
  459. // pl: Index of PEC layer. If there is none just send -1 //
  460. // x: Array containing the size parameters of the layers [0..L-1] //
  461. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  462. // nmax: Maximum number of multipolar expansion terms to be used for the //
  463. // calculations. Only use it if you know what you are doing, otherwise //
  464. // set this parameter to -1 and the function will calculate it. //
  465. // //
  466. // Output parameters: //
  467. // an, bn: Complex scattering amplitudes //
  468. // //
  469. // Return value: //
  470. // Number of multipolar expansion terms used for the calculations //
  471. //**********************************************************************************//
  472. template <typename FloatType>
  473. void MultiLayerMie<FloatType>::calcScattCoeffs() {
  474. isScaCoeffsCalc_ = false;
  475. const std::vector<FloatType>& x = size_param_;
  476. const std::vector<std::complex<FloatType> >& m = refractive_index_;
  477. const int& pl = PEC_layer_position_;
  478. const int L = refractive_index_.size();
  479. //************************************************************************//
  480. // Calculate the index of the first layer. It can be either 0 (default) //
  481. // or the index of the outermost PEC layer. In the latter case all layers //
  482. // below the PEC are discarded. //
  483. // ***********************************************************************//
  484. int fl = (pl > 0) ? pl : 0;
  485. if (nmax_preset_ <= 0) calcNmax(fl);
  486. else nmax_ = nmax_preset_;
  487. std::complex<FloatType> z1, z2;
  488. //**************************************************************************//
  489. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  490. // means that index = layer number - 1 or index = n - 1. The only exception //
  491. // are the arrays for representing D1, D3 and Q because they need a value //
  492. // for the index 0 (zero), hence it is important to consider this shift //
  493. // between different arrays. The change was done to optimize memory usage. //
  494. //**************************************************************************//
  495. // Allocate memory to the arrays
  496. std::vector<std::complex<FloatType> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  497. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  498. std::vector<std::vector<std::complex<FloatType> > > Q(L), Ha(L), Hb(L);
  499. for (int l = 0; l < L; l++) {
  500. Q[l].resize(nmax_ + 1);
  501. Ha[l].resize(nmax_);
  502. Hb[l].resize(nmax_);
  503. }
  504. an_.resize(nmax_);
  505. bn_.resize(nmax_);
  506. std::vector<std::complex<FloatType> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  507. //*************************************************//
  508. // Calculate D1 and D3 for z1 in the first layer //
  509. //*************************************************//
  510. if (fl == pl) { // PEC layer
  511. for (int n = 0; n <= nmax_; n++) {
  512. D1_mlxl[n] = std::complex<FloatType>(0.0, - 1.0);
  513. D3_mlxl[n] = std::complex<FloatType>(0.0, 1.0);
  514. }
  515. } else { // Regular layer
  516. z1 = x[fl]* m[fl];
  517. // Calculate D1 and D3
  518. calcD1D3(z1, D1_mlxl, D3_mlxl);
  519. }
  520. //******************************************************************//
  521. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  522. //******************************************************************//
  523. for (int n = 0; n < nmax_; n++) {
  524. Ha[fl][n] = D1_mlxl[n + 1];
  525. Hb[fl][n] = D1_mlxl[n + 1];
  526. }
  527. //*****************************************************//
  528. // Iteration from the second layer to the last one (L) //
  529. //*****************************************************//
  530. std::complex<FloatType> Temp, Num, Denom;
  531. std::complex<FloatType> G1, G2;
  532. for (int l = fl + 1; l < L; l++) {
  533. //************************************************************//
  534. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  535. //************************************************************//
  536. z1 = x[l]*m[l];
  537. z2 = x[l - 1]*m[l];
  538. //Calculate D1 and D3 for z1
  539. calcD1D3(z1, D1_mlxl, D3_mlxl);
  540. //Calculate D1 and D3 for z2
  541. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  542. //*************************************************//
  543. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  544. //*************************************************//
  545. // Upward recurrence for Q - equations (19a) and (19b)
  546. Num = std::complex<FloatType>(nmm::exp(-2.0*(z1.imag() - z2.imag())), 0.0)
  547. *std::complex<FloatType>(nmm::cos(-2.0*z2.real()) - nmm::exp(-2.0*z2.imag()), nmm::sin(-2.0*z2.real()));
  548. Denom = std::complex<FloatType>(nmm::cos(-2.0*z1.real()) - nmm::exp(-2.0*z1.imag()), nmm::sin(-2.0*z1.real()));
  549. Q[l][0] = Num/Denom;
  550. for (int n = 1; n <= nmax_; n++) {
  551. Num = (z1*D1_mlxl[n] + FloatType(n))*(FloatType(n) - z1*D3_mlxl[n - 1]);
  552. Denom = (z2*D1_mlxlM1[n] + FloatType(n))*(FloatType(n) - z2*D3_mlxlM1[n - 1]);
  553. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  554. }
  555. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  556. for (int n = 1; n <= nmax_; n++) {
  557. //Ha
  558. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  559. G1 = -D1_mlxlM1[n];
  560. G2 = -D3_mlxlM1[n];
  561. } else {
  562. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  563. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  564. } // end of if PEC
  565. Temp = Q[l][n]*G1;
  566. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  567. Denom = G2 - Temp;
  568. Ha[l][n - 1] = Num/Denom;
  569. //Hb
  570. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  571. G1 = Hb[l - 1][n - 1];
  572. G2 = Hb[l - 1][n - 1];
  573. } else {
  574. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  575. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  576. } // end of if PEC
  577. Temp = Q[l][n]*G1;
  578. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  579. Denom = (G2- Temp);
  580. Hb[l][n - 1] = (Num/ Denom);
  581. } // end of for Ha and Hb terms
  582. } // end of for layers iteration
  583. //**************************************//
  584. //Calculate Psi and Zeta for XL //
  585. //**************************************//
  586. // Calculate PsiXL and ZetaXL
  587. calcPsiZeta(std::complex<FloatType>(x[L - 1],0.0), PsiXL, ZetaXL);
  588. //*********************************************************************//
  589. // Finally, we calculate the scattering coefficients (an and bn) and //
  590. // the angular functions (Pi and Tau). Note that for these arrays the //
  591. // first layer is 0 (zero), in future versions all arrays will follow //
  592. // this convention to save memory. (13 Nov, 2014) //
  593. //*********************************************************************//
  594. for (int n = 0; n < nmax_; n++) {
  595. //********************************************************************//
  596. //Expressions for calculating an and bn coefficients are not valid if //
  597. //there is only one PEC layer (ie, for a simple PEC sphere). //
  598. //********************************************************************//
  599. if (pl < (L - 1)) {
  600. 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]);
  601. 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]);
  602. } else {
  603. an_[n] = calc_an(n + 1, x[L - 1], std::complex<FloatType>(0.0, 0.0), std::complex<FloatType>(1.0, 0.0), PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  604. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  605. }
  606. if (nmm::isnan(an_[n].real()) || nmm::isnan(an_[n].imag()) ||
  607. nmm::isnan(bn_[n].real()) || nmm::isnan(bn_[n].imag())
  608. ) {
  609. // TODO somehow notify Python users about it
  610. std::cout << "nmax value was chaned due to unexpected error. New values is "<< n
  611. << " (was "<<nmax_<<")"<<std::endl;
  612. nmax_ = n;
  613. break;
  614. }
  615. } // end of for an and bn terms
  616. isScaCoeffsCalc_ = true;
  617. } // end of MultiLayerMie::calcScattCoeffs()
  618. //**********************************************************************************//
  619. // This function calculates the actual scattering parameters and amplitudes //
  620. // //
  621. // Input parameters: //
  622. // L: Number of layers //
  623. // pl: Index of PEC layer. If there is none just send -1 //
  624. // x: Array containing the size parameters of the layers [0..L-1] //
  625. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  626. // nTheta: Number of scattering angles //
  627. // Theta: Array containing all the scattering angles where the scattering //
  628. // amplitudes will be calculated //
  629. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  630. // calculations. Only use it if you know what you are doing, otherwise //
  631. // set this parameter to -1 and the function will calculate it //
  632. // //
  633. // Output parameters: //
  634. // Qext: Efficiency factor for extinction //
  635. // Qsca: Efficiency factor for scattering //
  636. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  637. // Qbk: Efficiency factor for backscattering //
  638. // Qpr: Efficiency factor for the radiation pressure //
  639. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  640. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  641. // S1, S2: Complex scattering amplitudes //
  642. // //
  643. // Return value: //
  644. // Number of multipolar expansion terms used for the calculations //
  645. //**********************************************************************************//
  646. template <typename FloatType>
  647. void MultiLayerMie<FloatType>::RunMieCalculation() {
  648. if (size_param_.size() != refractive_index_.size())
  649. throw std::invalid_argument("Each size parameter should have only one index!");
  650. if (size_param_.size() == 0)
  651. throw std::invalid_argument("Initialize model first!");
  652. const std::vector<FloatType>& x = size_param_;
  653. //MarkUncalculated();
  654. // Calculate scattering coefficients
  655. if (!isScaCoeffsCalc_) calcScattCoeffs();
  656. // Initialize the scattering parameters
  657. Qext_ = 0.0;
  658. Qsca_ = 0.0;
  659. Qabs_ = 0.0;
  660. Qbk_ = 0.0;
  661. Qpr_ = 0.0;
  662. asymmetry_factor_ = 0.0;
  663. albedo_ = 0.0;
  664. // Initialize the scattering amplitudes
  665. std::vector<std::complex<FloatType> > tmp1(theta_.size(),std::complex<FloatType>(0.0, 0.0));
  666. S1_.swap(tmp1);
  667. S2_ = S1_;
  668. // Precalculate cos(theta) - gives about 5% speed up.
  669. std::vector<FloatType> costheta(theta_.size(), 0.0);
  670. for (unsigned int t = 0; t < theta_.size(); t++) {
  671. costheta[t] = nmm::cos(theta_[t]);
  672. }
  673. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  674. std::complex<FloatType> Qbktmp(0.0, 0.0);
  675. std::vector< std::complex<FloatType> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  676. // By using downward recurrence we avoid loss of precision due to float rounding errors
  677. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  678. // http://en.wikipedia.org/wiki/Loss_of_significance
  679. for (int n = nmax_ - 2; n >= 0; n--) {
  680. // for (int n = 0; n < nmax_; n++) {
  681. const int n1 = n + 1;
  682. if (mode_n_ == Modes::kAll) {
  683. // Equation (27)
  684. Qext_ += (n1 + n1 + 1.0) * (an_[n].real() + bn_[n].real());
  685. // Equation (28)
  686. Qsca_ += (n1 + n1 + 1.0) * (an_[n].real() * an_[n].real() + an_[n].imag() * an_[n].imag()
  687. + bn_[n].real() * bn_[n].real() + bn_[n].imag() * bn_[n].imag());
  688. // std::cout<<"n ="<< n1 << " ext:"<<Qext_ <<" sca:"<<Qsca_<<std::endl;
  689. // Equation (29)
  690. Qpr_ += ((n1 * (n1 + 2.0) / (n1 + 1.0)) * ((an_[n] * std::conj(an_[n1]) + bn_[n] * std::conj(bn_[n1])).real())
  691. + ((n1 + n1 + 1.0) / (n1 * (n1 + 1.0))) * (an_[n] * std::conj(bn_[n])).real());
  692. // Equation (33)
  693. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (an_[n] - bn_[n]);
  694. // Calculate the scattering amplitudes (S1 and S2) Equations (25a) - (25b)
  695. for (unsigned int t = 0; t < theta_.size(); t++) {
  696. calcPiTau(costheta[t], Pi, Tau);
  697. S1_[t] += calc_S1(n1, an_[n], bn_[n], Pi[n], Tau[n]);
  698. S2_[t] += calc_S2(n1, an_[n], bn_[n], Pi[n], Tau[n]);
  699. }
  700. continue;
  701. }
  702. if (n1 == mode_n_) {
  703. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  704. Qext_ += (n1 + n1 + 1.0) * (an_[n].real());
  705. Qsca_ += (n1 + n1 + 1.0) * (an_[n].real() * an_[n].real() + an_[n].imag() * an_[n].imag());
  706. Qpr_ += std::nan("");
  707. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (an_[n]);
  708. for (unsigned int t = 0; t < theta_.size(); t++) {
  709. calcPiTau(costheta[t], Pi, Tau);
  710. S1_[t] += calc_S1(n1, an_[n], static_cast<std::complex<FloatType>>(0), Pi[n], Tau[n]);
  711. S2_[t] += calc_S2(n1, an_[n], static_cast<std::complex<FloatType>>(0), Pi[n], Tau[n]);
  712. }
  713. }
  714. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  715. Qext_ += (n1 + n1 + 1.0) * (bn_[n].real());
  716. Qsca_ += (n1 + n1 + 1.0) * (bn_[n].real() * bn_[n].real() + bn_[n].imag() * bn_[n].imag());
  717. Qpr_ += std::nan("");
  718. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (bn_[n]);
  719. for (unsigned int t = 0; t < theta_.size(); t++) {
  720. calcPiTau(costheta[t], Pi, Tau);
  721. S1_[t] += calc_S1(n1, static_cast<std::complex<FloatType>>(0), bn_[n], Pi[n], Tau[n]);
  722. S2_[t] += calc_S2(n1, static_cast<std::complex<FloatType>>(0), bn_[n], Pi[n], Tau[n]);
  723. }
  724. }
  725. }
  726. }
  727. FloatType x2 = pow2(x.back());
  728. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  729. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  730. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  731. Qabs_ = Qext_ - Qsca_; // Equation (30)
  732. albedo_ = Qsca_/Qext_; // Equation (31)
  733. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  734. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  735. isMieCalculated_ = true;
  736. }
  737. //**********************************************************************************//
  738. // This function calculates the expansion coefficients inside the particle, //
  739. // required to calculate the near-field parameters. //
  740. // //
  741. // Input parameters: //
  742. // L: Number of layers //
  743. // pl: Index of PEC layer. If there is none just send -1 //
  744. // x: Array containing the size parameters of the layers [0..L-1] //
  745. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  746. // nmax: Maximum number of multipolar expansion terms to be used for the //
  747. // calculations. Only use it if you know what you are doing, otherwise //
  748. // set this parameter to -1 and the function will calculate it. //
  749. // //
  750. // Output parameters: //
  751. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  752. // //
  753. // Return value: //
  754. // Number of multipolar expansion terms used for the calculations //
  755. //**********************************************************************************//
  756. template <typename FloatType>
  757. void MultiLayerMie<FloatType>::calcExpanCoeffs() {
  758. if (!isScaCoeffsCalc_)
  759. throw std::invalid_argument("(calcExpanCoeffs) You should calculate external coefficients first!");
  760. isExpCoeffsCalc_ = false;
  761. std::complex<FloatType> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  762. const int L = refractive_index_.size();
  763. aln_.resize(L + 1);
  764. bln_.resize(L + 1);
  765. cln_.resize(L + 1);
  766. dln_.resize(L + 1);
  767. for (int l = 0; l <= L; l++) {
  768. aln_[l].resize(nmax_);
  769. bln_[l].resize(nmax_);
  770. cln_[l].resize(nmax_);
  771. dln_[l].resize(nmax_);
  772. }
  773. // Yang, paragraph under eq. A3
  774. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  775. for (int n = 0; n < nmax_; n++) {
  776. aln_[L][n] = an_[n];
  777. bln_[L][n] = bn_[n];
  778. cln_[L][n] = c_one;
  779. dln_[L][n] = c_one;
  780. }
  781. std::vector<std::complex<FloatType> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  782. std::vector<std::complex<FloatType> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  783. std::complex<FloatType> denomZeta, denomPsi, T1, T2, T3, T4;
  784. auto& m = refractive_index_;
  785. std::vector< std::complex<FloatType> > m1(L);
  786. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  787. m1[L - 1] = std::complex<FloatType> (1.0, 0.0);
  788. std::complex<FloatType> z, z1;
  789. for (int l = L - 1; l >= 0; l--) {
  790. if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!!
  791. for (int n = 0; n < nmax_; n++) {
  792. // aln
  793. aln_[l][n] = c_zero;
  794. // bln
  795. bln_[l][n] = c_zero;
  796. // cln
  797. cln_[l][n] = c_zero;
  798. // dln
  799. dln_[l][n] = c_zero;
  800. }
  801. } else { // Regular material, just do the calculation
  802. z = size_param_[l]*m[l];
  803. z1 = size_param_[l]*m1[l];
  804. calcD1D3(z, D1z, D3z);
  805. calcD1D3(z1, D1z1, D3z1);
  806. calcPsiZeta(z, Psiz, Zetaz);
  807. calcPsiZeta(z1, Psiz1, Zetaz1);
  808. for (int n = 0; n < nmax_; n++) {
  809. int n1 = n + 1;
  810. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  811. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  812. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  813. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  814. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  815. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  816. // aln
  817. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  818. // bln
  819. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  820. // cln
  821. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  822. // dln
  823. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  824. } // end of all n
  825. } // end PEC condition
  826. } // end of all l
  827. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  828. for (int n = 0; n < nmax_; ++n) {
  829. if (cabs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  830. else {
  831. //throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  832. std::cout<< std::setprecision(100)
  833. << "Warning: Potentially unstable calculation of aln[0]["
  834. << n << "] = "<< aln_[0][n] <<std::endl;
  835. aln_[0][n] = 0.0;
  836. }
  837. if (cabs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  838. else {
  839. //throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  840. std::cout<< std::setprecision(100)
  841. << "Warning: Potentially unstable calculation of bln[0]["
  842. << n << "] = "<< bln_[0][n] <<std::endl;
  843. bln_[0][n] = 0.0;
  844. }
  845. }
  846. isExpCoeffsCalc_ = true;
  847. } // end of void MultiLayerMie::calcExpanCoeffs()
  848. //**********************************************************************************//
  849. // This function calculates the electric (E) and magnetic (H) fields inside and //
  850. // around the particle. //
  851. // //
  852. // Input parameters (coordinates of the point): //
  853. // Rho: Radial distance //
  854. // Phi: Azimuthal angle //
  855. // Theta: Polar angle //
  856. // mode_n: mode order. //
  857. // -1 - use all modes (all_) //
  858. // 1 - use dipole mode only //
  859. // 2 - use quadrupole mode only //
  860. // ... //
  861. // mode_type: only used when mode_n != -1 //
  862. // 0 - electric only //
  863. // 1 - magnetic only //
  864. // //
  865. // //
  866. // Output parameters: //
  867. // E, H: Complex electric and magnetic fields //
  868. //**********************************************************************************//
  869. template <typename FloatType>
  870. void MultiLayerMie<FloatType>::calcFieldByComponents(const FloatType Rho,
  871. const FloatType Theta, const FloatType Phi,
  872. std::vector<std::complex<FloatType> >& E,
  873. std::vector<std::complex<FloatType> >& H) {
  874. std::complex<FloatType> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  875. // Vector containing precomputed integer powers of i to avoid computation
  876. std::vector<std::complex<FloatType> > ipow = {c_one, c_i, -c_one, -c_i};
  877. std::vector<std::complex<FloatType> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  878. std::vector<std::complex<FloatType> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  879. std::vector<std::complex<FloatType> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  880. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  881. int l = 0; // Layer number
  882. std::complex<FloatType> ml;
  883. // Initialize E and H
  884. for (int i = 0; i < 3; i++) {
  885. E[i] = c_zero;
  886. H[i] = c_zero;
  887. }
  888. if (Rho > size_param_.back()) {
  889. l = size_param_.size();
  890. ml = c_one;
  891. } else {
  892. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  893. if (Rho <= size_param_[i]) {
  894. l = i;
  895. }
  896. }
  897. ml = refractive_index_[l];
  898. }
  899. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  900. calcD1D3(Rho*ml, D1n, D3n);
  901. // Calculate Ricatti-Bessel functions
  902. calcPsiZeta(Rho*ml, Psi, Zeta);
  903. // Calculate angular functions Pi and Tau
  904. calcPiTau(nmm::cos(Theta), Pi, Tau);
  905. // for (int n = nmax_ - 2; n >= 0; n--) {
  906. for (int n = 0; n < nmax_-1; n++) {
  907. int n1 = n + 1;
  908. auto rn = static_cast<FloatType>(n1);
  909. // using BH 4.12 and 4.50
  910. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  911. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  912. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  913. std::complex<FloatType> En = ipow[n1 % 4]
  914. *static_cast<FloatType>((rn + rn + 1.0)/(rn*rn + rn));
  915. for (int i = 0; i < 3; i++) {
  916. auto Ediff = En*( cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  917. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  918. auto Hdiff = En*( -dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  919. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  920. if (nmm::isnan(Ediff.real()) || nmm::isnan(Ediff.imag()) ||
  921. nmm::isnan(Hdiff.real()) || nmm::isnan(Hdiff.imag())
  922. ) break;
  923. if (mode_n_ == Modes::kAll) {
  924. // electric field E [V m - 1] = EF*E0
  925. E[i] += Ediff;
  926. H[i] += Hdiff;
  927. continue;
  928. }
  929. if (n1 == mode_n_) {
  930. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  931. E[i] += En*( -c_i*dln_[l][n]*N1e1n[i]
  932. + c_i*aln_[l][n]*N3e1n[i]);
  933. H[i] += En*(-dln_[l][n]*M1e1n[i]
  934. +aln_[l][n]*M3e1n[i]);
  935. //std::cout << mode_n_;
  936. }
  937. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  938. E[i] += En*( cln_[l][n]*M1o1n[i]
  939. - bln_[l][n]*M3o1n[i]);
  940. H[i] += En*( -c_i*cln_[l][n]*N1o1n[i]
  941. + c_i*bln_[l][n]*N3o1n[i]);
  942. //std::cout << mode_n_;
  943. }
  944. //std::cout << std::endl;
  945. }
  946. //throw std::invalid_argument("Error! Unexpected mode for field evaluation!\n mode_n="+std::to_string(mode_n)+", mode_type="+std::to_string(mode_type)+"\n=====*****=====");
  947. }
  948. } // end of for all n
  949. // magnetic field
  950. std::complex<FloatType> hffact = ml/static_cast<FloatType>(cc_*mu_);
  951. for (int i = 0; i < 3; i++) {
  952. H[i] = hffact*H[i];
  953. }
  954. } // end of MultiLayerMie::calcFieldByComponents(...)
  955. //**********************************************************************************//
  956. // This function calculates complex electric and magnetic field in the surroundings //
  957. // and inside the particle. //
  958. // //
  959. // Input parameters: //
  960. // L: Number of layers //
  961. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  962. // x: Array containing the size parameters of the layers [0..L-1] //
  963. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  964. // nmax: Maximum number of multipolar expansion terms to be used for the //
  965. // calculations. Only use it if you know what you are doing, otherwise //
  966. // set this parameter to 0 (zero) and the function will calculate it. //
  967. // ncoord: Number of coordinate points //
  968. // Coords: Array containing all coordinates where the complex electric and //
  969. // magnetic fields will be calculated //
  970. // mode_n: mode order. //
  971. // -1 - use all modes (all_) //
  972. // 1 - use dipole mode only //
  973. // 2 - use quadrupole mode only //
  974. // ... //
  975. // mode_type: only used when mode_n != -1 //
  976. // 0 - electric only //
  977. // 1 - magnetic only //
  978. // //
  979. // Output parameters: //
  980. // E, H: Complex electric and magnetic field at the provided coordinates //
  981. // //
  982. // Return value: //
  983. // Number of multipolar expansion terms used for the calculations //
  984. //**********************************************************************************//
  985. template <typename FloatType>
  986. void MultiLayerMie<FloatType>::RunFieldCalculation() {
  987. FloatType Rho, Theta, Phi;
  988. // Calculate scattering coefficients an_ and bn_
  989. calcScattCoeffs();
  990. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  991. calcExpanCoeffs();
  992. long total_points = coords_[0].size();
  993. E_.resize(total_points);
  994. H_.resize(total_points);
  995. Es_.resize(total_points);
  996. Hs_.resize(total_points);
  997. for (auto& f : E_) f.resize(3);
  998. for (auto& f : H_) f.resize(3);
  999. for (auto& f : Es_) f.resize(3);
  1000. for (auto& f : Hs_) f.resize(3);
  1001. for (int point = 0; point < total_points; point++) {
  1002. const FloatType& Xp = coords_[0][point];
  1003. const FloatType& Yp = coords_[1][point];
  1004. const FloatType& Zp = coords_[2][point];
  1005. // Convert to spherical coordinates
  1006. Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1007. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1008. Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0;
  1009. // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc.
  1010. Phi = nmm::atan2(Yp,Xp);
  1011. // Avoid convergence problems due to Rho too small
  1012. if (Rho < 1e-5) Rho = 1e-5;
  1013. // std::cout << "Xp: "<<Xp<< " Yp: "<<Yp<< " Zp: "<<Zp<<std::endl;
  1014. // std::cout << " Rho: "<<Rho<<" Theta: "<<Theta<<" Phi:"<<Phi<<std::endl<<std::endl;
  1015. //*******************************************************//
  1016. // external scattering field = incident + scattered //
  1017. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1018. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1019. //*******************************************************//
  1020. // This array contains the fields in spherical coordinates
  1021. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  1022. // Do the actual calculation of electric and magnetic field
  1023. calcFieldByComponents(Rho, Theta, Phi, Es, Hs);
  1024. for (int sph_coord = 0; sph_coord<3; ++sph_coord) {
  1025. Es_[point][sph_coord] = Es[sph_coord];
  1026. Hs_[point][sph_coord] = Hs[sph_coord];
  1027. }
  1028. { //Now, convert the fields back to cartesian coordinates
  1029. using nmm::sin;
  1030. using nmm::cos;
  1031. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1032. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1033. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1034. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1035. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1036. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1037. }
  1038. } // end of for all field coordinates
  1039. } // end of MultiLayerMie::RunFieldCalculation()
  1040. } // end of namespace nmie
  1041. #endif // SRC_NMIE_IMPL_HPP_