nmie-wrapper.cc 56 KB

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