#ifndef SRC_NMIE_NEARFIELD_HPP_ #define SRC_NMIE_NEARFIELD_HPP_ //**********************************************************************************// // Copyright (C) 2009-2021 Ovidio Pena // // Copyright (C) 2013-2021 Konstantin Ladutenko // // // // This file is part of scattnlay // // // // This program is free software: you can redistribute it and/or modify // // it under the terms of the GNU General Public License as published by // // the Free Software Foundation, either version 3 of the License, or // // (at your option) any later version. // // // // This program is distributed in the hope that it will be useful, // // but WITHOUT ANY WARRANTY; without even the implied warranty of // // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the // // GNU General Public License for more details. // // // // The only additional remark is that we expect that all publications // // describing work using this software, or all commercial products // // using it, cite at least one of the following references: // // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by // // a multilayered sphere," Computer Physics Communications, // // vol. 180, Nov. 2009, pp. 2348-2354. // // [2] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie // // calculation of electromagnetic near-field for a multilayered // // sphere," Computer Physics Communications, vol. 214, May 2017, // // pp. 225-230. // // // // You should have received a copy of the GNU General Public License // // along with this program. If not, see . // //**********************************************************************************// //**********************************************************************************// // This class implements the algorithm for a multilayered sphere described by: // // [1] W. Yang, "Improved recursive algorithm for light scattering by a // // multilayered sphere,” Applied Optics, vol. 42, Mar. 2003, pp. 1710-1720. // // // // You can find the description of all the used equations in: // // [2] O. Pena and U. Pal, "Scattering of electromagnetic radiation by // // a multilayered sphere," Computer Physics Communications, // // vol. 180, Nov. 2009, pp. 2348-2354. // // [3] K. Ladutenko, U. Pal, A. Rivera, and O. Pena-Rodriguez, "Mie // // calculation of electromagnetic near-field for a multilayered // // sphere," Computer Physics Communications, vol. 214, May 2017, // // pp. 225-230. // // // // Hereinafter all equations numbers refer to [2] // //**********************************************************************************// #include #include #include #include //#include "nmie.hpp" namespace nmie { //class implementation //**********************************************************************************// // This function calculates the expansion coefficients inside the particle, // // required to calculate the near-field parameters. // // // // Input parameters: // // L: Number of layers // // pl: Index of PEC layer. If there is none just send -1 // // x: Array containing the size parameters of the layers [0..L-1] // // m: Array containing the relative refractive indexes of the layers [0..L-1] // // nmax: Maximum number of multipolar expansion terms to be used for the // // calculations. Only use it if you know what you are doing, otherwise // // set this parameter to -1 and the function will calculate it. // // // // Output parameters: // // aln, bln, cln, dln: Complex scattering amplitudes inside the particle // // // // Return value: // // Number of multipolar expansion terms used for the calculations // //**********************************************************************************// template void MultiLayerMie::calcExpanCoeffs() { if (!isScaCoeffsCalc_) throw std::invalid_argument("(calcExpanCoeffs) You should calculate external coefficients first!"); isExpCoeffsCalc_ = false; std::complex c_one(1.0, 0.0), c_zero(0.0, 0.0); const int L = refractive_index_.size(); aln_.resize(L + 1); bln_.resize(L + 1); cln_.resize(L + 1); dln_.resize(L + 1); for (int l = 0; l <= L; l++) { aln_[l].resize(nmax_); bln_[l].resize(nmax_); cln_[l].resize(nmax_); dln_[l].resize(nmax_); } // Yang, paragraph under eq. A3 // a^(L + 1)_n = a_n, d^(L + 1) = 1 ... for (int n = 0; n < nmax_; n++) { aln_[L][n] = an_[n]; bln_[L][n] = bn_[n]; cln_[L][n] = c_one; dln_[L][n] = c_one; } std::vector > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1); std::vector > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1); std::complex denomZeta, denomPsi, T1, T2, T3, T4; auto &m = refractive_index_; std::vector< std::complex > m1(L); for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1]; m1[L - 1] = std::complex (1.0, 0.0); std::complex z, z1; for (int l = L - 1; l >= 0; l--) { if (l <= PEC_layer_position_) { // We are inside a PEC. All coefficients must be zero!!! for (int n = 0; n < nmax_; n++) { // aln aln_[l][n] = c_zero; // bln bln_[l][n] = c_zero; // cln cln_[l][n] = c_zero; // dln dln_[l][n] = c_zero; } } else { // Regular material, just do the calculation z = size_param_[l]*m[l]; z1 = size_param_[l]*m1[l]; calcD1D3(z, D1z, D3z); calcD1D3(z1, D1z1, D3z1); calcPsiZeta(z, Psiz, Zetaz); calcPsiZeta(z1, Psiz1, Zetaz1); for (int n = 0; n < nmax_; n++) { int n1 = n + 1; denomZeta = Zetaz[n1]*(D1z[n1] - D3z[n1]); denomPsi = Psiz[n1]*(D1z[n1] - D3z[n1]); T1 = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1]; T2 = (bln_[l + 1][n]*Zetaz1[n1] - cln_[l + 1][n]*Psiz1[n1])*m[l]/m1[l]; T3 = (dln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - aln_[l + 1][n]*D3z1[n1]*Zetaz1[n1])*m[l]/m1[l]; T4 = cln_[l + 1][n]*D1z1[n1]*Psiz1[n1] - bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1]; // aln aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta; // bln bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta; // cln cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi; // dln dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi; } // end of all n } // end PEC condition } // end of all l // Check the result and change aln_[0][n] and aln_[0][n] for exact zero for (int n = 0; n < nmax_; ++n) { if (cabs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0; else { //throw std::invalid_argument("Unstable calculation of aln_[0][n]!"); std::cout<< std::setprecision(100) << "Warning: Potentially unstable calculation of aln[0][" << n << "] = "<< aln_[0][n] < hffact = ml/static_cast(cc_*mu_); for (int i = 0; i < 3; i++) { H[i] = hffact*H[i]; } } // end of MultiLayerMie::calcFieldByComponents(...) //**********************************************************************************// // This function calculates complex electric and magnetic field in the surroundings // // and inside the particle. // // // // Input parameters: // // L: Number of layers // // pl: Index of PEC layer. If there is none just send 0 (zero) // // x: Array containing the size parameters of the layers [0..L-1] // // m: Array containing the relative refractive indexes of the layers [0..L-1] // // nmax: Maximum number of multipolar expansion terms to be used for the // // calculations. Only use it if you know what you are doing, otherwise // // set this parameter to 0 (zero) and the function will calculate it. // // ncoord: Number of coordinate points // // Coords: Array containing all coordinates where the complex electric and // // magnetic fields will be calculated // // mode_n: mode order. // // -1 - use all modes (all_) // // 1 - use dipole mode only // // 2 - use quadrupole mode only // // ... // // mode_type: only used when mode_n != -1 // // 0 - electric only // // 1 - magnetic only // // // // Output parameters: // // E, H: Complex electric and magnetic field at the provided coordinates // // // // Return value: // // Number of multipolar expansion terms used for the calculations // //**********************************************************************************// template void MultiLayerMie::RunFieldCalculation() { // Calculate scattering coefficients an_ and bn_ calcScattCoeffs(); // Calculate expansion coefficients aln_, bln_, cln_, and dln_ calcExpanCoeffs(); Es_.clear(); Hs_.clear(); coords_polar_.clear(); long total_points = coords_[0].size(); for (int point = 0; point < total_points; point++) { const FloatType &Xp = coords_[0][point]; const FloatType &Yp = coords_[1][point]; const FloatType &Zp = coords_[2][point]; // Convert to spherical coordinates auto Rho = nmm::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp)); // If Rho=0 then Theta is undefined. Just set it to zero to avoid problems auto Theta = (Rho > 0.0) ? nmm::acos(Zp/Rho) : 0.0; // std::atan2 should take care of any special cases, e.g. Xp=Yp=0, etc. auto Phi = nmm::atan2(Yp,Xp); coords_polar_.push_back({Rho, Theta, Phi}); // Avoid convergence problems due to Rho too small if (Rho < 1e-5) Rho = 1e-5; //*******************************************************// // external scattering field = incident + scattered // // BH p.92 (4.37), 94 (4.45), 95 (4.50) // // assume: medium is non-absorbing; refim = 0; Uabs = 0 // //*******************************************************// // This array contains the fields in spherical coordinates std::vector > Es(3), Hs(3); // Do the actual calculation of electric and magnetic field std::vector > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1); std::vector Pi(nmax_), Tau(nmax_); std::complex ml; GetIndexAtRadius(Rho, ml); // Calculate logarithmic derivative of the Ricatti-Bessel functions calcD1D3(Rho*ml, D1n, D3n); // Calculate Ricatti-Bessel functions calcPsiZeta(Rho*ml, Psi, Zeta); // Calculate angular functions Pi and Tau calcPiTau(nmm::cos(Theta), Pi, Tau); calcFieldByComponents(Rho, Theta, Phi, Psi, D1n, Zeta, D3n, Pi, Tau, Es, Hs); Es_.push_back(Es); Hs_.push_back(Hs); } // end of for all field coordinates convertFieldsFromSphericalToCartesian(); } // end of MultiLayerMie::RunFieldCalculation() template double eval_delta(const unsigned int steps, const double from_value, const double to_value) { auto delta = std::abs(from_value - to_value); if (steps < 2) return delta; delta /= static_cast(steps-1); // We have a limited double precision evaluation of special functions, typically it is 1e-10. if ( (2.*delta)/std::abs(from_value+to_value) < 1e-9) throw std::invalid_argument("Error! The step is too fine, not supported!"); return delta; } // ml - refractive index // l - Layer number template template void MultiLayerMie::GetIndexAtRadius(const evalType Rho, std::complex &ml, unsigned int &l) { l = 0; if (Rho > size_param_.back()) { l = size_param_.size(); ml = std::complex(1.0, 0.0); } else { for (int i = size_param_.size() - 1; i >= 0 ; i--) { if (Rho <= size_param_[i]) { l = i; } } ml = ConvertComplex(refractive_index_[l]); } } template template void MultiLayerMie::GetIndexAtRadius(const evalType Rho, std::complex &ml) { unsigned int l; GetIndexAtRadius(Rho, ml, l); } template void MultiLayerMie::calcMieSeriesNeededToConverge(const FloatType Rho) { auto required_near_field_nmax = calcNmax(Rho); SetMaxTerms(required_near_field_nmax); // Calculate scattering coefficients an_ and bn_ calcScattCoeffs(); // We might be limited with available machine precision available_maximal_nmax_ = nmax_; // Calculate expansion coefficients aln_, bln_, cln_, and dln_ calcExpanCoeffs(); } template void MultiLayerMie::calcRadialOnlyDependantFunctions(const double from_Rho, const double to_Rho, const bool isIgnoreAvailableNmax, std::vector > > &Psi, std::vector > > &D1n, std::vector > > &Zeta, std::vector > > &D3n) { auto radius_points = Psi.size(); std::vector > > PsiZeta(radius_points); double delta_Rho = eval_delta(radius_points, from_Rho, to_Rho); for (unsigned int j=0; j < radius_points; j++) { auto Rho = static_cast(from_Rho + j*delta_Rho); // if (Rho < 1e-5) Rho = 1e-5; // TODO do we need this?. int near_field_nmax = calcNmax(Rho); // Skip if not enough terms in Mie series (i.e. required near field nmax > available terms ) if (near_field_nmax > available_maximal_nmax_ && !isIgnoreAvailableNmax) continue; if (near_field_nmax > available_maximal_nmax_) near_field_nmax = available_maximal_nmax_; Psi[j].resize(near_field_nmax + 1); D1n[j].resize(near_field_nmax + 1); Zeta[j].resize(near_field_nmax + 1); D3n[j].resize(near_field_nmax + 1); PsiZeta[j].resize(near_field_nmax + 1); std::complex ml; GetIndexAtRadius(Rho, ml); auto z = Rho*ml; evalDownwardD1(z, D1n[j]); evalUpwardPsi(z, D1n[j], Psi[j]); evalUpwardD3 (z, D1n[j], D3n[j], PsiZeta[j]); for (unsigned int k = 0; k < Zeta[j].size(); k++) { Zeta[j][k] = PsiZeta[j][k]/Psi[j][k]; } } } // input parameters: // outer_arc_points: will be increased to the nearest power of 2. template void MultiLayerMie::RunFieldCalculationPolar(const int outer_arc_points, const int radius_points, const double from_Rho, const double to_Rho, const double from_Theta, const double to_Theta, const double from_Phi, const double to_Phi, const bool isIgnoreAvailableNmax) { if (from_Rho > to_Rho || from_Theta > to_Theta || from_Phi > to_Phi || outer_arc_points < 1 || radius_points < 1 || from_Rho < 0.) throw std::invalid_argument("Error! Invalid argument for RunFieldCalculationPolar() !"); auto nmax_old = nmax_; int theta_points = 0, phi_points = 0; if (to_Theta-from_Theta > to_Phi-from_Phi) { theta_points = outer_arc_points; phi_points = static_cast((to_Phi-from_Phi)/(to_Theta-from_Theta) * outer_arc_points); } else { phi_points = outer_arc_points; theta_points = static_cast((to_Theta-from_Theta)/(to_Phi-from_Phi) * outer_arc_points); } if (theta_points == 0) theta_points = 1; if (phi_points == 0) phi_points = 1; calcMieSeriesNeededToConverge(to_Rho); std::vector > Pi(theta_points), Tau(theta_points); calcPiTauAllTheta(from_Theta, to_Theta, Pi, Tau); std::vector > > Psi(radius_points), D1n(radius_points), Zeta(radius_points), D3n(radius_points), PsiZeta(radius_points); calcRadialOnlyDependantFunctions(from_Rho, to_Rho, isIgnoreAvailableNmax, Psi, D1n, Zeta, D3n); double delta_Rho = eval_delta(radius_points, from_Rho, to_Rho); double delta_Theta = eval_delta(theta_points, from_Theta, to_Theta); double delta_Phi = eval_delta(phi_points, from_Phi, to_Phi); Es_.clear(); Hs_.clear(); coords_polar_.clear(); for (int j=0; j < radius_points; j++) { auto Rho = from_Rho + j * delta_Rho; std::vector< std::complex > Psi_dp = ConvertComplexVector(Psi[j]); std::vector< std::complex > Zeta_dp = ConvertComplexVector(Zeta[j]); std::vector< std::complex > D1n_dp = ConvertComplexVector(D1n[j]); std::vector< std::complex > D3n_dp = ConvertComplexVector(D3n[j]); for (int i = 0; i < theta_points; i++) { auto Theta = from_Theta + i * delta_Theta; std::vector Pi_dp = ConvertVector(Pi[i]); std::vector Tau_dp = ConvertVector(Tau[i]); for (int k = 0; k < phi_points; k++) { auto Phi = from_Phi + k * delta_Phi; coords_polar_.push_back({Rho, Theta, Phi}); // This array contains the fields in spherical coordinates // std::vector > Es(3), Hs(3); std::vector > Es(3), Hs(3); calcFieldByComponents( Rho, Theta, Phi, Psi_dp, D1n_dp, Zeta_dp, D3n_dp, Pi_dp, Tau_dp, Es, Hs // Rho, Theta, Phi, // Psi[j], D1n[j], Zeta[j], D3n[j], // Pi[i], Tau[i], Es, Hs ); Es_.push_back(ConvertComplexVector(Es)); Hs_.push_back(ConvertComplexVector(Hs)); } } } convertFieldsFromSphericalToCartesian(); nmax_ = nmax_old; } // Python interface } // end of namespace nmie #endif // SRC_NMIE_NEARFIELD_HPP_