123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317 |
- #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 ScattCoeffs(const unsigned int L, const int pl, std::vector<double>& x, std::vector<std::complex<double> >& m, const int nmax, std::vector<std::complex<double> >& an, std::vector<std::complex<double> >& bn) {
- if (x.size() != L || m.size() != L)
- throw std::invalid_argument("Declared number of layers do not fit x and m!");
- try {
- MultiLayerMie ml_mie;
- ml_mie.SetAllLayersSize(x);
- ml_mie.SetAllLayersIndex(m);
- ml_mie.SetPECLayer(pl);
- ml_mie.SetMaxTerms(nmax);
- ml_mie.calcScattCoeffs();
- an = ml_mie.GetAn();
- bn = ml_mie.GetBn();
- return ml_mie.GetMaxTerms();
- } 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, 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 ml_mie;
- ml_mie.SetAllLayersSize(x);
- ml_mie.SetAllLayersIndex(m);
- ml_mie.SetAngles(Theta);
- ml_mie.SetPECLayer(pl);
- ml_mie.SetMaxTerms(nmax);
- ml_mie.RunMieCalculation();
- *Qext = ml_mie.GetQext();
- *Qsca = ml_mie.GetQsca();
- *Qabs = ml_mie.GetQabs();
- *Qbk = ml_mie.GetQbk();
- *Qpr = ml_mie.GetQpr();
- *g = ml_mie.GetAsymmetryFactor();
- *Albedo = ml_mie.GetAlbedo();
- S1 = ml_mie.GetS1();
- S2 = ml_mie.GetS2();
- return ml_mie.GetMaxTerms();
- } 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 ml_mie;
- ml_mie.SetPECLayer(pl);
- ml_mie.SetAllLayersSize(x);
- ml_mie.SetAllLayersIndex(m);
- ml_mie.SetFieldCoords({Xp_vec, Yp_vec, Zp_vec});
- ml_mie.RunFieldCalculation();
- E = ml_mie.GetFieldE();
- H = ml_mie.GetFieldH();
- return ml_mie.GetMaxTerms();
- } 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::SetAllLayersSize(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::SetAllLayersIndex(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;
- }
-
-
-
- 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]) > 1.0e15) {
- 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 std::complex<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::calcScattCoeffs() {
- 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;
-
- calcScattCoeffs();
- 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::calcExpanCoeffs() {
- 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--) {
- if (l <= PEC_layer_position_) {
- for (int n = 0; n < nmax_; n++) {
-
- aln_[l][n] = c_zero;
-
- bln_[l][n] = c_zero;
-
- cln_[l][n] = c_zero;
-
- dln_[l][n] = c_zero;
- }
- } else {
- 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_[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) {
- if (std::abs(aln_[0][n]) < 1e-10) aln_[0][n] = 0.0;
- else {
-
- printf("Warning: Potentially unstable calculation of aln (aln[0][%i] = %g, %gi)\n", n, aln_[0][n].real(), aln_[0][n].imag());
- aln_[0][n] = 0.0;
- }
- if (std::abs(bln_[0][n]) < 1e-10) bln_[0][n] = 0.0;
- else {
-
- printf("Warning: Potentially unstable calculation of bln (bln[0][%i] = %g, %gi)\n", n, bln_[0][n].real(), bln_[0][n].imag());
- bln_[0][n] = 0.0;
- }
- }
- isExpCoeffsCalc_ = true;
- }
-
-
-
-
-
-
-
-
-
-
-
-
- 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_);
- 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*ml, D1n, D3n);
-
- calcPsiZeta(Rho*ml, 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*ml, Theta, Phi, Psi[n1], D1n[n1], Pi[n], Tau[n], rn, M1o1n, M1e1n, N1o1n, N1e1n);
- calcSpherHarm(Rho*ml, 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]);
- }
- }
-
- std::complex<double> hffact = ml/(cc_*mu_);
- for (int i = 0; i < 3; i++) {
- H[i] = hffact*H[i];
- }
- }
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- void MultiLayerMie::RunFieldCalculation() {
- double Rho, Theta, Phi;
-
- calcScattCoeffs();
-
- calcExpanCoeffs();
- 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];
- }
- }
- }
- }
|