nmie.cc 59 KB

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