nmie-wrapper.cc 80 KB

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