nmie-nearfield.hpp 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558
  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. std::complex<FloatType> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  81. const int L = refractive_index_.size();
  82. aln_.resize(L + 1);
  83. bln_.resize(L + 1);
  84. cln_.resize(L + 1);
  85. dln_.resize(L + 1);
  86. for (int l = 0; l <= L; l++) {
  87. aln_[l].resize(nmax_);
  88. bln_[l].resize(nmax_);
  89. cln_[l].resize(nmax_);
  90. dln_[l].resize(nmax_);
  91. }
  92. // Yang, paragraph under eq. A3
  93. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  94. for (int n = 0; n < nmax_; n++) {
  95. aln_[L][n] = an_[n];
  96. bln_[L][n] = bn_[n];
  97. cln_[L][n] = c_one;
  98. dln_[L][n] = c_one;
  99. }
  100. std::vector<std::complex<FloatType> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  101. std::vector<std::complex<FloatType> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  102. std::complex<FloatType> denomZeta, denomPsi, T1, T2, T3, T4;
  103. auto &m = refractive_index_;
  104. std::vector< std::complex<FloatType> > m1(L);
  105. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  106. m1[L - 1] = std::complex<FloatType> (1.0, 0.0);
  107. std::complex<FloatType> z, z1;
  108. for (int l = L - 1; l >= 0; l--) {
  109. if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!!
  110. for (int n = 0; n < nmax_; n++) {
  111. // aln
  112. aln_[l][n] = c_zero;
  113. // bln
  114. bln_[l][n] = c_zero;
  115. // cln
  116. cln_[l][n] = c_zero;
  117. // dln
  118. dln_[l][n] = c_zero;
  119. }
  120. } else { // Regular material, just do the calculation
  121. z = size_param_[l]*m[l];
  122. z1 = size_param_[l]*m1[l];
  123. calcD1D3(z, D1z, D3z);
  124. calcD1D3(z1, D1z1, D3z1);
  125. calcPsiZeta(z, Psiz, Zetaz);
  126. calcPsiZeta(z1, Psiz1, Zetaz1);
  127. for (int n = 0; n < nmax_; n++) {
  128. int n1 = n + 1;
  129. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  130. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  131. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  132. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  133. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  134. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  135. // aln
  136. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  137. // bln
  138. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  139. // cln
  140. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  141. // dln
  142. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  143. } // end of all n
  144. } // end PEC condition
  145. } // end of all l
  146. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  147. for (int n = 0; n < nmax_; ++n) {
  148. if (cabs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  149. else {
  150. //throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  151. std::cout<< std::setprecision(100)
  152. << "Warning: Potentially unstable calculation of aln[0]["
  153. << n << "] = "<< aln_[0][n] <<std::endl;
  154. aln_[0][n] = 0.0;
  155. }
  156. if (cabs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  157. else {
  158. //throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  159. std::cout<< std::setprecision(100)
  160. << "Warning: Potentially unstable calculation of bln[0]["
  161. << n << "] = "<< bln_[0][n] <<std::endl;
  162. bln_[0][n] = 0.0;
  163. }
  164. }
  165. isExpCoeffsCalc_ = true;
  166. } // end of void MultiLayerMie::calcExpanCoeffs()
  167. template <typename FloatType>
  168. void MultiLayerMie<FloatType>::convertFieldsFromSphericalToCartesian() {
  169. long total_points = coords_polar_.size();
  170. E_.clear(); H_.clear();
  171. for (int point=0; point < total_points; point++) {
  172. auto Theta = coords_polar_[point][1];
  173. auto Phi = coords_polar_[point][2];
  174. auto Es = Es_[point];
  175. auto Hs = Hs_[point];
  176. using nmm::sin;
  177. using nmm::cos;
  178. E_.push_back({ sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2],
  179. sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2],
  180. cos(Theta)*Es[0] - sin(Theta)*Es[1]});
  181. H_.push_back({ sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2],
  182. sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2],
  183. cos(Theta)*Hs[0] - sin(Theta)*Hs[1]});
  184. }
  185. } // end of void MultiLayerMie::convertFieldsFromSphericalToCartesian()
  186. //**********************************************************************************//
  187. // This function calculates the electric (E) and magnetic (H) fields inside and //
  188. // around the particle. //
  189. //
  190. // Main troubles of near-field evaluations originate from special functions
  191. // evaluation, so we expect that nmax needed for the convergence is the size
  192. // of Psi vector.
  193. // //
  194. // Input parameters (coordinates of the point): //
  195. // Rho: Radial distance //
  196. // Phi: Azimuthal angle //
  197. // Theta: Polar angle //
  198. // mode_n: mode order. //
  199. // -1 - use all modes (all_) //
  200. // 1 - use dipole mode only //
  201. // 2 - use quadrupole mode only //
  202. // ... //
  203. // mode_type: only used when mode_n != -1 //
  204. // 0 - electric only //
  205. // 1 - magnetic only //
  206. // //
  207. // //
  208. // Output parameters: //
  209. // E, H: Complex electric and magnetic fields //
  210. //**********************************************************************************//
  211. template <typename FloatType> template <typename evalType>
  212. void MultiLayerMie<FloatType>::calcFieldByComponents(const evalType Rho,
  213. const evalType Theta, const evalType Phi,
  214. const std::vector<std::complex<evalType> > &Psi,
  215. const std::vector<std::complex<evalType> > &D1n,
  216. const std::vector<std::complex<evalType> > &Zeta,
  217. const std::vector<std::complex<evalType> > &D3n,
  218. const std::vector<evalType> &Pi,
  219. const std::vector<evalType> &Tau,
  220. std::vector<std::complex<evalType> > &E,
  221. std::vector<std::complex<evalType> > &H) {
  222. auto nmax = Psi.size() - 1;
  223. std::complex<evalType> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  224. // Vector containing precomputed integer powers of i to avoid computation
  225. std::vector<std::complex<evalType> > ipow = {c_one, c_i, -c_one, -c_i};
  226. std::vector<std::complex<evalType> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  227. std::vector<std::complex<evalType> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  228. std::complex<evalType> ml;
  229. // Initialize E and H
  230. for (int i = 0; i < 3; i++) {
  231. E[i] = c_zero;
  232. H[i] = c_zero;
  233. }
  234. unsigned int l;
  235. GetIndexAtRadius(Rho, ml, l);
  236. for (unsigned int n = 0; n < nmax; n++) {
  237. int n1 = n + 1;
  238. auto rn = static_cast<evalType>(n1);
  239. // using BH 4.12 and 4.50
  240. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  241. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  242. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  243. std::complex<evalType> En = ipow[n1 % 4]
  244. *static_cast<evalType>((rn + rn + 1.0)/(rn*rn + rn));
  245. std::complex<evalType> Ediff, Hdiff;
  246. for (int i = 0; i < 3; i++) {
  247. Ediff = En*( cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  248. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  249. Hdiff = En*( -dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  250. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  251. if (nmm::isnan(Ediff.real()) || nmm::isnan(Ediff.imag()) ||
  252. nmm::isnan(Hdiff.real()) || nmm::isnan(Hdiff.imag())
  253. ) {
  254. std::cout << "Unexpected truncation during near-field evaluation at n = "<< n
  255. << " (of total nmax = "<<nmax<<")!!!"<<std::endl;
  256. break;
  257. }
  258. if (mode_n_ == Modes::kAll) {
  259. // electric field E [V m - 1] = EF*E0
  260. E[i] += Ediff;
  261. H[i] += Hdiff;
  262. continue;
  263. }
  264. if (n1 == mode_n_) {
  265. if (mode_type_ == Modes::kElectric || mode_type_ == Modes::kAll) {
  266. E[i] += En*( -c_i*dln_[l][n]*N1e1n[i]
  267. + c_i*aln_[l][n]*N3e1n[i]);
  268. H[i] += En*(-dln_[l][n]*M1e1n[i]
  269. +aln_[l][n]*M3e1n[i]);
  270. //std::cout << mode_n_;
  271. }
  272. if (mode_type_ == Modes::kMagnetic || mode_type_ == Modes::kAll) {
  273. E[i] += En*( cln_[l][n]*M1o1n[i]
  274. - bln_[l][n]*M3o1n[i]);
  275. H[i] += En*( -c_i*cln_[l][n]*N1o1n[i]
  276. + c_i*bln_[l][n]*N3o1n[i]);
  277. //std::cout << mode_n_;
  278. }
  279. //std::cout << std::endl;
  280. }
  281. //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=====*****=====");
  282. }
  283. if (nmm::isnan(Ediff.real()) || nmm::isnan(Ediff.imag()) ||
  284. nmm::isnan(Hdiff.real()) || nmm::isnan(Hdiff.imag())
  285. ) break;
  286. } // end of for all n
  287. // magnetic field
  288. std::complex<evalType> hffact = ml/static_cast<evalType>(cc_*mu_);
  289. for (int i = 0; i < 3; i++) {
  290. H[i] = hffact*H[i];
  291. }
  292. } // end of MultiLayerMie::calcFieldByComponents(...)
  293. //**********************************************************************************//
  294. // This function calculates complex electric and magnetic field in the surroundings //
  295. // and inside the particle. //
  296. // //
  297. // Input parameters: //
  298. // L: Number of layers //
  299. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  300. // x: Array containing the size parameters of the layers [0..L-1] //
  301. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  302. // nmax: Maximum number of multipolar expansion terms to be used for the //
  303. // calculations. Only use it if you know what you are doing, otherwise //
  304. // set this parameter to 0 (zero) and the function will calculate it. //
  305. // ncoord: Number of coordinate points //
  306. // Coords: Array containing all coordinates where the complex electric and //
  307. // magnetic fields will be calculated //
  308. // mode_n: mode order. //
  309. // -1 - use all modes (all_) //
  310. // 1 - use dipole mode only //
  311. // 2 - use quadrupole mode only //
  312. // ... //
  313. // mode_type: only used when mode_n != -1 //
  314. // 0 - electric only //
  315. // 1 - magnetic only //
  316. // //
  317. // Output parameters: //
  318. // E, H: Complex electric and magnetic field at the provided coordinates //
  319. // //
  320. // Return value: //
  321. // Number of multipolar expansion terms used for the calculations //
  322. //**********************************************************************************//
  323. template <typename FloatType>
  324. void MultiLayerMie<FloatType>::RunFieldCalculation() {
  325. // Calculate scattering coefficients an_ and bn_
  326. calcScattCoeffs();
  327. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  328. calcExpanCoeffs();
  329. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  330. long total_points = coords_[0].size();
  331. for (int point = 0; point < total_points; point++) {
  332. const FloatType &Xp = coords_[0][point];
  333. const FloatType &Yp = coords_[1][point];
  334. const FloatType &Zp = coords_[2][point];
  335. // Convert to spherical coordinates
  336. auto Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  337. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  338. auto Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0;
  339. // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc.
  340. auto Phi = nmm::atan2(Yp,Xp);
  341. coords_polar_.push_back({Rho, Theta, Phi});
  342. // Avoid convergence problems due to Rho too small
  343. if (Rho < 1e-5) Rho = 1e-5;
  344. //*******************************************************//
  345. // external scattering field = incident + scattered //
  346. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  347. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  348. //*******************************************************//
  349. // This array contains the fields in spherical coordinates
  350. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  351. // Do the actual calculation of electric and magnetic field
  352. std::vector<std::complex<FloatType> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  353. std::vector<FloatType> Pi(nmax_), Tau(nmax_);
  354. std::complex<FloatType> ml;
  355. GetIndexAtRadius(Rho, ml);
  356. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  357. calcD1D3(Rho*ml, D1n, D3n);
  358. // Calculate Ricatti-Bessel functions
  359. calcPsiZeta(Rho*ml, Psi, Zeta);
  360. // Calculate angular functions Pi and Tau
  361. calcPiTau(nmm::cos(Theta), Pi, Tau);
  362. calcFieldByComponents(Rho, Theta, Phi, Psi, D1n, Zeta, D3n, Pi, Tau, Es, Hs);
  363. Es_.push_back(Es);
  364. Hs_.push_back(Hs);
  365. } // end of for all field coordinates
  366. convertFieldsFromSphericalToCartesian();
  367. } // end of MultiLayerMie::RunFieldCalculation()
  368. template <typename FloatType>
  369. double eval_delta(const unsigned int steps, const double from_value, const double to_value) {
  370. auto delta = std::abs(from_value - to_value);
  371. if (steps < 2) return delta;
  372. delta /= static_cast<double>(steps-1);
  373. // We have a limited double precision evaluation of special functions, typically it is 1e-10.
  374. if ( (2.*delta)/std::abs(from_value+to_value) < 1e-9)
  375. throw std::invalid_argument("Error! The step is too fine, not supported!");
  376. return delta;
  377. }
  378. // ml - refractive index
  379. // l - Layer number
  380. template <typename FloatType> template <typename evalType>
  381. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  382. std::complex<evalType> &ml,
  383. unsigned int &l) {
  384. l = 0;
  385. if (Rho > size_param_.back()) {
  386. l = size_param_.size();
  387. ml = std::complex<FloatType>(1.0, 0.0);
  388. } else {
  389. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  390. if (Rho <= size_param_[i]) {
  391. l = i;
  392. }
  393. }
  394. ml = refractive_index_[l];
  395. }
  396. }
  397. template <typename FloatType> template <typename evalType>
  398. void MultiLayerMie<FloatType>::GetIndexAtRadius(const evalType Rho,
  399. std::complex<evalType> &ml) {
  400. unsigned int l;
  401. GetIndexAtRadius(Rho, ml, l);
  402. }
  403. template <typename FloatType>
  404. void MultiLayerMie<FloatType>::calcMieSeriesNeededToConverge(const FloatType Rho) {
  405. auto required_near_field_nmax = calcNmax(Rho);
  406. SetMaxTerms(required_near_field_nmax);
  407. // Calculate scattering coefficients an_ and bn_
  408. calcScattCoeffs();
  409. // We might be limited with available machine precision
  410. available_maximal_nmax_ = nmax_;
  411. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  412. calcExpanCoeffs();
  413. }
  414. template <typename FloatType>
  415. void MultiLayerMie<FloatType>::calcRadialOnlyDependantFunctions(const double from_Rho, const double to_Rho,
  416. const bool isIgnoreAvailableNmax,
  417. std::vector<std::vector<std::complex<FloatType> > > &Psi,
  418. std::vector<std::vector<std::complex<FloatType> > > &D1n,
  419. std::vector<std::vector<std::complex<FloatType> > > &Zeta,
  420. std::vector<std::vector<std::complex<FloatType> > > &D3n) {
  421. auto radius_points = Psi.size();
  422. std::vector<std::vector<std::complex<FloatType> > > PsiZeta(radius_points);
  423. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  424. for (unsigned int j=0; j < radius_points; j++) {
  425. auto Rho = static_cast<FloatType>(from_Rho + j*delta_Rho);
  426. // if (Rho < 1e-5) Rho = 1e-5; // TODO do we need this?.
  427. int near_field_nmax = calcNmax(Rho);
  428. // Skip if not enough terms in Mie series (i.e. required near field nmax > available terms )
  429. if (near_field_nmax > available_maximal_nmax_ && !isIgnoreAvailableNmax) continue;
  430. if (near_field_nmax > available_maximal_nmax_) near_field_nmax = available_maximal_nmax_;
  431. Psi[j].resize(near_field_nmax + 1); D1n[j].resize(near_field_nmax + 1);
  432. Zeta[j].resize(near_field_nmax + 1); D3n[j].resize(near_field_nmax + 1);
  433. PsiZeta[j].resize(near_field_nmax + 1);
  434. std::complex<FloatType> ml;
  435. GetIndexAtRadius(Rho, ml);
  436. auto z = Rho*ml;
  437. evalDownwardD1(z, D1n[j]);
  438. evalUpwardPsi(z, D1n[j], Psi[j]);
  439. evalUpwardD3 (z, D1n[j], D3n[j], PsiZeta[j]);
  440. for (unsigned int k = 0; k < Zeta[j].size(); k++) {
  441. Zeta[j][k] = PsiZeta[j][k]/Psi[j][k];
  442. }
  443. }
  444. }
  445. // input parameters:
  446. // outer_arc_points: will be increased to the nearest power of 2.
  447. template <typename FloatType>
  448. void MultiLayerMie<FloatType>::RunFieldCalculationPolar(const int outer_arc_points,
  449. const int radius_points,
  450. const double from_Rho, const double to_Rho,
  451. const double from_Theta, const double to_Theta,
  452. const double from_Phi, const double to_Phi,
  453. const bool isIgnoreAvailableNmax) {
  454. if (from_Rho > to_Rho || from_Theta > to_Theta || from_Phi > to_Phi
  455. || outer_arc_points < 1 || radius_points < 1
  456. || from_Rho < 0.)
  457. throw std::invalid_argument("Error! Invalid argument for RunFieldCalculationPolar() !");
  458. int theta_points = 0, phi_points = 0;
  459. if (to_Theta-from_Theta > to_Phi-from_Phi) {
  460. theta_points = outer_arc_points;
  461. phi_points = static_cast<int>((to_Phi-from_Phi)/(to_Theta-from_Theta) * outer_arc_points);
  462. } else {
  463. phi_points = outer_arc_points;
  464. theta_points = static_cast<int>((to_Theta-from_Theta)/(to_Phi-from_Phi) * outer_arc_points);
  465. }
  466. if (theta_points == 0) theta_points = 1;
  467. if (phi_points == 0) phi_points = 1;
  468. calcMieSeriesNeededToConverge(to_Rho);
  469. std::vector<std::vector<FloatType> > Pi(theta_points), Tau(theta_points);
  470. calcPiTauAllTheta(from_Theta, to_Theta, Pi, Tau);
  471. std::vector<std::vector<std::complex<FloatType> > > Psi(radius_points), D1n(radius_points),
  472. Zeta(radius_points), D3n(radius_points), PsiZeta(radius_points);
  473. calcRadialOnlyDependantFunctions(from_Rho, to_Rho, isIgnoreAvailableNmax,
  474. Psi, D1n, Zeta, D3n);
  475. double delta_Rho = eval_delta<double>(radius_points, from_Rho, to_Rho);
  476. double delta_Theta = eval_delta<double>(theta_points, from_Theta, to_Theta);
  477. double delta_Phi = eval_delta<double>(phi_points, from_Phi, to_Phi);
  478. Es_.clear(); Hs_.clear(); coords_polar_.clear();
  479. for (int j=0; j < radius_points; j++) {
  480. auto Rho = static_cast<FloatType>(from_Rho + j * delta_Rho);
  481. for (int i = 0; i < outer_arc_points; i++) {
  482. auto Theta = static_cast<FloatType>(from_Theta + i * delta_Theta);
  483. for (int k = 0; k < outer_arc_points; k++) {
  484. auto Phi = static_cast<FloatType>(from_Phi + k * delta_Phi);
  485. coords_polar_.push_back({Rho, Theta, Phi});
  486. // This array contains the fields in spherical coordinates
  487. std::vector<std::complex<FloatType> > Es(3), Hs(3);
  488. calcFieldByComponents(Rho, Theta, Phi,
  489. Psi[j], D1n[j], Zeta[j], D3n[j],
  490. Pi[i], Tau[i], Es, Hs);
  491. Es_.push_back(Es);
  492. Hs_.push_back(Hs);
  493. }
  494. }
  495. }
  496. convertFieldsFromSphericalToCartesian();
  497. }
  498. // Python interface
  499. } // end of namespace nmie
  500. #endif // SRC_NMIE_NEARFIELD_HPP_