nmie-wrapper.cc 68 KB

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