nmie-wrapper.cc 54 KB

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