nmie-core.cpp 79 KB

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