nmie-impl.cc 55 KB

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