nmie-impl.hpp 58 KB

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