nmie-wrapper.cc 42 KB

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