nmie.cc 84 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625
  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 "bessel.h"
  40. #include "nmie.h"
  41. #include <array>
  42. #include <algorithm>
  43. #include <cstdio>
  44. #include <cstdlib>
  45. #include <stdexcept>
  46. #include <vector>
  47. namespace nmie {
  48. //helpers
  49. template<class T> inline T pow2(const T value) {return value*value;}
  50. int round(double x) {
  51. return x >= 0 ? (int)(x + 0.5):(int)(x - 0.5);
  52. }
  53. //**********************************************************************************//
  54. // This function emulates a C call to calculate the actual scattering parameters //
  55. // and amplitudes. //
  56. // //
  57. // Input parameters: //
  58. // L: Number of layers //
  59. // pl: Index of PEC layer. If there is none just send -1 //
  60. // x: Array containing the size parameters of the layers [0..L-1] //
  61. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  62. // nTheta: Number of scattering angles //
  63. // Theta: Array containing all the scattering angles where the scattering //
  64. // amplitudes will be calculated //
  65. // nmax: Maximum number of multipolar expansion terms to be used for the //
  66. // calculations. Only use it if you know what you are doing, otherwise //
  67. // set this parameter to -1 and the function will calculate it //
  68. // //
  69. // Output parameters: //
  70. // Qext: Efficiency factor for extinction //
  71. // Qsca: Efficiency factor for scattering //
  72. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  73. // Qbk: Efficiency factor for backscattering //
  74. // Qpr: Efficiency factor for the radiation pressure //
  75. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  76. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  77. // S1, S2: Complex scattering amplitudes //
  78. // //
  79. // Return value: //
  80. // Number of multipolar expansion terms used for the calculations //
  81. //**********************************************************************************//
  82. int nMie(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  83. if (x.size() != L || m.size() != L)
  84. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  85. if (Theta.size() != nTheta)
  86. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  87. try {
  88. MultiLayerMie multi_layer_mie;
  89. multi_layer_mie.SetLayersSize(x);
  90. multi_layer_mie.SetLayersIndex(m);
  91. multi_layer_mie.SetAngles(Theta);
  92. multi_layer_mie.SetPECLayer(pl);
  93. multi_layer_mie.SetMaxTerms(nmax);
  94. multi_layer_mie.RunMieCalculation();
  95. *Qext = multi_layer_mie.GetQext();
  96. *Qsca = multi_layer_mie.GetQsca();
  97. *Qabs = multi_layer_mie.GetQabs();
  98. *Qbk = multi_layer_mie.GetQbk();
  99. *Qpr = multi_layer_mie.GetQpr();
  100. *g = multi_layer_mie.GetAsymmetryFactor();
  101. *Albedo = multi_layer_mie.GetAlbedo();
  102. S1 = multi_layer_mie.GetS1();
  103. S2 = multi_layer_mie.GetS2();
  104. } catch(const std::invalid_argument& ia) {
  105. // Will catch if multi_layer_mie fails or other errors.
  106. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  107. throw std::invalid_argument(ia);
  108. return -1;
  109. }
  110. return 0;
  111. }
  112. //**********************************************************************************//
  113. // This function is just a wrapper to call the full 'nMie' function with fewer //
  114. // parameters, it is here mainly for compatibility with older versions of the //
  115. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  116. // any limit for the maximum number of terms. //
  117. // //
  118. // Input parameters: //
  119. // L: Number of layers //
  120. // x: Array containing the size parameters of the layers [0..L-1] //
  121. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  122. // nTheta: Number of scattering angles //
  123. // Theta: Array containing all the scattering angles where the scattering //
  124. // amplitudes will be calculated //
  125. // //
  126. // Output parameters: //
  127. // Qext: Efficiency factor for extinction //
  128. // Qsca: Efficiency factor for scattering //
  129. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  130. // Qbk: Efficiency factor for backscattering //
  131. // Qpr: Efficiency factor for the radiation pressure //
  132. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  133. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  134. // S1, S2: Complex scattering amplitudes //
  135. // //
  136. // Return value: //
  137. // Number of multipolar expansion terms used for the calculations //
  138. //**********************************************************************************//
  139. int nMie(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  140. return nmie::nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  141. }
  142. //**********************************************************************************//
  143. // This function is just a wrapper to call the full 'nMie' function with fewer //
  144. // parameters, it is useful if you want to include a PEC layer but not a limit //
  145. // for the maximum number of terms. //
  146. // //
  147. // Input parameters: //
  148. // L: Number of layers //
  149. // pl: Index of PEC layer. If there is none just send -1 //
  150. // x: Array containing the size parameters of the layers [0..L-1] //
  151. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  152. // nTheta: Number of scattering angles //
  153. // Theta: Array containing all the scattering angles where the scattering //
  154. // amplitudes will be calculated //
  155. // //
  156. // Output parameters: //
  157. // Qext: Efficiency factor for extinction //
  158. // Qsca: Efficiency factor for scattering //
  159. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  160. // Qbk: Efficiency factor for backscattering //
  161. // Qpr: Efficiency factor for the radiation pressure //
  162. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  163. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  164. // S1, S2: Complex scattering amplitudes //
  165. // //
  166. // Return value: //
  167. // Number of multipolar expansion terms used for the calculations //
  168. //**********************************************************************************//
  169. int nMie(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  170. return nmie::nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  171. }
  172. //**********************************************************************************//
  173. // This function is just a wrapper to call the full 'nMie' function with fewer //
  174. // parameters, it is useful if you want to include a limit for the maximum number //
  175. // of terms but not a PEC layer. //
  176. // //
  177. // Input parameters: //
  178. // L: Number of layers //
  179. // x: Array containing the size parameters of the layers [0..L-1] //
  180. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  181. // nTheta: Number of scattering angles //
  182. // Theta: Array containing all the scattering angles where the scattering //
  183. // amplitudes will be calculated //
  184. // nmax: Maximum number of multipolar expansion terms to be used for the //
  185. // calculations. Only use it if you know what you are doing, otherwise //
  186. // set this parameter to -1 and the function will calculate it //
  187. // //
  188. // Output parameters: //
  189. // Qext: Efficiency factor for extinction //
  190. // Qsca: Efficiency factor for scattering //
  191. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  192. // Qbk: Efficiency factor for backscattering //
  193. // Qpr: Efficiency factor for the radiation pressure //
  194. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  195. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  196. // S1, S2: Complex scattering amplitudes //
  197. // //
  198. // Return value: //
  199. // Number of multipolar expansion terms used for the calculations //
  200. //**********************************************************************************//
  201. int nMie(const unsigned int L, std::vector<double>& x, std::vector<std::complex<double> >& m, const unsigned 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) {
  202. return nmie::nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  203. }
  204. //**********************************************************************************//
  205. // This function emulates a C call to calculate complex electric and magnetic field //
  206. // in the surroundings and inside (TODO) the particle. //
  207. // //
  208. // Input parameters: //
  209. // L: Number of layers //
  210. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  211. // x: Array containing the size parameters of the layers [0..L-1] //
  212. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  213. // nmax: Maximum number of multipolar expansion terms to be used for the //
  214. // calculations. Only use it if you know what you are doing, otherwise //
  215. // set this parameter to 0 (zero) and the function will calculate it. //
  216. // ncoord: Number of coordinate points //
  217. // Coords: Array containing all coordinates where the complex electric and //
  218. // magnetic fields will be calculated //
  219. // //
  220. // Output parameters: //
  221. // E, H: Complex electric and magnetic field at the provided coordinates //
  222. // //
  223. // Return value: //
  224. // Number of multipolar expansion terms used for the calculations //
  225. //**********************************************************************************//
  226. int nField(const unsigned int L, const int pl, const std::vector<double>& x, const std::vector<std::complex<double> >& m, const int nmax, const unsigned 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) {
  227. if (x.size() != L || m.size() != L)
  228. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  229. if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
  230. || E.size() != ncoord || H.size() != ncoord)
  231. throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
  232. for (auto f:E)
  233. if (f.size() != 3)
  234. throw std::invalid_argument("Field E is not 3D!");
  235. for (auto f:H)
  236. if (f.size() != 3)
  237. throw std::invalid_argument("Field H is not 3D!");
  238. try {
  239. MultiLayerMie multi_layer_mie;
  240. //multi_layer_mie.SetPECLayer(pl); // TODO add PEC layer to field plotting
  241. multi_layer_mie.SetLayersSize(x);
  242. multi_layer_mie.SetLayersIndex(m);
  243. multi_layer_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
  244. multi_layer_mie.RunFieldCalculation();
  245. E = multi_layer_mie.GetFieldE();
  246. H = multi_layer_mie.GetFieldH();
  247. } catch(const std::invalid_argument& ia) {
  248. // Will catch if multi_layer_mie fails or other errors.
  249. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  250. throw std::invalid_argument(ia);
  251. return - 1;
  252. }
  253. return 0;
  254. }
  255. // ********************************************************************** //
  256. // Returns previously calculated Qext //
  257. // ********************************************************************** //
  258. double MultiLayerMie::GetQext() {
  259. if (!isMieCalculated_)
  260. throw std::invalid_argument("You should run calculations before result request!");
  261. return Qext_;
  262. }
  263. // ********************************************************************** //
  264. // Returns previously calculated Qabs //
  265. // ********************************************************************** //
  266. double MultiLayerMie::GetQabs() {
  267. if (!isMieCalculated_)
  268. throw std::invalid_argument("You should run calculations before result request!");
  269. return Qabs_;
  270. }
  271. // ********************************************************************** //
  272. // Returns previously calculated Qsca //
  273. // ********************************************************************** //
  274. double MultiLayerMie::GetQsca() {
  275. if (!isMieCalculated_)
  276. throw std::invalid_argument("You should run calculations before result request!");
  277. return Qsca_;
  278. }
  279. // ********************************************************************** //
  280. // Returns previously calculated Qbk //
  281. // ********************************************************************** //
  282. double MultiLayerMie::GetQbk() {
  283. if (!isMieCalculated_)
  284. throw std::invalid_argument("You should run calculations before result request!");
  285. return Qbk_;
  286. }
  287. // ********************************************************************** //
  288. // Returns previously calculated Qpr //
  289. // ********************************************************************** //
  290. double MultiLayerMie::GetQpr() {
  291. if (!isMieCalculated_)
  292. throw std::invalid_argument("You should run calculations before result request!");
  293. return Qpr_;
  294. }
  295. // ********************************************************************** //
  296. // Returns previously calculated assymetry factor //
  297. // ********************************************************************** //
  298. double MultiLayerMie::GetAsymmetryFactor() {
  299. if (!isMieCalculated_)
  300. throw std::invalid_argument("You should run calculations before result request!");
  301. return asymmetry_factor_;
  302. }
  303. // ********************************************************************** //
  304. // Returns previously calculated Albedo //
  305. // ********************************************************************** //
  306. double MultiLayerMie::GetAlbedo() {
  307. if (!isMieCalculated_)
  308. throw std::invalid_argument("You should run calculations before result request!");
  309. return albedo_;
  310. }
  311. // ********************************************************************** //
  312. // Returns previously calculated S1 //
  313. // ********************************************************************** //
  314. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  315. if (!isMieCalculated_)
  316. throw std::invalid_argument("You should run calculations before result request!");
  317. return S1_;
  318. }
  319. // ********************************************************************** //
  320. // Returns previously calculated S2 //
  321. // ********************************************************************** //
  322. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  323. if (!isMieCalculated_)
  324. throw std::invalid_argument("You should run calculations before result request!");
  325. return S2_;
  326. }
  327. // ********************************************************************** //
  328. // Modify scattering (theta) angles //
  329. // ********************************************************************** //
  330. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  331. isExpCoeffsCalc_ = false;
  332. isScaCoeffsCalc_ = false;
  333. isMieCalculated_ = false;
  334. theta_ = angles;
  335. }
  336. // ********************************************************************** //
  337. // Modify size of all layers //
  338. // ********************************************************************** //
  339. void MultiLayerMie::SetLayersSize(const std::vector<double>& layer_size) {
  340. isExpCoeffsCalc_ = false;
  341. isScaCoeffsCalc_ = false;
  342. isMieCalculated_ = false;
  343. size_param_.clear();
  344. double prev_layer_size = 0.0;
  345. for (auto curr_layer_size : layer_size) {
  346. if (curr_layer_size <= 0.0)
  347. throw std::invalid_argument("Size parameter should be positive!");
  348. if (prev_layer_size > curr_layer_size)
  349. throw std::invalid_argument
  350. ("Size parameter for next layer should be larger than the previous one!");
  351. prev_layer_size = curr_layer_size;
  352. size_param_.push_back(curr_layer_size);
  353. }
  354. }
  355. // ********************************************************************** //
  356. // Modify refractive index of all layers //
  357. // ********************************************************************** //
  358. void MultiLayerMie::SetLayersIndex(const std::vector< std::complex<double> >& index) {
  359. isExpCoeffsCalc_ = false;
  360. isScaCoeffsCalc_ = false;
  361. isMieCalculated_ = false;
  362. refractive_index_ = index;
  363. }
  364. // ********************************************************************** //
  365. // Modify coordinates for field calculation //
  366. // ********************************************************************** //
  367. void MultiLayerMie::SetFieldCoords(const std::vector< std::vector<double> >& coords) {
  368. if (coords.size() != 3)
  369. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  370. if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
  371. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  372. coords_ = coords;
  373. }
  374. // ********************************************************************** //
  375. // ********************************************************************** //
  376. // ********************************************************************** //
  377. void MultiLayerMie::SetPECLayer(int layer_position) {
  378. isExpCoeffsCalc_ = false;
  379. isScaCoeffsCalc_ = false;
  380. isMieCalculated_ = false;
  381. if (layer_position < 0 && layer_position != -1)
  382. throw std::invalid_argument("Error! Layers are numbered from 0!");
  383. PEC_layer_position_ = layer_position;
  384. }
  385. // ********************************************************************** //
  386. // Set maximun number of terms to be used //
  387. // ********************************************************************** //
  388. void MultiLayerMie::SetMaxTerms(int nmax) {
  389. isExpCoeffsCalc_ = false;
  390. isScaCoeffsCalc_ = false;
  391. isMieCalculated_ = false;
  392. nmax_preset_ = nmax;
  393. }
  394. // ********************************************************************** //
  395. // ********************************************************************** //
  396. // ********************************************************************** //
  397. double MultiLayerMie::GetSizeParameter() {
  398. if (size_param_.size() > 0)
  399. return size_param_.back();
  400. else
  401. return 0;
  402. }
  403. // ********************************************************************** //
  404. // Clear layer information //
  405. // ********************************************************************** //
  406. void MultiLayerMie::ClearLayers() {
  407. isExpCoeffsCalc_ = false;
  408. isScaCoeffsCalc_ = false;
  409. isMieCalculated_ = false;
  410. size_param_.clear();
  411. refractive_index_.clear();
  412. }
  413. // ********************************************************************** //
  414. // ********************************************************************** //
  415. // ********************************************************************** //
  416. // Computational core
  417. // ********************************************************************** //
  418. // ********************************************************************** //
  419. // ********************************************************************** //
  420. // ********************************************************************** //
  421. // Calculate calcNstop - equation (17) //
  422. // ********************************************************************** //
  423. void MultiLayerMie::calcNstop() {
  424. const double& xL = size_param_.back();
  425. if (xL <= 8) {
  426. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  427. } else if (xL <= 4200) {
  428. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  429. } else {
  430. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  431. }
  432. }
  433. // ********************************************************************** //
  434. // Maximum number of terms required for the calculation //
  435. // ********************************************************************** //
  436. void MultiLayerMie::calcNmax(unsigned int first_layer) {
  437. int ri, riM1;
  438. const std::vector<double>& x = size_param_;
  439. const std::vector<std::complex<double> >& m = refractive_index_;
  440. calcNstop(); // Set initial nmax_ value
  441. for (unsigned int i = first_layer; i < x.size(); i++) {
  442. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  443. ri = round(std::abs(x[i]*m[i]));
  444. else
  445. ri = 0;
  446. nmax_ = std::max(nmax_, ri);
  447. // first layer is pec, if pec is present
  448. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  449. riM1 = round(std::abs(x[i - 1]* m[i]));
  450. else
  451. riM1 = 0;
  452. nmax_ = std::max(nmax_, riM1);
  453. }
  454. nmax_ += 15; // Final nmax_ value
  455. }
  456. // ********************************************************************** //
  457. // Calculate an - equation (5) //
  458. // ********************************************************************** //
  459. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  460. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  461. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  462. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  463. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  464. return Num/Denom;
  465. }
  466. // ********************************************************************** //
  467. // Calculate bn - equation (6) //
  468. // ********************************************************************** //
  469. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  470. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  471. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  472. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  473. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  474. return Num/Denom;
  475. }
  476. // ********************************************************************** //
  477. // Calculate an and bn for bulk sphere size x and index m //
  478. // equation (4.56) and (4.57) BH //
  479. // ********************************************************************** //
  480. void MultiLayerMie::calc_an_bn_bulk(std::vector<std::complex<double> >& an,
  481. std::vector<std::complex<double> >& bn,
  482. double x, std::complex<double> m) {
  483. //printf("==========\n m = %g,%g, x= %g\n", std::real(m), std::imag(m), x);
  484. std::vector<std::complex<double> > PsiX(nmax_ + 1), ZetaX(nmax_ + 1);
  485. std::vector<std::complex<double> > PsiMX(nmax_ + 1), ZetaMX(nmax_ + 1);
  486. // First, calculate the Riccati-Bessel functions
  487. calcPsiZeta(x, PsiX, ZetaX);
  488. calcPsiZeta(m*x, PsiMX, ZetaMX);
  489. std::vector<std::complex<double> > D1X(nmax_ + 1), D3X(nmax_ + 1);
  490. std::vector<std::complex<double> > D1MX(nmax_ + 1), D3MX(nmax_ + 1);
  491. // Calculate the logarithmic derivatives
  492. calcD1D3(x, D1X, D3X);
  493. calcD1D3(m*x, D1MX, D3MX);
  494. std::vector<std::complex<double> > dPsiX(nmax_ + 1), dZetaX(nmax_ + 1);
  495. std::vector<std::complex<double> > dPsiMX(nmax_ + 1);
  496. for (int i = 0; i < nmax_ + 1; ++i) {
  497. dPsiX[i] = D1X[i]*PsiX[i];
  498. dPsiMX[i] = D1MX[i]*PsiMX[i];
  499. //dZetaX[i] = D3X[i]*ZetaX[i];
  500. }
  501. bessel::calcZeta(nmax_, x, ZetaX, dZetaX);
  502. an.resize(nmax_);
  503. bn.resize(nmax_);
  504. for (int i = 0; i < nmax_; i++) {
  505. int n = i+1;
  506. std::complex<double> Num = m*PsiMX[n]*dPsiX[n] - PsiX[n]*dPsiMX[n];
  507. std::complex<double> Denom = m*PsiMX[n]*dZetaX[n] - ZetaX[n]*dPsiMX[n];
  508. an[i] = Num/Denom;
  509. Num = PsiMX[n]*dPsiX[n] - m*PsiX[n]*dPsiMX[n];
  510. Denom = PsiMX[n]*dZetaX[n] - m*ZetaX[n]*dPsiMX[n];
  511. bn[i] = Num/Denom;
  512. }
  513. // printf("dPsiX\n");
  514. // for (auto a: dPsiX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  515. // printf("\ndPsiMX\n");
  516. // for (auto a: dPsiMX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  517. // printf("\nPsiX\n");
  518. // for (auto a: PsiX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  519. // printf("\nPsiMX\n");
  520. // for (auto a: PsiMX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  521. // printf("\nZetaX\n");
  522. // for (auto a: ZetaX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  523. // printf("\ndZetaX\n");
  524. // for (auto a: dZetaX) printf("%10.5er%+10.5ei ",std::real(a), std::imag(a));
  525. // printf("\nsize param: %g\n", x);
  526. }
  527. // ********************************************************************** //
  528. // Calculates S1 - equation (25a) //
  529. // ********************************************************************** //
  530. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  531. double Pi, double Tau) {
  532. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  533. }
  534. // ********************************************************************** //
  535. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  536. // Pi and Tau) //
  537. // ********************************************************************** //
  538. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  539. double Pi, double Tau) {
  540. return calc_S1(n, an, bn, Tau, Pi);
  541. }
  542. //**********************************************************************************//
  543. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  544. // functions (D1 and D3) for a complex argument (z). //
  545. // Equations (16a), (16b) and (18a) - (18d) //
  546. // //
  547. // Input parameters: //
  548. // z: Complex argument to evaluate D1 and D3 //
  549. // nmax_: Maximum number of terms to calculate D1 and D3 //
  550. // //
  551. // Output parameters: //
  552. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  553. //**********************************************************************************//
  554. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  555. std::vector<std::complex<double> >& D1,
  556. std::vector<std::complex<double> >& D3) {
  557. // // Downward recurrence for D1 - equations (16a) and (16b)
  558. // D1[nmax_] = std::complex<double>(0.0, 0.0);
  559. // const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  560. // for (int n = nmax_; n > 0; n--) {
  561. // D1[n - 1] = static_cast<double>(n)*zinv - 1.0/(D1[n] + static_cast<double>(n)*zinv);
  562. // }
  563. // if (std::abs(D1[0]) > 100000.0)
  564. // throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  565. // // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  566. // PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  567. // *std::exp(-2.0*z.imag()));
  568. // D3[0] = std::complex<double>(0.0, 1.0);
  569. // for (int n = 1; n <= nmax_; n++) {
  570. // PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  571. // *(static_cast<double>(n)*zinv - D3[n - 1]);
  572. // D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  573. // }
  574. std::vector<std::complex<double> > ZetaZ(nmax_+1), dZetaZ(nmax_ + 1);
  575. bessel::calcZeta(nmax_, z, ZetaZ, dZetaZ );
  576. for (int n = 0; n < nmax_+1; ++n) {
  577. D3[n]=dZetaZ[n]/ZetaZ[n];
  578. }
  579. std::vector<std::complex<double> > PsiZ(nmax_+1), dPsiZ(nmax_ + 1);
  580. bessel::calcPsi(nmax_, z, PsiZ, dPsiZ );
  581. for (int n = 0; n < nmax_+1; ++n) {
  582. D1[n]=dPsiZ[n]/PsiZ[n];
  583. }
  584. }
  585. //**********************************************************************************//
  586. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  587. // complex argument (z). //
  588. // Equations (20a) - (21b) //
  589. // //
  590. // Input parameters: //
  591. // z: Complex argument to evaluate Psi and Zeta //
  592. // nmax: Maximum number of terms to calculate Psi and Zeta //
  593. // //
  594. // Output parameters: //
  595. // Psi, Zeta: Riccati-Bessel functions //
  596. //**********************************************************************************//
  597. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  598. std::vector<std::complex<double> >& Psi,
  599. std::vector<std::complex<double> >& Zeta) {
  600. // std::complex<double> c_i(0.0, 1.0);
  601. // std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
  602. // // First, calculate the logarithmic derivatives
  603. // calcD1D3(z, D1, D3);
  604. // // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  605. // Psi[0] = std::sin(z);
  606. // Zeta[0] = std::sin(z) - c_i*std::cos(z);
  607. // for (int n = 1; n <= nmax_; n++) {
  608. // Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  609. // Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  610. // }
  611. std::vector<std::complex<double> > dZetaZ(nmax_ + 1);
  612. bessel::calcZeta(nmax_, z, Zeta, dZetaZ);
  613. std::vector<std::complex<double> > dPsiZ(nmax_ + 1);
  614. bessel::calcPsi(nmax_, z, Psi, dPsiZ );
  615. }
  616. //**********************************************************************************//
  617. // This function calculates the spherical Bessel (jn) and Hankel (h1n) functions //
  618. // and their derivatives for a given complex value z. See pag. 87 B&H. //
  619. // //
  620. // Input parameters: //
  621. // z: Complex argument to evaluate jn and h1n //
  622. // nmax_: Maximum number of terms to calculate jn and h1n //
  623. // //
  624. // Output parameters: //
  625. // jn, h1n: Spherical Bessel and Hankel functions //
  626. // jnp, h1np: Derivatives of the spherical Bessel and Hankel functions //
  627. // //
  628. // What we actually calculate are the Ricatti-Bessel fucntions and then simply //
  629. // evaluate the spherical Bessel and Hankel functions and their derivatives //
  630. // using the relations: //
  631. // //
  632. // j[n] = Psi[n]/z //
  633. // j'[n] = j[n-1]-(n+1)*jn[n])/z //
  634. // h1[n] = Zeta[n]/z //
  635. // h1'[n] = h1[n-1]-(n+1)*h1n[n]/z //
  636. // //
  637. //**********************************************************************************//
  638. void MultiLayerMie::sbesjh(std::complex<double> z,
  639. std::vector<std::complex<double> >& jn, std::vector<std::complex<double> >& jnp,
  640. std::vector<std::complex<double> >& h1n, std::vector<std::complex<double> >& h1np) {
  641. // std::vector<std::complex<double> > Psi(nmax_ + 1), Zeta(nmax_ + 1);
  642. // // First, calculate the Riccati-Bessel functions
  643. // calcPsiZeta(z, Psi, Zeta);
  644. // // Now, calculate Spherical Bessel and Hankel functions and their derivatives
  645. // for (int n = 0; n <= nmax_; n++) {
  646. // jn[n] = Psi[n]/z;
  647. // h1n[n] = Zeta[n]/z;
  648. // if (n == 0) {
  649. // jnp[0] = -Psi[1]/z - jn[0]/z;
  650. // h1np[0] = -Zeta[1]/z - h1n[0]/z;
  651. // } else {
  652. // jnp[n] = jn[n - 1] - static_cast<double>(n + 1)*jn[n]/z;
  653. // h1np[n] = h1n[n - 1] - static_cast<double>(n + 1)*h1n[n]/z;
  654. // }
  655. // }
  656. std::vector< std::complex<double> > yn, ynp;
  657. int nm;
  658. bessel::csphjy (nmax_, z, nm, jn, jnp, yn, ynp );
  659. auto c_i = std::complex<double>(0.0,1.0);
  660. h1n.resize(nmax_+1);
  661. h1np.resize(nmax_+1);
  662. for (int i = 0; i < nmax_+1; ++i) {
  663. h1n[i] = jn[i] + c_i*yn[i];
  664. h1np[i] = jnp[i] + c_i*ynp[i];
  665. }
  666. }
  667. //**********************************************************************************//
  668. // This function calculates Pi and Tau for a given value of cos(Theta). //
  669. // Equations (26a) - (26c) //
  670. // //
  671. // Input parameters: //
  672. // nmax_: Maximum number of terms to calculate Pi and Tau //
  673. // nTheta: Number of scattering angles //
  674. // Theta: Array containing all the scattering angles where the scattering //
  675. // amplitudes will be calculated //
  676. // //
  677. // Output parameters: //
  678. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  679. //**********************************************************************************//
  680. void MultiLayerMie::calcPiTau(const double& costheta,
  681. std::vector<double>& Pi, std::vector<double>& Tau) {
  682. int i;
  683. //****************************************************//
  684. // Equations (26a) - (26c) //
  685. //****************************************************//
  686. // Initialize Pi and Tau
  687. Pi[0] = 1.0; // n=1
  688. Tau[0] = costheta;
  689. // Calculate the actual values
  690. if (nmax_ > 1) {
  691. Pi[1] = 3*costheta*Pi[0]; //n=2
  692. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  693. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  694. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  695. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  696. }
  697. }
  698. } // end of MultiLayerMie::calcPiTau(...)
  699. //**********************************************************************************//
  700. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  701. // required to calculate the near-field parameters. //
  702. // //
  703. // Input parameters: //
  704. // Rho: Radial distance //
  705. // Phi: Azimuthal angle //
  706. // Theta: Polar angle //
  707. // zn: Either the spherical Bessel or Hankel function of first kind //
  708. // dzn: Derivative of zn //
  709. // Pi, Tau: Angular functions Pi and Tau //
  710. // n: Order of vector spherical harmonics //
  711. // //
  712. // Output parameters: //
  713. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  714. //**********************************************************************************//
  715. void MultiLayerMie::calcSpherHarm(const double Rho, const double Theta, const double Phi,
  716. const std::complex<double>& zn, const std::complex<double>& dzn,
  717. const double& Pi, const double& Tau, const double& n,
  718. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  719. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  720. // using eq 4.50 in BH
  721. std::complex<double> c_zero(0.0, 0.0);
  722. std::complex<double> deriv = Rho*dzn + zn;
  723. using std::sin;
  724. using std::cos;
  725. Mo1n[0] = c_zero;
  726. Mo1n[1] = cos(Phi)*Pi*zn;
  727. Mo1n[2] = -sin(Phi)*Tau*zn;
  728. Me1n[0] = c_zero;
  729. Me1n[1] = -sin(Phi)*Pi*zn;
  730. Me1n[2] = -cos(Phi)*Tau*zn;
  731. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  732. No1n[1] = sin(Phi)*Tau*deriv/Rho;
  733. No1n[2] = cos(Phi)*Pi*deriv/Rho;
  734. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*zn/Rho;
  735. Ne1n[1] = cos(Phi)*Tau*deriv/Rho;
  736. Ne1n[2] = -sin(Phi)*Pi*deriv/Rho;
  737. } // end of MultiLayerMie::calcSpherHarm(...)
  738. //**********************************************************************************//
  739. // This function calculates the scattering coefficients required to calculate //
  740. // both the near- and far-field parameters. //
  741. // //
  742. // Input parameters: //
  743. // L: Number of layers //
  744. // pl: Index of PEC layer. If there is none just send -1 //
  745. // x: Array containing the size parameters of the layers [0..L-1] //
  746. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  747. // nmax: Maximum number of multipolar expansion terms to be used for the //
  748. // calculations. Only use it if you know what you are doing, otherwise //
  749. // set this parameter to -1 and the function will calculate it. //
  750. // //
  751. // Output parameters: //
  752. // an, bn: Complex scattering amplitudes //
  753. // //
  754. // Return value: //
  755. // Number of multipolar expansion terms used for the calculations //
  756. //**********************************************************************************//
  757. void MultiLayerMie::ScattCoeffs() {
  758. isScaCoeffsCalc_ = false;
  759. const std::vector<double>& x = size_param_;
  760. const std::vector<std::complex<double> >& m = refractive_index_;
  761. const int& pl = PEC_layer_position_;
  762. const int L = refractive_index_.size();
  763. //************************************************************************//
  764. // Calculate the index of the first layer. It can be either 0 (default) //
  765. // or the index of the outermost PEC layer. In the latter case all layers //
  766. // below the PEC are discarded. //
  767. // ***********************************************************************//
  768. int fl = (pl > 0) ? pl : 0;
  769. if (nmax_preset_ <= 0) calcNmax(fl);
  770. else nmax_ = nmax_preset_;
  771. std::complex<double> z1, z2;
  772. //**************************************************************************//
  773. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  774. // means that index = layer number - 1 or index = n - 1. The only exception //
  775. // are the arrays for representing D1, D3 and Q because they need a value //
  776. // for the index 0 (zero), hence it is important to consider this shift //
  777. // between different arrays. The change was done to optimize memory usage. //
  778. //**************************************************************************//
  779. // Allocate memory to the arrays
  780. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  781. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  782. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  783. for (int l = 0; l < L; l++) {
  784. Q[l].resize(nmax_ + 1);
  785. Ha[l].resize(nmax_);
  786. Hb[l].resize(nmax_);
  787. }
  788. an_.resize(nmax_);
  789. bn_.resize(nmax_);
  790. PsiZeta_.resize(nmax_ + 1);
  791. std::vector<std::complex<double> > 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 Psi and Zeta for XL //
  870. //**************************************//
  871. // Calculate PsiXL and ZetaXL
  872. calcPsiZeta(x[L - 1], PsiXL, ZetaXL);
  873. //*********************************************************************//
  874. // Finally, we calculate the scattering coefficients (an and bn) and //
  875. // the angular functions (Pi and Tau). Note that for these arrays the //
  876. // first layer is 0 (zero), in future versions all arrays will follow //
  877. // this convention to save memory. (13 Nov, 2014) //
  878. //*********************************************************************//
  879. for (int n = 0; n < nmax_; n++) {
  880. //********************************************************************//
  881. //Expressions for calculating an and bn coefficients are not valid if //
  882. //there is only one PEC layer (ie, for a simple PEC sphere). //
  883. //********************************************************************//
  884. if (pl < (L - 1)) {
  885. 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]);
  886. 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]);
  887. } else {
  888. 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]);
  889. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  890. }
  891. } // end of for an and bn terms
  892. isScaCoeffsCalc_ = true;
  893. } // end of MultiLayerMie::ScattCoeffs(...)
  894. //**********************************************************************************//
  895. // This function calculates the actual scattering parameters and amplitudes //
  896. // //
  897. // Input parameters: //
  898. // L: Number of layers //
  899. // pl: Index of PEC layer. If there is none just send -1 //
  900. // x: Array containing the size parameters of the layers [0..L-1] //
  901. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  902. // nTheta: Number of scattering angles //
  903. // Theta: Array containing all the scattering angles where the scattering //
  904. // amplitudes will be calculated //
  905. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  906. // calculations. Only use it if you know what you are doing, otherwise //
  907. // set this parameter to -1 and the function will calculate it //
  908. // //
  909. // Output parameters: //
  910. // Qext: Efficiency factor for extinction //
  911. // Qsca: Efficiency factor for scattering //
  912. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  913. // Qbk: Efficiency factor for backscattering //
  914. // Qpr: Efficiency factor for the radiation pressure //
  915. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  916. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  917. // S1, S2: Complex scattering amplitudes //
  918. // //
  919. // Return value: //
  920. // Number of multipolar expansion terms used for the calculations //
  921. //**********************************************************************************//
  922. void MultiLayerMie::RunMieCalculation() {
  923. if (size_param_.size() != refractive_index_.size())
  924. throw std::invalid_argument("Each size parameter should have only one index!");
  925. if (size_param_.size() == 0)
  926. throw std::invalid_argument("Initialize model first!");
  927. const std::vector<double>& x = size_param_;
  928. isExpCoeffsCalc_ = false;
  929. isScaCoeffsCalc_ = false;
  930. isMieCalculated_ = false;
  931. // Calculate scattering coefficients
  932. ScattCoeffs();
  933. if (!isScaCoeffsCalc_) // TODO seems to be unreachable
  934. throw std::invalid_argument("Calculation of scattering coefficients failed!");
  935. // Initialize the scattering parameters
  936. Qext_ = 0.0;
  937. Qsca_ = 0.0;
  938. Qabs_ = 0.0;
  939. Qbk_ = 0.0;
  940. Qpr_ = 0.0;
  941. asymmetry_factor_ = 0.0;
  942. albedo_ = 0.0;
  943. // Initialize the scattering amplitudes
  944. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  945. S1_.swap(tmp1);
  946. S2_ = S1_;
  947. std::vector<double> Pi(nmax_), Tau(nmax_);
  948. std::complex<double> Qbktmp(0.0, 0.0);
  949. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  950. // By using downward recurrence we avoid loss of precision due to float rounding errors
  951. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  952. // http://en.wikipedia.org/wiki/Loss_of_significance
  953. for (int i = nmax_ - 2; i >= 0; i--) {
  954. const int n = i + 1;
  955. // Equation (27)
  956. Qext_ += (n + n + 1.0)*(an_[i].real() + bn_[i].real());
  957. // Equation (28)
  958. Qsca_ += (n + n + 1.0)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  959. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  960. // Equation (29)
  961. Qpr_ += ((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  962. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  963. // Equation (33)
  964. Qbktmp += (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  965. // Calculate the scattering amplitudes (S1 and S2) //
  966. // Equations (25a) - (25b) //
  967. for (unsigned int t = 0; t < theta_.size(); t++) {
  968. calcPiTau(std::cos(theta_[t]), Pi, Tau);
  969. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
  970. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
  971. }
  972. }
  973. double x2 = pow2(x.back());
  974. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  975. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  976. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  977. Qabs_ = Qext_ - Qsca_; // Equation (30)
  978. albedo_ = Qsca_/Qext_; // Equation (31)
  979. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  980. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  981. isMieCalculated_ = true;
  982. }
  983. //**********************************************************************************//
  984. // This function calculates the expansion coefficients inside the particle, //
  985. // required to calculate the near-field parameters. //
  986. // //
  987. // Input parameters: //
  988. // L: Number of layers //
  989. // pl: Index of PEC layer. If there is none just send -1 //
  990. // x: Array containing the size parameters of the layers [0..L-1] //
  991. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  992. // nmax: Maximum number of multipolar expansion terms to be used for the //
  993. // calculations. Only use it if you know what you are doing, otherwise //
  994. // set this parameter to -1 and the function will calculate it. //
  995. // //
  996. // Output parameters: //
  997. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  998. // //
  999. // Return value: //
  1000. // Number of multipolar expansion terms used for the calculations //
  1001. //**********************************************************************************//
  1002. void MultiLayerMie::ExpanCoeffs() {
  1003. if (!isScaCoeffsCalc_)
  1004. throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
  1005. isExpCoeffsCalc_ = false;
  1006. std::complex<double> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  1007. const int L = refractive_index_.size();
  1008. aln_.resize(L + 1);
  1009. bln_.resize(L + 1);
  1010. cln_.resize(L + 1);
  1011. dln_.resize(L + 1);
  1012. for (int l = 0; l <= L; l++) {
  1013. aln_[l].resize(nmax_);
  1014. bln_[l].resize(nmax_);
  1015. cln_[l].resize(nmax_);
  1016. dln_[l].resize(nmax_);
  1017. }
  1018. // Yang, paragraph under eq. A3
  1019. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  1020. for (int i = 0; i < nmax_; ++i) {
  1021. aln_[L][i] = an_[i];
  1022. bln_[L][i] = bn_[i];
  1023. cln_[L][i] = c_one;
  1024. dln_[L][i] = c_one;
  1025. }
  1026. std::vector<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  1027. std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  1028. std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
  1029. auto& m = refractive_index_;
  1030. std::vector< std::complex<double> > m1(L);
  1031. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  1032. m1[L - 1] = std::complex<double> (1.0, 0.0);
  1033. std::complex<double> z, z1;
  1034. for (int l = L - 1; l >= 0; l--) {
  1035. z = size_param_[l]*m[l];
  1036. z1 = size_param_[l]*m1[l];
  1037. calcD1D3(z, D1z, D3z);
  1038. calcD1D3(z1, D1z1, D3z1);
  1039. calcPsiZeta(z, Psiz, Zetaz);
  1040. calcPsiZeta(z1, Psiz1, Zetaz1);
  1041. for (int n = 0; n < nmax_; n++) {
  1042. int n1 = n + 1;
  1043. denomZeta = m1[l]*Zetaz[n1]*(D1z[n1] - D3z[n1]);
  1044. denomPsi = m1[l]*Psiz[n1]*(D1z[n1] - D3z[n1]);
  1045. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  1046. T2 = bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1];
  1047. T3 = D1z1[n1]*dln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*aln_[l + 1][n]*Zetaz1[n1];
  1048. T4 = D1z1[n1]*cln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*bln_[l + 1][n]*Zetaz1[n1];
  1049. // aln
  1050. aln_[l][n] = (D1z[n1]*m1[l]*T1 + m[l]*T3)/denomZeta;
  1051. // bln
  1052. bln_[l][n] = (D1z[n1]*m[l]*T2 + m1[l]*T4)/denomZeta;
  1053. // cln
  1054. cln_[l][n] = (D3z[n1]*m[l]*T2 + m1[l]*T4)/denomPsi;
  1055. // dln
  1056. dln_[l][n] = (D3z[n1]*m1[l]*T1 + m[l]*T3)/denomPsi;
  1057. } // end of all n
  1058. } // end of all l
  1059. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  1060. for (int n = 0; n < nmax_; ++n) {
  1061. printf("n=%d, aln_=%g,%g, bln_=%g,%g \n", n, real(aln_[0][n]), imag(aln_[0][n]),
  1062. real(bln_[0][n]), imag(bln_[0][n]));
  1063. if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  1064. else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  1065. if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  1066. else throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  1067. }
  1068. isExpCoeffsCalc_ = true;
  1069. } // end of void MultiLayerMie::ExpanCoeffs()
  1070. //**********************************************************************************//
  1071. // This function calculates the expansion coefficients inside the particle, //
  1072. // required to calculate the near-field parameters. //
  1073. // //
  1074. // Input parameters: //
  1075. // L: Number of layers //
  1076. // pl: Index of PEC layer. If there is none just send -1 //
  1077. // x: Array containing the size parameters of the layers [0..L-1] //
  1078. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1079. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1080. // calculations. Only use it if you know what you are doing, otherwise //
  1081. // set this parameter to -1 and the function will calculate it. //
  1082. // //
  1083. // Output parameters: //
  1084. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  1085. // //
  1086. // Return value: //
  1087. // Number of multipolar expansion terms used for the calculations //
  1088. //**********************************************************************************//
  1089. void MultiLayerMie::ExpanCoeffsV2() {
  1090. if (!isScaCoeffsCalc_)
  1091. throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
  1092. isExpCoeffsCalc_ = false;
  1093. std::complex<double> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  1094. const int L = refractive_index_.size();
  1095. aln_.resize(L + 1);
  1096. bln_.resize(L + 1);
  1097. cln_.resize(L + 1);
  1098. dln_.resize(L + 1);
  1099. for (int l = 0; l <= L; l++) {
  1100. aln_[l].resize(nmax_);
  1101. bln_[l].resize(nmax_);
  1102. cln_[l].resize(nmax_);
  1103. dln_[l].resize(nmax_);
  1104. }
  1105. // Yang, paragraph under eq. A3
  1106. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  1107. for (int n = 0; n < nmax_; n++) {
  1108. aln_[L][n] = an_[n];
  1109. bln_[L][n] = bn_[n];
  1110. cln_[L][n] = c_one;
  1111. dln_[L][n] = c_one;
  1112. }
  1113. std::vector<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  1114. std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  1115. std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
  1116. std::vector<std::vector<std::complex<double> > > a(2);
  1117. a[0].resize(3);
  1118. a[1].resize(3);
  1119. auto& m = refractive_index_;
  1120. std::vector< std::complex<double> > m1(L);
  1121. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  1122. m1[L - 1] = std::complex<double> (1.0, 0.0);
  1123. std::complex<double> z, z1;
  1124. for (int l = L - 1; l >= 0; l--) {
  1125. z = size_param_[l]*m[l];
  1126. z1 = size_param_[l]*m1[l];
  1127. calcD1D3(z, D1z, D3z);
  1128. calcD1D3(z1, D1z1, D3z1);
  1129. calcPsiZeta(z, Psiz, Zetaz);
  1130. calcPsiZeta(z1, Psiz1, Zetaz1);
  1131. for (int n = 0; n < nmax_; n++) {
  1132. int n1 = n + 1;
  1133. a[0][0] = m1[l]*D3z[n1]*Zetaz[n1];
  1134. a[0][1] = -m1[l]*D1z[n1]*Psiz[n1];
  1135. a[0][2] = aln_[l + 1][n]*m[l]*D3z1[n1]*Zetaz1[n1];
  1136. a[0][2] -= dln_[l + 1][n]*m[l]*D1z1[n1]*Psiz1[n1];
  1137. a[1][0] = Zetaz[n1];
  1138. a[1][1] = -Psiz[n1];
  1139. a[1][2] = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  1140. // aln
  1141. aln_[l][n] = (a[0][2]*a[1][1] - a[0][1]*a[1][2])/(a[0][0]*a[1][1] - a[0][1]*a[1][0]);
  1142. // dln
  1143. dln_[l][n] = (a[0][2]*a[1][0] - a[0][0]*a[1][2])/(a[0][1]*a[1][0] - a[0][0]*a[0][1]);
  1144. /*for (int i = 0; i < 2; i++) {
  1145. for (int j = 0; j < 3; j++) {
  1146. printf("a[%i, %i] = %g,%g ", i, j, real(a[i][j]), imag(a[i][j]));
  1147. }
  1148. printf("\n");
  1149. }
  1150. printf("aln_[%i, %i] = %g,%g; dln_[%i, %i] = %g,%g\n\n", l, n, real(aln_[l][n]), imag(aln_[l][n]), l, n, real(dln_[l][n]), imag(dln_[l][n]));*/
  1151. a[0][0] = D3z[n1]*Zetaz[n1];
  1152. a[0][1] = -D1z[n1]*Psiz[n1];
  1153. a[0][2] = bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  1154. a[0][2] -= cln_[l + 1][n]*D1z1[n1]*Psiz1[n1];
  1155. a[1][0] = m1[l]*Zetaz[n1];
  1156. a[1][1] = -m1[l]*Psiz[n1];
  1157. a[1][2] = bln_[l + 1][n]*m[l]*Zetaz1[n1] - cln_[l + 1][n]*m[l]*Psiz1[n1];
  1158. // bln
  1159. bln_[l][n] = (a[0][2]*a[1][1] - a[0][1]*a[1][2])/(a[0][0]*a[1][1] - a[0][1]*a[1][0]);
  1160. // cln
  1161. cln_[l][n] = (a[0][2]*a[1][0] - a[0][0]*a[1][2])/(a[0][1]*a[1][0] - a[0][0]*a[0][1]);
  1162. /*denomZeta = m1[l]*Zetaz[n1]*(D1z[n1] - D3z[n1]);
  1163. denomPsi = m1[l]*Psiz[n1]*(D1z[n1] - D3z[n1]);
  1164. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  1165. T2 = bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1];
  1166. T3 = D1z1[n1]*dln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*aln_[l + 1][n]*Zetaz1[n1];
  1167. T4 = D1z1[n1]*cln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*bln_[l + 1][n]*Zetaz1[n1];
  1168. // aln
  1169. aln_[l][n] = (D1z[n1]*m1[l]*T1 + m[l]*T3)/denomZeta;
  1170. // bln
  1171. bln_[l][n] = (D1z[n1]*m[l]*T2 + m1[l]*T4)/denomZeta;
  1172. // cln
  1173. cln_[l][n] = (D3z[n1]*m[l]*T2 + m1[l]*T4)/denomPsi;
  1174. // dln
  1175. dln_[l][n] = (D3z[n1]*m1[l]*T1 + m[l]*T3)/denomPsi;*/
  1176. printf("aln_[%i, %i] = %g,%g; bln_[%i, %i] = %g,%g; cln_[%i, %i] = %g,%g; dln_[%i, %i] = %g,%g\n", l, n, real(aln_[l][n]), imag(aln_[l][n]), l, n, real(bln_[l][n]), imag(bln_[l][n]), l, n, real(cln_[l][n]), imag(cln_[l][n]), l, n, real(dln_[l][n]), imag(dln_[l][n]));
  1177. } // end of all n
  1178. } // end of all l
  1179. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  1180. for (int n = 0; n < nmax_; ++n) {
  1181. // printf("n=%d, aln_=%g,%g, bln_=%g,%g \n", n, real(aln_[0][n]), imag(aln_[0][n]),
  1182. // real(bln_[0][n]), imag(bln_[0][n]));
  1183. if (std::abs(aln_[0][n]) < 1e-1) aln_[0][n] = 0.0;
  1184. else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  1185. if (std::abs(bln_[0][n]) < 1e-1) bln_[0][n] = 0.0;
  1186. else throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  1187. }
  1188. isExpCoeffsCalc_ = true;
  1189. } // end of void MultiLayerMie::ExpanCoeffs()
  1190. // ********************************************************************** //
  1191. // external scattering field = incident + scattered //
  1192. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1193. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1194. // ********************************************************************** //
  1195. void MultiLayerMie::fieldExt(const double Rho, const double Theta, const double Phi,
  1196. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1197. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1198. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1199. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1200. std::vector<std::complex<double> > Ei(3, c_zero), Hi(3, c_zero), Es(3, c_zero), Hs(3, c_zero);
  1201. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1202. std::vector<double> Pi(nmax_), Tau(nmax_);
  1203. // Calculate spherical Bessel and Hankel functions
  1204. sbesjh(Rho, jn, jnp, h1n, h1np);
  1205. // Calculate angular functions Pi and Tau
  1206. calcPiTau(std::cos(Theta), Pi, Tau);
  1207. for (int n = 0; n < nmax_; n++) {
  1208. int n1 = n + 1;
  1209. double rn = static_cast<double>(n1);
  1210. // using BH 4.12 and 4.50
  1211. calcSpherHarm(Rho, Theta, Phi, h1n[n1], h1np[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1212. // scattered field: BH p.94 (4.45)
  1213. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1214. for (int i = 0; i < 3; i++) {
  1215. Es[i] = Es[i] + En*(c_i*an_[n]*N3e1n[i] - bn_[n]*M3o1n[i]);
  1216. Hs[i] = Hs[i] + En*(c_i*bn_[n]*N3o1n[i] + an_[n]*M3e1n[i]);
  1217. }
  1218. }
  1219. // incident E field: BH p.89 (4.21); cf. p.92 (4.37), p.93 (4.38)
  1220. // basis unit vectors = er, etheta, ephi
  1221. std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
  1222. {
  1223. using std::sin;
  1224. using std::cos;
  1225. Ei[0] = eifac*sin(Theta)*cos(Phi);
  1226. Ei[1] = eifac*cos(Theta)*cos(Phi);
  1227. Ei[2] = -eifac*sin(Phi);
  1228. }
  1229. // magnetic field
  1230. double hffact = 1.0/(cc_*mu_);
  1231. for (int i = 0; i < 3; i++) {
  1232. Hs[i] = hffact*Hs[i];
  1233. }
  1234. // incident H field: BH p.26 (2.43), p.89 (4.21)
  1235. std::complex<double> hffacta = hffact;
  1236. std::complex<double> hifac = eifac*hffacta;
  1237. {
  1238. using std::sin;
  1239. using std::cos;
  1240. Hi[0] = hifac*sin(Theta)*sin(Phi);
  1241. Hi[1] = hifac*cos(Theta)*sin(Phi);
  1242. Hi[2] = hifac*cos(Phi);
  1243. }
  1244. for (int i = 0; i < 3; i++) {
  1245. // electric field E [V m - 1] = EF*E0
  1246. E[i] = Ei[i] + Es[i];
  1247. H[i] = Hi[i] + Hs[i];
  1248. }
  1249. } // end of MultiLayerMie::fieldExt(...)
  1250. //**********************************************************************************//
  1251. // This function calculates the electric (E) and magnetic (H) fields inside and //
  1252. // around the particle. //
  1253. // //
  1254. // Input parameters (coordinates of the point): //
  1255. // Rho: Radial distance //
  1256. // Phi: Azimuthal angle //
  1257. // Theta: Polar angle //
  1258. // //
  1259. // Output parameters: //
  1260. // E, H: Complex electric and magnetic fields //
  1261. //**********************************************************************************//
  1262. void MultiLayerMie::calcField(const double Rho, const double Theta, const double Phi,
  1263. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1264. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1265. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1266. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1267. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  1268. std::vector<std::complex<double> > jn(nmax_ + 1), jnp(nmax_ + 1), h1n(nmax_ + 1), h1np(nmax_ + 1);
  1269. std::vector<double> Pi(nmax_), Tau(nmax_);
  1270. std::vector<std::complex<double> > Ei(3), Hi(3);
  1271. int l = 0; // Layer number
  1272. std::complex<double> ml;
  1273. // Initialize E and H
  1274. for (int i = 0; i < 3; i++) {
  1275. E[i] = c_zero;
  1276. H[i] = c_zero;
  1277. }
  1278. if (Rho > size_param_.back()) {
  1279. l = size_param_.size();
  1280. ml = c_one;
  1281. } else {
  1282. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  1283. if (Rho <= size_param_[i]) {
  1284. l = i;
  1285. }
  1286. }
  1287. ml = refractive_index_[l];
  1288. }
  1289. // Calculate spherical Bessel and Hankel functions and their derivatives
  1290. sbesjh(Rho*ml, jn, jnp, h1n, h1np);
  1291. // Calculate angular functions Pi and Tau
  1292. calcPiTau(std::cos(Theta), Pi, Tau);
  1293. printf("Thetd = %g, cos(Theta) = %g\n", Theta, std::cos(Theta));
  1294. printf("Pi:\n");
  1295. for (auto p : Pi) printf("%11.4e\n",p);
  1296. printf("Tau:\n");
  1297. for (auto p : Tau) printf("%11.4e\n",p);
  1298. for (int n = nmax_ - 2; n >= 0; n--) {
  1299. int n1 = n + 1;
  1300. double rn = static_cast<double>(n1);
  1301. // using BH 4.12 and 4.50
  1302. calcSpherHarm(Rho, Theta, Phi, jn[n1], jnp[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1303. calcSpherHarm(Rho, Theta, Phi, h1n[n1], h1np[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1304. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1305. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1306. for (int i = 0; i < 3; i++) {
  1307. // electric field E [V m - 1] = EF*E0
  1308. E[i] += En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  1309. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  1310. H[i] += En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  1311. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  1312. Ei[i] += En*(M1o1n[i] - c_i*N1e1n[i]);
  1313. Hi[i] += En*(-M1e1n[i] - c_i*N1o1n[i]);
  1314. }
  1315. } // end of for all n
  1316. //printf("rho = %10.5e; phi = %10.5eº; theta = %10.5eº; x[%i] = %10.5e; m[%i] = %10.5er%+10.5ei\n", Rho, Phi*180./PI_, Theta*180./PI_, l, size_param_[l], l, std::real(ml), std::imag(ml));
  1317. // magnetic field
  1318. double hffact = 1.0/(cc_*mu_);
  1319. for (int i = 0; i < 3; i++) {
  1320. H[i] = hffact*H[i];
  1321. Hi[i] *= hffact;
  1322. printf("E[%i] = %10.5er%+10.5ei; H[%i] = %10.5er%+10.5ei\n", i, std::real(E[i]), std::imag(E[i]), i, std::real(H[i]), std::imag(H[i]));
  1323. }
  1324. } // end of MultiLayerMie::calcField(...)
  1325. //**********************************************************************************//
  1326. // This function calculates complex electric and magnetic field in the surroundings //
  1327. // and inside the particle. //
  1328. // //
  1329. // Input parameters: //
  1330. // L: Number of layers //
  1331. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1332. // x: Array containing the size parameters of the layers [0..L-1] //
  1333. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1334. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1335. // calculations. Only use it if you know what you are doing, otherwise //
  1336. // set this parameter to 0 (zero) and the function will calculate it. //
  1337. // ncoord: Number of coordinate points //
  1338. // Coords: Array containing all coordinates where the complex electric and //
  1339. // magnetic fields will be calculated //
  1340. // //
  1341. // Output parameters: //
  1342. // E, H: Complex electric and magnetic field at the provided coordinates //
  1343. // //
  1344. // Return value: //
  1345. // Number of multipolar expansion terms used for the calculations //
  1346. //**********************************************************************************//
  1347. void MultiLayerMie::RunFieldCalculation() {
  1348. double Rho, Theta, Phi;
  1349. // Calculate scattering coefficients an_ and bn_
  1350. ScattCoeffs();
  1351. // std::vector<std::complex<double> > an1(nmax_), bn1(nmax_);
  1352. // calc_an_bn_bulk(an1, bn1, size_param_.back(), refractive_index_.back());
  1353. // for (int n = 0; n < nmax_; n++) {
  1354. // printf("an_[%i] = %10.5er%+10.5ei; an_bulk_[%i] = %10.5er%+10.5ei\n", n, std::real(an_[n]), std::imag(an_[n]), n, std::real(an1[n]), std::imag(an1[n]));
  1355. // printf("bn_[%i] = %10.5er%+10.5ei; bn_bulk_[%i] = %10.5er%+10.5ei\n", n, std::real(bn_[n]), std::imag(bn_[n]), n, std::real(bn1[n]), std::imag(bn1[n]));
  1356. // }
  1357. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  1358. ExpanCoeffs();
  1359. // for (int i = 0; i < nmax_; ++i) {
  1360. // printf("cln_[%i] = %10.5er%+10.5ei; dln_[%i] = %10.5er%+10.5ei\n", i, std::real(cln_[0][i]), std::imag(cln_[0][i]),
  1361. // i, std::real(dln_[0][i]), std::imag(dln_[0][i]));
  1362. // }
  1363. long total_points = coords_[0].size();
  1364. E_.resize(total_points);
  1365. H_.resize(total_points);
  1366. for (auto& f : E_) f.resize(3);
  1367. for (auto& f : H_) f.resize(3);
  1368. for (int point = 0; point < total_points; point++) {
  1369. const double& Xp = coords_[0][point];
  1370. const double& Yp = coords_[1][point];
  1371. const double& Zp = coords_[2][point];
  1372. // Convert to spherical coordinates
  1373. Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1374. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1375. Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1376. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1377. if (Xp == 0.0)
  1378. Phi = (Yp != 0.0) ? std::asin(Yp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1379. else
  1380. Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
  1381. // Avoid convergence problems due to Rho too small
  1382. if (Rho < 1e-5) Rho = 1e-5;
  1383. //printf("X = %g; Y = %g; Z = %g; pho = %g; phi = %g; theta = %g\n", Xp, Yp, Zp, Rho, Phi*180./PI_, Theta*180./PI_);
  1384. //*******************************************************//
  1385. // external scattering field = incident + scattered //
  1386. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1387. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1388. //*******************************************************//
  1389. // This array contains the fields in spherical coordinates
  1390. std::vector<std::complex<double> > Es(3), Hs(3);
  1391. // Firstly the easiest case: the field outside the particle
  1392. // if (Rho >= GetSizeParameter()) {
  1393. // fieldExt(Rho, Theta, Phi, Es, Hs);
  1394. // } else {
  1395. calcField(Rho, Theta, Phi, Es, Hs); //Should work fine both: inside and outside the particle
  1396. //}
  1397. { //Now, convert the fields back to cartesian coordinates
  1398. using std::sin;
  1399. using std::cos;
  1400. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1401. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1402. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1403. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1404. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1405. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1406. }
  1407. } // end of for all field coordinates
  1408. } // end of MultiLayerMie::RunFieldCalculation()
  1409. } // end of namespace nmie