12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562 |
- #include "bessel.h"
- #include "nmie.h"
- #include <array>
- #include <algorithm>
- #include <cstdio>
- #include <cstdlib>
- #include <stdexcept>
- #include <vector>
- namespace nmie {
-
- template<class T> inline T pow2(const T value) {return value*value;}
- int round(double x) {
- return x >= 0 ? (int)(x + 0.5):(int)(x - 0.5);
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- 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) {
- if (x.size() != L || m.size() != L)
- throw std::invalid_argument("Declared number of layers do not fit x and m!");
- if (Theta.size() != nTheta)
- throw std::invalid_argument("Declared number of sample for Theta is not correct!");
- try {
- MultiLayerMie multi_layer_mie;
- multi_layer_mie.SetLayersSize(x);
- multi_layer_mie.SetLayersIndex(m);
- multi_layer_mie.SetAngles(Theta);
- multi_layer_mie.SetPECLayer(pl);
- multi_layer_mie.SetMaxTerms(nmax);
- multi_layer_mie.RunMieCalculation();
- *Qext = multi_layer_mie.GetQext();
- *Qsca = multi_layer_mie.GetQsca();
- *Qabs = multi_layer_mie.GetQabs();
- *Qbk = multi_layer_mie.GetQbk();
- *Qpr = multi_layer_mie.GetQpr();
- *g = multi_layer_mie.GetAsymmetryFactor();
- *Albedo = multi_layer_mie.GetAlbedo();
- S1 = multi_layer_mie.GetS1();
- S2 = multi_layer_mie.GetS2();
- } catch(const std::invalid_argument& ia) {
-
- std::cerr << "Invalid argument: " << ia.what() << std::endl;
- throw std::invalid_argument(ia);
- return -1;
- }
- return 0;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- 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) {
- return nmie::nMie(L, -1, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- 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) {
- return nmie::nMie(L, pl, x, m, nTheta, Theta, -1, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- 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) {
- return nmie::nMie(L, -1, x, m, nTheta, Theta, nmax, Qext, Qsca, Qabs, Qbk, Qpr, g, Albedo, S1, S2);
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- 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) {
- if (x.size() != L || m.size() != L)
- throw std::invalid_argument("Declared number of layers do not fit x and m!");
- if (Xp_vec.size() != ncoord || Yp_vec.size() != ncoord || Zp_vec.size() != ncoord
- || E.size() != ncoord || H.size() != ncoord)
- throw std::invalid_argument("Declared number of coords do not fit Xp, Yp, Zp, E, or H!");
- for (auto f:E)
- if (f.size() != 3)
- throw std::invalid_argument("Field E is not 3D!");
- for (auto f:H)
- if (f.size() != 3)
- throw std::invalid_argument("Field H is not 3D!");
- try {
- MultiLayerMie multi_layer_mie;
-
- multi_layer_mie.SetLayersSize(x);
- multi_layer_mie.SetLayersIndex(m);
- multi_layer_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
- multi_layer_mie.RunFieldCalculation();
- E = multi_layer_mie.GetFieldE();
- H = multi_layer_mie.GetFieldH();
- } catch(const std::invalid_argument& ia) {
-
- std::cerr << "Invalid argument: " << ia.what() << std::endl;
- throw std::invalid_argument(ia);
- return - 1;
- }
- return 0;
- }
-
-
-
- double MultiLayerMie::GetQext() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return Qext_;
- }
-
-
-
- double MultiLayerMie::GetQabs() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return Qabs_;
- }
-
-
-
- double MultiLayerMie::GetQsca() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return Qsca_;
- }
-
-
-
- double MultiLayerMie::GetQbk() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return Qbk_;
- }
-
-
-
- double MultiLayerMie::GetQpr() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return Qpr_;
- }
-
-
-
- double MultiLayerMie::GetAsymmetryFactor() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return asymmetry_factor_;
- }
-
-
-
- double MultiLayerMie::GetAlbedo() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return albedo_;
- }
-
-
-
- std::vector<std::complex<double> > MultiLayerMie::GetS1() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return S1_;
- }
-
-
-
- std::vector<std::complex<double> > MultiLayerMie::GetS2() {
- if (!isMieCalculated_)
- throw std::invalid_argument("You should run calculations before result request!");
- return S2_;
- }
-
-
-
- void MultiLayerMie::SetAngles(const std::vector<double>& angles) {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- theta_ = angles;
- }
-
-
-
- void MultiLayerMie::SetLayersSize(const std::vector<double>& layer_size) {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- size_param_.clear();
- double prev_layer_size = 0.0;
- for (auto curr_layer_size : layer_size) {
- if (curr_layer_size <= 0.0)
- throw std::invalid_argument("Size parameter should be positive!");
- if (prev_layer_size > curr_layer_size)
- throw std::invalid_argument
- ("Size parameter for next layer should be larger than the previous one!");
- prev_layer_size = curr_layer_size;
- size_param_.push_back(curr_layer_size);
- }
- }
-
-
-
- void MultiLayerMie::SetLayersIndex(const std::vector< std::complex<double> >& index) {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- refractive_index_ = index;
- }
-
-
-
- void MultiLayerMie::SetFieldCoords(const std::vector< std::vector<double> >& coords) {
- if (coords.size() != 3)
- throw std::invalid_argument("Error! Wrong dimension of field monitor points!");
- if (coords[0].size() != coords[1].size() || coords[0].size() != coords[2].size())
- throw std::invalid_argument("Error! Missing coordinates for field monitor points!");
- coords_ = coords;
- }
-
-
-
- void MultiLayerMie::SetPECLayer(int layer_position) {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- if (layer_position < 0 && layer_position != -1)
- throw std::invalid_argument("Error! Layers are numbered from 0!");
- PEC_layer_position_ = layer_position;
- }
-
-
-
- void MultiLayerMie::SetMaxTerms(int nmax) {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- nmax_preset_ = nmax;
- }
-
-
-
- double MultiLayerMie::GetSizeParameter() {
- if (size_param_.size() > 0)
- return size_param_.back();
- else
- return 0;
- }
-
-
-
- void MultiLayerMie::ClearLayers() {
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
- size_param_.clear();
- refractive_index_.clear();
- }
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcNstop() {
- const double& xL = size_param_.back();
- if (xL <= 8) {
- nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 1);
- } else if (xL <= 4200) {
- nmax_ = round(xL + 4.05*pow(xL, 1.0/3.0) + 2);
- } else {
- nmax_ = round(xL + 4.0*pow(xL, 1.0/3.0) + 2);
- }
- }
-
-
-
- void MultiLayerMie::calcNmax(unsigned int first_layer) {
- int ri, riM1;
- const std::vector<double>& x = size_param_;
- const std::vector<std::complex<double> >& m = refractive_index_;
- calcNstop();
- for (unsigned int i = first_layer; i < x.size(); i++) {
- if (static_cast<int>(i) > PEC_layer_position_)
- ri = round(std::abs(x[i]*m[i]));
- else
- ri = 0;
- nmax_ = std::max(nmax_, ri);
-
- if ((i > first_layer) && (static_cast<int>(i - 1) > PEC_layer_position_))
- riM1 = round(std::abs(x[i - 1]* m[i]));
- else
- riM1 = 0;
- nmax_ = std::max(nmax_, riM1);
- }
- nmax_ += 15;
- }
-
-
-
- std::complex<double> MultiLayerMie::calc_an(int n, double XL, std::complex<double> Ha, std::complex<double> mL,
- std::complex<double> PsiXL, std::complex<double> ZetaXL,
- std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
- std::complex<double> Num = (Ha/mL + n/XL)*PsiXL - PsiXLM1;
- std::complex<double> Denom = (Ha/mL + n/XL)*ZetaXL - ZetaXLM1;
- return Num/Denom;
- }
-
-
-
- std::complex<double> MultiLayerMie::calc_bn(int n, double XL, std::complex<double> Hb, std::complex<double> mL,
- std::complex<double> PsiXL, std::complex<double> ZetaXL,
- std::complex<double> PsiXLM1, std::complex<double> ZetaXLM1) {
- std::complex<double> Num = (mL*Hb + n/XL)*PsiXL - PsiXLM1;
- std::complex<double> Denom = (mL*Hb + n/XL)*ZetaXL - ZetaXLM1;
- return Num/Denom;
- }
-
-
-
-
- void MultiLayerMie::calc_an_bn_bulk(std::vector<std::complex<double> >& an,
- std::vector<std::complex<double> >& bn,
- double x, std::complex<double> m) {
-
- std::vector<std::complex<double> > PsiX(nmax_ + 1), ZetaX(nmax_ + 1);
- std::vector<std::complex<double> > PsiMX(nmax_ + 1), ZetaMX(nmax_ + 1);
-
- calcPsiZeta(x, PsiX, ZetaX);
- calcPsiZeta(m*x, PsiMX, ZetaMX);
- std::vector<std::complex<double> > D1X(nmax_ + 1), D3X(nmax_ + 1);
- std::vector<std::complex<double> > D1MX(nmax_ + 1), D3MX(nmax_ + 1);
-
- calcD1D3(x, D1X, D3X);
- calcD1D3(m*x, D1MX, D3MX);
- std::vector<std::complex<double> > dPsiX(nmax_ + 1), dZetaX(nmax_ + 1);
- std::vector<std::complex<double> > dPsiMX(nmax_ + 1);
- for (int i = 0; i < nmax_ + 1; ++i) {
- dPsiX[i] = D1X[i]*PsiX[i];
- dPsiMX[i] = D1MX[i]*PsiMX[i];
-
- }
- bessel::calcZeta(nmax_, x, ZetaX, dZetaX);
- an.resize(nmax_);
- bn.resize(nmax_);
- for (int i = 0; i < nmax_; i++) {
- int n = i+1;
- std::complex<double> Num = m*PsiMX[n]*dPsiX[n] - PsiX[n]*dPsiMX[n];
- std::complex<double> Denom = m*PsiMX[n]*dZetaX[n] - ZetaX[n]*dPsiMX[n];
- an[i] = Num/Denom;
- Num = PsiMX[n]*dPsiX[n] - m*PsiX[n]*dPsiMX[n];
- Denom = PsiMX[n]*dZetaX[n] - m*ZetaX[n]*dPsiMX[n];
- bn[i] = Num/Denom;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
- }
-
-
-
- std::complex<double> MultiLayerMie::calc_S1(int n, std::complex<double> an, std::complex<double> bn,
- double Pi, double Tau) {
- return double(n + n + 1)*(Pi*an + Tau*bn)/double(n*n + n);
- }
-
-
-
-
- std::complex<double> MultiLayerMie::calc_S2(int n, std::complex<double> an, std::complex<double> bn,
- double Pi, double Tau) {
- return calc_S1(n, an, bn, Tau, Pi);
- }
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcD1D3(const std::complex<double> z,
- std::vector<std::complex<double> >& D1,
- std::vector<std::complex<double> >& D3) {
-
- D1[nmax_] = std::complex<double>(0.0, 0.0);
- const std::complex<double> zinv = std::complex<double>(1.0, 0.0)/z;
- for (int n = nmax_; n > 0; n--) {
- D1[n - 1] = static_cast<double>(n)*zinv - 1.0/(D1[n] + static_cast<double>(n)*zinv);
- }
- if (std::abs(D1[0]) > 100000.0)
- throw std::invalid_argument("Unstable D1! Please, try to change input parameters!\n");
-
- PsiZeta_[0] = 0.5*(1.0 - std::complex<double>(std::cos(2.0*z.real()), std::sin(2.0*z.real()))
- *std::exp(-2.0*z.imag()));
- D3[0] = std::complex<double>(0.0, 1.0);
- for (int n = 1; n <= nmax_; n++) {
- PsiZeta_[n] = PsiZeta_[n - 1]*(static_cast<double>(n)*zinv - D1[n - 1])
- *(static_cast<double>(n)*zinv - D3[n - 1]);
- D3[n] = D1[n] + std::complex<double>(0.0, 1.0)/PsiZeta_[n];
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcPsiZeta(std::complex<double> z,
- std::vector<std::complex<double> >& Psi,
- std::vector<std::complex<double> >& Zeta) {
- std::complex<double> c_i(0.0, 1.0);
- std::vector<std::complex<double> > D1(nmax_ + 1), D3(nmax_ + 1);
-
- calcD1D3(z, D1, D3);
-
- Psi[0] = std::sin(z);
- Zeta[0] = std::sin(z) - c_i*std::cos(z);
- for (int n = 1; n <= nmax_; n++) {
- Psi[n] = Psi[n - 1]*(static_cast<double>(n)/z - D1[n - 1]);
- Zeta[n] = Zeta[n - 1]*(static_cast<double>(n)/z - D3[n - 1]);
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcPiTau(const double& costheta,
- std::vector<double>& Pi, std::vector<double>& Tau) {
- int i;
-
-
-
-
- Pi[0] = 1.0;
- Tau[0] = costheta;
-
- if (nmax_ > 1) {
- Pi[1] = 3*costheta*Pi[0];
- Tau[1] = 2*costheta*Pi[1] - 3*Pi[0];
- for (i = 2; i < nmax_; i++) {
- Pi[i] = ((i + i + 1)*costheta*Pi[i - 1] - (i + 1)*Pi[i - 2])/i;
- Tau[i] = (i + 1)*costheta*Pi[i] - (i + 2)*Pi[i - 1];
- }
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcSpherHarm(const double Rho, const double Theta, const double Phi,
- const std::complex<double>& rn, const std::complex<double>& Dn,
- const double& Pi, const double& Tau, const double& n,
- std::vector<std::complex<double> >& Mo1n, std::vector<std::complex<double> >& Me1n,
- std::vector<std::complex<double> >& No1n, std::vector<std::complex<double> >& Ne1n) {
-
- std::complex<double> c_zero(0.0, 0.0);
- using std::sin;
- using std::cos;
- Mo1n[0] = c_zero;
- Mo1n[1] = cos(Phi)*Pi*rn/Rho;
- Mo1n[2] = -sin(Phi)*Tau*rn/Rho;
- Me1n[0] = c_zero;
- Me1n[1] = -sin(Phi)*Pi*rn/Rho;
- Me1n[2] = -cos(Phi)*Tau*rn/Rho;
- No1n[0] = sin(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
- No1n[1] = sin(Phi)*Tau*Dn*rn/Rho;
- No1n[2] = cos(Phi)*Pi*Dn*rn/Rho;
- Ne1n[0] = cos(Phi)*(n*n + n)*sin(Theta)*Pi*rn/Rho/Rho;
- Ne1n[1] = cos(Phi)*Tau*Dn*rn/Rho;
- Ne1n[2] = -sin(Phi)*Pi*Dn*rn/Rho;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::ScattCoeffs() {
- isScaCoeffsCalc_ = false;
- const std::vector<double>& x = size_param_;
- const std::vector<std::complex<double> >& m = refractive_index_;
- const int& pl = PEC_layer_position_;
- const int L = refractive_index_.size();
-
-
-
-
-
- int fl = (pl > 0) ? pl : 0;
- if (nmax_preset_ <= 0) calcNmax(fl);
- else nmax_ = nmax_preset_;
- std::complex<double> z1, z2;
-
-
-
-
-
-
-
-
- std::vector<std::complex<double> > D1_mlxl(nmax_ + 1), D1_mlxlM1(nmax_ + 1),
- D3_mlxl(nmax_ + 1), D3_mlxlM1(nmax_ + 1);
- std::vector<std::vector<std::complex<double> > > Q(L), Ha(L), Hb(L);
- for (int l = 0; l < L; l++) {
- Q[l].resize(nmax_ + 1);
- Ha[l].resize(nmax_);
- Hb[l].resize(nmax_);
- }
- an_.resize(nmax_);
- bn_.resize(nmax_);
- PsiZeta_.resize(nmax_ + 1);
- std::vector<std::complex<double> > PsiXL(nmax_ + 1), ZetaXL(nmax_ + 1);
-
-
-
- if (fl == pl) {
- for (int n = 0; n <= nmax_; n++) {
- D1_mlxl[n] = std::complex<double>(0.0, - 1.0);
- D3_mlxl[n] = std::complex<double>(0.0, 1.0);
- }
- } else {
- z1 = x[fl]* m[fl];
-
- calcD1D3(z1, D1_mlxl, D3_mlxl);
- }
-
-
-
- for (int n = 0; n < nmax_; n++) {
- Ha[fl][n] = D1_mlxl[n + 1];
- Hb[fl][n] = D1_mlxl[n + 1];
- }
-
-
-
- std::complex<double> Temp, Num, Denom;
- std::complex<double> G1, G2;
- for (int l = fl + 1; l < L; l++) {
-
-
-
- z1 = x[l]*m[l];
- z2 = x[l - 1]*m[l];
-
- calcD1D3(z1, D1_mlxl, D3_mlxl);
-
- calcD1D3(z2, D1_mlxlM1, D3_mlxlM1);
-
-
-
-
- Num = std::exp(-2.0*(z1.imag() - z2.imag()))
- *std::complex<double>(std::cos(-2.0*z2.real()) - std::exp(-2.0*z2.imag()), std::sin(-2.0*z2.real()));
- Denom = std::complex<double>(std::cos(-2.0*z1.real()) - std::exp(-2.0*z1.imag()), std::sin(-2.0*z1.real()));
- Q[l][0] = Num/Denom;
- for (int n = 1; n <= nmax_; n++) {
- Num = (z1*D1_mlxl[n] + double(n))*(double(n) - z1*D3_mlxl[n - 1]);
- Denom = (z2*D1_mlxlM1[n] + double(n))*(double(n) - z2*D3_mlxlM1[n - 1]);
- Q[l][n] = ((pow2(x[l - 1]/x[l])* Q[l][n - 1])*Num)/Denom;
- }
-
- for (int n = 1; n <= nmax_; n++) {
-
- if ((l - 1) == pl) {
- G1 = -D1_mlxlM1[n];
- G2 = -D3_mlxlM1[n];
- } else {
- G1 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D1_mlxlM1[n]);
- G2 = (m[l]*Ha[l - 1][n - 1]) - (m[l - 1]*D3_mlxlM1[n]);
- }
- Temp = Q[l][n]*G1;
- Num = (G2*D1_mlxl[n]) - (Temp*D3_mlxl[n]);
- Denom = G2 - Temp;
- Ha[l][n - 1] = Num/Denom;
-
- if ((l - 1) == pl) {
- G1 = Hb[l - 1][n - 1];
- G2 = Hb[l - 1][n - 1];
- } else {
- G1 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D1_mlxlM1[n]);
- G2 = (m[l - 1]*Hb[l - 1][n - 1]) - (m[l]*D3_mlxlM1[n]);
- }
- Temp = Q[l][n]*G1;
- Num = (G2*D1_mlxl[n]) - (Temp* D3_mlxl[n]);
- Denom = (G2- Temp);
- Hb[l][n - 1] = (Num/ Denom);
- }
- }
-
-
-
-
- calcPsiZeta(x[L - 1], PsiXL, ZetaXL);
-
-
-
-
-
-
- for (int n = 0; n < nmax_; n++) {
-
-
-
-
- if (pl < (L - 1)) {
- 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]);
- 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]);
- } else {
- 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]);
- bn_[n] = PsiXL[n + 1]/ZetaXL[n + 1];
- }
- }
- isScaCoeffsCalc_ = true;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::RunMieCalculation() {
- if (size_param_.size() != refractive_index_.size())
- throw std::invalid_argument("Each size parameter should have only one index!");
- if (size_param_.size() == 0)
- throw std::invalid_argument("Initialize model first!");
- const std::vector<double>& x = size_param_;
- isExpCoeffsCalc_ = false;
- isScaCoeffsCalc_ = false;
- isMieCalculated_ = false;
-
- ScattCoeffs();
- if (!isScaCoeffsCalc_)
- throw std::invalid_argument("Calculation of scattering coefficients failed!");
-
- Qext_ = 0.0;
- Qsca_ = 0.0;
- Qabs_ = 0.0;
- Qbk_ = 0.0;
- Qpr_ = 0.0;
- asymmetry_factor_ = 0.0;
- albedo_ = 0.0;
-
- std::vector<std::complex<double> > tmp1(theta_.size(),std::complex<double>(0.0, 0.0));
- S1_.swap(tmp1);
- S2_ = S1_;
- std::vector<double> Pi(nmax_), Tau(nmax_);
- std::complex<double> Qbktmp(0.0, 0.0);
- std::vector< std::complex<double> > Qbktmp_ch(nmax_ - 1, Qbktmp);
-
-
-
- for (int i = nmax_ - 2; i >= 0; i--) {
- const int n = i + 1;
-
- Qext_ += (n + n + 1.0)*(an_[i].real() + bn_[i].real());
-
- Qsca_ += (n + n + 1.0)*(an_[i].real()*an_[i].real() + an_[i].imag()*an_[i].imag()
- + bn_[i].real()*bn_[i].real() + bn_[i].imag()*bn_[i].imag());
-
- Qpr_ += ((n*(n + 2)/(n + 1))*((an_[i]*std::conj(an_[n]) + bn_[i]*std::conj(bn_[n])).real())
- + ((double)(n + n + 1)/(n*(n + 1)))*(an_[i]*std::conj(bn_[i])).real());
-
- Qbktmp += (double)(n + n + 1)*(1 - 2*(n % 2))*(an_[i]- bn_[i]);
-
-
- for (unsigned int t = 0; t < theta_.size(); t++) {
- calcPiTau(std::cos(theta_[t]), Pi, Tau);
- S1_[t] += calc_S1(n, an_[i], bn_[i], Pi[i], Tau[i]);
- S2_[t] += calc_S2(n, an_[i], bn_[i], Pi[i], Tau[i]);
- }
- }
- double x2 = pow2(x.back());
- Qext_ = 2.0*(Qext_)/x2;
- Qsca_ = 2.0*(Qsca_)/x2;
- Qpr_ = Qext_ - 4.0*(Qpr_)/x2;
- Qabs_ = Qext_ - Qsca_;
- albedo_ = Qsca_/Qext_;
- asymmetry_factor_ = (Qext_ - Qpr_)/Qsca_;
- Qbk_ = (Qbktmp.real()*Qbktmp.real() + Qbktmp.imag()*Qbktmp.imag())/x2;
- isMieCalculated_ = true;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::ExpanCoeffs() {
- if (!isScaCoeffsCalc_)
- throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
- isExpCoeffsCalc_ = false;
- std::complex<double> 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_);
- }
-
-
- 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<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
- std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
- std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
- auto& m = refractive_index_;
- std::vector< std::complex<double> > m1(L);
- for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
- m1[L - 1] = std::complex<double> (1.0, 0.0);
- std::complex<double> z, z1;
- for (int l = L - 1; l >= 0; l--) {
- 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 = (D1z1[n1]*dln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*aln_[l + 1][n]*Zetaz1[n1])*m[l]/m1[l];
- T4 = D1z1[n1]*cln_[l + 1][n]*Psiz1[n1] - D3z1[n1]*bln_[l + 1][n]*Zetaz1[n1];
-
- aln_[l][n] = (D1z[n1]*T1 + T3)/denomZeta;
-
- bln_[l][n] = (D1z[n1]*T2 + T4)/denomZeta;
-
- cln_[l][n] = (D3z[n1]*T2 + T4)/denomPsi;
-
- dln_[l][n] = (D3z[n1]*T1 + T3)/denomPsi;
-
- }
- }
-
- for (int n = 0; n < nmax_; ++n) {
- aln_[0][n] = 0.0;
- bln_[0][n] = 0.0;
-
-
-
-
-
-
- }
- isExpCoeffsCalc_ = true;
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::ExpanCoeffsV2() {
- if (!isScaCoeffsCalc_)
- throw std::invalid_argument("(ExpanCoeffs) You should calculate external coefficients first!");
- isExpCoeffsCalc_ = false;
- std::complex<double> 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_);
- }
-
-
- 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<std::complex<double> > D1z(nmax_ + 1), D1z1(nmax_ + 1), D3z(nmax_ + 1), D3z1(nmax_ + 1);
- std::vector<std::complex<double> > Psiz(nmax_ + 1), Psiz1(nmax_ + 1), Zetaz(nmax_ + 1), Zetaz1(nmax_ + 1);
- std::complex<double> denomZeta, denomPsi, T1, T2, T3, T4;
- std::vector<std::vector<std::complex<double> > > a(2);
- a[0].resize(3);
- a[1].resize(3);
- auto& m = refractive_index_;
- std::vector< std::complex<double> > m1(L);
- for (int l = 0; l < L - 1; l++) m1[l] = m[l + 1];
- m1[L - 1] = std::complex<double> (1.0, 0.0);
- std::complex<double> z, z1;
- for (int l = L - 1; l >= 0; l--) {
- 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;
- a[0][0] = m1[l]*D3z[n1]*Zetaz[n1];
- a[0][1] = -m1[l]*D1z[n1]*Psiz[n1];
- a[0][2] = aln_[l + 1][n]*m[l]*D3z1[n1]*Zetaz1[n1];
- a[0][2] -= dln_[l + 1][n]*m[l]*D1z1[n1]*Psiz1[n1];
- a[1][0] = Zetaz[n1];
- a[1][1] = -Psiz[n1];
- a[1][2] = aln_[l + 1][n]*Zetaz1[n1] - dln_[l + 1][n]*Psiz1[n1];
-
- aln_[l][n] = (a[0][2]*a[1][1] - a[0][1]*a[1][2])/(a[0][0]*a[1][1] - a[0][1]*a[1][0]);
-
- dln_[l][n] = (a[0][2]*a[1][0] - a[0][0]*a[1][2])/(a[0][1]*a[1][0] - a[0][0]*a[0][1]);
-
- a[0][0] = D3z[n1]*Zetaz[n1];
- a[0][1] = -D1z[n1]*Psiz[n1];
- a[0][2] = bln_[l + 1][n]*D3z1[n1]*Zetaz1[n1];
- a[0][2] -= cln_[l + 1][n]*D1z1[n1]*Psiz1[n1];
- a[1][0] = m1[l]*Zetaz[n1];
- a[1][1] = -m1[l]*Psiz[n1];
- a[1][2] = bln_[l + 1][n]*m[l]*Zetaz1[n1] - cln_[l + 1][n]*m[l]*Psiz1[n1];
-
- bln_[l][n] = (a[0][2]*a[1][1] - a[0][1]*a[1][2])/(a[0][0]*a[1][1] - a[0][1]*a[1][0]);
-
- cln_[l][n] = (a[0][2]*a[1][0] - a[0][0]*a[1][2])/(a[0][1]*a[1][0] - a[0][0]*a[0][1]);
-
- }
- }
-
- for (int n = 0; n < nmax_; ++n) {
-
-
- if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
- else throw std::invalid_argument("Unstable calculation of aln_[0][n]!");
- if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
- else throw std::invalid_argument("Unstable calculation of bln_[0][n]!");
- }
- isExpCoeffsCalc_ = true;
- }
-
-
-
-
-
- void MultiLayerMie::fieldExt(const double Rho, const double Theta, const double Phi,
- std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
- std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
- std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i};
- std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
- std::vector<std::complex<double> > Ei(3, c_zero), Hi(3, c_zero), Es(3, c_zero), Hs(3, c_zero);
- std::vector<std::complex<double> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
- std::vector<double> Pi(nmax_), Tau(nmax_);
-
- if (Rho < size_param_.back()) {
- for (int i = 0; i < 3; i++) {
- E[i] = c_zero;
- H[i] = c_zero;
- }
- return;
- }
-
- calcD1D3(Rho, D1n, D3n);
-
- calcPsiZeta(Rho, Psi, Zeta);
-
- calcPiTau(std::cos(Theta), Pi, Tau);
- for (int n = 0; n < nmax_; n++) {
- int n1 = n + 1;
- double rn = static_cast<double>(n1);
-
- calcSpherHarm(Rho, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
-
- std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
- for (int i = 0; i < 3; i++) {
- Es[i] = Es[i] + En*(c_i*an_[n]*N3e1n[i] - bn_[n]*M3o1n[i]);
- Hs[i] = Hs[i] + En*(c_i*bn_[n]*N3o1n[i] + an_[n]*M3e1n[i]);
- }
- }
-
-
- std::complex<double> eifac = std::exp(std::complex<double>(0.0, Rho*std::cos(Theta)));
- {
- using std::sin;
- using std::cos;
- Ei[0] = eifac*sin(Theta)*cos(Phi);
- Ei[1] = eifac*cos(Theta)*cos(Phi);
- Ei[2] = -eifac*sin(Phi);
- }
-
- double hffact = 1.0/(cc_*mu_);
- for (int i = 0; i < 3; i++) {
- Hs[i] = hffact*Hs[i];
- }
-
- std::complex<double> hffacta = hffact;
- std::complex<double> hifac = eifac*hffacta;
- {
- using std::sin;
- using std::cos;
- Hi[0] = hifac*sin(Theta)*sin(Phi);
- Hi[1] = hifac*cos(Theta)*sin(Phi);
- Hi[2] = hifac*cos(Phi);
- }
- for (int i = 0; i < 3; i++) {
-
- E[i] = Ei[i] + Es[i];
- H[i] = Hi[i] + Hs[i];
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::calcField(const double Rho, const double Theta, const double Phi,
- std::vector<std::complex<double> >& E, std::vector<std::complex<double> >& H) {
- std::complex<double> c_zero(0.0, 0.0), c_i(0.0, 1.0), c_one(1.0, 0.0);
- std::vector<std::complex<double> > ipow = {c_one, c_i, -c_one, -c_i};
- std::vector<std::complex<double> > M3o1n(3), M3e1n(3), N3o1n(3), N3e1n(3);
- std::vector<std::complex<double> > M1o1n(3), M1e1n(3), N1o1n(3), N1e1n(3);
- std::vector<std::complex<double> > Psi(nmax_ + 1), D1n(nmax_ + 1), Zeta(nmax_ + 1), D3n(nmax_ + 1);
- std::vector<double> Pi(nmax_), Tau(nmax_);
- std::vector<std::complex<double> > Ei(3), Hi(3);
- int l = 0;
- std::complex<double> ml;
-
- for (int i = 0; i < 3; i++) {
- E[i] = c_zero;
- H[i] = c_zero;
- }
- if (Rho > size_param_.back()) {
- l = size_param_.size();
- ml = c_one;
- } else {
- for (int i = size_param_.size() - 1; i >= 0 ; i--) {
- if (Rho <= size_param_[i]) {
- l = i;
-
- }
- }
- ml = refractive_index_[l];
- }
-
- calcD1D3(Rho, D1n, D3n);
-
- calcPsiZeta(Rho, Psi, Zeta);
-
- calcPiTau(std::cos(Theta), Pi, Tau);
- for (int n = nmax_ - 2; n >= 0; n--) {
- int n1 = n + 1;
- double rn = static_cast<double>(n1);
-
- calcSpherHarm(Rho, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
- calcSpherHarm(Rho, Theta, Phi, Zeta[n1], D3n[n1], Pi[n], Tau[n], rn, M3o1n, M3e1n, N3o1n, N3e1n);
-
-
-
-
- std::complex<double> En = ipow[n1 % 4]*(rn + rn + 1.0)/(rn*rn + rn);
- for (int i = 0; i < 3; i++) {
-
- E[i] += En*(cln_[l][n]*M1o1n[i] - c_i*dln_[l][n]*N1e1n[i]
- + c_i*aln_[l][n]*N3e1n[i] - bln_[l][n]*M3o1n[i]);
- H[i] += En*(-dln_[l][n]*M1e1n[i] - c_i*cln_[l][n]*N1o1n[i]
- + c_i*bln_[l][n]*N3o1n[i] + aln_[l][n]*M3e1n[i]);
- Ei[i] += En*(M1o1n[i] - c_i*N1e1n[i]);
- Hi[i] += En*(-M1e1n[i] - c_i*N1o1n[i]);
- }
- }
-
-
- double hffact = 1.0/(cc_*mu_);
- for (int i = 0; i < 3; i++) {
- H[i] = hffact*H[i];
- Hi[i] *= hffact;
-
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::RunFieldCalculation() {
- double Rho, Theta, Phi;
-
- ScattCoeffs();
-
-
-
-
-
-
-
- ExpanCoeffs();
-
-
-
-
-
- long total_points = coords_[0].size();
- E_.resize(total_points);
- H_.resize(total_points);
- for (auto& f : E_) f.resize(3);
- for (auto& f : H_) f.resize(3);
- for (int point = 0; point < total_points; point++) {
- const double& Xp = coords_[0][point];
- const double& Yp = coords_[1][point];
- const double& Zp = coords_[2][point];
-
- Rho = std::sqrt(pow2(Xp) + pow2(Yp) + pow2(Zp));
-
- Theta = (Rho > 0.0) ? std::acos(Zp/Rho) : 0.0;
-
- if (Xp == 0.0)
- Phi = (Yp != 0.0) ? std::asin(Yp/std::sqrt(pow2(Xp) + pow2(Yp))) : 0.0;
- else
- Phi = std::acos(Xp/std::sqrt(pow2(Xp) + pow2(Yp)));
-
- if (Rho < 1e-5) Rho = 1e-5;
-
-
-
-
-
-
-
- std::vector<std::complex<double> > Es(3), Hs(3);
-
-
-
-
- calcField(Rho, Theta, Phi, Es, Hs);
-
- {
- using std::sin;
- using std::cos;
- E_[point][0] = sin(Theta)*cos(Phi)*Es[0] + cos(Theta)*cos(Phi)*Es[1] - sin(Phi)*Es[2];
- E_[point][1] = sin(Theta)*sin(Phi)*Es[0] + cos(Theta)*sin(Phi)*Es[1] + cos(Phi)*Es[2];
- E_[point][2] = cos(Theta)*Es[0] - sin(Theta)*Es[1];
- H_[point][0] = sin(Theta)*cos(Phi)*Hs[0] + cos(Theta)*cos(Phi)*Hs[1] - sin(Phi)*Hs[2];
- H_[point][1] = sin(Theta)*sin(Phi)*Hs[0] + cos(Theta)*sin(Phi)*Hs[1] + cos(Phi)*Hs[2];
- H_[point][2] = cos(Theta)*Hs[0] - sin(Theta)*Hs[1];
- }
- }
- }
- }
|