nmie.cc 56 KB

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