nmie-wrapper.cc 47 KB

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