nmie-impl.hpp 57 KB

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