nmie-wrapper.cc 71 KB

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