nmie-wrapper.cc 63 KB

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