nmie-impl.hpp 58 KB

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