nmie.cc 70 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319
  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 scattering coefficients //
  54. // required to calculate both the near- and far-field parameters. //
  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. // nmax: Maximum number of multipolar expansion terms to be used for the //
  62. // calculations. Only use it if you know what you are doing, otherwise //
  63. // set this parameter to -1 and the function will calculate it. //
  64. // //
  65. // Output parameters: //
  66. // an, bn: Complex scattering amplitudes //
  67. // //
  68. // Return value: //
  69. // Number of multipolar expansion terms used for the calculations //
  70. //**********************************************************************************//
  71. int ScattCoeffs(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const int nmax, std::vector<std::complex<double> >& an, std::vector<std::complex<double> >& bn) {
  72. if (x.size() != L || m.size() != L)
  73. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  74. try {
  75. MultiLayerMie ml_mie;
  76. ml_mie.SetAllLayersSize(x);
  77. ml_mie.SetAllLayersIndex(m);
  78. ml_mie.SetPECLayer(pl);
  79. ml_mie.SetMaxTerms(nmax);
  80. ml_mie.calcScattCoeffs();
  81. an = ml_mie.GetAn();
  82. bn = ml_mie.GetBn();
  83. return ml_mie.GetMaxTerms();
  84. } catch(const std::invalid_argument& ia) {
  85. // Will catch if ml_mie fails or other errors.
  86. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  87. throw std::invalid_argument(ia);
  88. return -1;
  89. }
  90. return 0;
  91. }
  92. //**********************************************************************************//
  93. // This function emulates a C call to calculate the actual scattering parameters //
  94. // and amplitudes. //
  95. // //
  96. // Input parameters: //
  97. // L: Number of layers //
  98. // pl: Index of PEC layer. If there is none just send -1 //
  99. // x: Array containing the size parameters of the layers [0..L-1] //
  100. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  101. // nTheta: Number of scattering angles //
  102. // Theta: Array containing all the scattering angles where the scattering //
  103. // amplitudes will be calculated //
  104. // nmax: Maximum number of multipolar expansion terms to be used for the //
  105. // calculations. Only use it if you know what you are doing, otherwise //
  106. // set this parameter to -1 and the function will calculate it //
  107. // //
  108. // Output parameters: //
  109. // Qext: Efficiency factor for extinction //
  110. // Qsca: Efficiency factor for scattering //
  111. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  112. // Qbk: Efficiency factor for backscattering //
  113. // Qpr: Efficiency factor for the radiation pressure //
  114. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  115. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  116. // S1, S2: Complex scattering amplitudes //
  117. // //
  118. // Return value: //
  119. // Number of multipolar expansion terms used for the calculations //
  120. //**********************************************************************************//
  121. 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) {
  122. if (x.size() != L || m.size() != L)
  123. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  124. if (Theta.size() != nTheta)
  125. throw std::invalid_argument("Declared number of sample for Theta is not correct!");
  126. try {
  127. MultiLayerMie ml_mie;
  128. ml_mie.SetAllLayersSize(x);
  129. ml_mie.SetAllLayersIndex(m);
  130. ml_mie.SetAngles(Theta);
  131. ml_mie.SetPECLayer(pl);
  132. ml_mie.SetMaxTerms(nmax);
  133. ml_mie.RunMieCalculation();
  134. *Qext = ml_mie.GetQext();
  135. *Qsca = ml_mie.GetQsca();
  136. *Qabs = ml_mie.GetQabs();
  137. *Qbk = ml_mie.GetQbk();
  138. *Qpr = ml_mie.GetQpr();
  139. *g = ml_mie.GetAsymmetryFactor();
  140. *Albedo = ml_mie.GetAlbedo();
  141. S1 = ml_mie.GetS1();
  142. S2 = ml_mie.GetS2();
  143. return ml_mie.GetMaxTerms();
  144. } catch(const std::invalid_argument& ia) {
  145. // Will catch if ml_mie fails or other errors.
  146. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  147. throw std::invalid_argument(ia);
  148. return -1;
  149. }
  150. return 0;
  151. }
  152. //**********************************************************************************//
  153. // This function is just a wrapper to call the full 'nMie' function with fewer //
  154. // parameters, it is here mainly for compatibility with older versions of the //
  155. // program. Also, you can use it if you neither have a PEC layer nor want to define //
  156. // any limit for the maximum number of terms. //
  157. // //
  158. // Input parameters: //
  159. // L: Number of layers //
  160. // x: Array containing the size parameters of the layers [0..L-1] //
  161. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  162. // nTheta: Number of scattering angles //
  163. // Theta: Array containing all the scattering angles where the scattering //
  164. // amplitudes will be calculated //
  165. // //
  166. // Output parameters: //
  167. // Qext: Efficiency factor for extinction //
  168. // Qsca: Efficiency factor for scattering //
  169. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  170. // Qbk: Efficiency factor for backscattering //
  171. // Qpr: Efficiency factor for the radiation pressure //
  172. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  173. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  174. // S1, S2: Complex scattering amplitudes //
  175. // //
  176. // Return value: //
  177. // Number of multipolar expansion terms used for the calculations //
  178. //**********************************************************************************//
  179. 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) {
  180. return nmie::nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  181. }
  182. //**********************************************************************************//
  183. // This function is just a wrapper to call the full 'nMie' function with fewer //
  184. // parameters, it is useful if you want to include a PEC layer but not a limit //
  185. // for the maximum number of terms. //
  186. // //
  187. // Input parameters: //
  188. // L: Number of layers //
  189. // pl: Index of PEC layer. If there is none just send -1 //
  190. // x: Array containing the size parameters of the layers [0..L-1] //
  191. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  192. // nTheta: Number of scattering angles //
  193. // Theta: Array containing all the scattering angles where the scattering //
  194. // amplitudes will be calculated //
  195. // //
  196. // Output parameters: //
  197. // Qext: Efficiency factor for extinction //
  198. // Qsca: Efficiency factor for scattering //
  199. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  200. // Qbk: Efficiency factor for backscattering //
  201. // Qpr: Efficiency factor for the radiation pressure //
  202. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  203. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  204. // S1, S2: Complex scattering amplitudes //
  205. // //
  206. // Return value: //
  207. // Number of multipolar expansion terms used for the calculations //
  208. //**********************************************************************************//
  209. 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) {
  210. return nmie::nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  211. }
  212. //**********************************************************************************//
  213. // This function is just a wrapper to call the full 'nMie' function with fewer //
  214. // parameters, it is useful if you want to include a limit for the maximum number //
  215. // of terms but not a PEC layer. //
  216. // //
  217. // Input parameters: //
  218. // L: Number of layers //
  219. // x: Array containing the size parameters of the layers [0..L-1] //
  220. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  221. // nTheta: Number of scattering angles //
  222. // Theta: Array containing all the scattering angles where the scattering //
  223. // amplitudes will be calculated //
  224. // nmax: Maximum number of multipolar expansion terms to be used for the //
  225. // calculations. Only use it if you know what you are doing, otherwise //
  226. // set this parameter to -1 and the function will calculate it //
  227. // //
  228. // Output parameters: //
  229. // Qext: Efficiency factor for extinction //
  230. // Qsca: Efficiency factor for scattering //
  231. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  232. // Qbk: Efficiency factor for backscattering //
  233. // Qpr: Efficiency factor for the radiation pressure //
  234. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  235. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  236. // S1, S2: Complex scattering amplitudes //
  237. // //
  238. // Return value: //
  239. // Number of multipolar expansion terms used for the calculations //
  240. //**********************************************************************************//
  241. 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) {
  242. return nmie::nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
  243. }
  244. //**********************************************************************************//
  245. // This function emulates a C call to calculate complex electric and magnetic field //
  246. // in the surroundings and inside the particle. //
  247. // //
  248. // Input parameters: //
  249. // L: Number of layers //
  250. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  251. // x: Array containing the size parameters of the layers [0..L-1] //
  252. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  253. // nmax: Maximum number of multipolar expansion terms to be used for the //
  254. // calculations. Only use it if you know what you are doing, otherwise //
  255. // set this parameter to 0 (zero) and the function will calculate it. //
  256. // ncoord: Number of coordinate points //
  257. // Coords: Array containing all coordinates where the complex electric and //
  258. // magnetic fields will be calculated //
  259. // //
  260. // Output parameters: //
  261. // E, H: Complex electric and magnetic field at the provided coordinates //
  262. // //
  263. // Return value: //
  264. // Number of multipolar expansion terms used for the calculations //
  265. //**********************************************************************************//
  266. 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) {
  267. if (x.size() != L || m.size() != L)
  268. throw std::invalid_argument("Declared number of layers do not fit x and m!");
  269. if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
  270. || E.size() != ncoord || H.size() != ncoord)
  271. throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
  272. for (auto f:E)
  273. if (f.size() != 3)
  274. throw std::invalid_argument("Field E is not 3D!");
  275. for (auto f:H)
  276. if (f.size() != 3)
  277. throw std::invalid_argument("Field H is not 3D!");
  278. try {
  279. MultiLayerMie ml_mie;
  280. ml_mie.SetPECLayer(pl);
  281. ml_mie.SetAllLayersSize(x);
  282. ml_mie.SetAllLayersIndex(m);
  283. ml_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
  284. ml_mie.RunFieldCalculation();
  285. E = ml_mie.GetFieldE();
  286. H = ml_mie.GetFieldH();
  287. return ml_mie.GetMaxTerms();
  288. } catch(const std::invalid_argument& ia) {
  289. // Will catch if ml_mie fails or other errors.
  290. std::cerr << "Invalid argument: " << ia.what() << std::endl;
  291. throw std::invalid_argument(ia);
  292. return - 1;
  293. }
  294. return 0;
  295. }
  296. // ********************************************************************** //
  297. // Returns previously calculated Qext //
  298. // ********************************************************************** //
  299. double MultiLayerMie::GetQext() {
  300. if (!isMieCalculated_)
  301. throw std::invalid_argument("You should run calculations before result request!");
  302. return Qext_;
  303. }
  304. // ********************************************************************** //
  305. // Returns previously calculated Qabs //
  306. // ********************************************************************** //
  307. double MultiLayerMie::GetQabs() {
  308. if (!isMieCalculated_)
  309. throw std::invalid_argument("You should run calculations before result request!");
  310. return Qabs_;
  311. }
  312. // ********************************************************************** //
  313. // Returns previously calculated Qsca //
  314. // ********************************************************************** //
  315. double MultiLayerMie::GetQsca() {
  316. if (!isMieCalculated_)
  317. throw std::invalid_argument("You should run calculations before result request!");
  318. return Qsca_;
  319. }
  320. // ********************************************************************** //
  321. // Returns previously calculated Qbk //
  322. // ********************************************************************** //
  323. double MultiLayerMie::GetQbk() {
  324. if (!isMieCalculated_)
  325. throw std::invalid_argument("You should run calculations before result request!");
  326. return Qbk_;
  327. }
  328. // ********************************************************************** //
  329. // Returns previously calculated Qpr //
  330. // ********************************************************************** //
  331. double MultiLayerMie::GetQpr() {
  332. if (!isMieCalculated_)
  333. throw std::invalid_argument("You should run calculations before result request!");
  334. return Qpr_;
  335. }
  336. // ********************************************************************** //
  337. // Returns previously calculated assymetry factor //
  338. // ********************************************************************** //
  339. double MultiLayerMie::GetAsymmetryFactor() {
  340. if (!isMieCalculated_)
  341. throw std::invalid_argument("You should run calculations before result request!");
  342. return asymmetry_factor_;
  343. }
  344. // ********************************************************************** //
  345. // Returns previously calculated Albedo //
  346. // ********************************************************************** //
  347. double MultiLayerMie::GetAlbedo() {
  348. if (!isMieCalculated_)
  349. throw std::invalid_argument("You should run calculations before result request!");
  350. return albedo_;
  351. }
  352. // ********************************************************************** //
  353. // Returns previously calculated S1 //
  354. // ********************************************************************** //
  355. std::vector<std::complex<double> > MultiLayerMie::GetS1() {
  356. if (!isMieCalculated_)
  357. throw std::invalid_argument("You should run calculations before result request!");
  358. return S1_;
  359. }
  360. // ********************************************************************** //
  361. // Returns previously calculated S2 //
  362. // ********************************************************************** //
  363. std::vector<std::complex<double> > MultiLayerMie::GetS2() {
  364. if (!isMieCalculated_)
  365. throw std::invalid_argument("You should run calculations before result request!");
  366. return S2_;
  367. }
  368. // ********************************************************************** //
  369. // Modify scattering (theta) angles //
  370. // ********************************************************************** //
  371. void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
  372. isExpCoeffsCalc_ = false;
  373. isScaCoeffsCalc_ = false;
  374. isMieCalculated_ = false;
  375. theta_ = angles;
  376. }
  377. // ********************************************************************** //
  378. // Modify size of all layers //
  379. // ********************************************************************** //
  380. void MultiLayerMie::SetAllLayersSize(const std::vector<double>& layer_size) {
  381. isExpCoeffsCalc_ = false;
  382. isScaCoeffsCalc_ = false;
  383. isMieCalculated_ = false;
  384. size_param_.clear();
  385. double prev_layer_size = 0.0;
  386. for (auto curr_layer_size : layer_size) {
  387. if (curr_layer_size <= 0.0)
  388. throw std::invalid_argument("Size parameter should be positive!");
  389. if (prev_layer_size > curr_layer_size)
  390. throw std::invalid_argument
  391. ("Size parameter for next layer should be larger than the previous one!");
  392. prev_layer_size = curr_layer_size;
  393. size_param_.push_back(curr_layer_size);
  394. }
  395. }
  396. // ********************************************************************** //
  397. // Modify refractive index of all layers //
  398. // ********************************************************************** //
  399. void MultiLayerMie::SetAllLayersIndex(const std::vector< std::complex<double> >& index) {
  400. isExpCoeffsCalc_ = false;
  401. isScaCoeffsCalc_ = false;
  402. isMieCalculated_ = false;
  403. refractive_index_ = index;
  404. }
  405. // ********************************************************************** //
  406. // Modify coordinates for field calculation //
  407. // ********************************************************************** //
  408. void MultiLayerMie::SetFieldCoords(const std::vector< std::vector<double> >& coords) {
  409. if (coords.size() != 3)
  410. throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
  411. if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
  412. throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
  413. coords_ = coords;
  414. }
  415. // ********************************************************************** //
  416. // ********************************************************************** //
  417. // ********************************************************************** //
  418. void MultiLayerMie::SetPECLayer(int layer_position) {
  419. isExpCoeffsCalc_ = false;
  420. isScaCoeffsCalc_ = false;
  421. isMieCalculated_ = false;
  422. if (layer_position < 0 && layer_position != -1)
  423. throw std::invalid_argument("Error! Layers are numbered from 0!");
  424. PEC_layer_position_ = layer_position;
  425. }
  426. // ********************************************************************** //
  427. // Set maximun number of terms to be used //
  428. // ********************************************************************** //
  429. void MultiLayerMie::SetMaxTerms(int nmax) {
  430. isExpCoeffsCalc_ = false;
  431. isScaCoeffsCalc_ = false;
  432. isMieCalculated_ = false;
  433. nmax_preset_ = nmax;
  434. }
  435. // ********************************************************************** //
  436. // ********************************************************************** //
  437. // ********************************************************************** //
  438. double MultiLayerMie::GetSizeParameter() {
  439. if (size_param_.size() > 0)
  440. return size_param_.back();
  441. else
  442. return 0;
  443. }
  444. // ********************************************************************** //
  445. // Clear layer information //
  446. // ********************************************************************** //
  447. void MultiLayerMie::ClearLayers() {
  448. isExpCoeffsCalc_ = false;
  449. isScaCoeffsCalc_ = false;
  450. isMieCalculated_ = false;
  451. size_param_.clear();
  452. refractive_index_.clear();
  453. }
  454. // ********************************************************************** //
  455. // ********************************************************************** //
  456. // ********************************************************************** //
  457. // Computational core
  458. // ********************************************************************** //
  459. // ********************************************************************** //
  460. // ********************************************************************** //
  461. // ********************************************************************** //
  462. // Calculate calcNstop - equation (17) //
  463. // ********************************************************************** //
  464. void MultiLayerMie::calcNstop() {
  465. const double& xL = size_param_.back();
  466. if (xL <= 8) {
  467. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
  468. } else if (xL <= 4200) {
  469. nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
  470. } else {
  471. nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
  472. }
  473. }
  474. // ********************************************************************** //
  475. // Maximum number of terms required for the calculation //
  476. // ********************************************************************** //
  477. void MultiLayerMie::calcNmax(unsigned int first_layer) {
  478. int ri, riM1;
  479. const std::vector<double>& x = size_param_;
  480. const std::vector<std::complex<double> >& m = refractive_index_;
  481. calcNstop(); // Set initial nmax_ value
  482. for (unsigned int i = first_layer; i < x.size(); i++) {
  483. if (static_cast<int>(i) > PEC_layer_position_) // static_cast used to avoid warning
  484. ri = round(std::abs(x[i]*m[i]));
  485. else
  486. ri = 0;
  487. nmax_ = std::max(nmax_, ri);
  488. // first layer is pec, if pec is present
  489. if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
  490. riM1 = round(std::abs(x[i - 1]* m[i]));
  491. else
  492. riM1 = 0;
  493. nmax_ = std::max(nmax_, riM1);
  494. }
  495. nmax_ += 15; // Final nmax_ value
  496. }
  497. // ********************************************************************** //
  498. // Calculate an - equation (5) //
  499. // ********************************************************************** //
  500. std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
  501. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  502. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  503. std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
  504. std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
  505. return Num/Denom;
  506. }
  507. // ********************************************************************** //
  508. // Calculate bn - equation (6) //
  509. // ********************************************************************** //
  510. std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
  511. std::complex<double> PsiXL, std::complex<double> ZetaXL,
  512. std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
  513. std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
  514. std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
  515. return Num/Denom;
  516. }
  517. // ********************************************************************** //
  518. // Calculates S1 - equation (25a) //
  519. // ********************************************************************** //
  520. std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
  521. double Pi, double Tau) {
  522. return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
  523. }
  524. // ********************************************************************** //
  525. // Calculates S2 - equation (25b) (it's the same as (25a), just switches //
  526. // Pi and Tau) //
  527. // ********************************************************************** //
  528. std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
  529. double Pi, double Tau) {
  530. return calc_S1(n, an, bn, Tau, Pi);
  531. }
  532. //**********************************************************************************//
  533. // This function calculates the logarithmic derivatives of the Riccati-Bessel //
  534. // functions (D1 and D3) for a complex argument (z). //
  535. // Equations (16a), (16b) and (18a) - (18d) //
  536. // //
  537. // Input parameters: //
  538. // z: Complex argument to evaluate D1 and D3 //
  539. // nmax_: Maximum number of terms to calculate D1 and D3 //
  540. // //
  541. // Output parameters: //
  542. // D1, D3: Logarithmic derivatives of the Riccati-Bessel functions //
  543. //**********************************************************************************//
  544. void MultiLayerMie::calcD1D3(const std::complex<double> z,
  545. std::vector<std::complex<double> >& D1,
  546. std::vector<std::complex<double> >& D3) {
  547. // Downward recurrence for D1 - equations (16a) and (16b)
  548. D1[nmax_] = std::complex<double>(0.0, 0.0);
  549. const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
  550. for (int n = nmax_; n > 0; n--) {
  551. D1[n - 1] = static_cast<double>(n)*zinv - 1.0/(D1[n] + static_cast<double>(n)*zinv);
  552. }
  553. if (std::abs(D1[0]) > 1.0e15) {
  554. throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
  555. //printf("Warning: Potentially unstable D1! Please, try to change input parameters!\n");
  556. }
  557. // Upward recurrence for PsiZeta and D3 - equations (18a) - (18d)
  558. PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
  559. *std::exp(-2.0*z.imag()));
  560. D3[0] = std::complex<double>(0.0, 1.0);
  561. for (int n = 1; n <= nmax_; n++) {
  562. PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
  563. *(static_cast<double>(n)*zinv - D3[n - 1]);
  564. D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
  565. }
  566. }
  567. //**********************************************************************************//
  568. // This function calculates the Riccati-Bessel functions (Psi and Zeta) for a //
  569. // complex argument (z). //
  570. // Equations (20a) - (21b) //
  571. // //
  572. // Input parameters: //
  573. // z: Complex argument to evaluate Psi and Zeta //
  574. // nmax: Maximum number of terms to calculate Psi and Zeta //
  575. // //
  576. // Output parameters: //
  577. // Psi, Zeta: Riccati-Bessel functions //
  578. //**********************************************************************************//
  579. void MultiLayerMie::calcPsiZeta(std::complex<double> z,
  580. std::vector<std::complex<double> >& Psi,
  581. std::vector<std::complex<double> >& Zeta) {
  582. std::complex<double> c_i(0.0, 1.0);
  583. std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
  584. // First, calculate the logarithmic derivatives
  585. calcD1D3(z, D1, D3);
  586. // Now, use the upward recurrence to calculate Psi and Zeta - equations (20a) - (21b)
  587. Psi[0] = std::sin(z);
  588. Zeta[0] = std::sin(z) - c_i*std::cos(z);
  589. for (int n = 1; n <= nmax_; n++) {
  590. Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
  591. Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
  592. }
  593. }
  594. //**********************************************************************************//
  595. // This function calculates Pi and Tau for a given value of cos(Theta). //
  596. // Equations (26a) - (26c) //
  597. // //
  598. // Input parameters: //
  599. // nmax_: Maximum number of terms to calculate Pi and Tau //
  600. // nTheta: Number of scattering angles //
  601. // Theta: Array containing all the scattering angles where the scattering //
  602. // amplitudes will be calculated //
  603. // //
  604. // Output parameters: //
  605. // Pi, Tau: Angular functions Pi and Tau, as defined in equations (26a) - (26c) //
  606. //**********************************************************************************//
  607. void MultiLayerMie::calcPiTau(const double& costheta,
  608. std::vector<double>& Pi, std::vector<double>& Tau) {
  609. int i;
  610. //****************************************************//
  611. // Equations (26a) - (26c) //
  612. //****************************************************//
  613. // Initialize Pi and Tau
  614. Pi[0] = 1.0; // n=1
  615. Tau[0] = costheta;
  616. // Calculate the actual values
  617. if (nmax_ > 1) {
  618. Pi[1] = 3*costheta*Pi[0]; //n=2
  619. Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
  620. for (i = 2; i < nmax_; i++) { //n=[3..nmax_]
  621. Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
  622. Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
  623. }
  624. }
  625. } // end of MultiLayerMie::calcPiTau(...)
  626. //**********************************************************************************//
  627. // This function calculates vector spherical harmonics (eq. 4.50, p. 95 BH), //
  628. // required to calculate the near-field parameters. //
  629. // //
  630. // Input parameters: //
  631. // Rho: Radial distance //
  632. // Phi: Azimuthal angle //
  633. // Theta: Polar angle //
  634. // rn: Either the spherical Ricatti-Bessel function of first or third kind //
  635. // Dn: Logarithmic derivative of rn //
  636. // Pi, Tau: Angular functions Pi and Tau //
  637. // n: Order of vector spherical harmonics //
  638. // //
  639. // Output parameters: //
  640. // Mo1n, Me1n, No1n, Ne1n: Complex vector spherical harmonics //
  641. //**********************************************************************************//
  642. void MultiLayerMie::calcSpherHarm(const std::complex<double> Rho, const double Theta, const double Phi,
  643. const std::complex<double>& rn, const std::complex<double>& Dn,
  644. const double& Pi, const double& Tau, const double& n,
  645. std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
  646. std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
  647. // using eq 4.50 in BH
  648. std::complex<double> c_zero(0.0, 0.0);
  649. using std::sin;
  650. using std::cos;
  651. Mo1n[0] = c_zero;
  652. Mo1n[1] = cos(Phi)*Pi*rn/Rho;
  653. Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
  654. Me1n[0] = c_zero;
  655. Me1n[1] = -sin(Phi)*Pi*rn/Rho;
  656. Me1n[2] = -cos(Phi)*Tau*rn/Rho;
  657. No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  658. No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
  659. No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
  660. Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
  661. Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
  662. Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
  663. } // end of MultiLayerMie::calcSpherHarm(...)
  664. //**********************************************************************************//
  665. // This function calculates the scattering coefficients required to calculate //
  666. // both the near- and far-field parameters. //
  667. // //
  668. // Input parameters: //
  669. // L: Number of layers //
  670. // pl: Index of PEC layer. If there is none just send -1 //
  671. // x: Array containing the size parameters of the layers [0..L-1] //
  672. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  673. // nmax: Maximum number of multipolar expansion terms to be used for the //
  674. // calculations. Only use it if you know what you are doing, otherwise //
  675. // set this parameter to -1 and the function will calculate it. //
  676. // //
  677. // Output parameters: //
  678. // an, bn: Complex scattering amplitudes //
  679. // //
  680. // Return value: //
  681. // Number of multipolar expansion terms used for the calculations //
  682. //**********************************************************************************//
  683. void MultiLayerMie::calcScattCoeffs() {
  684. isScaCoeffsCalc_ = false;
  685. const std::vector<double>& x = size_param_;
  686. const std::vector<std::complex<double> >& m = refractive_index_;
  687. const int& pl = PEC_layer_position_;
  688. const int L = refractive_index_.size();
  689. //************************************************************************//
  690. // Calculate the index of the first layer. It can be either 0 (default) //
  691. // or the index of the outermost PEC layer. In the latter case all layers //
  692. // below the PEC are discarded. //
  693. // ***********************************************************************//
  694. int fl = (pl > 0) ? pl : 0;
  695. if (nmax_preset_ <= 0) calcNmax(fl);
  696. else nmax_ = nmax_preset_;
  697. std::complex<double> z1, z2;
  698. //**************************************************************************//
  699. // Note that since Fri, Nov 14, 2014 all arrays start from 0 (zero), which //
  700. // means that index = layer number - 1 or index = n - 1. The only exception //
  701. // are the arrays for representing D1, D3 and Q because they need a value //
  702. // for the index 0 (zero), hence it is important to consider this shift //
  703. // between different arrays. The change was done to optimize memory usage. //
  704. //**************************************************************************//
  705. // Allocate memory to the arrays
  706. std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
  707. D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
  708. std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
  709. for (int l = 0; l < L; l++) {
  710. Q[l].resize(nmax_ + 1);
  711. Ha[l].resize(nmax_);
  712. Hb[l].resize(nmax_);
  713. }
  714. an_.resize(nmax_);
  715. bn_.resize(nmax_);
  716. PsiZeta_.resize(nmax_ + 1);
  717. std::vector<std::complex<double> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
  718. //*************************************************//
  719. // Calculate D1 and D3 for z1 in the first layer //
  720. //*************************************************//
  721. if (fl == pl) { // PEC layer
  722. for (int n = 0; n <= nmax_; n++) {
  723. D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
  724. D3_mlxl[n] = std::complex<double>(0.0, 1.0);
  725. }
  726. } else { // Regular layer
  727. z1 = x[fl]* m[fl];
  728. // Calculate D1 and D3
  729. calcD1D3(z1, D1_mlxl, D3_mlxl);
  730. }
  731. //******************************************************************//
  732. // Calculate Ha and Hb in the first layer - equations (7a) and (8a) //
  733. //******************************************************************//
  734. for (int n = 0; n < nmax_; n++) {
  735. Ha[fl][n] = D1_mlxl[n + 1];
  736. Hb[fl][n] = D1_mlxl[n + 1];
  737. }
  738. //*****************************************************//
  739. // Iteration from the second layer to the last one (L) //
  740. //*****************************************************//
  741. std::complex<double> Temp, Num, Denom;
  742. std::complex<double> G1, G2;
  743. for (int l = fl + 1; l < L; l++) {
  744. //************************************************************//
  745. //Calculate D1 and D3 for z1 and z2 in the layers fl + 1..L //
  746. //************************************************************//
  747. z1 = x[l]*m[l];
  748. z2 = x[l - 1]*m[l];
  749. //Calculate D1 and D3 for z1
  750. calcD1D3(z1, D1_mlxl, D3_mlxl);
  751. //Calculate D1 and D3 for z2
  752. calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
  753. //*************************************************//
  754. //Calculate Q, Ha and Hb in the layers fl + 1..L //
  755. //*************************************************//
  756. // Upward recurrence for Q - equations (19a) and (19b)
  757. Num = std::exp(-2.0*(z1.imag() - z2.imag()))
  758. *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
  759. Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
  760. Q[l][0] = Num/Denom;
  761. for (int n = 1; n <= nmax_; n++) {
  762. Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
  763. Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
  764. Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
  765. }
  766. // Upward recurrence for Ha and Hb - equations (7b), (8b) and (12) - (15)
  767. for (int n = 1; n <= nmax_; n++) {
  768. //Ha
  769. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  770. G1 = -D1_mlxlM1[n];
  771. G2 = -D3_mlxlM1[n];
  772. } else {
  773. G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
  774. G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
  775. } // end of if PEC
  776. Temp = Q[l][n]*G1;
  777. Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
  778. Denom = G2 - Temp;
  779. Ha[l][n - 1] = Num/Denom;
  780. //Hb
  781. if ((l - 1) == pl) { // The layer below the current one is a PEC layer
  782. G1 = Hb[l - 1][n - 1];
  783. G2 = Hb[l - 1][n - 1];
  784. } else {
  785. G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
  786. G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
  787. } // end of if PEC
  788. Temp = Q[l][n]*G1;
  789. Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
  790. Denom = (G2- Temp);
  791. Hb[l][n - 1] = (Num/ Denom);
  792. } // end of for Ha and Hb terms
  793. } // end of for layers iteration
  794. //**************************************//
  795. //Calculate Psi and Zeta for XL //
  796. //**************************************//
  797. // Calculate PsiXL and ZetaXL
  798. calcPsiZeta(x[L - 1], PsiXL, ZetaXL);
  799. //*********************************************************************//
  800. // Finally, we calculate the scattering coefficients (an and bn) and //
  801. // the angular functions (Pi and Tau). Note that for these arrays the //
  802. // first layer is 0 (zero), in future versions all arrays will follow //
  803. // this convention to save memory. (13 Nov, 2014) //
  804. //*********************************************************************//
  805. for (int n = 0; n < nmax_; n++) {
  806. //********************************************************************//
  807. //Expressions for calculating an and bn coefficients are not valid if //
  808. //there is only one PEC layer (ie, for a simple PEC sphere). //
  809. //********************************************************************//
  810. if (pl < (L - 1)) {
  811. 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]);
  812. 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]);
  813. } else {
  814. 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]);
  815. bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
  816. }
  817. } // end of for an and bn terms
  818. isScaCoeffsCalc_ = true;
  819. } // end of MultiLayerMie::calcScattCoeffs()
  820. //**********************************************************************************//
  821. // This function calculates the actual scattering parameters and amplitudes //
  822. // //
  823. // Input parameters: //
  824. // L: Number of layers //
  825. // pl: Index of PEC layer. If there is none just send -1 //
  826. // x: Array containing the size parameters of the layers [0..L-1] //
  827. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  828. // nTheta: Number of scattering angles //
  829. // Theta: Array containing all the scattering angles where the scattering //
  830. // amplitudes will be calculated //
  831. // nmax_: Maximum number of multipolar expansion terms to be used for the //
  832. // calculations. Only use it if you know what you are doing, otherwise //
  833. // set this parameter to -1 and the function will calculate it //
  834. // //
  835. // Output parameters: //
  836. // Qext: Efficiency factor for extinction //
  837. // Qsca: Efficiency factor for scattering //
  838. // Qabs: Efficiency factor for absorption (Qabs = Qext - Qsca) //
  839. // Qbk: Efficiency factor for backscattering //
  840. // Qpr: Efficiency factor for the radiation pressure //
  841. // g: Asymmetry factor (g = (Qext-Qpr)/Qsca) //
  842. // Albedo: Single scattering albedo (Albedo = Qsca/Qext) //
  843. // S1, S2: Complex scattering amplitudes //
  844. // //
  845. // Return value: //
  846. // Number of multipolar expansion terms used for the calculations //
  847. //**********************************************************************************//
  848. void MultiLayerMie::RunMieCalculation() {
  849. if (size_param_.size() != refractive_index_.size())
  850. throw std::invalid_argument("Each size parameter should have only one index!");
  851. if (size_param_.size() == 0)
  852. throw std::invalid_argument("Initialize model first!");
  853. const std::vector<double>& x = size_param_;
  854. isExpCoeffsCalc_ = false;
  855. isScaCoeffsCalc_ = false;
  856. isMieCalculated_ = false;
  857. // Calculate scattering coefficients
  858. calcScattCoeffs();
  859. // Initialize the scattering parameters
  860. Qext_ = 0.0;
  861. Qsca_ = 0.0;
  862. Qabs_ = 0.0;
  863. Qbk_ = 0.0;
  864. Qpr_ = 0.0;
  865. asymmetry_factor_ = 0.0;
  866. albedo_ = 0.0;
  867. // Initialize the scattering amplitudes
  868. std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
  869. S1_.swap(tmp1);
  870. S2_ = S1_;
  871. std::vector<double> Pi(nmax_), Tau(nmax_);
  872. std::complex<double> Qbktmp(0.0, 0.0);
  873. std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
  874. // By using downward recurrence we avoid loss of precision due to float rounding errors
  875. // See: https://docs.oracle.com/cd/E19957-01/806-3568/ncg_goldberg.html
  876. // http://en.wikipedia.org/wiki/Loss_of_significance
  877. for (int i = nmax_ - 2; i >= 0; i--) {
  878. const int n = i + 1;
  879. // Equation (27)
  880. Qext_ += (n + n + 1.0)*(an_[i].real() + bn_[i].real());
  881. // Equation (28)
  882. Qsca_ += (n + n + 1.0)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
  883. + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
  884. // Equation (29)
  885. Qpr_ += ((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
  886. + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
  887. // Equation (33)
  888. Qbktmp += (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
  889. // Calculate the scattering amplitudes (S1 and S2) //
  890. // Precalculate cos(theta) - gives about 5% speed up.
  891. std::vector<double> costheta(theta_.size(), 0.0);
  892. for (int t = 0; t < theta_.size(); t++) {
  893. costheta[t] = std::cos(theta_[t]);
  894. }
  895. // Equations (25a) - (25b) //
  896. for (unsigned int t = 0; t < theta_.size(); t++) {
  897. calcPiTau(costheta[t], Pi, Tau);
  898. S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
  899. S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
  900. }
  901. }
  902. double x2 = pow2(x.back());
  903. Qext_ = 2.0*(Qext_)/x2; // Equation (27)
  904. Qsca_ = 2.0*(Qsca_)/x2; // Equation (28)
  905. Qpr_ = Qext_ - 4.0*(Qpr_)/x2; // Equation (29)
  906. Qabs_ = Qext_ - Qsca_; // Equation (30)
  907. albedo_ = Qsca_/Qext_; // Equation (31)
  908. asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_; // Equation (32)
  909. Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2; // Equation (33)
  910. isMieCalculated_ = true;
  911. }
  912. //**********************************************************************************//
  913. // This function calculates the expansion coefficients inside the particle, //
  914. // required to calculate the near-field parameters. //
  915. // //
  916. // Input parameters: //
  917. // L: Number of layers //
  918. // pl: Index of PEC layer. If there is none just send -1 //
  919. // x: Array containing the size parameters of the layers [0..L-1] //
  920. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  921. // nmax: Maximum number of multipolar expansion terms to be used for the //
  922. // calculations. Only use it if you know what you are doing, otherwise //
  923. // set this parameter to -1 and the function will calculate it. //
  924. // //
  925. // Output parameters: //
  926. // aln, bln, cln, dln: Complex scattering amplitudes inside the particle //
  927. // //
  928. // Return value: //
  929. // Number of multipolar expansion terms used for the calculations //
  930. //**********************************************************************************//
  931. void MultiLayerMie::calcExpanCoeffs() {
  932. if (!isScaCoeffsCalc_)
  933. throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
  934. isExpCoeffsCalc_ = false;
  935. std::complex<double> c_one(1.0, 0.0), c_zero(0.0, 0.0);
  936. const int L = refractive_index_.size();
  937. aln_.resize(L + 1);
  938. bln_.resize(L + 1);
  939. cln_.resize(L + 1);
  940. dln_.resize(L + 1);
  941. for (int l = 0; l <= L; l++) {
  942. aln_[l].resize(nmax_);
  943. bln_[l].resize(nmax_);
  944. cln_[l].resize(nmax_);
  945. dln_[l].resize(nmax_);
  946. }
  947. // Yang, paragraph under eq. A3
  948. // a^(L + 1)_n = a_n, d^(L + 1) = 1 ...
  949. for (int n = 0; n < nmax_; n++) {
  950. aln_[L][n] = an_[n];
  951. bln_[L][n] = bn_[n];
  952. cln_[L][n] = c_one;
  953. dln_[L][n] = c_one;
  954. }
  955. std::vector<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
  956. std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
  957. std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
  958. auto& m = refractive_index_;
  959. std::vector< std::complex<double> > m1(L);
  960. for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
  961. m1[L - 1] = std::complex<double> (1.0, 0.0);
  962. std::complex<double> z, z1;
  963. for (int l = L - 1; l >= 0; l--) {
  964. if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!!
  965. for (int n = 0; n < nmax_; n++) {
  966. // aln
  967. aln_[l][n] = c_zero;
  968. // bln
  969. bln_[l][n] = c_zero;
  970. // cln
  971. cln_[l][n] = c_zero;
  972. // dln
  973. dln_[l][n] = c_zero;
  974. }
  975. } else { // Regular material, just do the calculation
  976. z = size_param_[l]*m[l];
  977. z1 = size_param_[l]*m1[l];
  978. calcD1D3(z, D1z, D3z);
  979. calcD1D3(z1, D1z1, D3z1);
  980. calcPsiZeta(z, Psiz, Zetaz);
  981. calcPsiZeta(z1, Psiz1, Zetaz1);
  982. for (int n = 0; n < nmax_; n++) {
  983. int n1 = n + 1;
  984. denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]);
  985. denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]);
  986. T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
  987. T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l];
  988. T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l];
  989. T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
  990. // aln
  991. aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
  992. // bln
  993. bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
  994. // cln
  995. cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
  996. // dln
  997. dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
  998. } // end of all n
  999. } // end PEC condition
  1000. } // end of all l
  1001. // Check the result and change aln_[0][n] and aln_[0][n] for exact zero
  1002. for (int n = 0; n < nmax_; ++n) {
  1003. if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
  1004. else {
  1005. //throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
  1006. printf("Warning: Potentially unstable calculation of aln (aln[0][%i] = %g, %gi)\n", n, aln_[0][n].real(), aln_[0][n].imag());
  1007. aln_[0][n] = 0.0;
  1008. }
  1009. if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
  1010. else {
  1011. //throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
  1012. printf("Warning: Potentially unstable calculation of bln (bln[0][%i] = %g, %gi) pl=%d\n", n, bln_[0][n].real(), bln_[0][n].imag(), PEC_layer_position_);
  1013. bln_[0][n] = 0.0;
  1014. }
  1015. }
  1016. isExpCoeffsCalc_ = true;
  1017. } // end of void MultiLayerMie::calcExpanCoeffs()
  1018. //**********************************************************************************//
  1019. // This function calculates the electric (E) and magnetic (H) fields inside and //
  1020. // around the particle. //
  1021. // //
  1022. // Input parameters (coordinates of the point): //
  1023. // Rho: Radial distance //
  1024. // Phi: Azimuthal angle //
  1025. // Theta: Polar angle //
  1026. // //
  1027. // Output parameters: //
  1028. // E, H: Complex electric and magnetic fields //
  1029. //**********************************************************************************//
  1030. void MultiLayerMie::calcField(const double Rho, const double Theta, const double Phi,
  1031. std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
  1032. std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
  1033. std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i}; // Vector containing precomputed integer powers of i to avoid computation
  1034. std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
  1035. std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
  1036. std::vector<std::complex<double> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
  1037. std::vector<double> Pi(nmax_), Tau(nmax_);
  1038. int l = 0; // Layer number
  1039. std::complex<double> ml;
  1040. // Initialize E and H
  1041. for (int i = 0; i < 3; i++) {
  1042. E[i] = c_zero;
  1043. H[i] = c_zero;
  1044. }
  1045. if (Rho > size_param_.back()) {
  1046. l = size_param_.size();
  1047. ml = c_one;
  1048. } else {
  1049. for (int i = size_param_.size() - 1; i >= 0 ; i--) {
  1050. if (Rho <= size_param_[i]) {
  1051. l = i;
  1052. }
  1053. }
  1054. ml = refractive_index_[l];
  1055. }
  1056. // Calculate logarithmic derivative of the Ricatti-Bessel functions
  1057. calcD1D3(Rho*ml, D1n, D3n);
  1058. // Calculate Ricatti-Bessel functions
  1059. calcPsiZeta(Rho*ml, Psi, Zeta);
  1060. // Calculate angular functions Pi and Tau
  1061. calcPiTau(std::cos(Theta), Pi, Tau);
  1062. for (int n = nmax_ - 2; n >= 0; n--) {
  1063. int n1 = n + 1;
  1064. double rn = static_cast<double>(n1);
  1065. // using BH 4.12 and 4.50
  1066. calcSpherHarm(Rho*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
  1067. calcSpherHarm(Rho*ml, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
  1068. // Total field in the lth layer: eqs. (1) and (2) in Yang, Appl. Opt., 42 (2003) 1710-1720
  1069. std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
  1070. for (int i = 0; i < 3; i++) {
  1071. // electric field E [V m - 1] = EF*E0
  1072. E[i] += En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
  1073. + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
  1074. H[i] += En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
  1075. + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
  1076. }
  1077. } // end of for all n
  1078. // magnetic field
  1079. std::complex<double> hffact = ml/(cc_*mu_);
  1080. for (int i = 0; i < 3; i++) {
  1081. H[i] = hffact*H[i];
  1082. }
  1083. } // end of MultiLayerMie::calcField(...)
  1084. //**********************************************************************************//
  1085. // This function calculates complex electric and magnetic field in the surroundings //
  1086. // and inside the particle. //
  1087. // //
  1088. // Input parameters: //
  1089. // L: Number of layers //
  1090. // pl: Index of PEC layer. If there is none just send 0 (zero) //
  1091. // x: Array containing the size parameters of the layers [0..L-1] //
  1092. // m: Array containing the relative refractive indexes of the layers [0..L-1] //
  1093. // nmax: Maximum number of multipolar expansion terms to be used for the //
  1094. // calculations. Only use it if you know what you are doing, otherwise //
  1095. // set this parameter to 0 (zero) and the function will calculate it. //
  1096. // ncoord: Number of coordinate points //
  1097. // Coords: Array containing all coordinates where the complex electric and //
  1098. // magnetic fields will be calculated //
  1099. // //
  1100. // Output parameters: //
  1101. // E, H: Complex electric and magnetic field at the provided coordinates //
  1102. // //
  1103. // Return value: //
  1104. // Number of multipolar expansion terms used for the calculations //
  1105. //**********************************************************************************//
  1106. void MultiLayerMie::RunFieldCalculation() {
  1107. double Rho, Theta, Phi;
  1108. // Calculate scattering coefficients an_ and bn_
  1109. calcScattCoeffs();
  1110. // Calculate expansion coefficients aln_, bln_, cln_, and dln_
  1111. calcExpanCoeffs();
  1112. long total_points = coords_[0].size();
  1113. E_.resize(total_points);
  1114. H_.resize(total_points);
  1115. for (auto& f : E_) f.resize(3);
  1116. for (auto& f : H_) f.resize(3);
  1117. for (int point = 0; point < total_points; point++) {
  1118. const double& Xp = coords_[0][point];
  1119. const double& Yp = coords_[1][point];
  1120. const double& Zp = coords_[2][point];
  1121. // Convert to spherical coordinates
  1122. Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
  1123. // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems
  1124. Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
  1125. // If Xp=Yp=0 then Phi is undefined. Just set it to zero to avoid problems
  1126. if (Xp == 0.0)
  1127. Phi = (Yp != 0.0) ? std::asin(Yp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
  1128. else
  1129. Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
  1130. // Avoid convergence problems due to Rho too small
  1131. if (Rho < 1e-5) Rho = 1e-5;
  1132. //*******************************************************//
  1133. // external scattering field = incident + scattered //
  1134. // BH p.92 (4.37), 94 (4.45), 95 (4.50) //
  1135. // assume: medium is non-absorbing; refim = 0; Uabs = 0 //
  1136. //*******************************************************//
  1137. // This array contains the fields in spherical coordinates
  1138. std::vector<std::complex<double> > Es(3), Hs(3);
  1139. // Do the actual calculation of electric and magnetic field
  1140. calcField(Rho, Theta, Phi, Es, Hs);
  1141. { //Now, convert the fields back to cartesian coordinates
  1142. using std::sin;
  1143. using std::cos;
  1144. E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
  1145. E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
  1146. E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
  1147. H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
  1148. H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
  1149. H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
  1150. }
  1151. } // end of for all field coordinates
  1152. } // end of MultiLayerMie::RunFieldCalculation()
  1153. } // end of namespace nmie