nmie-wrapper.cc 50 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089
  1. ///
  2. /// @file nmie-wrapper.cc
  3. /// @author Ladutenko Konstantin <kostyfisik at gmail (.) com>
  4. /// @date Tue Sep 3 00:38:27 2013
  5. /// @copyright 2013 Ladutenko Konstantin
  6. ///
  7. /// nmie-wrapper is free software: you can redistribute it and/or modify
  8. /// it under the terms of the GNU General Public License as published by
  9. /// the Free Software Foundation, either version 3 of the License, or
  10. /// (at your option) any later version.
  11. ///
  12. /// nmie-wrapper is distributed in the hope that it will be useful,
  13. /// but WITHOUT ANY WARRANTY; without even the implied warranty of
  14. /// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  15. /// GNU General Public License for more details.
  16. ///
  17. /// You should have received a copy of the GNU General Public License
  18. /// along with nmie-wrapper. If not, see <http://www.gnu.org/licenses/>.
  19. ///
  20. /// nmie-wrapper uses nmie.c from scattnlay by Ovidio Pena
  21. /// <ovidio@bytesfall.com> as a linked library. He has an additional condition to
  22. /// his library:
  23. // The only additional condition is that we expect that all publications //
  24. // describing work using this software , or all commercial products //
  25. // using it, cite the following reference: //
  26. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  27. // a multilayered sphere," Computer Physics Communications, //
  28. // vol. 180, Nov. 2009, pp. 2348-2354. //
  29. ///
  30. /// @brief Wrapper class around nMie function for ease of use
  31. ///
  32. #include "nmie-wrapper.h"
  33. //#include "nmie.h"
  34. #include <array>
  35. #include <algorithm>
  36. #include <cstdio>
  37. #include <cstdlib>
  38. #include <stdexcept>
  39. #include <vector>
  40. namespace nmie {
  41. //helpers
  42. template<class T> inline T pow2(const T value) {return value*value;}
  43. #define round(x) ((x) >= 0 ? (int)((x) + 0.5):(int)((x) - 0.5))
  44. // ********************************************************************** //
  45. // ********************************************************************** //
  46. // ********************************************************************** //
  47. //emulate C call.
  48. int nMie_wrapper(int L, std::vector<double> x, std::vector<std::complex<double> > m,
  49. int nTheta, std::vector<double> Theta,
  50. double *Qext, double *Qsca, double *Qabs, double *Qbk, double *Qpr, double *g, double *Albedo,
  51. std::vector<std::complex<double> >& S1, std::vector<std::complex<double> >& S2) {
  52. if (x.size() != L || m.size() != L)
  53. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  54. if (Theta.size() != nTheta)
  55. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  56. try {
  57. MultiLayerMie multi_layer_mie;
  58. multi_layer_mie.SetWidthSP(x);
  59. multi_layer_mie.SetIndexSP(m);
  60. multi_layer_mie.SetAngles(Theta);
  61. multi_layer_mie.RunMieCalculations();
  62. *Qext = multi_layer_mie.GetQext();
  63. *Qsca = multi_layer_mie.GetQsca();
  64. *Qabs = multi_layer_mie.GetQabs();
  65. *Qbk = multi_layer_mie.GetQbk();
  66. *Qpr = multi_layer_mie.GetQpr();
  67. *g = multi_layer_mie.GetAsymmetryFactor();
  68. *Albedo = multi_layer_mie.GetAlbedo();
  69. S1 = multi_layer_mie.GetS1();
  70. S2 = multi_layer_mie.GetS2();
  71. } catch( const std::invalid_argument& ia ) {
  72. // Will catch if multi_layer_mie fails or other errors.
  73. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  74. return -1;
  75. }
  76. return 0;
  77. }
  78. // ********************************************************************** //
  79. // ********************************************************************** //
  80. // ********************************************************************** //
  81. double MultiLayerMie::GetQext() {
  82. if (!isMieCalculated_)
  83. throw std::invalid_argument("You should run calculations before result reques!");
  84. return Qext_;
  85. }
  86. // ********************************************************************** //
  87. // ********************************************************************** //
  88. // ********************************************************************** //
  89. double MultiLayerMie::GetQabs() {
  90. if (!isMieCalculated_)
  91. throw std::invalid_argument("You should run calculations before result reques!");
  92. return Qabs_;
  93. }
  94. // ********************************************************************** //
  95. // ********************************************************************** //
  96. // ********************************************************************** //
  97. double MultiLayerMie::GetQsca() {
  98. if (!isMieCalculated_)
  99. throw std::invalid_argument("You should run calculations before result reques!");
  100. return Qsca_;
  101. }
  102. // ********************************************************************** //
  103. // ********************************************************************** //
  104. // ********************************************************************** //
  105. double MultiLayerMie::GetQbk() {
  106. if (!isMieCalculated_)
  107. throw std::invalid_argument("You should run calculations before result reques!");
  108. return Qbk_;
  109. }
  110. // ********************************************************************** //
  111. // ********************************************************************** //
  112. // ********************************************************************** //
  113. double MultiLayerMie::GetQpr() {
  114. if (!isMieCalculated_)
  115. throw std::invalid_argument("You should run calculations before result reques!");
  116. return Qpr_;
  117. }
  118. // ********************************************************************** //
  119. // ********************************************************************** //
  120. // ********************************************************************** //
  121. double MultiLayerMie::GetAsymmetryFactor() {
  122. if (!isMieCalculated_)
  123. throw std::invalid_argument("You should run calculations before result reques!");
  124. return asymmetry_factor_;
  125. }
  126. // ********************************************************************** //
  127. // ********************************************************************** //
  128. // ********************************************************************** //
  129. double MultiLayerMie::GetAlbedo() {
  130. if (!isMieCalculated_)
  131. throw std::invalid_argument("You should run calculations before result reques!");
  132. return albedo_;
  133. }
  134. // ********************************************************************** //
  135. // ********************************************************************** //
  136. // ********************************************************************** //
  137. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  138. return S1_;
  139. }
  140. // ********************************************************************** //
  141. // ********************************************************************** //
  142. // ********************************************************************** //
  143. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  144. return S2_;
  145. }
  146. // ********************************************************************** //
  147. // ********************************************************************** //
  148. // ********************************************************************** //
  149. void MultiLayerMie::SetAngles(std::vector<double> angles) {
  150. isMieCalculated_ = false;
  151. theta_.clear();
  152. for (auto value : angles) theta_.push_back(value);
  153. } // end of SetAngles()
  154. // ********************************************************************** //
  155. // ********************************************************************** //
  156. // ********************************************************************** //
  157. void MultiLayerMie::SetWidthSP(std::vector<double> size_parameter) {
  158. isMieCalculated_ = false;
  159. size_parameter_.clear();
  160. double prev_size_parameter = 0.0;
  161. for (auto layer_size_parameter : size_parameter) {
  162. if (layer_size_parameter <= 0.0)
  163. throw std::invalid_argument("Size parameter should be positive!");
  164. if (prev_size_parameter > layer_size_parameter)
  165. throw std::invalid_argument
  166. ("Size parameter for next layer should be larger than the previous one!");
  167. prev_size_parameter = layer_size_parameter;
  168. size_parameter_.push_back(layer_size_parameter);
  169. }
  170. }
  171. // end of void MultiLayerMie::SetWidthSP(...);
  172. // ********************************************************************** //
  173. // ********************************************************************** //
  174. // ********************************************************************** //
  175. void MultiLayerMie::SetIndexSP(std::vector< std::complex<double> > index) {
  176. isMieCalculated_ = false;
  177. index_.clear();
  178. for (auto value : index) index_.push_back(value);
  179. } // end of void MultiLayerMie::SetIndexSP(...);
  180. // ********************************************************************** //
  181. // ********************************************************************** //
  182. // ********************************************************************** //
  183. void MultiLayerMie::SetPEC(int layer_position) {
  184. if (layer_position < 0)
  185. throw std::invalid_argument("Error! Layers are numbered from 0!");
  186. PEC_layer_position_ = layer_position;
  187. }
  188. // ********************************************************************** //
  189. // ********************************************************************** //
  190. // ********************************************************************** //
  191. void MultiLayerMie::SetMaxTermsNumber(int nmax) {
  192. nmax_ = nmax;
  193. }
  194. // ********************************************************************** //
  195. // ********************************************************************** //
  196. // ********************************************************************** //
  197. void MultiLayerMie::GenerateSizeParameter() {
  198. size_parameter_.clear();
  199. double radius = 0.0;
  200. for (auto width : target_width_) {
  201. radius += width;
  202. size_parameter_.push_back(2*PI*radius / wavelength_);
  203. }
  204. for (auto width : coating_width_) {
  205. radius += width;
  206. size_parameter_.push_back(2*PI*radius / wavelength_);
  207. }
  208. total_radius_ = radius;
  209. } // end of void MultiLayerMie::GenerateSizeParameter();
  210. // ********************************************************************** //
  211. // ********************************************************************** //
  212. // ********************************************************************** //
  213. double MultiLayerMie::GetTotalRadius() {
  214. if (total_radius_ == 0) GenerateSizeParameter();
  215. return total_radius_;
  216. } // end of double MultiLayerMie::GetTotalRadius();
  217. // ********************************************************************** //
  218. // ********************************************************************** //
  219. // ********************************************************************** //
  220. std::vector< std::array<double,5> >
  221. MultiLayerMie::GetSpectra(double from_WL, double to_WL, int samples) {
  222. std::vector< std::array<double,5> > spectra;
  223. double step_WL = (to_WL - from_WL)/ static_cast<double>(samples);
  224. double wavelength_backup = wavelength_;
  225. long fails = 0;
  226. for (double WL = from_WL; WL < to_WL; WL += step_WL) {
  227. double Qext, Qsca, Qabs, Qbk;
  228. wavelength_ = WL;
  229. try {
  230. RunMieCalculations();
  231. } catch( const std::invalid_argument& ia ) {
  232. fails++;
  233. continue;
  234. }
  235. //printf("%3.1f ",WL);
  236. spectra.push_back({wavelength_, Qext, Qsca, Qabs, Qbk});
  237. } // end of for each WL in spectra
  238. printf("fails %li\n",fails);
  239. wavelength_ = wavelength_backup;
  240. return spectra;
  241. }
  242. // ********************************************************************** //
  243. // ********************************************************************** //
  244. // ********************************************************************** //
  245. // Calculate Nstop - equation (17)
  246. //
  247. void MultiLayerMie::Nstop() {
  248. const double& xL = size_parameter_.back();
  249. if (xL <= 8) {
  250. nmax_ = round(xL + 4*pow(xL, 1/3) + 1);
  251. } else if (xL <= 4200) {
  252. nmax_ = round(xL + 4.05*pow(xL, 1/3) + 2);
  253. } else {
  254. nmax_ = round(xL + 4*pow(xL, 1/3) + 2);
  255. }
  256. }
  257. // ********************************************************************** //
  258. // ********************************************************************** //
  259. // ********************************************************************** //
  260. void MultiLayerMie::Nmax(int first_layer) {
  261. int ri, riM1;
  262. const std::vector<double>& x = size_parameter_;
  263. const std::vector<std::complex<double> >& m = index_;
  264. const int& pl = PEC_layer_position_;
  265. Nstop(); // Set initial nmax_ value
  266. for (int i = first_layer; i < x.size(); i++) {
  267. if (i > PEC_layer_position_)
  268. ri = round(std::abs(x[i]*m[i]));
  269. else
  270. ri = 0;
  271. nmax_ = std::max(nmax_, ri);
  272. // first layer is pec, if pec is present
  273. if ((i > first_layer) && ((i - 1) > PEC_layer_position_))
  274. riM1 = round(std::abs(x[i - 1]* m[i]));
  275. // TODO Ovidio, should we check?
  276. // riM2 = round(std::abs(x[i]* m[i-1]))
  277. else
  278. riM1 = 0;
  279. nmax_ = std::max(nmax_, riM1);
  280. }
  281. nmax_ += 15; // Final nmax_ value
  282. }
  283. //**********************************************************************************//
  284. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  285. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  286. // //
  287. // Input parameters: //
  288. // z: Real argument to evaluate jn and h1n //
  289. // nmax_: Maximum number of terms to calculate jn and h1n //
  290. // //
  291. // Output parameters: //
  292. // jn, h1n: Spherical Bessel and Hankel functions //
  293. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  294. // //
  295. // The implementation follows the algorithm by I.J. Thompson and A.R. Barnett, //
  296. // Comp. Phys. Comm. 47 (1987) 245-257. //
  297. // //
  298. // Complex spherical Bessel functions from n=0..nmax_-1 for z in the upper half //
  299. // plane (Im(z) > -3). //
  300. // //
  301. // j[n] = j/n(z) Regular solution: j[0]=sin(z)/z //
  302. // j'[n] = d[j/n(z)]/dz //
  303. // h1[n] = h[0]/n(z) Irregular Hankel function: //
  304. // h1'[n] = d[h[0]/n(z)]/dz h1[0] = j0(z) + i*y0(z) //
  305. // = (sin(z)-i*cos(z))/z //
  306. // = -i*exp(i*z)/z //
  307. // Using complex CF1, and trigonometric forms for n=0 solutions. //
  308. //**********************************************************************************//
  309. void MultiLayerMie::sbesjh(std::complex<double> z, std::vector<std::complex<double> >& jn, std::vector<std::complex<double> >& jnp, std::vector<std::complex<double> >& h1n, std::vector<std::complex<double> >& h1np) {
  310. const int limit = 20000;
  311. const double accur = 1.0e-12;
  312. const double tm30 = 1e-30;
  313. int n;
  314. double absc;
  315. std::complex<double> zi, w;
  316. std::complex<double> pl, f, b, d, c, del, jn0, jndb, h1nldb, h1nbdb;
  317. absc = std::abs(std::real(z)) + std::abs(std::imag(z));
  318. if ((absc < accur) || (std::imag(z) < -3.0)) {
  319. throw std::invalid_argument("TODO add error description for condition if ((absc < accur) || (std::imag(z) < -3.0))");
  320. }
  321. zi = 1.0/z;
  322. w = zi + zi;
  323. pl = double(nmax_)*zi;
  324. f = pl + zi;
  325. b = f + f + zi;
  326. d = 0.0;
  327. c = f;
  328. for (n = 0; n < limit; n++) {
  329. d = b - d;
  330. c = b - 1.0/c;
  331. absc = std::abs(std::real(d)) + std::abs(std::imag(d));
  332. if (absc < tm30) {
  333. d = tm30;
  334. }
  335. absc = std::abs(std::real(c)) + std::abs(std::imag(c));
  336. if (absc < tm30) {
  337. c = tm30;
  338. }
  339. d = 1.0/d;
  340. del = d*c;
  341. f = f*del;
  342. b += w;
  343. absc = std::abs(std::real(del - 1.0)) + std::abs(std::imag(del - 1.0));
  344. if (absc < accur) {
  345. // We have obtained the desired accuracy
  346. break;
  347. }
  348. }
  349. if (absc > accur) {
  350. throw std::invalid_argument("We were not able to obtain the desired accuracy");
  351. }
  352. jn[nmax_ - 1] = tm30;
  353. jnp[nmax_ - 1] = f*jn[nmax_ - 1];
  354. // Downward recursion to n=0 (N.B. Coulomb Functions)
  355. for (n = nmax_ - 2; n >= 0; n--) {
  356. jn[n] = pl*jn[n + 1] + jnp[n + 1];
  357. jnp[n] = pl*jn[n] - jn[n + 1];
  358. pl = pl - zi;
  359. }
  360. // Calculate the n=0 Bessel Functions
  361. jn0 = zi*std::sin(z);
  362. h1n[0] = std::exp(std::complex<double>(0.0, 1.0)*z)*zi*(-std::complex<double>(0.0, 1.0));
  363. h1np[0] = h1n[0]*(std::complex<double>(0.0, 1.0) - zi);
  364. // Rescale j[n], j'[n], converting to spherical Bessel functions.
  365. // Recur h1[n], h1'[n] as spherical Bessel functions.
  366. w = 1.0/jn[0];
  367. pl = zi;
  368. for (n = 0; n < nmax_; n++) {
  369. jn[n] = jn0*(w*jn[n]);
  370. jnp[n] = jn0*(w*jnp[n]) - zi*jn[n];
  371. if (n != 0) {
  372. h1n[n] = (pl - zi)*h1n[n - 1] - h1np[n - 1];
  373. // check if hankel is increasing (upward stable)
  374. if (std::abs(h1n[n]) < std::abs(h1n[n - 1])) {
  375. jndb = z;
  376. h1nldb = h1n[n];
  377. h1nbdb = h1n[n - 1];
  378. }
  379. pl += zi;
  380. h1np[n] = -(pl*h1n[n]) + h1n[n - 1];
  381. }
  382. }
  383. }
  384. //**********************************************************************************//
  385. // This function calculates the spherical Bessel functions (bj and by) and the //
  386. // logarithmic derivative (bd) for a given complex value z. See pag. 87 B&H. //
  387. // //
  388. // Input parameters: //
  389. // z: Complex argument to evaluate bj, by and bd //
  390. // nmax_: Maximum number of terms to calculate bj, by and bd //
  391. // //
  392. // Output parameters: //
  393. // bj, by: Spherical Bessel functions //
  394. // bd: Logarithmic derivative //
  395. //**********************************************************************************//
  396. void MultiLayerMie::sphericalBessel(std::complex<double> z, std::vector<std::complex<double> >& bj, std::vector<std::complex<double> >& by, std::vector<std::complex<double> >& bd) {
  397. std::vector<std::complex<double> > jn, jnp, h1n, h1np;
  398. jn.resize(nmax_);
  399. jnp.resize(nmax_);
  400. h1n.resize(nmax_);
  401. h1np.resize(nmax_);
  402. sbesjh(z, jn, jnp, h1n, h1np);
  403. for (int n = 0; n < nmax_; n++) {
  404. bj[n] = jn[n];
  405. by[n] = (h1n[n] - jn[n])/std::complex<double>(0.0, 1.0);
  406. bd[n] = jnp[n]/jn[n] + 1.0/z;
  407. }
  408. }
  409. // ********************************************************************** //
  410. // ********************************************************************** //
  411. // ********************************************************************** //
  412. // Calculate an - equation (5)
  413. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  414. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  415. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  416. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  417. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  418. return Num/Denom;
  419. }
  420. // ********************************************************************** //
  421. // ********************************************************************** //
  422. // ********************************************************************** //
  423. // Calculate bn - equation (6)
  424. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  425. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  426. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  427. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  428. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  429. return Num/Denom;
  430. }
  431. // ********************************************************************** //
  432. // ********************************************************************** //
  433. // ********************************************************************** //
  434. // Calculates S1 - equation (25a)
  435. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  436. double Pi, double Tau) {
  437. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  438. }
  439. // ********************************************************************** //
  440. // ********************************************************************** //
  441. // ********************************************************************** //
  442. // Calculates S2 - equation (25b) (it's the same as (25a), just switches Pi and Tau)
  443. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  444. double Pi, double Tau) {
  445. return calc_S1(n, an, bn, Tau, Pi);
  446. }
  447. //**********************************************************************************//
  448. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  449. // real argument (x). //
  450. // Equations (20a) - (21b) //
  451. // //
  452. // Input parameters: //
  453. // x: Real argument to evaluate Psi and Zeta //
  454. // nmax: Maximum number of terms to calculate Psi and Zeta //
  455. // //
  456. // Output parameters: //
  457. // Psi, Zeta: Riccati-Bessel functions //
  458. //**********************************************************************************//
  459. void MultiLayerMie::calcPsiZeta(double x,
  460. std::vector<std::complex<double> > D1,
  461. std::vector<std::complex<double> > D3,
  462. std::vector<std::complex<double> >& Psi,
  463. std::vector<std::complex<double> >& Zeta) {
  464. int n;
  465. //Upward recurrence for Psi and Zeta - equations (20a) - (21b)
  466. Psi[0] = std::complex<double>(sin(x), 0);
  467. Zeta[0] = std::complex<double>(sin(x), -cos(x));
  468. for (n = 1; n <= nmax_; n++) {
  469. Psi[n] = Psi[n - 1]*(n/x - D1[n - 1]);
  470. Zeta[n] = Zeta[n - 1]*(n/x - D3[n - 1]);
  471. }
  472. }
  473. //**********************************************************************************//
  474. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  475. // functions (D1 and D3) for a complex argument (z). //
  476. // Equations (16a), (16b) and (18a) - (18d) //
  477. // //
  478. // Input parameters: //
  479. // z: Complex argument to evaluate D1 and D3 //
  480. // nmax_: Maximum number of terms to calculate D1 and D3 //
  481. // //
  482. // Output parameters: //
  483. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  484. //**********************************************************************************//
  485. void MultiLayerMie::calcD1D3(std::complex<double> z,
  486. std::vector<std::complex<double> >& D1,
  487. std::vector<std::complex<double> >& D3) {
  488. int n;
  489. std::vector<std::complex<double> > PsiZeta;
  490. PsiZeta.resize(nmax_ + 1);
  491. // Downward recurrence for D1 - equations (16a) and (16b)
  492. D1[nmax_] = std::complex<double>(0.0, 0.0);
  493. for (n = nmax_; n > 0; n--) {
  494. D1[n - 1] = double(n)/z - 1.0/(D1[n] + double(n)/z);
  495. }
  496. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  497. PsiZeta[0] = 0.5*(1.0 - std::complex<double>(cos(2.0*z.real()), sin(2.0*z.real()))*exp(-2.0*z.imag()));
  498. D3[0] = std::complex<double>(0.0, 1.0);
  499. for (n = 1; n <= nmax_; n++) {
  500. PsiZeta[n] = PsiZeta[n - 1]*(double(n)/z - D1[n - 1])*(double(n)/z- D3[n - 1]);
  501. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta[n];
  502. }
  503. }
  504. //**********************************************************************************//
  505. // This function calculates Pi and Tau for all values of Theta. //
  506. // Equations (26a) - (26c) //
  507. // //
  508. // Input parameters: //
  509. // nmax_: Maximum number of terms to calculate Pi and Tau //
  510. // nTheta: Number of scattering angles //
  511. // Theta: Array containing all the scattering angles where the scattering //
  512. // amplitudes will be calculated //
  513. // //
  514. // Output parameters: //
  515. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  516. //**********************************************************************************//
  517. void MultiLayerMie::calcPiTau(double Theta, std::vector<double>& Pi, std::vector<double>& Tau) {
  518. int n;
  519. //****************************************************//
  520. // Equations (26a) - (26c) //
  521. //****************************************************//
  522. for (n = 0; n < nmax_; n++) {
  523. if (n == 0) {
  524. // Initialize Pi and Tau
  525. Pi[n] = 1.0;
  526. Tau[n] = (n + 1)*cos(Theta);
  527. } else {
  528. // Calculate the actual values
  529. Pi[n] = ((n == 1) ? ((n + n + 1)*cos(Theta)*Pi[n - 1]/n)
  530. : (((n + n + 1)*cos(Theta)*Pi[n - 1] - (n + 1)*Pi[n - 2])/n));
  531. Tau[n] = (n + 1)*cos(Theta)*Pi[n] - (n + 2)*Pi[n - 1];
  532. }
  533. }
  534. }
  535. //**********************************************************************************//
  536. // This function calculates the scattering coefficients required to calculate //
  537. // both the near- and far-field parameters. //
  538. // //
  539. // Input parameters: //
  540. // L: Number of layers //
  541. // pl: Index of PEC layer. If there is none just send -1 //
  542. // x: Array containing the size parameters of the layers [0..L-1] //
  543. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  544. // nmax: Maximum number of multipolar expansion terms to be used for the //
  545. // calculations. Only use it if you know what you are doing, otherwise //
  546. // set this parameter to -1 and the function will calculate it. //
  547. // //
  548. // Output parameters: //
  549. // an, bn: Complex scattering amplitudes //
  550. // //
  551. // Return value: //
  552. // Number of multipolar expansion terms used for the calculations //
  553. //**********************************************************************************//
  554. void MultiLayerMie::ScattCoeffs(int L,
  555. std::vector<std::complex<double> >& an, std::vector<std::complex<double> >& bn) {
  556. //************************************************************************//
  557. // Calculate the index of the first layer. It can be either 0 (default) //
  558. // or the index of the outermost PEC layer. In the latter case all layers //
  559. // below the PEC are discarded. //
  560. //************************************************************************//
  561. const std::vector<double>& x = size_parameter_;
  562. const std::vector<std::complex<double> >& m = index_;
  563. const int& pl = PEC_layer_position_;
  564. // TODO, is it possible for PEC to have a zero index? If yes than is should be:
  565. // int fl = (pl > -1) ? pl : 0;
  566. int fl = (pl > 0) ? pl : 0;
  567. if (nmax_ <= 0) Nmax(fl);
  568. std::complex<double> z1, z2;
  569. std::complex<double> Num, Denom;
  570. std::complex<double> G1, G2;
  571. std::complex<double> Temp;
  572. int n, l;
  573. //**************************************************************************//
  574. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  575. // means that index = layer number - 1 or index = n - 1. The only exception //
  576. // are the arrays for representing D1, D3 and Q because they need a value //
  577. // for the index 0 (zero), hence it is important to consider this shift //
  578. // between different arrays. The change was done to optimize memory usage. //
  579. //**************************************************************************//
  580. // Allocate memory to the arrays
  581. std::vector<std::vector<std::complex<double> > > D1_mlxl, D1_mlxlM1;
  582. D1_mlxl.resize(L);
  583. D1_mlxlM1.resize(L);
  584. std::vector<std::vector<std::complex<double> > > D3_mlxl, D3_mlxlM1;
  585. D3_mlxl.resize(L);
  586. D3_mlxlM1.resize(L);
  587. std::vector<std::vector<std::complex<double> > > Q;
  588. Q.resize(L);
  589. std::vector<std::vector<std::complex<double> > > Ha, Hb;
  590. Ha.resize(L);
  591. Hb.resize(L);
  592. for (l = 0; l < L; l++) {
  593. D1_mlxl[l].resize(nmax_ + 1);
  594. D1_mlxlM1[l].resize(nmax_ + 1);
  595. D3_mlxl[l].resize(nmax_ + 1);
  596. D3_mlxlM1[l].resize(nmax_ + 1);
  597. Q[l].resize(nmax_ + 1);
  598. Ha[l].resize(nmax_);
  599. Hb[l].resize(nmax_);
  600. }
  601. an.resize(nmax_);
  602. bn.resize(nmax_);
  603. std::vector<std::complex<double> > D1XL, D3XL;
  604. D1XL.resize(nmax_ + 1);
  605. D3XL.resize(nmax_ + 1);
  606. std::vector<std::complex<double> > PsiXL, ZetaXL;
  607. PsiXL.resize(nmax_ + 1);
  608. ZetaXL.resize(nmax_ + 1);
  609. //*************************************************//
  610. // Calculate D1 and D3 for z1 in the first layer //
  611. //*************************************************//
  612. if (fl == pl) { // PEC layer
  613. for (n = 0; n <= nmax_; n++) {
  614. D1_mlxl[fl][n] = std::complex<double>(0.0, -1.0);
  615. D3_mlxl[fl][n] = std::complex<double>(0.0, 1.0);
  616. }
  617. } else { // Regular layer
  618. z1 = x[fl]* m[fl];
  619. // Calculate D1 and D3
  620. calcD1D3(z1, D1_mlxl[fl], D3_mlxl[fl]);
  621. }
  622. //******************************************************************//
  623. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  624. //******************************************************************//
  625. for (n = 0; n < nmax_; n++) {
  626. Ha[fl][n] = D1_mlxl[fl][n + 1];
  627. Hb[fl][n] = D1_mlxl[fl][n + 1];
  628. }
  629. //*****************************************************//
  630. // Iteration from the second layer to the last one (L) //
  631. //*****************************************************//
  632. for (l = fl + 1; l < L; l++) {
  633. //************************************************************//
  634. //Calculate D1 and D3 for z1 and z2 in the layers fl+1..L //
  635. //************************************************************//
  636. z1 = x[l]*m[l];
  637. z2 = x[l - 1]*m[l];
  638. //Calculate D1 and D3 for z1
  639. calcD1D3(z1, D1_mlxl[l], D3_mlxl[l]);
  640. //Calculate D1 and D3 for z2
  641. calcD1D3(z2, D1_mlxlM1[l], D3_mlxlM1[l]);
  642. //*********************************************//
  643. //Calculate Q, Ha and Hb in the layers fl+1..L //
  644. //*********************************************//
  645. // Upward recurrence for Q - equations (19a) and (19b)
  646. Num = exp(-2.0*(z1.imag() - z2.imag()))*std::complex<double>(cos(-2.0*z2.real()) - exp(-2.0*z2.imag()), sin(-2.0*z2.real()));
  647. Denom = std::complex<double>(cos(-2.0*z1.real()) - exp(-2.0*z1.imag()), sin(-2.0*z1.real()));
  648. Q[l][0] = Num/Denom;
  649. for (n = 1; n <= nmax_; n++) {
  650. Num = (z1*D1_mlxl[l][n] + double(n))*(double(n) - z1*D3_mlxl[l][n - 1]);
  651. Denom = (z2*D1_mlxlM1[l][n] + double(n))*(double(n) - z2*D3_mlxlM1[l][n - 1]);
  652. Q[l][n] = (((x[l - 1]*x[l - 1])/(x[l]*x[l])* Q[l][n - 1])*Num)/Denom;
  653. }
  654. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  655. for (n = 1; n <= nmax_; n++) {
  656. //Ha
  657. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  658. G1 = -D1_mlxlM1[l][n];
  659. G2 = -D3_mlxlM1[l][n];
  660. } else {
  661. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[l][n]);
  662. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[l][n]);
  663. }
  664. Temp = Q[l][n]*G1;
  665. Num = (G2*D1_mlxl[l][n]) - (Temp*D3_mlxl[l][n]);
  666. Denom = G2 - Temp;
  667. Ha[l][n - 1] = Num/Denom;
  668. //Hb
  669. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  670. G1 = Hb[l - 1][n - 1];
  671. G2 = Hb[l - 1][n - 1];
  672. } else {
  673. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[l][n]);
  674. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[l][n]);
  675. }
  676. Temp = Q[l][n]*G1;
  677. Num = (G2*D1_mlxl[l][n]) - (Temp* D3_mlxl[l][n]);
  678. Denom = (G2- Temp);
  679. Hb[l][n - 1] = (Num/ Denom);
  680. }
  681. }
  682. //**************************************//
  683. //Calculate D1, D3, Psi and Zeta for XL //
  684. //**************************************//
  685. // Calculate D1XL and D3XL
  686. calcD1D3(x[L - 1], D1XL, D3XL);
  687. // Calculate PsiXL and ZetaXL
  688. calcPsiZeta(x[L - 1], D1XL, D3XL, PsiXL, ZetaXL);
  689. //*********************************************************************//
  690. // Finally, we calculate the scattering coefficients (an and bn) and //
  691. // the angular functions (Pi and Tau). Note that for these arrays the //
  692. // first layer is 0 (zero), in future versions all arrays will follow //
  693. // this convention to save memory. (13 Nov, 2014) //
  694. //*********************************************************************//
  695. for (n = 0; n < nmax_; n++) {
  696. //********************************************************************//
  697. //Expressions for calculating an and bn coefficients are not valid if //
  698. //there is only one PEC layer (ie, for a simple PEC sphere). //
  699. //********************************************************************//
  700. if (pl < (L - 1)) {
  701. an[n] = calc_an(n + 1, x[L - 1], Ha[L - 1][n], m[L - 1], PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  702. bn[n] = calc_bn(n + 1, x[L - 1], Hb[L - 1][n], m[L - 1], PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  703. } else {
  704. an[n] = calc_an(n + 1, x[L - 1], std::complex<double>(0.0, 0.0), std::complex<double>(1.0, 0.0), PsiXL[n + 1], ZetaXL[n + 1], PsiXL[n], ZetaXL[n]);
  705. bn[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  706. }
  707. }
  708. }
  709. // ********************************************************************** //
  710. // ********************************************************************** //
  711. // ********************************************************************** //
  712. //**********************************************************************************//
  713. // This function calculates the actual scattering parameters and amplitudes //
  714. // //
  715. // Input parameters: //
  716. // L: Number of layers //
  717. // pl: Index of PEC layer. If there is none just send -1 //
  718. // x: Array containing the size parameters of the layers [0..L-1] //
  719. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  720. // nTheta: Number of scattering angles //
  721. // Theta: Array containing all the scattering angles where the scattering //
  722. // amplitudes will be calculated //
  723. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  724. // calculations. Only use it if you know what you are doing, otherwise //
  725. // set this parameter to -1 and the function will calculate it //
  726. // //
  727. // Output parameters: //
  728. // Qext: Efficiency factor for extinction //
  729. // Qsca: Efficiency factor for scattering //
  730. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  731. // Qbk: Efficiency factor for backscattering //
  732. // Qpr: Efficiency factor for the radiation pressure //
  733. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  734. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  735. // S1, S2: Complex scattering amplitudes //
  736. // //
  737. // Return value: //
  738. // Number of multipolar expansion terms used for the calculations //
  739. //**********************************************************************************//
  740. void MultiLayerMie::RunMieCalculations() {
  741. if (size_parameter_.size() != index_.size())
  742. throw std::invalid_argument("Each size parameter should have only one index!");
  743. int i, n, t;
  744. std::vector<std::complex<double> > an, bn;
  745. std::complex<double> Qbktmp;
  746. const std::vector<double>& x = size_parameter_;
  747. const std::vector<std::complex<double> >& m = index_;
  748. int L = index_.size();
  749. // Calculate scattering coefficients
  750. ScattCoeffs(L, an, bn);
  751. std::vector<double> Pi, Tau;
  752. Pi.resize(nmax_);
  753. Tau.resize(nmax_);
  754. double x2 = x[L - 1]*x[L - 1];
  755. // Initialize the scattering parameters
  756. Qext_ = 0;
  757. Qsca_ = 0;
  758. Qabs_ = 0;
  759. Qbk_ = 0;
  760. Qbktmp = std::complex<double>(0.0, 0.0);
  761. Qpr_ = 0;
  762. asymmetry_factor_ = 0;
  763. albedo_ = 0;
  764. // Initialize the scattering amplitudes
  765. int nTheta = theta_.size();
  766. S1_.resize(nTheta);
  767. S2_.resize(nTheta);
  768. for (t = 0; t < nTheta; t++) {
  769. S1_[t] = std::complex<double>(0.0, 0.0);
  770. S2_[t] = std::complex<double>(0.0, 0.0);
  771. }
  772. // By using downward recurrence we avoid loss of precision due to float rounding errors
  773. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  774. // http://en.wikipedia.org/wiki/Loss_of_significance
  775. for (i = nmax_ - 2; i >= 0; i--) {
  776. n = i + 1;
  777. // Equation (27)
  778. Qext_ += (n + n + 1)*(an[i].real() + bn[i].real());
  779. // Equation (28)
  780. Qsca_ += (n + n + 1)*(an[i].real()*an[i].real() + an[i].imag()*an[i].imag() + bn[i].real()*bn[i].real() + bn[i].imag()*bn[i].imag());
  781. // Equation (29) TODO We must check carefully this equation. If we
  782. // remove the typecast to double then the result changes. Which is
  783. // the correct one??? Ovidio (2014/12/10) With cast ratio will
  784. // give double, without cast (n + n + 1)/(n*(n + 1)) will be
  785. // rounded to integer. Tig (2015/02/24)
  786. Qpr_ += ((n*(n + 2)/(n + 1))*((an[i]*std::conj(an[n]) + bn[i]*std::conj(bn[n])).real()) + ((double)(n + n + 1)/(n*(n + 1)))*(an[i]*std::conj(bn[i])).real());
  787. // Equation (33)
  788. Qbktmp = Qbktmp + (double)(n + n + 1)*(1 - 2*(n % 2))*(an[i]- bn[i]);
  789. //****************************************************//
  790. // Calculate the scattering amplitudes (S1 and S2) //
  791. // Equations (25a) - (25b) //
  792. //****************************************************//
  793. for (t = 0; t < nTheta; t++) {
  794. calcPiTau(theta_[t], Pi, Tau);
  795. S1_[t] += calc_S1(n, an[i], bn[i], Pi[i], Tau[i]);
  796. S2_[t] += calc_S2(n, an[i], bn[i], Pi[i], Tau[i]);
  797. }
  798. }
  799. Qext_ = 2*(Qext_)/x2; // Equation (27)
  800. Qsca_ = 2*(Qsca_)/x2; // Equation (28)
  801. Qpr_ = Qext_ - 4*(Qpr_)/x2; // Equation (29)
  802. Qabs_ = Qext_ - Qsca_; // Equation (30)
  803. albedo_ = Qsca_ / Qext_; // Equation (31)
  804. asymmetry_factor_ = (Qext_ - Qpr_) / Qsca_; // Equation (32)
  805. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  806. isMieCalculated_ = true;
  807. //return nmax;
  808. }
  809. // ********************************************************************** //
  810. // ********************************************************************** //
  811. // ********************************************************************** //
  812. // external scattering field = incident + scattered
  813. // BH p.92 (4.37), 94 (4.45), 95 (4.50)
  814. // assume: medium is non-absorbing; refim = 0; Uabs = 0
  815. void MultiLayerMie::fieldExt(double Rho, double Phi, double Theta, std::vector<double> Pi, std::vector<double> Tau,
  816. std::vector<std::complex<double> > an, std::vector<std::complex<double> > bn,
  817. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  818. int i, n;
  819. double rn = 0.0;
  820. std::complex<double> zn, xxip, encap;
  821. std::vector<std::complex<double> > vm3o1n, vm3e1n, vn3o1n, vn3e1n;
  822. vm3o1n.resize(3);
  823. vm3e1n.resize(3);
  824. vn3o1n.resize(3);
  825. vn3e1n.resize(3);
  826. std::vector<std::complex<double> > Ei, Hi, Es, Hs;
  827. Ei.resize(3);
  828. Hi.resize(3);
  829. Es.resize(3);
  830. Hs.resize(3);
  831. for (i = 0; i < 3; i++) {
  832. Ei[i] = std::complex<double>(0.0, 0.0);
  833. Hi[i] = std::complex<double>(0.0, 0.0);
  834. Es[i] = std::complex<double>(0.0, 0.0);
  835. Hs[i] = std::complex<double>(0.0, 0.0);
  836. }
  837. std::vector<std::complex<double> > bj, by, bd;
  838. bj.resize(nmax_);
  839. by.resize(nmax_);
  840. bd.resize(nmax_);
  841. // Calculate spherical Bessel and Hankel functions
  842. sphericalBessel(Rho, bj, by, bd);
  843. for (n = 0; n < nmax_; n++) {
  844. rn = double(n + 1);
  845. zn = bj[n] + std::complex<double>(0.0, 1.0)*by[n];
  846. xxip = Rho*(bj[n] + std::complex<double>(0.0, 1.0)*by[n]) - rn*zn;
  847. vm3o1n[0] = std::complex<double>(0.0, 0.0);
  848. vm3o1n[1] = std::cos(Phi)*Pi[n]*zn;
  849. vm3o1n[2] = -(std::sin(Phi)*Tau[n]*zn);
  850. vm3e1n[0] = std::complex<double>(0.0, 0.0);
  851. vm3e1n[1] = -(std::sin(Phi)*Pi[n]*zn);
  852. vm3e1n[2] = -(std::cos(Phi)*Tau[n]*zn);
  853. vn3o1n[0] = std::sin(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  854. vn3o1n[1] = std::sin(Phi)*Tau[n]*xxip/Rho;
  855. vn3o1n[2] = std::cos(Phi)*Pi[n]*xxip/Rho;
  856. vn3e1n[0] = std::cos(Phi)*rn*(rn + 1.0)*std::sin(Theta)*Pi[n]*zn/Rho;
  857. vn3e1n[1] = std::cos(Phi)*Tau[n]*xxip/Rho;
  858. vn3e1n[2] = -(std::sin(Phi)*Pi[n]*xxip/Rho);
  859. // scattered field: BH p.94 (4.45)
  860. encap = std::pow(std::complex<double>(0.0, 1.0), rn)*(2.0*rn + 1.0)/(rn*(rn + 1.0));
  861. for (i = 0; i < 3; i++) {
  862. Es[i] = Es[i] + encap*(std::complex<double>(0.0, 1.0)*an[n]*vn3e1n[i] - bn[n]*vm3o1n[i]);
  863. Hs[i] = Hs[i] + encap*(std::complex<double>(0.0, 1.0)*bn[n]*vn3o1n[i] + an[n]*vm3e1n[i]);
  864. }
  865. }
  866. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  867. // basis unit vectors = er, etheta, ephi
  868. std::complex<double> eifac = std::exp(std::complex<double>(0.0, 1.0)*Rho*std::cos(Theta));
  869. Ei[0] = eifac*std::sin(Theta)*std::cos(Phi);
  870. Ei[1] = eifac*std::cos(Theta)*std::cos(Phi);
  871. Ei[2] = -(eifac*std::sin(Phi));
  872. // magnetic field
  873. double hffact = 1.0/(cc*mu);
  874. for (i = 0; i < 3; i++) {
  875. Hs[i] = hffact*Hs[i];
  876. }
  877. // incident H field: BH p.26 (2.43), p.89 (4.21)
  878. std::complex<double> hffacta = hffact;
  879. std::complex<double> hifac = eifac*hffacta;
  880. Hi[0] = hifac*std::sin(Theta)*std::sin(Phi);
  881. Hi[1] = hifac*std::cos(Theta)*std::sin(Phi);
  882. Hi[2] = hifac*std::cos(Phi);
  883. for (i = 0; i < 3; i++) {
  884. // electric field E [V m-1] = EF*E0
  885. E[i] = Ei[i] + Es[i];
  886. H[i] = Hi[i] + Hs[i];
  887. }
  888. }
  889. // ********************************************************************** //
  890. // ********************************************************************** //
  891. // ********************************************************************** //
  892. //**********************************************************************************//
  893. // This function calculates complex electric and magnetic field in the surroundings //
  894. // and inside (TODO) the particle. //
  895. // //
  896. // Input parameters: //
  897. // L: Number of layers //
  898. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  899. // x: Array containing the size parameters of the layers [0..L-1] //
  900. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  901. // nmax: Maximum number of multipolar expansion terms to be used for the //
  902. // calculations. Only use it if you know what you are doing, otherwise //
  903. // set this parameter to 0 (zero) and the function will calculate it. //
  904. // ncoord: Number of coordinate points //
  905. // Coords: Array containing all coordinates where the complex electric and //
  906. // magnetic fields will be calculated //
  907. // //
  908. // Output parameters: //
  909. // E, H: Complex electric and magnetic field at the provided coordinates //
  910. // //
  911. // Return value: //
  912. // Number of multipolar expansion terms used for the calculations //
  913. //**********************************************************************************//
  914. // int MultiLayerMie::nField(int L, int pl, std::vector<double> x, std::vector<std::complex<double> > m, int nmax,
  915. // int ncoord, std::vector<double> Xp, std::vector<double> Yp, std::vector<double> Zp,
  916. // std::vector<std::vector<std::complex<double> > >& E, std::vector<std::vector<std::complex<double> > >& H) {
  917. // int i, c;
  918. // double Rho, Phi, Theta;
  919. // std::vector<std::complex<double> > an, bn;
  920. // // This array contains the fields in spherical coordinates
  921. // std::vector<std::complex<double> > Es, Hs;
  922. // Es.resize(3);
  923. // Hs.resize(3);
  924. // // Calculate scattering coefficients
  925. // ScattCoeffs(L, pl, an, bn);
  926. // std::vector<double> Pi, Tau;
  927. // Pi.resize(nmax_);
  928. // Tau.resize(nmax_);
  929. // for (c = 0; c < ncoord; c++) {
  930. // // Convert to spherical coordinates
  931. // Rho = sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c] + Zp[c]*Zp[c]);
  932. // if (Rho < 1e-3) {
  933. // // Avoid convergence problems
  934. // Rho = 1e-3;
  935. // }
  936. // Phi = acos(Xp[c]/sqrt(Xp[c]*Xp[c] + Yp[c]*Yp[c]));
  937. // Theta = acos(Xp[c]/Rho);
  938. // calcPiTau(Theta, Pi, Tau);
  939. // //*******************************************************//
  940. // // external scattering field = incident + scattered //
  941. // // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  942. // // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  943. // //*******************************************************//
  944. // // Firstly the easiest case: the field outside the particle
  945. // if (Rho >= x[L - 1]) {
  946. // fieldExt(Rho, Phi, Theta, Pi, Tau, an, bn, Es, Hs);
  947. // } else {
  948. // // TODO, for now just set all the fields to zero
  949. // for (i = 0; i < 3; i++) {
  950. // Es[i] = std::complex<double>(0.0, 0.0);
  951. // Hs[i] = std::complex<double>(0.0, 0.0);
  952. // }
  953. // }
  954. // //Now, convert the fields back to cartesian coordinates
  955. // E[c][0] = std::sin(Theta)*std::cos(Phi)*Es[0] + std::cos(Theta)*std::cos(Phi)*Es[1] - std::sin(Phi)*Es[2];
  956. // E[c][1] = std::sin(Theta)*std::sin(Phi)*Es[0] + std::cos(Theta)*std::sin(Phi)*Es[1] + std::cos(Phi)*Es[2];
  957. // E[c][2] = std::cos(Theta)*Es[0] - std::sin(Theta)*Es[1];
  958. // H[c][0] = std::sin(Theta)*std::cos(Phi)*Hs[0] + std::cos(Theta)*std::cos(Phi)*Hs[1] - std::sin(Phi)*Hs[2];
  959. // H[c][1] = std::sin(Theta)*std::sin(Phi)*Hs[0] + std::cos(Theta)*std::sin(Phi)*Hs[1] + std::cos(Phi)*Hs[2];
  960. // H[c][2] = std::cos(Theta)*Hs[0] - std::sin(Theta)*Hs[1];
  961. // }
  962. // return nmax;
  963. // } // end of int nField()
  964. } // end of namespace nmie