nmie-wrapper.cc 50 KB

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