nmie-nearfield.hpp 35 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703
  1. #ifndef SRC_NMIE_NEARFIELD_HPP_
  2. #define SRC_NMIE_NEARFIELD_HPP_
  3. //**********************************************************************************//
  4. // Copyright (C) 2009-2021 Ovidio Pena <ovidio@bytesfall.com> //
  5. // Copyright (C) 2013-2021 Konstantin Ladutenko <kostyfisik@gmail.com> //
  6. // //
  7. // This file is part of scattnlay //
  8. // //
  9. // This program is free software: you can redistribute it and/or modify //
  10. // it under the terms of the GNU General Public License as published by //
  11. // the Free Software Foundation, either version 3 of the License, or //
  12. // (at your option) any later version. //
  13. // //
  14. // This program is distributed in the hope that it will be useful, //
  15. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  16. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  17. // GNU General Public License for more details. //
  18. // //
  19. // The only additional remark is that we expect that all publications //
  20. // describing work using this software, or all commercial products //
  21. // using it, cite at least one of the following references: //
  22. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  23. // a multilayered sphere," Computer Physics Communications, //
  24. // vol. 180, Nov. 2009, pp. 2348-2354. //
  25. // [2] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie //
  26. // calculation of electromagnetic near-field for a multilayered //
  27. // sphere," Computer Physics Communications, vol. 214, May 2017, //
  28. // pp. 225-230. //
  29. // //
  30. // You should have received a copy of the GNU General Public License //
  31. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  32. //**********************************************************************************//
  33. //**********************************************************************************//
  34. // This class implements the algorithm for a multilayered sphere described by: //
  35. // [1] W. Yang, "Improved recursive algorithm for light scattering by a //
  36. // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. //
  37. // //
  38. // You can find the description of all the used equations in: //
  39. // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  40. // a multilayered sphere," Computer Physics Communications, //
  41. // vol. 180, Nov. 2009, pp. 2348-2354. //
  42. // [3] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie //
  43. // calculation of electromagnetic near-field for a multilayered //
  44. // sphere," Computer Physics Communications, vol. 214, May 2017, //
  45. // pp. 225-230. //
  46. // //
  47. // Hereinafter all equations numbers refer to [2] //
  48. //**********************************************************************************//
  49. #include <iostream>
  50. #include <iomanip>
  51. #include <stdexcept>
  52. #include <vector>
  53. //#include "nmie.hpp"
  54. namespace nmie {
  55. //class implementation
  56. //**********************************************************************************//
  57. // This function calculates the expansion coefficients inside the particle, //
  58. // required to calculate the near-field parameters. //
  59. // //
  60. // Input parameters: //
  61. // L: Number of layers //
  62. // pl: Index of PEC layer. If there is none just send -1 //
  63. // x: Array containing the size parameters of the layers [0..L-1] //
  64. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  65. // nmax: Maximum number of multipolar expansion terms to be used for the //
  66. // calculations. Only use it if you know what you are doing, otherwise //
  67. // set this parameter to -1 and the function will calculate it. //
  68. // //
  69. // Output parameters: //
  70. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  71. // //
  72. // Return value: //
  73. // Number of multipolar expansion terms used for the calculations //
  74. //**********************************************************************************//
  75. template <typename FloatType>
  76. void MultiLayerMie<FloatType>::calcExpanCoeffs() {
  77. if (!isScaCoeffsCalc_)
  78. throw std::invalid_argument("(calcExpanCoeffs) You should calculate external coefficients first!");
  79. isExpCoeffsCalc_ = false;
  80. aln_.clear(); bln_.clear(); cln_.clear(); dln_.clear();
  81. std::complex<FloatType> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  82. const int L = refractive_index_.size();
  83. aln_.resize(L + 1);
  84. bln_.resize(L + 1);
  85. cln_.resize(L + 1);
  86. dln_.resize(L + 1);
  87. for (int l = 0; l <= L; l++) {
  88. aln_[l].resize(nmax_, static_cast<FloatType>(0.0));
  89. bln_[l].resize(nmax_, static_cast<FloatType>(0.0));
  90. cln_[l].resize(nmax_, static_cast<FloatType>(0.0));
  91. dln_[l].resize(nmax_, static_cast<FloatType>(0.0));
  92. }
  93. // Yang, paragraph under eq. A3
  94. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  95. for (int n = 0; n < nmax_; n++) {
  96. aln_[L][n] = an_[n];
  97. bln_[L][n] = bn_[n];
  98. cln_[L][n] = c_one;
  99. dln_[L][n] = c_one;
  100. }
  101. std::vector<std::complex<FloatType> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  102. std::vector<std::complex<FloatType> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  103. std::complex<FloatType> denomZeta, denomPsi, T1, T2, T3, T4;
  104. auto &m = refractive_index_;
  105. std::vector< std::complex<FloatType> > m1(L);
  106. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  107. m1[L - 1] = std::complex<FloatType> (1.0, 0.0);
  108. std::complex<FloatType> z, z1;
  109. for (int l = L - 1; l >= 0; l--) {
  110. if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!!
  111. for (int n = 0; n < nmax_; n++) {
  112. // aln
  113. aln_[l][n] = c_zero;
  114. // bln
  115. bln_[l][n] = c_zero;
  116. // cln
  117. cln_[l][n] = c_zero;
  118. // dln
  119. dln_[l][n] = c_zero;
  120. }
  121. } else { // Regular material, just do the calculation
  122. z = size_param_[l]*m[l];
  123. z1 = size_param_[l]*m1[l];
  124. calcD1D3(z, D1z, D3z);
  125. calcD1D3(z1, D1z1, D3z1);
  126. calcPsiZeta(z, Psiz, Zetaz);
  127. calcPsiZeta(z1, Psiz1, Zetaz1);
  128. for (int n = 0; n < nmax_; n++) {
  129. int n1 = n + 1;
  130. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  131. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  132. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  133. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  134. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  135. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  136. // aln
  137. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  138. // bln
  139. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  140. // cln
  141. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  142. // dln
  143. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  144. } // end of all n
  145. } // end PEC condition
  146. } // end of all l
  147. int print_precision = 16;
  148. #ifdef MULTI_PRECISION
  149. print_precision = MULTI_PRECISION;
  150. #endif
  151. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  152. int print_count = 0;
  153. for (int n = 0; n < nmax_; ++n) {
  154. if (cabs(aln_[0][n]) > 1e-10 && print_count < 2) {
  155. print_count++;
  156. std::cout<< std::setprecision(print_precision)
  157. << "Warning: Potentially unstable calculation of aln[0]["
  158. << n << "] = "<< aln_[0][n] <<std::endl;
  159. }
  160. if (cabs(bln_[0][n]) > 1e-10 && print_count < 2) {
  161. print_count++;
  162. std::cout<< std::setprecision(print_precision)
  163. << "Warning: Potentially unstable calculation of bln[0]["
  164. << n << "] = "<< bln_[0][n] <<std::endl;
  165. }
  166. aln_[0][n] = 0.0;
  167. bln_[0][n] = 0.0;
  168. }
  169. isExpCoeffsCalc_ = true;
  170. } // end of void MultiLayerMie::calcExpanCoeffs()
  171. template <typename FloatType>
  172. void MultiLayerMie<FloatType>::convertFieldsFromSphericalToCartesian() {
  173. long total_points = coords_polar_.size();
  174. E_.clear(); H_.clear();
  175. Eabs_.clear(); Habs_.clear();
  176. for (int point=0; point < total_points; point++) {
  177. auto Theta = coords_polar_[point][1];
  178. auto Phi = coords_polar_[point][2];
  179. auto Es = Es_[point];
  180. auto Hs = Hs_[point];
  181. using nmm::sin;
  182. using nmm::cos;
  183. E_.push_back({ sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2],
  184. sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2],
  185. cos(Theta)*Es[0] - sin(Theta)*Es[1]});
  186. H_.push_back({ sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2],
  187. sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2],
  188. cos(Theta)*Hs[0] - sin(Theta)*Hs[1]});
  189. Eabs_.push_back(vabs(E_.back()));
  190. Habs_.push_back(vabs(H_.back()));
  191. }
  192. } // end of void MultiLayerMie::convertFieldsFromSphericalToCartesian()
  193. //**********************************************************************************//
  194. // This function calculates the electric (E) and magnetic (H) fields inside and //
  195. // around the particle. //
  196. //
  197. // Main trouble of near-field evaluations is supposed to originate from special functions
  198. // evaluation, so we expect that nmax needed for the convergence is the size
  199. // of Psi vector.
  200. // //
  201. // Input parameters (coordinates of the point): //
  202. // Rho: Radial distance //
  203. // Phi: Azimuthal angle //
  204. // Theta: Polar angle //
  205. // mode_n: mode order. //
  206. // -1 - use all modes (all_) //
  207. // 1 - use dipole mode only //
  208. // 2 - use quadrupole mode only //
  209. // ... //
  210. // mode_type: only used when mode_n != -1 //
  211. // 0 - electric only //
  212. // 1 - magnetic only //
  213. // //
  214. // //
  215. // Output parameters: //
  216. // E, H: Complex electric and magnetic fields //
  217. //**********************************************************************************//
  218. template <typename FloatType> template <typename evalType>
  219. void MultiLayerMie<FloatType>::calcFieldByComponents(const evalType Rho,
  220. const evalType Theta, const evalType Phi,
  221. const std::vector<std::complex<evalType> > &Psi,
  222. const std::vector<std::complex<evalType> > &D1n,
  223. const std::vector<std::complex<evalType> > &Zeta,
  224. const std::vector<std::complex<evalType> > &D3n,
  225. const std::vector<evalType> &Pi,
  226. const std::vector<evalType> &Tau,
  227. std::vector<std::complex<evalType> > &E,
  228. std::vector<std::complex<evalType> > &H,
  229. std::vector<bool> &isConvergedE,
  230. std::vector<bool> &isConvergedH,
  231. bool isMarkUnconverged) {
  232. auto nmax = Psi.size() - 1;
  233. std::complex<evalType> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  234. // auto c_nan = ConvertComplex<FloatType>(std::complex<double>(std::nan(""), std::nan("")));
  235. // Vector containing precomputed integer powers of i to avoid computation
  236. std::vector<std::complex<evalType> > ipow = {c_one, c_i, -c_one, -c_i};
  237. std::vector<std::complex<evalType> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  238. std::vector<std::complex<evalType> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  239. std::complex<evalType> ml;
  240. // Initialize E and H
  241. for (int i = 0; i < 3; i++) {
  242. E[i] = c_zero;
  243. H[i] = c_zero;
  244. }
  245. const unsigned L = refractive_index_.size();
  246. for (int n = 0; n < nmax_; n++) {
  247. cln_[L][n] = c_zero;
  248. dln_[L][n] = c_zero;
  249. }
  250. unsigned int l;
  251. GetIndexAtRadius(Rho, ml, l);
  252. isConvergedE = {false, false, false}, isConvergedH = {false, false, false};
  253. // evalType E0 = 0, H0=0;
  254. for (unsigned int n = 0; n < nmax; n++) {
  255. int n1 = n + 1;
  256. auto rn = static_cast<evalType>(n1);
  257. // using BH 4.12 and 4.50
  258. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  259. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  260. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  261. std::complex<evalType> En = ipow[n1 % 4]
  262. *static_cast<evalType>((rn + rn + 1.0)/(rn*rn + rn));
  263. std::complex<evalType> Ediff, Hdiff;
  264. std::complex<FloatType> Ediff_ft, Hdiff_ft;
  265. auto aln = ConvertComplex<evalType>(aln_[l][n]);
  266. auto bln = ConvertComplex<evalType>(bln_[l][n]);
  267. auto cln = ConvertComplex<evalType>(cln_[l][n]);
  268. auto dln = ConvertComplex<evalType>(dln_[l][n]);
  269. for (int i = 0; i < 3; i++) {
  270. if (isConvergedE[i] && isConvergedH[i]) continue;
  271. Ediff = En*( cln*M1o1n[i] - c_i*dln*N1e1n[i]
  272. + c_i*aln*N3e1n[i] - bln*M3o1n[i]);
  273. Hdiff = En*( -dln*M1e1n[i] - c_i*cln*N1o1n[i]
  274. + c_i*bln*N3o1n[i] + aln*M3e1n[i]);
  275. Ediff_ft = ConvertComplex<FloatType>(Ediff);
  276. Hdiff_ft = ConvertComplex<FloatType>(Hdiff);
  277. if ( nmm::isnan(Ediff_ft.real()) || nmm::isnan(Ediff_ft.imag()) ||
  278. nmm::isnan(Hdiff_ft.real()) || nmm::isnan(Hdiff_ft.imag()) ) {
  279. std::cout << "Unexpected truncation during near-field evaluation at n = "<< n
  280. << " (of total nmax = "<<nmax<<")!!!"<<std::endl;
  281. break;
  282. }
  283. if (n!=0) {
  284. if (cabs(Ediff) == 0) isConvergedE[i] = true;
  285. if (cabs(Hdiff) == 0) isConvergedH[i] = true;
  286. if (cabs(E[i]) != 0.)
  287. if (cabs(Ediff)/cabs(E[i]) < nearfield_convergence_threshold_) isConvergedE[i] = true;
  288. if (cabs(H[i]) != 0.)
  289. if (cabs(Hdiff)/cabs(H[i]) < nearfield_convergence_threshold_) isConvergedH[i] = true;
  290. }
  291. if (isConvergedE[i]) Ediff = c_zero;
  292. if (isConvergedH[i]) Hdiff = c_zero;
  293. if ((!isConvergedH[i] || !isConvergedE[i]) && n==nmax-1 && GetFieldConvergence()) {
  294. std::cout<<"Econv:"<<cabs(Ediff)/cabs(E[i])<<" Hconv:"<<cabs(Hdiff)/cabs(H[i])<<std::endl;
  295. }
  296. if (mode_n_ == Modes::kAll) {
  297. // electric field E [V m - 1] = EF*E0
  298. E[i] += Ediff;
  299. H[i] += Hdiff;
  300. continue;
  301. }
  302. if (n1 == mode_n_) {
  303. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  304. E[i] += En*( -c_i*dln*N1e1n[i]
  305. + c_i*aln*N3e1n[i]);
  306. H[i] += En*(-dln*M1e1n[i]
  307. +aln*M3e1n[i]);
  308. //std::cout << mode_n_;
  309. }
  310. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  311. E[i] += En*( cln*M1o1n[i]
  312. - bln*M3o1n[i]);
  313. H[i] += En*( -c_i*cln*N1o1n[i]
  314. + c_i*bln*N3o1n[i]);
  315. //std::cout << mode_n_;
  316. }
  317. //std::cout << std::endl;
  318. }
  319. //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=====*****=====");
  320. }
  321. if (nmm::isnan(Ediff_ft.real()) || nmm::isnan(Ediff_ft.imag()) ||
  322. nmm::isnan(Hdiff_ft.real()) || nmm::isnan(Hdiff_ft.imag())
  323. ) break;
  324. } // end of for all n
  325. // Add the incident field
  326. if(l==L) {
  327. // using nmm::sin_t;
  328. // using nmm::cos_t;
  329. const auto z = Rho*cos_t(Theta);
  330. const auto Ex = std::complex<evalType>(cos_t(z), sin_t(z));
  331. E[0] += Ex*cos_t(Phi)*sin_t(Theta);
  332. E[1] += Ex*cos_t(Phi)*cos_t(Theta);
  333. E[2] += -Ex*sin_t(Phi);
  334. const auto Hy = Ex;
  335. H[0] += Hy*sin_t(Theta)*sin_t(Phi);
  336. H[1] += Hy*cos_t(Theta)*sin_t(Phi);
  337. H[2] += Hy*cos_t(Phi);
  338. }
  339. if( (!isConvergedE[0] || !isConvergedE[1] ||!isConvergedE[2] ||
  340. !isConvergedH[0] || !isConvergedH[1] ||!isConvergedH[2] ) && GetFieldConvergence()) {
  341. std::cout << "Field evaluation failed to converge an nmax = "<< nmax << std::endl;
  342. std::cout << "Near-field convergence threshold: "<<nearfield_convergence_threshold_<<std::endl;
  343. if (isMarkUnconverged) { //mark as NaN
  344. for(auto &ee :E) ee /= c_zero;
  345. for(auto &ee :H) ee /= c_zero;
  346. }
  347. }
  348. // magnetic field
  349. std::complex<evalType> hffact = ml/static_cast<evalType>(cc_*mu_);
  350. for (int i = 0; i < 3; i++) {
  351. H[i] = hffact*H[i];
  352. }
  353. } // end of MultiLayerMie::calcFieldByComponents(...)
  354. //**********************************************************************************//
  355. // This function calculates complex electric and magnetic field in the surroundings //
  356. // and inside the particle. //
  357. // //
  358. // Input parameters: //
  359. // L: Number of layers //
  360. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  361. // x: Array containing the size parameters of the layers [0..L-1] //
  362. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  363. // nmax: Maximum number of multipolar expansion terms to be used for the //
  364. // calculations. Only use it if you know what you are doing, otherwise //
  365. // set this parameter to 0 (zero) and the function will calculate it. //
  366. // ncoord: Number of coordinate points //
  367. // Coords: Array containing all coordinates where the complex electric and //
  368. // magnetic fields will be calculated //
  369. // mode_n: mode order. //
  370. // -1 - use all modes (all_) //
  371. // 1 - use dipole mode only //
  372. // 2 - use quadrupole mode only //
  373. // ... //
  374. // mode_type: only used when mode_n != -1 //
  375. // 0 - electric only //
  376. // 1 - magnetic only //
  377. // //
  378. // Output parameters: //
  379. // E, H: Complex electric and magnetic field at the provided coordinates //
  380. // //
  381. // Return value: //
  382. // Number of multipolar expansion terms used for the calculations //
  383. //**********************************************************************************//
  384. template <typename FloatType>
  385. void MultiLayerMie<FloatType>::RunFieldCalculation(bool isMarkUnconverged) {
  386. // Calculate scattering coefficients an_ and bn_
  387. calcScattCoeffs();
  388. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  389. calcExpanCoeffs();
  390. std::vector<bool> isConvergedE = {false, false, false}, isConvergedH = {false, false, false};
  391. isConvergedE_ = {true, true, true}, isConvergedH_ = {true, true, true};
  392. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  393. long total_points = coords_[0].size();
  394. for (int point = 0; point < total_points; point++) {
  395. const FloatType &Xp = coords_[0][point];
  396. const FloatType &Yp = coords_[1][point];
  397. const FloatType &Zp = coords_[2][point];
  398. // Convert to spherical coordinates
  399. auto Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  400. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  401. auto Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0;
  402. // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc.
  403. auto Phi = nmm::atan2(Yp,Xp);
  404. coords_polar_.push_back({Rho, Theta, Phi});
  405. // Avoid convergence problems due to Rho too small
  406. if (Rho < 1e-5) Rho = 1e-5;
  407. //*******************************************************//
  408. // external scattering field = incident + scattered //
  409. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  410. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  411. //*******************************************************//
  412. // This array contains the fields in spherical coordinates
  413. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  414. // Do the actual calculation of electric and magnetic field
  415. std::vector<std::complex<FloatType> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  416. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  417. std::complex<FloatType> ml;
  418. GetIndexAtRadius(Rho, ml);
  419. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  420. calcD1D3(Rho*ml, D1n, D3n);
  421. // Calculate Ricatti-Bessel functions
  422. calcPsiZeta(Rho*ml, Psi, Zeta);
  423. // Calculate angular functions Pi and Tau
  424. calcPiTau(nmm::cos(Theta), Pi, Tau);
  425. calcFieldByComponents(Rho, Theta, Phi, Psi, D1n, Zeta, D3n, Pi, Tau, Es, Hs,
  426. isConvergedE, isConvergedH, isMarkUnconverged);
  427. UpdateConvergenceStatus(isConvergedE, isConvergedH);
  428. Es_.push_back(Es);
  429. Hs_.push_back(Hs);
  430. } // end of for all field coordinates
  431. convertFieldsFromSphericalToCartesian();
  432. } // end of MultiLayerMie::RunFieldCalculation()
  433. // TODO do we really need this eval_delta()?
  434. template <typename FloatType>
  435. double eval_delta(const unsigned int steps, const double from_value, const double to_value) {
  436. auto delta = std::abs(from_value - to_value);
  437. if (steps < 2) return delta;
  438. delta /= static_cast<double>(steps-1);
  439. // We have a limited double precision evaluation of special functions, typically it is 1e-10.
  440. if ( (2.*delta)/std::abs(from_value+to_value) < 1e-9)
  441. throw std::invalid_argument("Error! The step is too fine, not supported!");
  442. return delta;
  443. }
  444. // ml - refractive index
  445. // l - Layer number
  446. template <typename FloatType> template <typename evalType>
  447. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  448. std::complex<evalType> &ml,
  449. unsigned int &l) {
  450. l = 0;
  451. if (Rho > size_param_.back()) {
  452. l = size_param_.size();
  453. ml = std::complex<evalType>(1.0, 0.0);
  454. } else {
  455. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  456. if (Rho <= size_param_[i]) {
  457. l = i;
  458. }
  459. }
  460. ml = ConvertComplex<evalType>(refractive_index_[l]);
  461. }
  462. }
  463. template <typename FloatType> template <typename evalType>
  464. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  465. std::complex<evalType> &ml) {
  466. unsigned int l;
  467. GetIndexAtRadius(Rho, ml, l);
  468. }
  469. template <typename FloatType>
  470. void MultiLayerMie<FloatType>::calcMieSeriesNeededToConverge(const FloatType Rho, int nmax_in) {
  471. if (nmax_in < 1) {
  472. auto required_near_field_nmax = calcNmax(Rho);
  473. SetMaxTerms(required_near_field_nmax);
  474. } else {
  475. SetMaxTerms(nmax_in);
  476. }
  477. // Calculate scattering coefficients an_ and bn_
  478. calcScattCoeffs();
  479. // We might be limited with available machine precision
  480. available_maximal_nmax_ = nmax_;
  481. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  482. calcExpanCoeffs();
  483. }
  484. template <typename FloatType>
  485. void MultiLayerMie<FloatType>::calcRadialOnlyDependantFunctions(const double from_Rho, const double to_Rho,
  486. std::vector<std::vector<std::complex<FloatType> > > &Psi,
  487. std::vector<std::vector<std::complex<FloatType> > > &D1n,
  488. std::vector<std::vector<std::complex<FloatType> > > &Zeta,
  489. std::vector<std::vector<std::complex<FloatType> > > &D3n,
  490. int nmax_in) {
  491. auto radius_points = Psi.size();
  492. std::vector<std::vector<std::complex<FloatType> > > PsiZeta(radius_points);
  493. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  494. for (unsigned int j=0; j < radius_points; j++) {
  495. auto Rho = static_cast<FloatType>(from_Rho + j*delta_Rho);
  496. // if (Rho < 1e-5) Rho = 1e-5; // TODO do we need this?.
  497. int near_field_nmax = nmax_in;
  498. if (nmax_in < 1) near_field_nmax = calcNmax(Rho);
  499. // Skip if not enough terms in Mie series (i.e. required near field nmax > available terms )
  500. if (near_field_nmax > available_maximal_nmax_) near_field_nmax = available_maximal_nmax_;
  501. Psi[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0)); D1n[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  502. Zeta[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0)); D3n[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  503. PsiZeta[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  504. std::complex<FloatType> ml;
  505. GetIndexAtRadius(Rho, ml);
  506. auto z = Rho*ml;
  507. evalDownwardD1(z, D1n[j]);
  508. evalUpwardPsi(z, D1n[j], Psi[j]);
  509. evalUpwardD3 (z, D1n[j], D3n[j], PsiZeta[j]);
  510. for (unsigned int k = 0; k < Zeta[j].size(); k++) {
  511. Zeta[j][k] = PsiZeta[j][k]/Psi[j][k];
  512. }
  513. }
  514. }
  515. // input parameters:
  516. // outer_arc_points: will be increased to the nearest power of 2.
  517. template <typename FloatType>
  518. void MultiLayerMie<FloatType>::RunFieldCalculationPolar(const int outer_arc_points,
  519. const int radius_points,
  520. const double from_Rho, const double to_Rho,
  521. const double from_Theta, const double to_Theta,
  522. const double from_Phi, const double to_Phi,
  523. const bool isMarkUnconverged,
  524. int nmax_in) {
  525. if (from_Rho > to_Rho || from_Theta > to_Theta || from_Phi > to_Phi
  526. || outer_arc_points < 1 || radius_points < 1
  527. || from_Rho < 0.)
  528. throw std::invalid_argument("Error! Invalid argument for RunFieldCalculationPolar() !");
  529. // auto nmax_old = nmax_;
  530. int theta_points = 0, phi_points = 0;
  531. if (to_Theta-from_Theta > to_Phi-from_Phi) {
  532. theta_points = outer_arc_points;
  533. phi_points = static_cast<int>((to_Phi-from_Phi)/(to_Theta-from_Theta) * outer_arc_points);
  534. } else {
  535. phi_points = outer_arc_points;
  536. theta_points = static_cast<int>((to_Theta-from_Theta)/(to_Phi-from_Phi) * outer_arc_points);
  537. }
  538. if (theta_points == 0) theta_points = 1;
  539. if (phi_points == 0) phi_points = 1;
  540. calcMieSeriesNeededToConverge(to_Rho, nmax_in);
  541. std::vector<std::vector<FloatType> > Pi(theta_points), Tau(theta_points);
  542. calcPiTauAllTheta(from_Theta, to_Theta, Pi, Tau);
  543. std::vector<std::vector<std::complex<FloatType> > > Psi(radius_points), D1n(radius_points),
  544. Zeta(radius_points), D3n(radius_points), PsiZeta(radius_points);
  545. calcRadialOnlyDependantFunctions(from_Rho, to_Rho,
  546. Psi, D1n, Zeta, D3n,
  547. nmax_in);
  548. // std::cout<<"Done evaluation of special functions."<<std::endl;
  549. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  550. double delta_Theta = eval_delta<double>(theta_points, from_Theta, to_Theta);
  551. double delta_Phi = eval_delta<double>(phi_points, from_Phi, to_Phi);
  552. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  553. std::vector<bool> isConvergedE = {false, false, false}, isConvergedH = {false, false, false};
  554. isConvergedE_ = {true, true, true}, isConvergedH_ = {true, true, true};
  555. for (int j=0; j < radius_points; j++) {
  556. auto Rho = from_Rho + j * delta_Rho;
  557. std::vector< std::complex<double> > Psi_dp = ConvertComplexVector<double>(Psi[j]);
  558. std::vector< std::complex<double> > Zeta_dp = ConvertComplexVector<double>(Zeta[j]);
  559. std::vector< std::complex<double> > D1n_dp = ConvertComplexVector<double>(D1n[j]);
  560. std::vector< std::complex<double> > D3n_dp = ConvertComplexVector<double>(D3n[j]);
  561. for (int i = 0; i < theta_points; i++) {
  562. auto Theta = from_Theta + i * delta_Theta;
  563. std::vector<double> Pi_dp = ConvertVector<double>(Pi[i]);
  564. std::vector<double> Tau_dp = ConvertVector<double>(Tau[i]);
  565. for (int k = 0; k < phi_points; k++) {
  566. auto Phi = from_Phi + k * delta_Phi;
  567. coords_polar_.push_back({Rho, Theta, Phi});
  568. std::vector<std::complex<double> > Es(3), Hs(3);
  569. calcFieldByComponents( Rho, Theta, Phi,
  570. Psi_dp, D1n_dp, Zeta_dp, D3n_dp, Pi_dp, Tau_dp,
  571. Es, Hs, isConvergedE, isConvergedH,
  572. isMarkUnconverged);
  573. UpdateConvergenceStatus(isConvergedE, isConvergedH);
  574. Es_.push_back(ConvertComplexVector<FloatType>(Es));
  575. Hs_.push_back(ConvertComplexVector<FloatType>(Hs));
  576. }
  577. }
  578. }
  579. convertFieldsFromSphericalToCartesian();
  580. }
  581. template <typename FloatType>
  582. void MultiLayerMie<FloatType>::UpdateConvergenceStatus(std::vector<bool> isConvergedE, std::vector<bool> isConvergedH) {
  583. for (int i = 0; i< 3; i++) isConvergedE_[i] = isConvergedE_[i] && isConvergedE[i];
  584. for (int i = 0; i< 3; i++) isConvergedH_[i] = isConvergedH_[i] && isConvergedH[i];
  585. }
  586. template <typename FloatType>
  587. bool MultiLayerMie<FloatType>::GetFieldConvergence () {
  588. bool convergence = true;
  589. for (auto conv:isConvergedE_) convergence = convergence && conv;
  590. for (auto conv:isConvergedH_) convergence = convergence && conv;
  591. return convergence;
  592. }
  593. template <typename FloatType>
  594. void MultiLayerMie<FloatType>::RunFieldCalculationCartesian(const int side_points,
  595. const double relative_side_length,
  596. const int plane_selected,
  597. const double at_x, const double at_y,
  598. const double at_z,
  599. const bool isMarkUnconverged,
  600. const int nmax_in) {
  601. SetMaxTerms(nmax_in);
  602. std::vector<FloatType> Xp(0), Yp(0), Zp(0);
  603. if (size_param_.size()<1) throw "Expect size_param_ to have at least one element before running a simulation";
  604. const FloatType max_coord_value = size_param_.back() * relative_side_length;
  605. const FloatType dx = max_coord_value*2/( (side_points<2 ? 2 : side_points) - 1.0);
  606. const FloatType dy = dx, dz = dx;
  607. auto push_coords = [&](
  608. const int nx, const int ny, const int nz,
  609. const FloatType xi, const FloatType yi, const FloatType zi) {
  610. for (int i = 0; i < nx; i++) {
  611. for (int j = 0; j < ny; j++) {
  612. for (int k = 0; k < nz; k++) {
  613. Xp.push_back(xi + static_cast<FloatType>(i) * dx);
  614. Yp.push_back(yi + static_cast<FloatType>(j) * dy);
  615. Zp.push_back(zi + static_cast<FloatType>(k) * dz);
  616. }
  617. }
  618. }
  619. };
  620. if (plane_selected == Planes::kEk ) {
  621. push_coords( side_points, 1, side_points,
  622. -max_coord_value+at_x, at_y, -max_coord_value+at_z);
  623. }
  624. if (plane_selected == Planes::kHk ) {
  625. push_coords( 1, side_points, side_points,
  626. at_x, -max_coord_value+at_y, -max_coord_value+at_z);
  627. }
  628. if (plane_selected == Planes::kEH) {
  629. push_coords( side_points, side_points, 1,
  630. -max_coord_value+at_x, -max_coord_value+at_y, at_z);
  631. }
  632. const unsigned int total_size = side_points*side_points;
  633. if (Xp.size() != total_size || Yp.size() != total_size || Zp.size() != total_size)
  634. throw std::invalid_argument("Error! Wrong dimension of field monitor points for cartesian grid!");
  635. SetFieldCoords({Xp, Yp, Zp});
  636. RunFieldCalculation(isMarkUnconverged);
  637. } // end of void MultiLayerMie<FloatType>::RunFieldCalculationCartesian(...)
  638. } // end of namespace nmie
  639. #endif // SRC_NMIE_NEARFIELD_HPP_