nmie-nearfield.hpp 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709
  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] << " which is expected to be exact zero!"<<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] << " which is expected to be exact zero!" <<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. std::vector< std::complex<evalType> > Ediff_prev = {{0.,0.},{0.,0.},{0.,0.}},
  255. Hdiff_prev = {{0.,0.},{0.,0.},{0.,0.}};
  256. for (unsigned int n = 0; n < nmax; n++) {
  257. if ( isConvergedE[0] && isConvergedE[1] && isConvergedE[2]
  258. && isConvergedH[0] && isConvergedH[1] && isConvergedH[2]) {
  259. std::cout<<"Near-field early convergence at nmax = "<<n+1<<std::endl;
  260. break;
  261. }
  262. int n1 = n + 1;
  263. auto rn = static_cast<evalType>(n1);
  264. // using BH 4.12 and 4.50
  265. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  266. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  267. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  268. std::complex<evalType> En = ipow[n1 % 4]
  269. *static_cast<evalType>((rn + rn + 1.0)/(rn*rn + rn));
  270. std::complex<evalType> Ediff, Hdiff;
  271. std::complex<FloatType> Ediff_ft, Hdiff_ft;
  272. auto aln = ConvertComplex<evalType>(aln_[l][n]);
  273. auto bln = ConvertComplex<evalType>(bln_[l][n]);
  274. auto cln = ConvertComplex<evalType>(cln_[l][n]);
  275. auto dln = ConvertComplex<evalType>(dln_[l][n]);
  276. for (int i = 0; i < 3; i++) {
  277. if (isConvergedE[i] && isConvergedH[i]) continue; // TODO is it safe?
  278. Ediff = En*( cln*M1o1n[i] - c_i*dln*N1e1n[i]
  279. + c_i*aln*N3e1n[i] - bln*M3o1n[i]);
  280. Hdiff = En*( -dln*M1e1n[i] - c_i*cln*N1o1n[i]
  281. + c_i*bln*N3o1n[i] + aln*M3e1n[i]);
  282. Ediff_ft = ConvertComplex<FloatType>(Ediff);
  283. Hdiff_ft = ConvertComplex<FloatType>(Hdiff);
  284. if ( nmm::isnan(Ediff_ft.real()) || nmm::isnan(Ediff_ft.imag()) ||
  285. nmm::isnan(Hdiff_ft.real()) || nmm::isnan(Hdiff_ft.imag()) ) {
  286. std::cout << "Unexpected truncation during near-field evaluation at n = "<< n
  287. << " (of total nmax = "<<nmax<<")!!!"<<std::endl;
  288. break;
  289. }
  290. if (n>0) {
  291. if (
  292. (cabs(Ediff_prev[i]) <= cabs(E[i]) * nearfield_convergence_threshold_)
  293. && (cabs(Ediff) <= cabs(E[i]) * nearfield_convergence_threshold_)
  294. ) isConvergedE[i] = true;
  295. if (
  296. (cabs(Hdiff_prev[i]) <= cabs(H[i]) * nearfield_convergence_threshold_)
  297. && (cabs(Hdiff) <= cabs(H[i]) * nearfield_convergence_threshold_)
  298. ) isConvergedH[i] = true;
  299. }
  300. Ediff_prev[i] = Ediff;
  301. Hdiff_prev[i] = Hdiff;
  302. if ((!isConvergedH[i] || !isConvergedE[i]) && n==nmax-1 && GetFieldConvergence()) {
  303. std::cout<<"Econv:"<<cabs(Ediff)/cabs(E[i])<<" Hconv:"<<cabs(Hdiff)/cabs(H[i])<<std::endl;
  304. }
  305. if (mode_n_ == Modes::kAll) {
  306. // electric field E [V m - 1] = EF*E0
  307. E[i] += Ediff;
  308. H[i] += Hdiff;
  309. continue;
  310. }
  311. if (n == 0) {
  312. }
  313. if (n1 == mode_n_) {
  314. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  315. E[i] += En*( -c_i*dln*N1e1n[i]
  316. + c_i*aln*N3e1n[i]);
  317. H[i] += En*(-dln*M1e1n[i]
  318. +aln*M3e1n[i]);
  319. //std::cout << mode_n_;
  320. }
  321. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  322. E[i] += En*( cln*M1o1n[i]
  323. - bln*M3o1n[i]);
  324. H[i] += En*( -c_i*cln*N1o1n[i]
  325. + c_i*bln*N3o1n[i]);
  326. //std::cout << mode_n_;
  327. }
  328. //std::cout << std::endl;
  329. }
  330. //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=====*****=====");
  331. }
  332. if (nmm::isnan(Ediff_ft.real()) || nmm::isnan(Ediff_ft.imag()) ||
  333. nmm::isnan(Hdiff_ft.real()) || nmm::isnan(Hdiff_ft.imag())
  334. ) break;
  335. } // end of for all n
  336. // Add the incident field
  337. if(l==L) {
  338. const auto z = Rho*cos_t(Theta);
  339. const auto Ex = std::complex<evalType>(cos_t(z), sin_t(z));
  340. E[0] += Ex*cos_t(Phi)*sin_t(Theta);
  341. E[1] += Ex*cos_t(Phi)*cos_t(Theta);
  342. E[2] += -Ex*sin_t(Phi);
  343. const auto Hy = Ex;
  344. H[0] += Hy*sin_t(Theta)*sin_t(Phi);
  345. H[1] += Hy*cos_t(Theta)*sin_t(Phi);
  346. H[2] += Hy*cos_t(Phi);
  347. }
  348. if( (!isConvergedE[0] || !isConvergedE[1] ||!isConvergedE[2] ||
  349. !isConvergedH[0] || !isConvergedH[1] ||!isConvergedH[2] ) && GetFieldConvergence()) {
  350. std::cout << "Field evaluation failed to converge an nmax = "<< nmax << std::endl;
  351. std::cout << "Near-field convergence threshold: "<<nearfield_convergence_threshold_<<std::endl;
  352. if (isMarkUnconverged) { //mark as NaN
  353. for(auto &ee :E) ee /= c_zero;
  354. for(auto &ee :H) ee /= c_zero;
  355. }
  356. }
  357. // magnetic field
  358. std::complex<evalType> hffact = ml/static_cast<evalType>(nmie::cc_*nmie::mu_);
  359. for (int i = 0; i < 3; i++) {
  360. H[i] = hffact*H[i];
  361. }
  362. } // end of MultiLayerMie::calcFieldByComponents(...)
  363. //**********************************************************************************//
  364. // This function calculates complex electric and magnetic field in the surroundings //
  365. // and inside the particle. //
  366. // //
  367. // Input parameters: //
  368. // L: Number of layers //
  369. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  370. // x: Array containing the size parameters of the layers [0..L-1] //
  371. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  372. // nmax: Maximum number of multipolar expansion terms to be used for the //
  373. // calculations. Only use it if you know what you are doing, otherwise //
  374. // set this parameter to 0 (zero) and the function will calculate it. //
  375. // ncoord: Number of coordinate points //
  376. // Coords: Array containing all coordinates where the complex electric and //
  377. // magnetic fields will be calculated //
  378. // mode_n: mode order. //
  379. // -1 - use all modes (all_) //
  380. // 1 - use dipole mode only //
  381. // 2 - use quadrupole mode only //
  382. // ... //
  383. // mode_type: only used when mode_n != -1 //
  384. // 0 - electric only //
  385. // 1 - magnetic only //
  386. // //
  387. // Output parameters: //
  388. // E, H: Complex electric and magnetic field at the provided coordinates //
  389. // //
  390. // Return value: //
  391. // Number of multipolar expansion terms used for the calculations //
  392. //**********************************************************************************//
  393. template <typename FloatType>
  394. void MultiLayerMie<FloatType>::RunFieldCalculation(bool isMarkUnconverged) {
  395. // Calculate scattering coefficients an_ and bn_
  396. calcScattCoeffs();
  397. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  398. calcExpanCoeffs();
  399. std::vector<bool> isConvergedE = {false, false, false}, isConvergedH = {false, false, false};
  400. isConvergedE_ = {true, true, true}, isConvergedH_ = {true, true, true};
  401. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  402. long total_points = coords_[0].size();
  403. for (int point = 0; point < total_points; point++) {
  404. const FloatType &Xp = coords_[0][point];
  405. const FloatType &Yp = coords_[1][point];
  406. const FloatType &Zp = coords_[2][point];
  407. // Convert to spherical coordinates
  408. auto Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  409. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  410. auto Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0;
  411. // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc.
  412. auto Phi = nmm::atan2(Yp,Xp);
  413. coords_polar_.push_back({Rho, Theta, Phi});
  414. // Avoid convergence problems due to Rho too small
  415. if (Rho < 1e-5) Rho = 1e-5;
  416. //*******************************************************//
  417. // external scattering field = incident + scattered //
  418. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  419. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  420. //*******************************************************//
  421. // This array contains the fields in spherical coordinates
  422. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  423. // Do the actual calculation of electric and magnetic field
  424. std::vector<std::complex<FloatType> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  425. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  426. std::complex<FloatType> ml;
  427. GetIndexAtRadius(Rho, ml);
  428. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  429. calcD1D3(Rho*ml, D1n, D3n);
  430. // Calculate Ricatti-Bessel functions
  431. calcPsiZeta(Rho*ml, Psi, Zeta);
  432. // Calculate angular functions Pi and Tau
  433. calcPiTau(nmm::cos(Theta), Pi, Tau);
  434. calcFieldByComponents(Rho, Theta, Phi, Psi, D1n, Zeta, D3n, Pi, Tau, Es, Hs,
  435. isConvergedE, isConvergedH, isMarkUnconverged);
  436. UpdateConvergenceStatus(isConvergedE, isConvergedH);
  437. Es_.push_back(Es);
  438. Hs_.push_back(Hs);
  439. } // end of for all field coordinates
  440. convertFieldsFromSphericalToCartesian();
  441. } // end of MultiLayerMie::RunFieldCalculation()
  442. // TODO do we really need this eval_delta()?
  443. template <typename FloatType>
  444. double eval_delta(const unsigned int steps, const double from_value, const double to_value) {
  445. auto delta = std::abs(from_value - to_value);
  446. if (steps < 2) return delta;
  447. delta /= static_cast<double>(steps-1);
  448. // We have a limited double precision evaluation of special functions, typically it is 1e-10.
  449. if ( (2.*delta)/std::abs(from_value+to_value) < 1e-9)
  450. throw std::invalid_argument("Error! The step is too fine, not supported!");
  451. return delta;
  452. }
  453. // ml - refractive index
  454. // l - Layer number
  455. template <typename FloatType> template <typename evalType>
  456. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  457. std::complex<evalType> &ml,
  458. unsigned int &l) {
  459. l = 0;
  460. if (Rho > size_param_.back()) {
  461. l = size_param_.size();
  462. ml = std::complex<evalType>(1.0, 0.0);
  463. } else {
  464. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  465. if (Rho <= size_param_[i]) {
  466. l = i;
  467. }
  468. }
  469. ml = ConvertComplex<evalType>(refractive_index_[l]);
  470. }
  471. }
  472. template <typename FloatType> template <typename evalType>
  473. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  474. std::complex<evalType> &ml) {
  475. unsigned int l;
  476. GetIndexAtRadius(Rho, ml, l);
  477. }
  478. template <typename FloatType>
  479. void MultiLayerMie<FloatType>::calcMieSeriesNeededToConverge(const FloatType Rho, int nmax_in) {
  480. if (nmax_in < 1) {
  481. auto required_near_field_nmax = calcNmax(Rho);
  482. SetMaxTerms(required_near_field_nmax);
  483. } else {
  484. SetMaxTerms(nmax_in);
  485. }
  486. // Calculate scattering coefficients an_ and bn_
  487. calcScattCoeffs();
  488. // We might be limited with available machine precision
  489. available_maximal_nmax_ = nmax_;
  490. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  491. calcExpanCoeffs();
  492. }
  493. template <typename FloatType>
  494. void MultiLayerMie<FloatType>::calcRadialOnlyDependantFunctions(const double from_Rho, const double to_Rho,
  495. std::vector<std::vector<std::complex<FloatType> > > &Psi,
  496. std::vector<std::vector<std::complex<FloatType> > > &D1n,
  497. std::vector<std::vector<std::complex<FloatType> > > &Zeta,
  498. std::vector<std::vector<std::complex<FloatType> > > &D3n,
  499. int nmax_in) {
  500. auto radius_points = Psi.size();
  501. std::vector<std::vector<std::complex<FloatType> > > PsiZeta(radius_points);
  502. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  503. for (unsigned int j=0; j < radius_points; j++) {
  504. auto Rho = static_cast<FloatType>(from_Rho + j*delta_Rho);
  505. // if (Rho < 1e-5) Rho = 1e-5; // TODO do we need this?.
  506. int near_field_nmax = nmax_in;
  507. if (nmax_in < 1) near_field_nmax = calcNmax(Rho);
  508. // Skip if not enough terms in Mie series (i.e. required near field nmax > available terms )
  509. if (near_field_nmax > available_maximal_nmax_) near_field_nmax = available_maximal_nmax_;
  510. Psi[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0)); D1n[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  511. Zeta[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0)); D3n[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  512. PsiZeta[j].resize(near_field_nmax + 1, static_cast<FloatType>(0.0));
  513. std::complex<FloatType> ml;
  514. GetIndexAtRadius(Rho, ml);
  515. auto z = Rho*ml;
  516. evalDownwardD1(z, D1n[j]);
  517. evalUpwardPsi(z, D1n[j], Psi[j]);
  518. evalUpwardD3 (z, D1n[j], D3n[j], PsiZeta[j]);
  519. for (unsigned int k = 0; k < Zeta[j].size(); k++) {
  520. Zeta[j][k] = PsiZeta[j][k]/Psi[j][k];
  521. }
  522. }
  523. }
  524. // input parameters:
  525. // outer_arc_points: will be increased to the nearest power of 2.
  526. template <typename FloatType>
  527. void MultiLayerMie<FloatType>::RunFieldCalculationPolar(const int outer_arc_points,
  528. const int radius_points,
  529. const double from_Rho, const double to_Rho,
  530. const double from_Theta, const double to_Theta,
  531. const double from_Phi, const double to_Phi,
  532. const bool isMarkUnconverged,
  533. int nmax_in) {
  534. if (from_Rho > to_Rho || from_Theta > to_Theta || from_Phi > to_Phi
  535. || outer_arc_points < 1 || radius_points < 1
  536. || from_Rho < 0.)
  537. throw std::invalid_argument("Error! Invalid argument for RunFieldCalculationPolar() !");
  538. // auto nmax_old = nmax_;
  539. int theta_points = 0, phi_points = 0;
  540. if (to_Theta-from_Theta > to_Phi-from_Phi) {
  541. theta_points = outer_arc_points;
  542. phi_points = static_cast<int>((to_Phi-from_Phi)/(to_Theta-from_Theta) * outer_arc_points);
  543. } else {
  544. phi_points = outer_arc_points;
  545. theta_points = static_cast<int>((to_Theta-from_Theta)/(to_Phi-from_Phi) * outer_arc_points);
  546. }
  547. if (theta_points == 0) theta_points = 1;
  548. if (phi_points == 0) phi_points = 1;
  549. calcMieSeriesNeededToConverge(to_Rho, nmax_in);
  550. std::vector<std::vector<FloatType> > Pi(theta_points), Tau(theta_points);
  551. calcPiTauAllTheta(from_Theta, to_Theta, Pi, Tau);
  552. std::vector<std::vector<std::complex<FloatType> > > Psi(radius_points), D1n(radius_points),
  553. Zeta(radius_points), D3n(radius_points), PsiZeta(radius_points);
  554. calcRadialOnlyDependantFunctions(from_Rho, to_Rho,
  555. Psi, D1n, Zeta, D3n,
  556. nmax_in);
  557. // std::cout<<"Done evaluation of special functions."<<std::endl;
  558. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  559. double delta_Theta = eval_delta<double>(theta_points, from_Theta, to_Theta);
  560. double delta_Phi = eval_delta<double>(phi_points, from_Phi, to_Phi);
  561. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  562. std::vector<bool> isConvergedE = {false, false, false}, isConvergedH = {false, false, false};
  563. isConvergedE_ = {true, true, true}, isConvergedH_ = {true, true, true};
  564. for (int j=0; j < radius_points; j++) {
  565. auto Rho = from_Rho + j * delta_Rho;
  566. std::vector< std::complex<double> > Psi_dp = ConvertComplexVector<double>(Psi[j]);
  567. std::vector< std::complex<double> > Zeta_dp = ConvertComplexVector<double>(Zeta[j]);
  568. std::vector< std::complex<double> > D1n_dp = ConvertComplexVector<double>(D1n[j]);
  569. std::vector< std::complex<double> > D3n_dp = ConvertComplexVector<double>(D3n[j]);
  570. for (int i = 0; i < theta_points; i++) {
  571. auto Theta = from_Theta + i * delta_Theta;
  572. std::vector<double> Pi_dp = ConvertVector<double>(Pi[i]);
  573. std::vector<double> Tau_dp = ConvertVector<double>(Tau[i]);
  574. for (int k = 0; k < phi_points; k++) {
  575. auto Phi = from_Phi + k * delta_Phi;
  576. coords_polar_.push_back({Rho, Theta, Phi});
  577. std::vector<std::complex<double> > Es(3), Hs(3);
  578. calcFieldByComponents( Rho, Theta, Phi,
  579. Psi_dp, D1n_dp, Zeta_dp, D3n_dp, Pi_dp, Tau_dp,
  580. Es, Hs, isConvergedE, isConvergedH,
  581. isMarkUnconverged);
  582. UpdateConvergenceStatus(isConvergedE, isConvergedH);
  583. Es_.push_back(ConvertComplexVector<FloatType>(Es));
  584. Hs_.push_back(ConvertComplexVector<FloatType>(Hs));
  585. }
  586. }
  587. }
  588. convertFieldsFromSphericalToCartesian();
  589. }
  590. template <typename FloatType>
  591. void MultiLayerMie<FloatType>::UpdateConvergenceStatus(std::vector<bool> isConvergedE, std::vector<bool> isConvergedH) {
  592. for (int i = 0; i< 3; i++) isConvergedE_[i] = isConvergedE_[i] && isConvergedE[i];
  593. for (int i = 0; i< 3; i++) isConvergedH_[i] = isConvergedH_[i] && isConvergedH[i];
  594. }
  595. template <typename FloatType>
  596. bool MultiLayerMie<FloatType>::GetFieldConvergence () {
  597. bool convergence = true;
  598. for (auto conv:isConvergedE_) convergence = convergence && conv;
  599. for (auto conv:isConvergedH_) convergence = convergence && conv;
  600. return convergence;
  601. }
  602. template <typename FloatType>
  603. void MultiLayerMie<FloatType>::RunFieldCalculationCartesian(const int first_side_points,
  604. const int second_side_points,
  605. const double relative_side_length,
  606. const int plane_selected,
  607. const double at_x, const double at_y,
  608. const double at_z,
  609. const bool isMarkUnconverged,
  610. const int nmax_in) {
  611. SetMaxTerms(nmax_in);
  612. std::vector<FloatType> Xp(0), Yp(0), Zp(0);
  613. if (size_param_.size()<1) throw "Expect size_param_ to have at least one element before running a simulation";
  614. const FloatType total_R = size_param_.back();
  615. const FloatType second_side_max_coord_value = total_R * relative_side_length;
  616. // TODO add test if side_1_points <= 1 or side_2_points <= 1
  617. const FloatType space_step = second_side_max_coord_value*2/( (second_side_points<2 ? 2 : second_side_points) - 1.0);
  618. auto push_coords = [&](const int nx, const int ny, const int nz) {
  619. const FloatType xi = at_x*total_R - space_step*(nx-1)/2;
  620. const FloatType yi = at_y*total_R - space_step*(ny-1)/2;
  621. const FloatType zi = at_z*total_R - space_step*(nz-1)/2;
  622. for (int i = 0; i < nx; i++) {
  623. for (int j = 0; j < ny; j++) {
  624. for (int k = 0; k < nz; k++) {
  625. Xp.push_back(xi + static_cast<FloatType>(i) * space_step);
  626. Yp.push_back(yi + static_cast<FloatType>(j) * space_step);
  627. Zp.push_back(zi + static_cast<FloatType>(k) * space_step);
  628. }
  629. }
  630. }
  631. };
  632. // TODO add test to check that side_2_points is for z-axis
  633. if (plane_selected == Planes::kEk) push_coords(first_side_points, 1, second_side_points);
  634. if (plane_selected == Planes::kHk) push_coords(1, first_side_points, second_side_points);
  635. if (plane_selected == Planes::kEH) push_coords(first_side_points, second_side_points, 1);
  636. const unsigned int total_size = first_side_points*second_side_points;
  637. if (Xp.size() != total_size || Yp.size() != total_size || Zp.size() != total_size)
  638. throw std::invalid_argument("Error! Wrong dimension of field monitor points for cartesian grid!");
  639. SetFieldCoords({Xp, Yp, Zp});
  640. RunFieldCalculation(isMarkUnconverged);
  641. } // end of void MultiLayerMie<FloatType>::RunFieldCalculationCartesian(...)
  642. } // end of namespace nmie
  643. #endif // SRC_NMIE_NEARFIELD_HPP_