nmie-impl.hpp 59 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192
  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. // ********************************************************************** //
  238. // Calculate calcNstop - equation (17) //
  239. // ********************************************************************** //
  240. template <typename FloatType>
  241. void MultiLayerMie<FloatType>::calcNstop() {
  242. const FloatType& xL = size_param_.back();
  243. if (xL <= 8) {
  244. nmax_ = newround(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  245. } else if (xL <= 4200) {
  246. nmax_ = newround(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  247. } else {
  248. nmax_ = newround(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  249. }
  250. }
  251. // ********************************************************************** //
  252. // Maximum number of terms required for the calculation //
  253. // ********************************************************************** //
  254. template <typename FloatType>
  255. void MultiLayerMie<FloatType>::calcNmax(unsigned int first_layer) {
  256. int ri, riM1;
  257. const std::vector<FloatType>& x = size_param_;
  258. const std::vector<std::complex<FloatType> >& m = refractive_index_;
  259. calcNstop(); // Set initial nmax_ value
  260. for (unsigned int i = first_layer; i < x.size(); i++) {
  261. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  262. ri = newround(cabs(x[i]*m[i]));
  263. else
  264. ri = 0;
  265. nmax_ = std::max(nmax_, ri);
  266. // first layer is pec, if pec is present
  267. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  268. riM1 = newround(cabs(x[i - 1]* m[i]));
  269. else
  270. riM1 = 0;
  271. nmax_ = std::max(nmax_, riM1);
  272. }
  273. nmax_ += 15; // Final nmax_ value
  274. // nmax_ *= nmax_;
  275. // printf("using nmax %i\n", nmax_);
  276. }
  277. // ********************************************************************** //
  278. // Calculate an - equation (5) //
  279. // ********************************************************************** //
  280. template <typename FloatType>
  281. std::complex<FloatType> MultiLayerMie<FloatType>::
  282. calc_an(int n, FloatType XL, std::complex<FloatType> Ha, std::complex<FloatType> mL,
  283. std::complex<FloatType> PsiXL, std::complex<FloatType> ZetaXL,
  284. std::complex<FloatType> PsiXLM1, std::complex<FloatType> ZetaXLM1) {
  285. std::complex<FloatType> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  286. std::complex<FloatType> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  287. // std::cout<< std::setprecision(100)
  288. // << "Ql " << PsiXL
  289. // << std::endl;
  290. return Num/Denom;
  291. }
  292. // ********************************************************************** //
  293. // Calculate bn - equation (6) //
  294. // ********************************************************************** //
  295. template <typename FloatType>
  296. std::complex<FloatType> MultiLayerMie<FloatType>::calc_bn(int n, FloatType XL, std::complex<FloatType> Hb, std::complex<FloatType> mL,
  297. std::complex<FloatType> PsiXL, std::complex<FloatType> ZetaXL,
  298. std::complex<FloatType> PsiXLM1, std::complex<FloatType> ZetaXLM1) {
  299. std::complex<FloatType> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  300. std::complex<FloatType> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  301. return Num/Denom;
  302. }
  303. // ********************************************************************** //
  304. // Calculates S1 - equation (25a) //
  305. // ********************************************************************** //
  306. template <typename FloatType>
  307. std::complex<FloatType> MultiLayerMie<FloatType>::calc_S1(int n, std::complex<FloatType> an, std::complex<FloatType> bn,
  308. FloatType Pi, FloatType Tau) {
  309. return FloatType(n + n + 1)*(Pi*an + Tau*bn)/FloatType(n*n + n);
  310. }
  311. // ********************************************************************** //
  312. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  313. // Pi and Tau) //
  314. // ********************************************************************** //
  315. template <typename FloatType>
  316. std::complex<FloatType> MultiLayerMie<FloatType>::calc_S2(int n, std::complex<FloatType> an, std::complex<FloatType> bn,
  317. FloatType Pi, FloatType Tau) {
  318. return calc_S1(n, an, bn, Tau, Pi);
  319. }
  320. //**********************************************************************************//
  321. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  322. // functions (D1 and D3) for a complex argument (z). //
  323. // Equations (16a), (16b) and (18a) - (18d) //
  324. // //
  325. // Input parameters: //
  326. // z: Complex argument to evaluate D1 and D3 //
  327. // nmax_: Maximum number of terms to calculate D1 and D3 //
  328. // //
  329. // Output parameters: //
  330. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  331. //**********************************************************************************//
  332. template <typename FloatType>
  333. void MultiLayerMie<FloatType>::calcD1D3(const std::complex<FloatType> z,
  334. std::vector<std::complex<FloatType> >& D1,
  335. std::vector<std::complex<FloatType> >& D3) {
  336. evalDownwardD1(z, D1);
  337. int lnmx = evalKapteynNumberOfLostSignificantDigits(nmax_, z);
  338. std::vector<std::complex<FloatType> > r;
  339. if (lnmx < 4) {
  340. r.resize(nmax_+1);
  341. evalForwardR(z, r);
  342. } else {
  343. int valid_digits = 6;
  344. int nstar = getNStar(nmax_, z, valid_digits);
  345. r.resize(nstar);
  346. evalBackwardR(z,r);
  347. }
  348. convertRtoD1(z, r, D1);
  349. // TODO: Do we need this check?
  350. // if (cabs(D1[0]) > 1.0e15) {
  351. // throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  352. // //printf("Warning: Potentially unstable D1! Please, try to change input parameters!\n");
  353. // }
  354. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  355. PsiZeta_[0] = static_cast<FloatType>(0.5)*(static_cast<FloatType>(1.0) - std::complex<FloatType>(nmm::cos(2.0*z.real()), nmm::sin(2.0*z.real()))
  356. *static_cast<FloatType>(nmm::exp(-2.0*z.imag())));
  357. D3[0] = std::complex<FloatType>(0.0, 1.0);
  358. const std::complex<FloatType> zinv = std::complex<FloatType>(1.0, 0.0)/z;
  359. for (int n = 1; n <= nmax_; n++) {
  360. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<FloatType>(n)*zinv - D1[n - 1])
  361. *(static_cast<FloatType>(n)*zinv - D3[n - 1]);
  362. D3[n] = D1[n] + std::complex<FloatType>(0.0, 1.0)/PsiZeta_[n];
  363. }
  364. }
  365. //**********************************************************************************//
  366. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  367. // complex argument (z). //
  368. // Equations (20a) - (21b) //
  369. // //
  370. // Input parameters: //
  371. // z: Complex argument to evaluate Psi and Zeta //
  372. // nmax: Maximum number of terms to calculate Psi and Zeta //
  373. // //
  374. // Output parameters: //
  375. // Psi, Zeta: Riccati-Bessel functions //
  376. //**********************************************************************************//
  377. template <typename FloatType>
  378. void MultiLayerMie<FloatType>::calcPsiZeta(std::complex<FloatType> z,
  379. std::vector<std::complex<FloatType> >& Psi,
  380. std::vector<std::complex<FloatType> >& Zeta) {
  381. std::complex<FloatType> c_i(0.0, 1.0);
  382. std::vector<std::complex<FloatType> > D1(nmax_ + 1), D3(nmax_ + 1);
  383. // First, calculate the logarithmic derivatives
  384. calcD1D3(z, D1, D3);
  385. // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  386. Psi[0] = std::sin(z);
  387. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  388. for (int n = 1; n <= nmax_; n++) {
  389. Psi[n] = Psi[n - 1]*(std::complex<FloatType>(n,0.0)/z - D1[n - 1]);
  390. Zeta[n] = Zeta[n - 1]*(std::complex<FloatType>(n,0.0)/z - D3[n - 1]);
  391. }
  392. }
  393. //**********************************************************************************//
  394. // This function calculates Pi and Tau for a given value of cos(Theta). //
  395. // Equations (26a) - (26c) //
  396. // //
  397. // Input parameters: //
  398. // nmax_: Maximum number of terms to calculate Pi and Tau //
  399. // nTheta: Number of scattering angles //
  400. // Theta: Array containing all the scattering angles where the scattering //
  401. // amplitudes will be calculated //
  402. // //
  403. // Output parameters: //
  404. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  405. //**********************************************************************************//
  406. template <typename FloatType>
  407. void MultiLayerMie<FloatType>::calcPiTau(const FloatType& costheta,
  408. std::vector<FloatType>& Pi, std::vector<FloatType>& Tau) {
  409. int i;
  410. //****************************************************//
  411. // Equations (26a) - (26c) //
  412. //****************************************************//
  413. // Initialize Pi and Tau
  414. Pi[0] = 1.0; // n=1
  415. Tau[0] = costheta;
  416. // Calculate the actual values
  417. if (nmax_ > 1) {
  418. Pi[1] = 3*costheta*Pi[0]; //n=2
  419. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  420. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  421. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  422. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  423. }
  424. }
  425. } // end of MultiLayerMie::calcPiTau(...)
  426. //**********************************************************************************//
  427. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  428. // required to calculate the near-field parameters. //
  429. // //
  430. // Input parameters: //
  431. // Rho: Radial distance //
  432. // Phi: Azimuthal angle //
  433. // Theta: Polar angle //
  434. // rn: Either the spherical Ricatti-Bessel function of first or third kind //
  435. // Dn: Logarithmic derivative of rn //
  436. // Pi, Tau: Angular functions Pi and Tau //
  437. // n: Order of vector spherical harmonics //
  438. // //
  439. // Output parameters: //
  440. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  441. //**********************************************************************************//
  442. template <typename FloatType>
  443. void MultiLayerMie<FloatType>::calcSpherHarm(const std::complex<FloatType> Rho, const FloatType Theta, const FloatType Phi,
  444. const std::complex<FloatType>& rn, const std::complex<FloatType>& Dn,
  445. const FloatType& Pi, const FloatType& Tau, const FloatType& n,
  446. std::vector<std::complex<FloatType> >& Mo1n, std::vector<std::complex<FloatType> >& Me1n,
  447. std::vector<std::complex<FloatType> >& No1n, std::vector<std::complex<FloatType> >& Ne1n) {
  448. // using eq 4.50 in BH
  449. std::complex<FloatType> c_zero(0.0, 0.0);
  450. using nmm::sin;
  451. using nmm::cos;
  452. Mo1n[0] = c_zero;
  453. Mo1n[1] = cos(Phi)*Pi*rn/Rho;
  454. Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
  455. Me1n[0] = c_zero;
  456. Me1n[1] = -sin(Phi)*Pi*rn/Rho;
  457. Me1n[2] = -cos(Phi)*Tau*rn/Rho;
  458. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  459. No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
  460. No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
  461. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  462. Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
  463. Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
  464. } // end of MultiLayerMie::calcSpherHarm(...)
  465. //**********************************************************************************//
  466. // This function calculates the scattering coefficients required to calculate //
  467. // both the near- and far-field parameters. //
  468. // //
  469. // Input parameters: //
  470. // L: Number of layers //
  471. // pl: Index of PEC layer. If there is none just send -1 //
  472. // x: Array containing the size parameters of the layers [0..L-1] //
  473. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  474. // nmax: Maximum number of multipolar expansion terms to be used for the //
  475. // calculations. Only use it if you know what you are doing, otherwise //
  476. // set this parameter to -1 and the function will calculate it. //
  477. // //
  478. // Output parameters: //
  479. // an, bn: Complex scattering amplitudes //
  480. // //
  481. // Return value: //
  482. // Number of multipolar expansion terms used for the calculations //
  483. //**********************************************************************************//
  484. template <typename FloatType>
  485. void MultiLayerMie<FloatType>::calcScattCoeffs() {
  486. isScaCoeffsCalc_ = false;
  487. const std::vector<FloatType>& x = size_param_;
  488. const std::vector<std::complex<FloatType> >& m = refractive_index_;
  489. const int& pl = PEC_layer_position_;
  490. const int L = refractive_index_.size();
  491. //************************************************************************//
  492. // Calculate the index of the first layer. It can be either 0 (default) //
  493. // or the index of the outermost PEC layer. In the latter case all layers //
  494. // below the PEC are discarded. //
  495. // ***********************************************************************//
  496. int fl = (pl > 0) ? pl : 0;
  497. if (nmax_preset_ <= 0) calcNmax(fl);
  498. else nmax_ = nmax_preset_;
  499. std::complex<FloatType> z1, z2;
  500. //**************************************************************************//
  501. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  502. // means that index = layer number - 1 or index = n - 1. The only exception //
  503. // are the arrays for representing D1, D3 and Q because they need a value //
  504. // for the index 0 (zero), hence it is important to consider this shift //
  505. // between different arrays. The change was done to optimize memory usage. //
  506. //**************************************************************************//
  507. // Allocate memory to the arrays
  508. std::vector<std::complex<FloatType> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  509. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  510. std::vector<std::vector<std::complex<FloatType> > > Q(L), Ha(L), Hb(L);
  511. for (int l = 0; l < L; l++) {
  512. Q[l].resize(nmax_ + 1);
  513. Ha[l].resize(nmax_);
  514. Hb[l].resize(nmax_);
  515. }
  516. an_.resize(nmax_);
  517. bn_.resize(nmax_);
  518. PsiZeta_.resize(nmax_ + 1);
  519. std::vector<std::complex<FloatType> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  520. //*************************************************//
  521. // Calculate D1 and D3 for z1 in the first layer //
  522. //*************************************************//
  523. if (fl == pl) { // PEC layer
  524. for (int n = 0; n <= nmax_; n++) {
  525. D1_mlxl[n] = std::complex<FloatType>(0.0, - 1.0);
  526. D3_mlxl[n] = std::complex<FloatType>(0.0, 1.0);
  527. }
  528. } else { // Regular layer
  529. z1 = x[fl]* m[fl];
  530. // Calculate D1 and D3
  531. calcD1D3(z1, D1_mlxl, D3_mlxl);
  532. }
  533. //******************************************************************//
  534. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  535. //******************************************************************//
  536. for (int n = 0; n < nmax_; n++) {
  537. Ha[fl][n] = D1_mlxl[n + 1];
  538. Hb[fl][n] = D1_mlxl[n + 1];
  539. }
  540. //*****************************************************//
  541. // Iteration from the second layer to the last one (L) //
  542. //*****************************************************//
  543. std::complex<FloatType> Temp, Num, Denom;
  544. std::complex<FloatType> G1, G2;
  545. for (int l = fl + 1; l < L; l++) {
  546. //************************************************************//
  547. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  548. //************************************************************//
  549. z1 = x[l]*m[l];
  550. z2 = x[l - 1]*m[l];
  551. //Calculate D1 and D3 for z1
  552. calcD1D3(z1, D1_mlxl, D3_mlxl);
  553. //Calculate D1 and D3 for z2
  554. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  555. //*************************************************//
  556. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  557. //*************************************************//
  558. // Upward recurrence for Q - equations (19a) and (19b)
  559. Num = std::complex<FloatType>(nmm::exp(-2.0*(z1.imag() - z2.imag())), 0.0)
  560. *std::complex<FloatType>(nmm::cos(-2.0*z2.real()) - nmm::exp(-2.0*z2.imag()), nmm::sin(-2.0*z2.real()));
  561. Denom = std::complex<FloatType>(nmm::cos(-2.0*z1.real()) - nmm::exp(-2.0*z1.imag()), nmm::sin(-2.0*z1.real()));
  562. Q[l][0] = Num/Denom;
  563. for (int n = 1; n <= nmax_; n++) {
  564. Num = (z1*D1_mlxl[n] + FloatType(n))*(FloatType(n) - z1*D3_mlxl[n - 1]);
  565. Denom = (z2*D1_mlxlM1[n] + FloatType(n))*(FloatType(n) - z2*D3_mlxlM1[n - 1]);
  566. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  567. }
  568. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  569. for (int n = 1; n <= nmax_; n++) {
  570. //Ha
  571. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  572. G1 = -D1_mlxlM1[n];
  573. G2 = -D3_mlxlM1[n];
  574. } else {
  575. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  576. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  577. } // end of if PEC
  578. Temp = Q[l][n]*G1;
  579. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  580. Denom = G2 - Temp;
  581. Ha[l][n - 1] = Num/Denom;
  582. //Hb
  583. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  584. G1 = Hb[l - 1][n - 1];
  585. G2 = Hb[l - 1][n - 1];
  586. } else {
  587. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  588. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  589. } // end of if PEC
  590. Temp = Q[l][n]*G1;
  591. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  592. Denom = (G2- Temp);
  593. Hb[l][n - 1] = (Num/ Denom);
  594. } // end of for Ha and Hb terms
  595. } // end of for layers iteration
  596. //**************************************//
  597. //Calculate Psi and Zeta for XL //
  598. //**************************************//
  599. // Calculate PsiXL and ZetaXL
  600. calcPsiZeta(std::complex<FloatType>(x[L - 1],0.0), PsiXL, ZetaXL);
  601. //*********************************************************************//
  602. // Finally, we calculate the scattering coefficients (an and bn) and //
  603. // the angular functions (Pi and Tau). Note that for these arrays the //
  604. // first layer is 0 (zero), in future versions all arrays will follow //
  605. // this convention to save memory. (13 Nov, 2014) //
  606. //*********************************************************************//
  607. for (int n = 0; n < nmax_; n++) {
  608. //********************************************************************//
  609. //Expressions for calculating an and bn coefficients are not valid if //
  610. //there is only one PEC layer (ie, for a simple PEC sphere). //
  611. //********************************************************************//
  612. if (pl < (L - 1)) {
  613. 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]);
  614. 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]);
  615. } else {
  616. 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]);
  617. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  618. }
  619. } // end of for an and bn terms
  620. isScaCoeffsCalc_ = true;
  621. } // end of MultiLayerMie::calcScattCoeffs()
  622. //**********************************************************************************//
  623. // This function calculates the actual scattering parameters and amplitudes //
  624. // //
  625. // Input parameters: //
  626. // L: Number of layers //
  627. // pl: Index of PEC layer. If there is none just send -1 //
  628. // x: Array containing the size parameters of the layers [0..L-1] //
  629. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  630. // nTheta: Number of scattering angles //
  631. // Theta: Array containing all the scattering angles where the scattering //
  632. // amplitudes will be calculated //
  633. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  634. // calculations. Only use it if you know what you are doing, otherwise //
  635. // set this parameter to -1 and the function will calculate it //
  636. // //
  637. // Output parameters: //
  638. // Qext: Efficiency factor for extinction //
  639. // Qsca: Efficiency factor for scattering //
  640. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  641. // Qbk: Efficiency factor for backscattering //
  642. // Qpr: Efficiency factor for the radiation pressure //
  643. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  644. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  645. // S1, S2: Complex scattering amplitudes //
  646. // //
  647. // Return value: //
  648. // Number of multipolar expansion terms used for the calculations //
  649. //**********************************************************************************//
  650. template <typename FloatType>
  651. void MultiLayerMie<FloatType>::RunMieCalculation() {
  652. if (size_param_.size() != refractive_index_.size())
  653. throw std::invalid_argument("Each size parameter should have only one index!");
  654. if (size_param_.size() == 0)
  655. throw std::invalid_argument("Initialize model first!");
  656. const std::vector<FloatType>& x = size_param_;
  657. //MarkUncalculated();
  658. // Calculate scattering coefficients
  659. if (!isScaCoeffsCalc_) calcScattCoeffs();
  660. // Initialize the scattering parameters
  661. Qext_ = 0.0;
  662. Qsca_ = 0.0;
  663. Qabs_ = 0.0;
  664. Qbk_ = 0.0;
  665. Qpr_ = 0.0;
  666. asymmetry_factor_ = 0.0;
  667. albedo_ = 0.0;
  668. // Initialize the scattering amplitudes
  669. std::vector<std::complex<FloatType> > tmp1(theta_.size(),std::complex<FloatType>(0.0, 0.0));
  670. S1_.swap(tmp1);
  671. S2_ = S1_;
  672. // Precalculate cos(theta) - gives about 5% speed up.
  673. std::vector<FloatType> costheta(theta_.size(), 0.0);
  674. for (unsigned int t = 0; t < theta_.size(); t++) {
  675. costheta[t] = nmm::cos(theta_[t]);
  676. }
  677. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  678. std::complex<FloatType> Qbktmp(0.0, 0.0);
  679. std::vector< std::complex<FloatType> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  680. // By using downward recurrence we avoid loss of precision due to float rounding errors
  681. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  682. // http://en.wikipedia.org/wiki/Loss_of_significance
  683. for (int n = nmax_ - 2; n >= 0; n--) {
  684. const int n1 = n + 1;
  685. if (mode_n_ == Modes::kAll) {
  686. // Equation (27)
  687. Qext_ += (n1 + n1 + 1.0) * (an_[n].real() + bn_[n].real());
  688. // Equation (28)
  689. Qsca_ += (n1 + n1 + 1.0) * (an_[n].real() * an_[n].real() + an_[n].imag() * an_[n].imag()
  690. + bn_[n].real() * bn_[n].real() + bn_[n].imag() * bn_[n].imag());
  691. // Equation (29)
  692. Qpr_ += ((n1 * (n1 + 2.0) / (n1 + 1.0)) * ((an_[n] * std::conj(an_[n1]) + bn_[n] * std::conj(bn_[n1])).real())
  693. + ((n1 + n1 + 1.0) / (n1 * (n1 + 1.0))) * (an_[n] * std::conj(bn_[n])).real());
  694. // Equation (33)
  695. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (an_[n] - bn_[n]);
  696. // Calculate the scattering amplitudes (S1 and S2) Equations (25a) - (25b)
  697. for (unsigned int t = 0; t < theta_.size(); t++) {
  698. calcPiTau(costheta[t], Pi, Tau);
  699. S1_[t] += calc_S1(n1, an_[n], bn_[n], Pi[n], Tau[n]);
  700. S2_[t] += calc_S2(n1, an_[n], bn_[n], Pi[n], Tau[n]);
  701. }
  702. continue;
  703. }
  704. if (n1 == mode_n_) {
  705. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  706. Qext_ += (n1 + n1 + 1.0) * (an_[n].real());
  707. Qsca_ += (n1 + n1 + 1.0) * (an_[n].real() * an_[n].real() + an_[n].imag() * an_[n].imag());
  708. Qpr_ += std::nan("");
  709. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (an_[n]);
  710. for (unsigned int t = 0; t < theta_.size(); t++) {
  711. calcPiTau(costheta[t], Pi, Tau);
  712. S1_[t] += calc_S1(n1, an_[n], static_cast<std::complex<FloatType>>(0), Pi[n], Tau[n]);
  713. S2_[t] += calc_S2(n1, an_[n], static_cast<std::complex<FloatType>>(0), Pi[n], Tau[n]);
  714. }
  715. }
  716. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  717. Qext_ += (n1 + n1 + 1.0) * (bn_[n].real());
  718. Qsca_ += (n1 + n1 + 1.0) * (bn_[n].real() * bn_[n].real() + bn_[n].imag() * bn_[n].imag());
  719. Qpr_ += std::nan("");
  720. Qbktmp += (FloatType) (n1 + n1 + 1.0) * (1.0 - 2.0 * (n1 % 2)) * (bn_[n]);
  721. for (unsigned int t = 0; t < theta_.size(); t++) {
  722. calcPiTau(costheta[t], Pi, Tau);
  723. S1_[t] += calc_S1(n1, static_cast<std::complex<FloatType>>(0), bn_[n], Pi[n], Tau[n]);
  724. S2_[t] += calc_S2(n1, static_cast<std::complex<FloatType>>(0), bn_[n], Pi[n], Tau[n]);
  725. }
  726. }
  727. }
  728. }
  729. FloatType x2 = pow2(x.back());
  730. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  731. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  732. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  733. Qabs_ = Qext_ - Qsca_; // Equation (30)
  734. albedo_ = Qsca_/Qext_; // Equation (31)
  735. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  736. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  737. isMieCalculated_ = true;
  738. }
  739. //**********************************************************************************//
  740. // This function calculates the expansion coefficients inside the particle, //
  741. // required to calculate the near-field parameters. //
  742. // //
  743. // Input parameters: //
  744. // L: Number of layers //
  745. // pl: Index of PEC layer. If there is none just send -1 //
  746. // x: Array containing the size parameters of the layers [0..L-1] //
  747. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  748. // nmax: Maximum number of multipolar expansion terms to be used for the //
  749. // calculations. Only use it if you know what you are doing, otherwise //
  750. // set this parameter to -1 and the function will calculate it. //
  751. // //
  752. // Output parameters: //
  753. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  754. // //
  755. // Return value: //
  756. // Number of multipolar expansion terms used for the calculations //
  757. //**********************************************************************************//
  758. template <typename FloatType>
  759. void MultiLayerMie<FloatType>::calcExpanCoeffs() {
  760. if (!isScaCoeffsCalc_)
  761. throw std::invalid_argument("(calcExpanCoeffs) You should calculate external coefficients first!");
  762. isExpCoeffsCalc_ = false;
  763. std::complex<FloatType> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  764. const int L = refractive_index_.size();
  765. aln_.resize(L + 1);
  766. bln_.resize(L + 1);
  767. cln_.resize(L + 1);
  768. dln_.resize(L + 1);
  769. for (int l = 0; l <= L; l++) {
  770. aln_[l].resize(nmax_);
  771. bln_[l].resize(nmax_);
  772. cln_[l].resize(nmax_);
  773. dln_[l].resize(nmax_);
  774. }
  775. // Yang, paragraph under eq. A3
  776. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  777. for (int n = 0; n < nmax_; n++) {
  778. aln_[L][n] = an_[n];
  779. bln_[L][n] = bn_[n];
  780. cln_[L][n] = c_one;
  781. dln_[L][n] = c_one;
  782. }
  783. std::vector<std::complex<FloatType> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  784. std::vector<std::complex<FloatType> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  785. std::complex<FloatType> denomZeta, denomPsi, T1, T2, T3, T4;
  786. auto& m = refractive_index_;
  787. std::vector< std::complex<FloatType> > m1(L);
  788. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  789. m1[L - 1] = std::complex<FloatType> (1.0, 0.0);
  790. std::complex<FloatType> z, z1;
  791. for (int l = L - 1; l >= 0; l--) {
  792. if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!!
  793. for (int n = 0; n < nmax_; n++) {
  794. // aln
  795. aln_[l][n] = c_zero;
  796. // bln
  797. bln_[l][n] = c_zero;
  798. // cln
  799. cln_[l][n] = c_zero;
  800. // dln
  801. dln_[l][n] = c_zero;
  802. }
  803. } else { // Regular material, just do the calculation
  804. z = size_param_[l]*m[l];
  805. z1 = size_param_[l]*m1[l];
  806. calcD1D3(z, D1z, D3z);
  807. calcD1D3(z1, D1z1, D3z1);
  808. calcPsiZeta(z, Psiz, Zetaz);
  809. calcPsiZeta(z1, Psiz1, Zetaz1);
  810. for (int n = 0; n < nmax_; n++) {
  811. int n1 = n + 1;
  812. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  813. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  814. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  815. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  816. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  817. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  818. // aln
  819. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  820. // bln
  821. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  822. // cln
  823. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  824. // dln
  825. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  826. } // end of all n
  827. } // end PEC condition
  828. } // end of all l
  829. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  830. for (int n = 0; n < nmax_; ++n) {
  831. if (cabs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  832. else {
  833. //throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  834. std::cout<< std::setprecision(100)
  835. << "Warning: Potentially unstable calculation of aln[0]["
  836. << n << "] = "<< aln_[0][n] <<std::endl;
  837. aln_[0][n] = 0.0;
  838. }
  839. if (cabs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  840. else {
  841. //throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  842. std::cout<< std::setprecision(100)
  843. << "Warning: Potentially unstable calculation of bln[0]["
  844. << n << "] = "<< bln_[0][n] <<std::endl;
  845. bln_[0][n] = 0.0;
  846. }
  847. }
  848. isExpCoeffsCalc_ = true;
  849. } // end of void MultiLayerMie::calcExpanCoeffs()
  850. //**********************************************************************************//
  851. // This function calculates the electric (E) and magnetic (H) fields inside and //
  852. // around the particle. //
  853. // //
  854. // Input parameters (coordinates of the point): //
  855. // Rho: Radial distance //
  856. // Phi: Azimuthal angle //
  857. // Theta: Polar angle //
  858. // mode_n: mode order. //
  859. // -1 - use all modes (all_) //
  860. // 1 - use dipole mode only //
  861. // 2 - use quadrupole mode only //
  862. // ... //
  863. // mode_type: only used when mode_n != -1 //
  864. // 0 - electric only //
  865. // 1 - magnetic only //
  866. // //
  867. // //
  868. // Output parameters: //
  869. // E, H: Complex electric and magnetic fields //
  870. //**********************************************************************************//
  871. template <typename FloatType>
  872. void MultiLayerMie<FloatType>::calcFieldByComponents(const FloatType Rho,
  873. const FloatType Theta, const FloatType Phi,
  874. std::vector<std::complex<FloatType> >& E,
  875. std::vector<std::complex<FloatType> >& H) {
  876. std::complex<FloatType> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  877. // Vector containing precomputed integer powers of i to avoid computation
  878. std::vector<std::complex<FloatType> > ipow = {c_one, c_i, -c_one, -c_i};
  879. std::vector<std::complex<FloatType> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  880. std::vector<std::complex<FloatType> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  881. std::vector<std::complex<FloatType> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  882. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  883. int l = 0; // Layer number
  884. std::complex<FloatType> ml;
  885. // Initialize E and H
  886. for (int i = 0; i < 3; i++) {
  887. E[i] = c_zero;
  888. H[i] = c_zero;
  889. }
  890. if (Rho > size_param_.back()) {
  891. l = size_param_.size();
  892. ml = c_one;
  893. } else {
  894. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  895. if (Rho <= size_param_[i]) {
  896. l = i;
  897. }
  898. }
  899. ml = refractive_index_[l];
  900. }
  901. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  902. calcD1D3(Rho*ml, D1n, D3n);
  903. // Calculate Ricatti-Bessel functions
  904. calcPsiZeta(Rho*ml, Psi, Zeta);
  905. // Calculate angular functions Pi and Tau
  906. calcPiTau(nmm::cos(Theta), Pi, Tau);
  907. for (int n = nmax_ - 2; n >= 0; n--) {
  908. int n1 = n + 1;
  909. auto rn = static_cast<FloatType>(n1);
  910. // using BH 4.12 and 4.50
  911. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  912. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  913. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  914. std::complex<FloatType> En = ipow[n1 % 4]
  915. *static_cast<FloatType>((rn + rn + 1.0)/(rn*rn + rn));
  916. for (int i = 0; i < 3; i++) {
  917. if (mode_n_ == Modes::kAll) {
  918. // electric field E [V m - 1] = EF*E0
  919. E[i] += En*( cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  920. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  921. H[i] += En*( -dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  922. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  923. continue;
  924. }
  925. if (n1 == mode_n_) {
  926. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  927. E[i] += En*( -c_i*dln_[l][n]*N1e1n[i]
  928. + c_i*aln_[l][n]*N3e1n[i]);
  929. H[i] += En*(-dln_[l][n]*M1e1n[i]
  930. +aln_[l][n]*M3e1n[i]);
  931. //std::cout << mode_n_;
  932. }
  933. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  934. E[i] += En*( cln_[l][n]*M1o1n[i]
  935. - bln_[l][n]*M3o1n[i]);
  936. H[i] += En*( -c_i*cln_[l][n]*N1o1n[i]
  937. + c_i*bln_[l][n]*N3o1n[i]);
  938. //std::cout << mode_n_;
  939. }
  940. //std::cout << std::endl;
  941. }
  942. //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=====*****=====");
  943. }
  944. } // end of for all n
  945. // magnetic field
  946. std::complex<FloatType> hffact = ml/static_cast<FloatType>(cc_*mu_);
  947. for (int i = 0; i < 3; i++) {
  948. H[i] = hffact*H[i];
  949. }
  950. } // end of MultiLayerMie::calcFieldByComponents(...)
  951. //**********************************************************************************//
  952. // This function calculates complex electric and magnetic field in the surroundings //
  953. // and inside the particle. //
  954. // //
  955. // Input parameters: //
  956. // L: Number of layers //
  957. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  958. // x: Array containing the size parameters of the layers [0..L-1] //
  959. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  960. // nmax: Maximum number of multipolar expansion terms to be used for the //
  961. // calculations. Only use it if you know what you are doing, otherwise //
  962. // set this parameter to 0 (zero) and the function will calculate it. //
  963. // ncoord: Number of coordinate points //
  964. // Coords: Array containing all coordinates where the complex electric and //
  965. // magnetic fields will be calculated //
  966. // mode_n: mode order. //
  967. // -1 - use all modes (all_) //
  968. // 1 - use dipole mode only //
  969. // 2 - use quadrupole mode only //
  970. // ... //
  971. // mode_type: only used when mode_n != -1 //
  972. // 0 - electric only //
  973. // 1 - magnetic only //
  974. // //
  975. // Output parameters: //
  976. // E, H: Complex electric and magnetic field at the provided coordinates //
  977. // //
  978. // Return value: //
  979. // Number of multipolar expansion terms used for the calculations //
  980. //**********************************************************************************//
  981. template <typename FloatType>
  982. void MultiLayerMie<FloatType>::RunFieldCalculation() {
  983. FloatType Rho, Theta, Phi;
  984. // Calculate scattering coefficients an_ and bn_
  985. calcScattCoeffs();
  986. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  987. calcExpanCoeffs();
  988. long total_points = coords_[0].size();
  989. E_.resize(total_points);
  990. H_.resize(total_points);
  991. Es_.resize(total_points);
  992. Hs_.resize(total_points);
  993. for (auto& f : E_) f.resize(3);
  994. for (auto& f : H_) f.resize(3);
  995. for (auto& f : Es_) f.resize(3);
  996. for (auto& f : Hs_) f.resize(3);
  997. for (int point = 0; point < total_points; point++) {
  998. const FloatType& Xp = coords_[0][point];
  999. const FloatType& Yp = coords_[1][point];
  1000. const FloatType& Zp = coords_[2][point];
  1001. // Convert to spherical coordinates
  1002. Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1003. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1004. Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0;
  1005. // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc.
  1006. Phi = nmm::atan2(Yp,Xp);
  1007. // Avoid convergence problems due to Rho too small
  1008. if (Rho < 1e-5) Rho = 1e-5;
  1009. // std::cout << "Xp: "<<Xp<< " Yp: "<<Yp<< " Zp: "<<Zp<<std::endl;
  1010. // std::cout << " Rho: "<<Rho<<" Theta: "<<Theta<<" Phi:"<<Phi<<std::endl<<std::endl;
  1011. //*******************************************************//
  1012. // external scattering field = incident + scattered //
  1013. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1014. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1015. //*******************************************************//
  1016. // This array contains the fields in spherical coordinates
  1017. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  1018. // Do the actual calculation of electric and magnetic field
  1019. calcFieldByComponents(Rho, Theta, Phi, Es, Hs);
  1020. for (int sph_coord = 0; sph_coord<3; ++sph_coord) {
  1021. Es_[point][sph_coord] = Es[sph_coord];
  1022. Hs_[point][sph_coord] = Hs[sph_coord];
  1023. }
  1024. { //Now, convert the fields back to cartesian coordinates
  1025. using nmm::sin;
  1026. using nmm::cos;
  1027. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1028. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1029. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1030. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1031. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1032. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1033. }
  1034. } // end of for all field coordinates
  1035. } // end of MultiLayerMie::RunFieldCalculation()
  1036. } // end of namespace nmie
  1037. #endif // SRC_NMIE_IMPL_HPP_