37 #include <Faddeeva/Faddeeva.hh>
46 return Complex(0, 2) * (Constant::inv_sqrt_pi - z *
w);
88 Eigen::Ref<Eigen::VectorXcd> F,
89 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
93 const Numeric& magnetic_magnitude,
94 const Numeric& doppler_constant,
102 switch (lineshape_type) {
104 case LineShape::Type::SDVP:
114 case LineShape::Type::VP:
125 case LineShape::Type::DP:
135 case LineShape::Type::LP:
137 F, dF,
data,
f_grid, zeeman_df, magnetic_magnitude, line.
F0(), X);
141 switch (mirroring_type) {
142 case Absorption::MirroringType::None:
145 case Absorption::MirroringType::Lorentz: {
147 Eigen::VectorXcd Fm(F.size());
161 case Absorption::MirroringType::SameAsLineShape: {
163 Eigen::VectorXcd Fm(F.size());
165 switch (lineshape_type) {
166 case LineShape::Type::DP:
176 case LineShape::Type::LP:
186 case LineShape::Type::VP:
198 case LineShape::Type::SDVP:
220 case Absorption::NormalizationType::None:
222 case Absorption::NormalizationType::VVH:
225 case Absorption::NormalizationType::VVW:
235 Eigen::Ref<Eigen::VectorXcd> F,
236 Eigen::Ref<Eigen::MatrixXcd> dF,
238 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
240 const Numeric& magnetic_magnitude,
244 const Index& line_ind,
249 constexpr
Complex cpi(0, Constant::pi);
250 constexpr
Complex iz(0.0, 1.0);
253 auto nppd = derivatives_data_position.
nelem();
256 const Numeric F0 = F0_noshift + lso.
D0 + zeeman_df * magnetic_magnitude + lso.
DV;
259 auto z =
data.col(0);
264 F.noalias() = z.cwiseInverse();
267 dw.noalias() = -Constant::pi * F.array().square().matrix();
269 for (
auto iq = 0; iq < nppd; iq++) {
270 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
272 if (deriv == JacPropMatType::Temperature)
275 dF.col(iq).noalias() = -iz *
dw;
276 else if ((deriv == JacPropMatType::LineCenter or
279 dF.col(iq).noalias() = iz *
dw;
282 dF.col(iq).noalias() =
dw;
285 dF.col(iq).noalias() = iz *
dw;
287 dF.col(iq).noalias() = iz * zeeman_df *
dw;
288 else if (deriv == JacPropMatType::VMR) {
292 dF.col(iq).setZero();
299 Eigen::Ref<Eigen::VectorXcd> F,
300 Eigen::Ref<Eigen::MatrixXcd> dF,
302 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
304 const Numeric& magnetic_magnitude,
309 const Index& line_ind,
315 constexpr
Complex iz(0.0, 1.0);
318 auto nppd = derivatives_data_position.
nelem();
321 const Numeric F0 = F0_noshift + zeeman_df * magnetic_magnitude + lso.
D0 + lso.
DV;
323 const Numeric invGD = 1.0 / GD;
324 const Numeric dGD_dT = dGD_div_F0_dT *
F0 - GD_div_F0 * (dT.
D0 + dT.
DV);
327 const Numeric fac = Constant::inv_sqrt_pi * invGD;
330 auto z =
data.col(0);
337 F.noalias() =
fac * z.unaryExpr(&
w);
340 dw.noalias() = 2 * (
Complex(0,
fac * Constant::inv_sqrt_pi) -
341 z.cwiseProduct(F).array())
344 for (
auto iq = 0; iq < nppd; iq++) {
345 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
348 dF.col(iq).noalias() = invGD *
dw;
349 else if (deriv == JacPropMatType::Temperature)
350 dF.col(iq).noalias() =
352 F * dGD_dT * invGD -
dw.cwiseProduct(z) * dGD_dT * invGD;
353 else if ((deriv == JacPropMatType::LineCenter or
356 dF.col(iq).noalias() = -F /
F0 -
dw * invGD -
dw.cwiseProduct(z) /
F0;
359 dF.col(iq).noalias() = iz *
dw * invGD;
362 dF.col(iq).noalias() = -
dw * invGD;
364 dF.col(iq).noalias() =
dw * (-zeeman_df * invGD);
365 else if (deriv == JacPropMatType::VMR and
367 dF.col(iq).noalias() =
370 dF.col(iq).setZero();
376 Eigen::Ref<Eigen::VectorXcd> F,
377 Eigen::Ref<Eigen::MatrixXcd> dF,
379 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
381 const Numeric& magnetic_magnitude,
385 const Index& line_ind,
388 const Numeric& dGD_div_F0_dT) {
389 auto nppd = derivatives_data_position.
nelem();
392 const Numeric F0 = F0_noshift + zeeman_df * magnetic_magnitude;
393 const Numeric invGD = 1.0 / (GD_div_F0 *
F0);
396 auto x =
data.col(0);
397 auto mx2 =
data.col(1);
399 x.noalias() = (
f_grid.array() -
F0).matrix() * invGD;
400 mx2.noalias() = -
data.col(0).cwiseAbs2();
401 F.noalias() = invGD * Constant::inv_sqrt_pi * mx2.array().exp().matrix();
403 for (
auto iq = 0; iq < nppd; iq++) {
404 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
407 dF.col(iq).noalias() = -2 * invGD * F.cwiseProduct(
x);
408 else if (deriv == JacPropMatType::Temperature)
409 dF.col(iq).noalias() =
410 -dGD_div_F0_dT *
F0 * invGD * (2.0 * F.cwiseProduct(mx2) + F);
411 else if (deriv == JacPropMatType::LineCenter and
413 dF.col(iq).noalias() = -F.cwiseProduct(mx2) * (2 /
F0) +
414 F.cwiseProduct(
x) * 2 * (invGD - 1 /
F0);
415 else if (deriv == JacPropMatType::VMR)
416 dF.col(iq).setZero();
421 Eigen::Ref<Eigen::VectorXcd> F,
422 Eigen::Ref<Eigen::MatrixXcd> dF,
423 const Eigen::Ref<Eigen::VectorXcd> Fm,
424 const Eigen::Ref<Eigen::MatrixXcd> dFm,
426 const bool with_mirroring,
428 const Index& line_ind,
433 auto nppd = derivatives_data_position.
nelem();
438 if (with_mirroring) {
439 dF.noalias() += dFm *
conj(
LM);
441 for (
auto iq = 0; iq < nppd; iq++) {
442 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
444 if (deriv == JacPropMatType::Temperature) {
446 dF.col(iq).noalias() += F * c + Fm *
conj(c);
449 dF.col(iq).noalias() = F + Fm;
452 dF.col(iq).noalias() =
Complex(0, -1) * (F - Fm);
453 else if (deriv == JacPropMatType::VMR and
456 dF.col(iq).noalias() += F * c + Fm *
conj(c);
460 for (
auto iq = 0; iq < nppd; iq++) {
461 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
463 if (deriv == JacPropMatType::Temperature)
464 dF.col(iq).noalias() += F *
Complex(dT.
G, -dT.
Y);
467 dF.col(iq).noalias() = F;
470 dF.col(iq).noalias() =
Complex(0, -1) * F;
471 else if (deriv == JacPropMatType::VMR and
473 dF.col(iq).noalias() += F *
Complex(dVMR.
G, -dVMR.
Y);
478 if (with_mirroring) F.noalias() += Fm *
std::conj(
LM);
482 Eigen::Ref<Eigen::VectorXcd> F,
483 Eigen::Ref<Eigen::MatrixXcd> dF,
484 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
488 const Index& line_ind,
491 auto nf =
f_grid.size(), nppd = derivatives_data_position.
nelem();
495 (Constant::h) / (2.0 * Constant::k * T) /
496 std::sinh((Constant::h *
F0) / (2.0 * Constant::k * T)) * invF0;
503 (2.0 * std::tanh(
F0 * Constant::h / (2.0 * Constant::k * T)))) /
504 (Constant::k * T * T);
508 for (
auto iv = 0; iv < nf; iv++) {
512 for (
auto iq = 0; iq < nppd; iq++) {
513 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
516 if (deriv == JacPropMatType::Temperature)
517 dF(iv, iq) += dmafac_dT_div_fun * F[iv];
518 else if (deriv == JacPropMatType::LineCenter and
521 (-invF0 - Constant::h / (2.0 * Constant::k * T *
522 std::tanh(
F0 * Constant::h /
523 (2.0 * Constant::k * T)))) *
526 dF(iv, iq) += (2.0 /
f_grid[iv]) * F[iv];
532 Eigen::Ref<Eigen::VectorXcd> F,
533 Eigen::Ref<Eigen::MatrixXcd> dF,
535 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
539 const Index& line_ind,
542 auto nppd = derivatives_data_position.
nelem();
545 const Numeric kT = 2.0 * Constant::k * T;
546 const Numeric c1 = Constant::h / kT;
549 const Numeric tanh_f0part = std::tanh(c1 *
F0);
553 auto tanh_f =
data.col(0);
554 auto ftanh_f =
data.col(1);
556 tanh_f.noalias() = (c1 *
f_grid.array()).tanh().matrix();
557 ftanh_f.noalias() =
f_grid.cwiseProduct(tanh_f) / denom;
558 F.array() *= ftanh_f.array();
560 dF.array().colwise() *= ftanh_f.array();
561 for (
auto iq = 0; iq < nppd; iq++) {
562 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
564 if (deriv == JacPropMatType::Temperature)
565 dF.col(iq).noalias() +=
567 ((denom -
F0 / tanh_f0part) * F - denom * ftanh_f.cwiseProduct(F) +
568 f_grid.cwiseProduct(F).cwiseQuotient(tanh_f));
570 dF.col(iq).noalias() +=
572 c1 * (F.cwiseQuotient(tanh_f) - F.cwiseProduct(tanh_f));
573 else if (deriv == JacPropMatType::LineCenter and
575 dF.col(iq).noalias() +=
576 (-1.0 /
F0 + c1 * tanh_f0part - c1 / tanh_f0part) * F;
581 Eigen::Ref<Eigen::VectorXcd> F,
582 Eigen::Ref<Eigen::MatrixXcd> dF,
583 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
586 const Index& line_ind,
589 auto nf =
f_grid.size(), nppd = derivatives_data_position.
nelem();
593 const Numeric invF02 = invF0 * invF0;
595 for (
auto iv = 0; iv < nf; iv++) {
602 for (
auto iq = 0; iq < nppd; iq++) {
603 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
609 if (deriv == JacPropMatType::LineCenter and
611 dF(iv, iq) -= 2.0 * invF0 * F[iv];
613 dF(iv, iq) += 2.0 /
f_grid[iv] * F[iv];
629 return S0 * K1 * K2 * QT0 / QT;
633 Eigen::Ref<Eigen::VectorXcd> F,
634 Eigen::Ref<Eigen::MatrixXcd> dF,
635 Eigen::Ref<Eigen::VectorXcd>
N,
636 Eigen::Ref<Eigen::MatrixXcd> dN,
644 const Index& line_ind,
648 auto nppd = derivatives_data_position.
nelem();
655 const Numeric invQT = 1.0 / QT;
656 const Numeric S = line.
I0() * isotopic_ratio * QT0 * invQT * K1 * K2;
660 for (
auto iq = 0; iq < nppd; iq++) {
661 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
663 if (deriv == JacPropMatType::Temperature)
664 dF.col(iq).noalias() +=
668 else if (deriv == JacPropMatType::LineStrength and
670 dF.col(iq).noalias() = F / line.
I0();
671 else if (deriv == JacPropMatType::LineCenter and
673 dF.col(iq).noalias() +=
685 Eigen::Ref<Eigen::VectorXcd> F,
686 Eigen::Ref<Eigen::MatrixXcd> dF,
687 Eigen::Ref<Eigen::VectorXcd>
N,
688 Eigen::Ref<Eigen::MatrixXcd> dN,
700 const Index& line_ind,
704 auto nppd = derivatives_data_position.
nelem();
716 const Numeric invQT = 1.0 / QT;
717 const Numeric QT_ratio = QT0 * invQT;
719 const Numeric dS_dS0_abs = isotopic_ratio * QT_ratio * K1 * K2 * K3;
720 const Numeric S_abs = line.
I0() * dS_dS0_abs;
721 const Numeric dS_dS0_src = isotopic_ratio * QT_ratio * K1 * K2 * K4;
722 const Numeric S_src = line.
I0() * dS_dS0_src;
724 dN.noalias() = dF * (S_src - S_abs);
726 for (
auto iq = 0; iq < nppd; iq++) {
727 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
729 if (deriv == JacPropMatType::Temperature) {
745 dN.col(iq).noalias() += F * (dS_dT_src - dS_dT_abs);
746 dF.col(iq).noalias() += F * dS_dT_abs;
747 }
else if (deriv == JacPropMatType::LineStrength and
749 dF.col(iq).noalias() = F * dS_dS0_abs;
750 dN.col(iq).noalias() = F * (dS_dS0_src - dS_dS0_abs);
751 }
else if (deriv == JacPropMatType::LineCenter and
763 dN.col(iq).noalias() += F * (dS_dF0_src - dS_dF0_abs);
764 dF.col(iq).noalias() += F * dS_dF0_abs;
765 }
else if (deriv == JacPropMatType::NLTE) {
770 dN.col(iq).noalias() = -F * dS_dTl_abs;
771 dF.col(iq).noalias() = -dN.col(iq);
777 dN.col(iq).noalias() = F * (dS_dTu_src - dS_dTu_abs);
778 dF.col(iq).noalias() = F * dS_dTu_abs;
783 N.noalias() = F * (S_src - S_abs);
788 Eigen::Ref<Eigen::MatrixXcd> dF,
790 const Index& line_ind,
796 auto nppd = derivatives_data_position.
nelem();
798 for (
auto iq = 0; iq < nppd; iq++) {
800 derivatives_data[derivatives_data_position[iq]];
808 return std::sqrt(Constant::doppler_broadening_const_squared * T / mass);
819 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
824 const bool need_cutoff = (fmax > fmin);
828 Index i_f_max = nf - 1;
831 while (start_cutoff < nf and fmin >
f_grid[start_cutoff])
833 while (i_f_max >= start_cutoff and fmax <
f_grid[i_f_max])
837 nelem_cutoff = i_f_max - start_cutoff + 1;
845 Eigen::Ref<Eigen::VectorXcd> F,
846 Eigen::Ref<Eigen::MatrixXcd> dF,
847 Eigen::Ref<Eigen::VectorXcd>
N,
848 Eigen::Ref<Eigen::MatrixXcd> dN,
857 const Index& line_ind,
861 auto nppd = derivatives_data_position.
nelem();
865 constexpr
Numeric c1 = Constant::h / (4 * Constant::pi);
879 const Numeric k = c3 * (r1 *
x - r2) * (A21 / c2);
882 const Numeric e = c3 * r2 * A21;
885 const Numeric exp_T = std::exp(Constant::h / Constant::k *
F0 / T),
886 b = c2 / (exp_T - 1);
889 const Numeric ratio = e / b - k;
892 dN.noalias() = dF * ratio;
896 for (
auto iq = 0; iq < nppd; iq++) {
897 const auto& deriv = derivatives_data[derivatives_data_position[iq]];
899 if (deriv == JacPropMatType::Temperature)
900 dN.col(iq).noalias() +=
901 F * e * Constant::h *
F0 * exp_T / (c2 * Constant::k * T * T);
902 else if (deriv == JacPropMatType::LineCenter and
905 Constant::h * exp_T / (c2 * Constant::k * T) - 3.0 * b /
F0;
906 Numeric de_df0 = c1 * r2 * A21;
907 Numeric dk_df0 = c1 * (r1 *
x - r2) * (A21 / c2) - 3.0 * k /
F0;
909 dN.col(iq).noalias() += F * (e * done_over_b_df0 + de_df0 / b - dk_df0);
910 dF.col(iq).noalias() += F * dk_df0;
911 }
else if (deriv == JacPropMatType::NLTE) {
913 const Numeric dk_dr2 = -c3 * A21 / c2, de_dr2 = c3 * A21,
914 dratio_dr2 = de_dr2 / b - dk_dr2;
915 dN.col(iq).noalias() = F * dratio_dr2;
916 dF.col(iq).noalias() = F * dk_dr2;
918 const Numeric dk_dr1 = c3 *
x * A21 / c2;
919 dF.col(iq).noalias() = F * dk_dr1;
925 N.noalias() = F * ratio;
932 Eigen::Ref<Eigen::MatrixXcd> dF,
933 const Eigen::Ref<const Eigen::VectorXd>
f_grid,
935 const Numeric& magnetic_magnitude_si,
940 const Index& line_ind,
943 const Numeric& dGD_div_F0_dT_si,
946 using Constant::inv_sqrt_pi;
951 using Constant::sqrt_ln_2;
952 using Constant::sqrt_pi;
961 freq2kaycm(F0_noshift_si + zeeman_df_si * magnetic_magnitude_si);
962 const Numeric GamD = GD_div_F0_si * sg0 / sqrt_ln_2;
963 const auto lso =
si2cgs(lso_si);
964 const auto dT =
si2cgs(dT_si);
965 const auto dV =
si2cgs(dVMR_si);
968 const Numeric cte = sqrt_ln_2 / GamD;
972 const Complex c0(lso.G0, -lso.D0);
973 const Complex c2(lso.G2, -lso.D2);
974 const Complex c0t = (1 - lso.ETA) * (c0 - 1.5 * c2) + lso.FVC;
975 const Complex c2t = (1 - lso.ETA) * c2;
980 for (
auto iv = 0; iv <
f_grid.size(); iv++) {
986 Complex Z1, Z2, Zb, W1, W2, Wb;
993 Aterm = sqrt_pi * cte * W1;
996 Bterm = sqrt_pi * cte * ((1 -
pow2(Z1)) * W1 + Z1 * inv_sqrt_pi);
998 Bterm = cte * (sqrt_pi * W1 + 0.5 / Z1 - 0.75 /
pow3(Z1));
1004 Z2 = sqrtXY + sqrtY;
1009 Aterm = sqrt_pi * cte * (W1 - W2);
1010 Bterm = (-1 + sqrt_pi / (2 * sqrtY) * (1 -
pow2(Z1)) * W1 -
1011 sqrt_pi / (2 * sqrtY) * (1 -
pow2(Z2)) * W2) /
1021 if (
abs(sqrtX) <= 4e3) {
1025 Aterm = (2 * sqrt_pi / c2t) * (inv_sqrt_pi - sqrtX * Wb);
1028 (-1 + 2 * sqrt_pi * (1 - X - 2 * Y) * (inv_sqrt_pi - sqrtX * Wb) +
1029 2 * sqrt_pi * sqrtXY * W1);
1031 Aterm = (1 / c2t) * (1 / X - 1.5 /
pow2(X));
1032 Bterm = (1 / c2t) * (-1 + (1 - X - 2 * Y) * (1 / X - 1.5 /
pow2(X)) +
1033 2 * sqrt_pi * sqrtXY * W1);
1036 Z1 = sqrtXY - sqrtY;
1037 Z2 = Z1 + 2 * sqrtY;
1043 Aterm = sqrt_pi * cte * (W1 - W2);
1044 Bterm = (-1 + sqrt_pi / (2 * sqrtY) * (1 -
pow2(Z1)) * W1 -
1045 sqrt_pi / (2 * sqrtY) * (1 -
pow2(Z2)) * W2) /
1049 F(iv) = Aterm / (pi * (((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm +
1050 Bterm * c2 * lso.ETA + 1));
1052 for (
auto iq = 0; iq < derivatives_data_position.
nelem(); iq++) {
1054 derivatives_data[derivatives_data_position[iq]];
1057 Complex dc0 = 0, dc2 = 0, dc0t = 0, dc2t = 0;
1058 if (rt == JacPropMatType::Temperature) {
1059 dcte = (-(dGD_div_F0_dT_si * sg0 / Constant::sqrt_ln_2) / GamD) *
1063 dc0t = (-(c0 - 1.5 * c2) * dT.ETA + dT.FVC) +
1064 ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1065 dc2t = (-c2 * dT.ETA) + ((1 - lso.ETA) * dc2);
1066 }
else if (rt == JacPropMatType::VMR and
1070 dc0t = (-(c0 - 1.5 * c2) * dV.ETA + dV.FVC) +
1071 ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1072 dc2t = (-c2 * dV.ETA) + ((1 - lso.ETA) * dc2);
1077 dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1078 dc2t = ((1 - lso.ETA) * dc2);
1083 dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1084 dc2t = ((1 - lso.ETA) * dc2);
1089 dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1094 dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1100 dc0t = (-c0 + 1.5 * c2);
1104 const Complex dY = (-2 * dcte / cte - 2 * dc2t / c2t) * Y;
1109 else if (rt == JacPropMatType::LineCenter and
1120 if (
abs(c2t) == 0) {
1124 else if (rt == JacPropMatType::LineCenter and
1133 const Complex dW1 = iz * dZ1 *
dw(Z1, W1);
1135 dAterm = sqrt_pi * (W1 * dcte + cte * dW1);
1139 -(sqrt_pi * ((
pow2(Z1) - 1) * dW1 + 2 * W1 * Z1 * dZ1) - dZ1) *
1141 (sqrt_pi * (
pow2(Z1) - 1) * W1 - Z1) * dcte;
1144 ((sqrt_pi * W1 *
pow3(Z1) + 0.5 *
pow2(Z1) - 0.75) * Z1 * dcte +
1145 (sqrt_pi *
pow4(Z1) * dW1 - 0.5 *
pow2(Z1) * dZ1 + 2.25 * dZ1) *
1148 }
else if (
abs(X) <= 3e-8 *
abs(Y)) {
1152 else if (rt == JacPropMatType::LineCenter and
1161 const Complex dZ2 = dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1164 const Complex dW1 = iz * dZ1 *
dw(Z1, W1);
1165 const Complex dW2 = iz * dZ2 *
dw(Z2, W2);
1167 dAterm = sqrt_pi * ((W1 - W2) * dcte + (dW1 - dW2) * cte);
1170 (((
pow2(Z1) - 1) * W1 - (
pow2(Z2) - 1) * W2) * dY +
1172 (-(
pow2(Z1) - 1) * dW1 + (
pow2(Z2) - 1) * dW2 -
1173 2 * W1 * Z1 * dZ1 + 2 * W2 * Z2 * dZ2) *
1177 (sqrt_pi * (
pow2(Z1) - 1) * W1 -
1178 sqrt_pi * (
pow2(Z2) - 1) * W2 + 2 * sqrtY) *
1180 (4 * Y * sqrtY *
pow2(c2t));
1181 }
else if (
abs(Y) <= 1e-15 *
abs(X)) {
1182 const Complex dZ1 = (dX + dY) / (2 * sqrtXY);
1183 const Complex dW1 = iz * dZ1 *
dw(Z1, W1);
1184 if (
abs(sqrtX) <= 4e3) {
1185 const Complex dZb = dX / (2 * sqrtX);
1186 const Complex dWb = iz * dZb *
dw(Zb, Wb);
1188 dAterm = (-sqrt_pi * (Wb * dX + 2 * X * dWb) * c2t +
1189 2 * (sqrt_pi * Wb * sqrtX - 1) * sqrtX * dc2t) /
1190 (sqrtX *
pow2(c2t));
1194 (2 * (sqrt_pi * Wb * sqrtX - 1) * (X + 2 * Y - 1) +
1195 2 * sqrt_pi * sqrtXY * W1 - 1) *
1198 ((sqrt_pi * Wb * sqrtX - 1) * (dX + 2 * dY) +
1199 sqrt_pi * sqrtXY * dW1) *
1201 sqrt_pi * (Wb * dX + 2 * X * dWb) * sqrtXY * (X + 2 * Y - 1) +
1202 sqrt_pi * (dX + dY) * W1 * sqrtX) *
1204 (sqrtXY * sqrtX *
pow2(c2t));
1206 dAterm = ((-X + 3.0) * c2t * dX - (X - 1.5) * X * dc2t) /
1209 dBterm = (((-2 * sqrt_pi * sqrtXY * W1 + 1) *
pow2(X) +
1210 (X - 1.5) * (X + 2 * Y - 1)) *
1212 ((X - 3.0) * sqrtXY * (X + 2 * Y - 1) * dX -
1213 (X - 1.5) * sqrtXY * (dX + 2 * dY) * X +
1214 2 * sqrt_pi * (X + Y) *
pow3(X) * dW1 +
1215 sqrt_pi * (dX + dY) * W1 *
pow3(X)) *
1220 const Complex dZ1 = -dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1222 const Complex dZ2 = dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1225 const Complex dW1 = iz * dZ1 *
dw(Z1, W1);
1226 const Complex dW2 = iz * dZ2 *
dw(Z2, W2);
1228 dAterm = sqrt_pi * ((W1 - W2) * dcte + (dW1 - dW2) * cte);
1231 (((
pow2(Z1) - 1) * W1 - (
pow2(Z2) - 1) * W2) * dY +
1233 (-(
pow2(Z1) - 1) * dW1 + (
pow2(Z2) - 1) * dW2 -
1234 2 * W1 * Z1 * dZ1 + 2 * W2 * Z2 * dZ2) *
1238 (sqrt_pi * (
pow2(Z1) - 1) * W1 -
1239 sqrt_pi * (
pow2(Z2) - 1) * W2 + 2 * sqrtY) *
1241 (4 * Y * sqrtY *
pow2(c2t));
1245 if (rt == JacPropMatType::Temperature) {
1248 }
else if (rt == JacPropMatType::VMR) {
1259 ((((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm + Bterm * c2 * lso.ETA +
1262 (((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * dAterm +
1263 ((c0 - 1.5 * c2) * dETA + (dc0 - 1.5 * dc2) * lso.ETA - dFVC) *
1265 Bterm * c2 * dETA + Bterm * lso.ETA * dc2 + c2 * lso.ETA * dBterm) *
1267 (pi *
pow2(((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm +
1268 Bterm * c2 * lso.ETA + 1));
1274 for (
auto iq = 0; iq < derivatives_data_position.
nelem(); iq++) {
1276 derivatives_data[derivatives_data_position[iq]];
1309 const bool no_negatives,
1313 const Index nj = derivatives_data_active.
nelem();
1324 Eigen::Matrix<Numeric, 1, 1> fc;
1325 auto& Fc = scratch.
Fc;
1326 auto& Nc = scratch.
Nc;
1327 auto& dFc = scratch.
dFc;
1328 auto& dNc = scratch.
dNc;
1329 auto& datac = scratch.
datac;
1352 if (band.
Cutoff() == Absorption::CutoffType::LineByLineOffset and i>0) {
1365 const auto f = f_full.middleRows(
start,
nelem);
1371 const auto dXdT = do_temperature ?
1375 const auto dXdVMR = do_vmr.test ?
1379 const Index nz = zeeman ?
1382 for (
Index iz=0; iz<nz; iz++) {
1392 case LineShape::Type::DP:
1393 set_doppler(F, dF,
data, f, dfdH, H, band.
F0(i), DC, band, i, derivatives_data, derivatives_data_active, dDCdT);
1394 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1395 set_doppler(Fc, dFc, datac, fc, dfdH, H, band.
F0(i), DC, band, i, derivatives_data, derivatives_data_active, dDCdT);
1398 case LineShape::Type::SDVP:
1399 set_htp(F, dF, f, dfdH, H, band.
F0(i), DC, X, band, i, derivatives_data, derivatives_data_active, dDCdT, dXdT, dXdVMR);
1400 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1401 set_htp(Fc, dFc, fc, dfdH, H, band.
F0(i), DC, X, band, i, derivatives_data, derivatives_data_active, dDCdT, dXdT, dXdVMR);
1403 case LineShape::Type::LP:
1404 set_lorentz(F, dF,
data, f, dfdH, H, band.
F0(i), X, band, i, derivatives_data, derivatives_data_active, dXdT, dXdVMR);
1405 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1406 set_lorentz(Fc, dFc, datac, fc, dfdH, H, band.
F0(i), X, band, i, derivatives_data, derivatives_data_active, dXdT, dXdVMR);
1408 case LineShape::Type::VP:
1409 set_voigt(F, dF,
data, f, dfdH, H, band.
F0(i), DC, X, band, i, derivatives_data, derivatives_data_active, dDCdT, dXdT, dXdVMR);
1410 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1411 set_voigt(Fc, dFc, datac, fc, dfdH, H, band.
F0(i), DC, X, band, i, derivatives_data, derivatives_data_active, dDCdT, dXdT, dXdVMR);
1416 if (band.
Cutoff() not_eq Absorption::CutoffType::None) {
1418 for (
Index ij = 0; ij < nj; ij++) {
1419 dF.col(ij).array() -= dFc[ij];
1424 const bool with_mirroring =
1425 band.
Mirroring() not_eq Absorption::MirroringType::None and
1428 case Absorption::MirroringType::None:
1431 case Absorption::MirroringType::Lorentz:
1432 set_lorentz(
N, dN,
data, f, -dfdH, H, -band.
F0(i),
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1433 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1434 set_lorentz(Nc, dNc, datac, fc, -dfdH, H, -band.
F0(i),
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1436 case Absorption::MirroringType::SameAsLineShape:
1438 case LineShape::Type::DP:
1439 set_doppler(
N, dN,
data, f, -dfdH, H, -band.
F0(i), -DC, band, i, derivatives_data, derivatives_data_active, -dDCdT);
1440 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1441 set_doppler(Nc, dNc, datac, fc, -dfdH, H, -band.
F0(i), -DC, band, i, derivatives_data, derivatives_data_active, -dDCdT);
1443 case LineShape::Type::LP:
1444 set_lorentz(
N, dN,
data, f, -dfdH, H, -band.
F0(i),
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1445 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1446 set_lorentz(Nc, dNc, datac, fc, -dfdH, H, -band.
F0(i),
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1448 case LineShape::Type::VP:
1449 set_voigt(
N, dN,
data, f, -dfdH, H, -band.
F0(i), -DC,
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, -dDCdT, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1450 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1451 set_voigt(Nc, dNc, datac, fc, -dfdH, H, -band.
F0(i), -DC,
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, -dDCdT, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1454 case LineShape::Type::SDVP:
1456 set_htp(
N, dN, f, -dfdH, H, -band.
F0(i), -DC,
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, -dDCdT, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1457 if (band.
Cutoff() not_eq Absorption::CutoffType::None)
1458 set_htp(Nc, dNc, fc, -dfdH, H, -band.
F0(i), -DC,
LineShape::mirroredOutput(X), band, i, derivatives_data, derivatives_data_active, -dDCdT, do_temperature ?
LineShape::mirroredOutput(dXdT) : empty_output, do_vmr.test ?
LineShape::mirroredOutput(dXdVMR) : empty_output);
1465 if (band.
Cutoff() not_eq Absorption::CutoffType::None and with_mirroring) {
1467 for (
Index ij = 0; ij < nj; ij++) {
1468 dN.col(ij).array() -= dNc[ij];
1474 apply_linemixing_scaling_and_mirroring(F, dF,
N, dN, X, with_mirroring, band, i, derivatives_data, derivatives_data_active, dXdT, dXdVMR);
1482 case Absorption::NormalizationType::None:
1484 case Absorption::NormalizationType::VVH:
1487 case Absorption::NormalizationType::VVW:
1488 apply_VVW_scaling(F, dF, f, band.
F0(i), band, i, derivatives_data, derivatives_data_active);
1498 case Absorption::PopulationType::ByHITRANRosenkranzRelmat:
1499 case Absorption::PopulationType::ByLTE:
1500 apply_linestrength_scaling_by_lte(F, dF,
N, dN, band.
Line(i), T, band.
T0(), isot_ratio, QT, QT0, band, i, derivatives_data, derivatives_data_active, dQTdT);
1502 case Absorption::PopulationType::ByNLTEVibrationalTemperatures: {
1504 apply_linestrength_scaling_by_vibrational_nlte(F, dF,
N, dN, band.
Line(i), T, band.
T0(), nlte_data.T_upp, nlte_data.T_low, nlte_data.E_upp, nlte_data.E_low, isot_ratio, QT, QT0, band, i, derivatives_data, derivatives_data_active, dQTdT);
1506 case Absorption::PopulationType::ByNLTEPopulationDistribution: {
1508 apply_linestrength_from_nlte_level_distributions(F, dF,
N, dN, nlte_data.r_low, nlte_data.r_upp, band.
g_low(i), band.
g_upp(i), band.
A(i), band.
F0(i), T, band, i, derivatives_data, derivatives_data_active);
1530 auto reset_zeroes = (sum.
F.array().real() < 0);
1532 sum.
N = reset_zeroes.select(
Complex(0, 0), sum.
N);
1533 for (
Index ij=0; ij<nj; ij++)
1534 sum.
dF.col(ij) = reset_zeroes.select(
Complex(0, 0), sum.
dF.col(ij));
1535 for (
Index ij=0; ij<nj; ij++)
1536 sum.
dN.col(ij) = reset_zeroes.select(
Complex(0, 0), sum.
dN.col(ij));
1537 sum.
F = reset_zeroes.select(
Complex(0, 0), sum.
F);