nmie.cc 67 KB

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