nmie.cc 77 KB

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