ARTS  2.4.0(git:4fb77825)
linefunctions.cc
Go to the documentation of this file.
1 /* Copyright (C) 2017
2  * Richard Larsson <ric.larsson@gmail.com>
3  *
4  * This program is free software; you can redistribute it and/or modify it
5  * under the terms of the GNU General Public License as published by the
6  * Free Software Foundation; either version 2, or (at your option) any
7  * later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program; if not, write to the Free Software
16  * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
17  * USA. */
18 
35 #include "linefunctions.h"
36 #include <Eigen/Core>
37 #include <Faddeeva/Faddeeva.hh>
38 #include "constants.h"
39 #include "linescaling.h"
40 
42 inline Complex w(Complex z) noexcept { return Faddeeva::w(z); }
43 
45 constexpr Complex dw(Complex z, Complex w) noexcept {
46  return Complex(0, 2) * (Constant::inv_sqrt_pi - z * w);
47 }
48 
50 constexpr Complex pCqSDHC_to_arts(Complex x) noexcept {
51  using Constant::c;
52  using Constant::pow2;
54 
55  return conj(hitran2arts_linestrength(x) / pow2(c));
56 }
57 
60  using Constant::c;
61  using Constant::pow4;
63 
65 }
66 
69  using Constant::c;
70  using Constant::pow4;
72 
73  return Complex(0, -1) *
75 }
76 
79  using Constant::c;
80  using Constant::pow4;
82 
84  pow4(c);
85 }
86 
88  Eigen::Ref<Eigen::VectorXcd> F,
89  const Eigen::Ref<const Eigen::VectorXd> f_grid,
90  const Absorption::SingleLine& line,
91  const Numeric& temperature,
92  const Numeric& zeeman_df,
93  const Numeric& magnetic_magnitude,
94  const Numeric& doppler_constant,
95  const LineShape::Output& X,
96  const LineShape::Type lineshape_type,
97  const Absorption::MirroringType mirroring_type,
98  const Absorption::NormalizationType norm_type)
99 {
100  Eigen::MatrixXcd dF(0, 0), data(F.size(), Linefunctions::ExpectedDataSize());
101 
102  switch (lineshape_type) {
104  case LineShape::Type::SDVP:
105  set_htp(F,
106  dF,
107  f_grid,
108  zeeman_df,
109  magnetic_magnitude,
110  line.F0(),
111  doppler_constant,
112  X);
113  break;
114  case LineShape::Type::VP:
115  set_voigt(F,
116  dF,
117  data,
118  f_grid,
119  zeeman_df,
120  magnetic_magnitude,
121  line.F0(),
122  doppler_constant,
123  X);
124  break;
125  case LineShape::Type::DP:
126  set_doppler(F,
127  dF,
128  data,
129  f_grid,
130  zeeman_df,
131  magnetic_magnitude,
132  line.F0(),
133  doppler_constant);
134  break;
135  case LineShape::Type::LP:
136  set_lorentz(
137  F, dF, data, f_grid, zeeman_df, magnetic_magnitude, line.F0(), X);
138  break;
139  }
140 
141  switch (mirroring_type) {
142  case Absorption::MirroringType::None:
144  break;
145  case Absorption::MirroringType::Lorentz: {
146  // Set the mirroring computational vectors and size them as needed
147  Eigen::VectorXcd Fm(F.size());
148 
149  set_lorentz(Fm,
150  dF,
151  data,
152  f_grid,
153  -zeeman_df,
154  magnetic_magnitude,
155  -line.F0(),
157 
158  // Apply mirroring; FIXME: Add conjugate?
159  F.noalias() += Fm;
160  } break;
161  case Absorption::MirroringType::SameAsLineShape: {
162  // Set the mirroring computational vectors and size them as needed
163  Eigen::VectorXcd Fm(F.size());
164 
165  switch (lineshape_type) {
166  case LineShape::Type::DP:
167  set_doppler(Fm,
168  dF,
169  data,
170  f_grid,
171  -zeeman_df,
172  magnetic_magnitude,
173  -line.F0(),
174  -doppler_constant);
175  break;
176  case LineShape::Type::LP:
177  set_lorentz(Fm,
178  dF,
179  data,
180  f_grid,
181  -zeeman_df,
182  magnetic_magnitude,
183  -line.F0(),
185  break;
186  case LineShape::Type::VP:
187  set_voigt(Fm,
188  dF,
189  data,
190  f_grid,
191  -zeeman_df,
192  magnetic_magnitude,
193  -line.F0(),
194  -doppler_constant,
196  break;
198  case LineShape::Type::SDVP:
199  // WARNING: This mirroring is not tested and it might require, e.g., FVC to be treated differently
200  set_htp(Fm,
201  dF,
202  f_grid,
203  -zeeman_df,
204  magnetic_magnitude,
205  -line.F0(),
206  -doppler_constant,
208  break;
209  }
210 
211  F.noalias() += Fm;
212  break;
213  }
214  }
215 
216  // Line normalization if necessary
217  // The user sets this by setting LSM LNT followed by a type
218  // that is internally interpreted to mean some kind of lineshape normalization
219  switch (norm_type) {
220  case Absorption::NormalizationType::None:
221  break;
222  case Absorption::NormalizationType::VVH:
223  apply_VVH_scaling(F, dF, data, f_grid, line.F0(), temperature);
224  break;
225  case Absorption::NormalizationType::VVW:
226  apply_VVW_scaling(F, dF, f_grid, line.F0());
227  break;
229  apply_rosenkranz_quadratic_scaling(F, dF, f_grid, line.F0(), temperature);
230  break;
231  }
232 }
233 
235  Eigen::Ref<Eigen::VectorXcd> F,
236  Eigen::Ref<Eigen::MatrixXcd> dF,
237  Eigen::Ref<Eigen::Matrix<Complex, Eigen::Dynamic, ExpectedDataSize()>> data,
238  const Eigen::Ref<const Eigen::VectorXd> f_grid,
239  const Numeric& zeeman_df,
240  const Numeric& magnetic_magnitude,
241  const Numeric& F0_noshift,
242  const LineShape::Output& lso,
243  const AbsorptionLines& band,
244  const Index& line_ind,
245  const ArrayOfRetrievalQuantity& derivatives_data,
246  const ArrayOfIndex& derivatives_data_position,
247  const LineShape::Output& dT,
248  const LineShape::Output& dVMR) {
249  constexpr Complex cpi(0, Constant::pi);
250  constexpr Complex iz(0.0, 1.0);
251 
252  // Size of the problem
253  auto nppd = derivatives_data_position.nelem();
254 
255  // The central frequency
256  const Numeric F0 = F0_noshift + lso.D0 + zeeman_df * magnetic_magnitude + lso.DV;
257 
258  // Naming req data blocks
259  auto z = data.col(0);
260  auto dw = data.col(1);
261 
262  z.noalias() =
263  (Constant::pi * Complex(lso.G0, F0) - cpi * f_grid.array()).matrix();
264  F.noalias() = z.cwiseInverse();
265 
266  if (nppd) {
267  dw.noalias() = -Constant::pi * F.array().square().matrix();
268 
269  for (auto iq = 0; iq < nppd; iq++) {
270  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
271 
272  if (deriv == JacPropMatType::Temperature)
273  dF.col(iq).noalias() = Complex(dT.G0, dT.D0 + dT.DV) * dw;
274  else if (is_frequency_parameter(deriv))
275  dF.col(iq).noalias() = -iz * dw;
276  else if ((deriv == JacPropMatType::LineCenter or
277  is_pressure_broadening_DV(deriv)) and
278  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
279  dF.col(iq).noalias() = iz * dw;
280  else if (is_pressure_broadening_G0(deriv) and
281  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
282  dF.col(iq).noalias() = dw;
283  else if (is_pressure_broadening_D0(deriv) and
284  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
285  dF.col(iq).noalias() = iz * dw;
286  else if (is_magnetic_parameter(deriv))
287  dF.col(iq).noalias() = iz * zeeman_df * dw;
288  else if (deriv == JacPropMatType::VMR) {
289  if (Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
290  dF.col(iq).noalias() = Complex(dVMR.G0, dVMR.D0 + dVMR.DV) * dw;
291  else
292  dF.col(iq).setZero();
293  }
294  }
295  }
296 }
297 
299  Eigen::Ref<Eigen::VectorXcd> F,
300  Eigen::Ref<Eigen::MatrixXcd> dF,
301  Eigen::Ref<Eigen::Matrix<Complex, Eigen::Dynamic, ExpectedDataSize()>> data,
302  const Eigen::Ref<const Eigen::VectorXd> f_grid,
303  const Numeric& zeeman_df,
304  const Numeric& magnetic_magnitude,
305  const Numeric& F0_noshift,
306  const Numeric& GD_div_F0,
307  const LineShape::Output& lso,
308  const AbsorptionLines& band,
309  const Index& line_ind,
310  const ArrayOfRetrievalQuantity& derivatives_data,
311  const ArrayOfIndex& derivatives_data_position,
312  const Numeric& dGD_div_F0_dT,
313  const LineShape::Output& dT,
314  const LineShape::Output& dVMR) {
315  constexpr Complex iz(0.0, 1.0);
316 
317  // Size of problem
318  auto nppd = derivatives_data_position.nelem();
319 
320  // Doppler broadening and line center
321  const Numeric F0 = F0_noshift + zeeman_df * magnetic_magnitude + lso.D0 + lso.DV;
322  const Numeric GD = GD_div_F0 * F0;
323  const Numeric invGD = 1.0 / GD;
324  const Numeric dGD_dT = dGD_div_F0_dT * F0 - GD_div_F0 * (dT.D0 + dT.DV);
325 
326  // constant normalization factor for Voigt
327  const Numeric fac = Constant::inv_sqrt_pi * invGD;
328 
329  // Naming req data blocks
330  auto z = data.col(0);
331  auto dw = data.col(1);
332 
333  // Frequency grid
334  z.noalias() = invGD * (Complex(-F0, lso.G0) + f_grid.array()).matrix();
335 
336  // Line shape
337  F.noalias() = fac * z.unaryExpr(&w);
338 
339  if (nppd) {
340  dw.noalias() = 2 * (Complex(0, fac * Constant::inv_sqrt_pi) -
341  z.cwiseProduct(F).array())
342  .matrix();
343 
344  for (auto iq = 0; iq < nppd; iq++) {
345  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
346 
347  if (is_frequency_parameter(deriv))
348  dF.col(iq).noalias() = invGD * dw;
349  else if (deriv == JacPropMatType::Temperature)
350  dF.col(iq).noalias() =
351  dw * Complex(-dT.D0 - dT.DV, dT.G0) * invGD -
352  F * dGD_dT * invGD - dw.cwiseProduct(z) * dGD_dT * invGD;
353  else if ((deriv == JacPropMatType::LineCenter or
354  is_pressure_broadening_DV(deriv)) and
355  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
356  dF.col(iq).noalias() = -F / F0 - dw * invGD - dw.cwiseProduct(z) / F0;
357  else if (is_pressure_broadening_G0(deriv) and
358  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
359  dF.col(iq).noalias() = iz * dw * invGD;
360  else if (is_pressure_broadening_D0(deriv) and
361  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
362  dF.col(iq).noalias() = -dw * invGD;
363  else if (is_magnetic_parameter(deriv))
364  dF.col(iq).noalias() = dw * (-zeeman_df * invGD);
365  else if (deriv == JacPropMatType::VMR and
366  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
367  dF.col(iq).noalias() =
368  dw * Complex(-dVMR.D0 - dVMR.DV, dVMR.G0) * invGD;
369  else
370  dF.col(iq).setZero();
371  }
372  }
373 }
374 
376  Eigen::Ref<Eigen::VectorXcd> F,
377  Eigen::Ref<Eigen::MatrixXcd> dF,
378  Eigen::Ref<Eigen::Matrix<Complex, Eigen::Dynamic, ExpectedDataSize()>> data,
379  const Eigen::Ref<const Eigen::VectorXd> f_grid,
380  const Numeric& zeeman_df,
381  const Numeric& magnetic_magnitude,
382  const Numeric& F0_noshift,
383  const Numeric& GD_div_F0,
384  const AbsorptionLines& band,
385  const Index& line_ind,
386  const ArrayOfRetrievalQuantity& derivatives_data,
387  const ArrayOfIndex& derivatives_data_position,
388  const Numeric& dGD_div_F0_dT) {
389  auto nppd = derivatives_data_position.nelem();
390 
391  // Doppler broadening and line center
392  const Numeric F0 = F0_noshift + zeeman_df * magnetic_magnitude;
393  const Numeric invGD = 1.0 / (GD_div_F0 * F0);
394 
395  // Naming data blocks
396  auto x = data.col(0);
397  auto mx2 = data.col(1);
398 
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();
402 
403  for (auto iq = 0; iq < nppd; iq++) {
404  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
405 
406  if (is_frequency_parameter(deriv))
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
412  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
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(); // Must reset incase other lineshapes are mixed in
417  }
418 }
419 
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,
425  const LineShape::Output& X,
426  const bool with_mirroring,
427  const AbsorptionLines& band,
428  const Index& line_ind,
429  const ArrayOfRetrievalQuantity& derivatives_data,
430  const ArrayOfIndex& derivatives_data_position,
431  const LineShape::Output& dT,
432  const LineShape::Output& dVMR) {
433  auto nppd = derivatives_data_position.nelem();
434 
435  const Complex LM = Complex(1.0 + X.G, -X.Y);
436 
437  dF *= LM;
438  if (with_mirroring) {
439  dF.noalias() += dFm * conj(LM);
440 
441  for (auto iq = 0; iq < nppd; iq++) {
442  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
443 
444  if (deriv == JacPropMatType::Temperature) {
445  const auto c = Complex(dT.G, -dT.Y);
446  dF.col(iq).noalias() += F * c + Fm * conj(c);
447  } else if (is_pressure_broadening_G(deriv) and
448  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
449  dF.col(iq).noalias() = F + Fm;
450  else if (is_pressure_broadening_Y(deriv) and
451  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
452  dF.col(iq).noalias() = Complex(0, -1) * (F - Fm);
453  else if (deriv == JacPropMatType::VMR and
454  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind)) {
455  const auto c = Complex(dVMR.G, -dVMR.Y);
456  dF.col(iq).noalias() += F * c + Fm * conj(c);
457  }
458  }
459  } else {
460  for (auto iq = 0; iq < nppd; iq++) {
461  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
462 
463  if (deriv == JacPropMatType::Temperature)
464  dF.col(iq).noalias() += F * Complex(dT.G, -dT.Y);
465  else if (is_pressure_broadening_G(deriv) and
466  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
467  dF.col(iq).noalias() = F;
468  else if (is_pressure_broadening_Y(deriv) and
469  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
470  dF.col(iq).noalias() = Complex(0, -1) * F;
471  else if (deriv == JacPropMatType::VMR and
472  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
473  dF.col(iq).noalias() += F * Complex(dVMR.G, -dVMR.Y);
474  }
475  }
476 
477  F *= LM;
478  if (with_mirroring) F.noalias() += Fm * std::conj(LM);
479 }
480 
482  Eigen::Ref<Eigen::VectorXcd> F,
483  Eigen::Ref<Eigen::MatrixXcd> dF,
484  const Eigen::Ref<const Eigen::VectorXd> f_grid,
485  const Numeric& F0,
486  const Numeric& T,
487  const AbsorptionLines& band,
488  const Index& line_ind,
489  const ArrayOfRetrievalQuantity& derivatives_data,
490  const ArrayOfIndex& derivatives_data_position) {
491  auto nf = f_grid.size(), nppd = derivatives_data_position.nelem();
492 
493  const Numeric invF0 = 1.0 / F0;
494  const Numeric mafac =
495  (Constant::h) / (2.0 * Constant::k * T) /
496  std::sinh((Constant::h * F0) / (2.0 * Constant::k * T)) * invF0;
497 
498  Numeric dmafac_dT_div_fun = 0;
499  if (do_temperature_jacobian(derivatives_data))
500  dmafac_dT_div_fun =
501  -(Constant::k * T -
502  F0 * Constant::h /
503  (2.0 * std::tanh(F0 * Constant::h / (2.0 * Constant::k * T)))) /
504  (Constant::k * T * T);
505 
506  Numeric fun;
507 
508  for (auto iv = 0; iv < nf; iv++) {
509  fun = mafac * (f_grid[iv] * f_grid[iv]);
510  F[iv] *= fun;
511 
512  for (auto iq = 0; iq < nppd; iq++) {
513  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
514 
515  dF(iv, iq) *= fun;
516  if (deriv == JacPropMatType::Temperature)
517  dF(iv, iq) += dmafac_dT_div_fun * F[iv];
518  else if (deriv == JacPropMatType::LineCenter and
519  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
520  dF(iv, iq) +=
521  (-invF0 - Constant::h / (2.0 * Constant::k * T *
522  std::tanh(F0 * Constant::h /
523  (2.0 * Constant::k * T)))) *
524  F[iv];
525  else if (is_frequency_parameter(deriv))
526  dF(iv, iq) += (2.0 / f_grid[iv]) * F[iv];
527  }
528  }
529 }
530 
532  Eigen::Ref<Eigen::VectorXcd> F,
533  Eigen::Ref<Eigen::MatrixXcd> dF,
534  Eigen::Ref<Eigen::Matrix<Complex, Eigen::Dynamic, ExpectedDataSize()>> data,
535  const Eigen::Ref<const Eigen::VectorXd> f_grid,
536  const Numeric& F0,
537  const Numeric& T,
538  const AbsorptionLines& band,
539  const Index& line_ind,
540  const ArrayOfRetrievalQuantity& derivatives_data,
541  const ArrayOfIndex& derivatives_data_position) {
542  auto nppd = derivatives_data_position.nelem();
543 
544  // 2kT is constant for the loop
545  const Numeric kT = 2.0 * Constant::k * T;
546  const Numeric c1 = Constant::h / kT;
547 
548  // denominator is constant for the loop
549  const Numeric tanh_f0part = std::tanh(c1 * F0);
550  const Numeric denom = F0 * tanh_f0part;
551 
552  // Name the data
553  auto tanh_f = data.col(0);
554  auto ftanh_f = data.col(1);
555 
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();
559 
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]];
563 
564  if (deriv == JacPropMatType::Temperature)
565  dF.col(iq).noalias() +=
566  -c1 / T *
567  ((denom - F0 / tanh_f0part) * F - denom * ftanh_f.cwiseProduct(F) +
568  f_grid.cwiseProduct(F).cwiseQuotient(tanh_f));
569  else if (is_frequency_parameter(deriv))
570  dF.col(iq).noalias() +=
571  F.cwiseQuotient(f_grid) +
572  c1 * (F.cwiseQuotient(tanh_f) - F.cwiseProduct(tanh_f));
573  else if (deriv == JacPropMatType::LineCenter and
574  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
575  dF.col(iq).noalias() +=
576  (-1.0 / F0 + c1 * tanh_f0part - c1 / tanh_f0part) * F;
577  }
578 }
579 
581  Eigen::Ref<Eigen::VectorXcd> F,
582  Eigen::Ref<Eigen::MatrixXcd> dF,
583  const Eigen::Ref<const Eigen::VectorXd> f_grid,
584  const Numeric& F0,
585  const AbsorptionLines& band,
586  const Index& line_ind,
587  const ArrayOfRetrievalQuantity& derivatives_data,
588  const ArrayOfIndex& derivatives_data_position) {
589  auto nf = f_grid.size(), nppd = derivatives_data_position.nelem();
590 
591  // denominator is constant for the loop
592  const Numeric invF0 = 1.0 / F0;
593  const Numeric invF02 = invF0 * invF0;
594 
595  for (auto iv = 0; iv < nf; iv++) {
596  // Set the factor
597  const Numeric fac = f_grid[iv] * f_grid[iv] * invF02;
598 
599  // Set the line shape
600  F[iv] *= fac;
601 
602  for (auto iq = 0; iq < nppd; iq++) {
603  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
604 
605  // The factor is applied to all partial derivatives
606  dF(iv, iq) *= fac;
607 
608  // These partial derivatives are special
609  if (deriv == JacPropMatType::LineCenter and
610  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
611  dF(iv, iq) -= 2.0 * invF0 * F[iv];
612  else if (is_frequency_parameter(deriv))
613  dF(iv, iq) += 2.0 / f_grid[iv] * F[iv];
614  }
615  }
616 }
617 
619  Numeric E0,
620  Numeric F0,
621  Numeric QT0,
622  Numeric T0,
623  Numeric QT,
624  Numeric T) {
625  const Numeric gamma = stimulated_emission(T, F0);
626  const Numeric gamma_ref = stimulated_emission(T0, F0);
627  const Numeric K1 = boltzman_ratio(T, T0, E0);
628  const Numeric K2 = stimulated_relative_emission(gamma, gamma_ref);
629  return S0 * K1 * K2 * QT0 / QT;
630 }
631 
633  Eigen::Ref<Eigen::VectorXcd> F,
634  Eigen::Ref<Eigen::MatrixXcd> dF,
635  Eigen::Ref<Eigen::VectorXcd> N,
636  Eigen::Ref<Eigen::MatrixXcd> dN,
637  const Absorption::SingleLine& line,
638  const Numeric& T,
639  const Numeric& T0,
640  const Numeric& isotopic_ratio,
641  const Numeric& QT,
642  const Numeric& QT0,
643  const AbsorptionLines& band,
644  const Index& line_ind,
645  const ArrayOfRetrievalQuantity& derivatives_data,
646  const ArrayOfIndex& derivatives_data_position,
647  const Numeric& dQT_dT) {
648  auto nppd = derivatives_data_position.nelem();
649 
650  const Numeric gamma = stimulated_emission(T, line.F0());
651  const Numeric gamma_ref = stimulated_emission(T0, line.F0());
652  const Numeric K1 = boltzman_ratio(T, T0, line.E0());
653  const Numeric K2 = stimulated_relative_emission(gamma, gamma_ref);
654 
655  const Numeric invQT = 1.0 / QT;
656  const Numeric S = line.I0() * isotopic_ratio * QT0 * invQT * K1 * K2;
657 
658  F *= S;
659  dF *= S;
660  for (auto iq = 0; iq < nppd; iq++) {
661  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
662 
663  if (deriv == JacPropMatType::Temperature)
664  dF.col(iq).noalias() +=
665  F * (dstimulated_relative_emission_dT(gamma, gamma_ref, line.F0(), T) /
666  K2 +
667  dboltzman_ratio_dT_div_boltzmann_ratio(T, line.E0()) - invQT * dQT_dT);
668  else if (deriv == JacPropMatType::LineStrength and
669  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
670  dF.col(iq).noalias() = F / line.I0(); //nb. overwrite
671  else if (deriv == JacPropMatType::LineCenter and
672  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind))
673  dF.col(iq).noalias() +=
674  F *
675  dstimulated_relative_emission_dF0(gamma, gamma_ref, T, T0) /
676  K2;
677  }
678 
679  // No NLTE variables
680  N.setZero();
681  dN.setZero();
682 }
683 
685  Eigen::Ref<Eigen::VectorXcd> F,
686  Eigen::Ref<Eigen::MatrixXcd> dF,
687  Eigen::Ref<Eigen::VectorXcd> N,
688  Eigen::Ref<Eigen::MatrixXcd> dN,
689  const Absorption::SingleLine& line,
690  const Numeric& T,
691  const Numeric& T0,
692  const Numeric& Tu,
693  const Numeric& Tl,
694  const Numeric& Evu,
695  const Numeric& Evl,
696  const Numeric& isotopic_ratio,
697  const Numeric& QT,
698  const Numeric& QT0,
699  const AbsorptionLines& band,
700  const Index& line_ind,
701  const ArrayOfRetrievalQuantity& derivatives_data,
702  const ArrayOfIndex& derivatives_data_position,
703  const Numeric& dQT_dT) {
704  auto nppd = derivatives_data_position.nelem();
705 
706  const Numeric gamma = stimulated_emission(T,line.F0());
707  const Numeric gamma_ref = stimulated_emission(T0,line.F0());
708  const Numeric r_low = boltzman_ratio(Tl,T,Evl);
709  const Numeric r_upp = boltzman_ratio(Tu,T,Evu);
710 
711  const Numeric K1 = boltzman_ratio(T,T0,line.E0());
712  const Numeric K2 = stimulated_relative_emission(gamma,gamma_ref);
713  const Numeric K3 = absorption_nlte_ratio(gamma,r_upp,r_low);
714  const Numeric K4 = r_upp;
715 
716  const Numeric invQT = 1.0 / QT;
717  const Numeric QT_ratio = QT0 * invQT;
718 
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;
723 
724  dN.noalias() = dF * (S_src - S_abs);
725  dF *= S_abs;
726  for (auto iq = 0; iq < nppd; iq++) {
727  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
728 
729  if (deriv == JacPropMatType::Temperature) {
730  const Numeric dS_dT_abs =
731  S_abs *
732  (dstimulated_relative_emission_dT(gamma, gamma_ref, line.F0(), T) /
733  K2 +
734  dboltzman_ratio_dT(K1, T, line.E0()) / K1 +
735  dabsorption_nlte_rate_dT(gamma, T, line.F0(), Evl, Evu, K4, r_low) /
736  K3 -
737  invQT * dQT_dT);
738  const Numeric dS_dT_src =
739  S_src *
740  (dstimulated_relative_emission_dT(gamma, gamma_ref, line.F0(), T) /
741  K2 +
742  dboltzman_ratio_dT(K1, T, line.E0()) / K1 -
743  dboltzman_ratio_dT(K4, T, Evu) / K4 - invQT * dQT_dT);
744 
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
748  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind)) {
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
752  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind)) {
753  const Numeric dS_dF0_abs =
754  S_abs *
755  (dstimulated_relative_emission_dF0(gamma, gamma_ref, T, T0) /
756  K2 +
757  dabsorption_nlte_rate_dF0(gamma, T, K4, r_low) / K3);
758  const Numeric dS_dF0_src =
759  S_src *
760  dstimulated_relative_emission_dF0(gamma, gamma_ref, T, T0) /
761  K2;
762 
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) {
766  if (Absorption::id_in_line_lower(band, deriv.QuantumIdentity(), line_ind)) {
767  const Numeric dS_dTl_abs =
768  S_abs * dabsorption_nlte_rate_dTl(gamma, T, Tl, Evl, r_low) / K3;
769 
770  dN.col(iq).noalias() = -F * dS_dTl_abs;
771  dF.col(iq).noalias() = -dN.col(iq);
772  } else if (Absorption::id_in_line_upper(band, deriv.QuantumIdentity(), line_ind)) {
773  const Numeric dS_dTu_abs =
774  S_abs * dabsorption_nlte_rate_dTu(gamma, T, Tu, Evu, K4) / K3;
775  const Numeric dS_dTu_src = S_src * dboltzman_ratio_dT(K4, Tu, Evu) / K4;
776 
777  dN.col(iq).noalias() = F * (dS_dTu_src - dS_dTu_abs);
778  dF.col(iq).noalias() = F * dS_dTu_abs;
779  }
780  }
781  }
782 
783  N.noalias() = F * (S_src - S_abs);
784  F *= S_abs;
785 }
786 
788  Eigen::Ref<Eigen::MatrixXcd> dF,
789  const AbsorptionLines& band,
790  const Index& line_ind,
791  const ArrayOfRetrievalQuantity& derivatives_data,
792  const ArrayOfIndex& derivatives_data_position,
793  const Numeric& T,
794  const Numeric& P,
795  const Vector& vmrs) {
796  auto nppd = derivatives_data_position.nelem();
797 
798  for (auto iq = 0; iq < nppd; iq++) {
799  const RetrievalQuantity& rt =
800  derivatives_data[derivatives_data_position[iq]];
801  if (is_lineshape_parameter(rt) and
802  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
803  dF.col(iq) *= band.ShapeParameter_dInternal(line_ind, T, P, vmrs, rt);
804  }
805 }
806 
808  return std::sqrt(Constant::doppler_broadening_const_squared * T / mass);
809 }
810 
812  const Numeric& dc) {
813  return dc / (2 * T);
814 }
815 
817  Index& start_cutoff,
818  Index& nelem_cutoff,
819  const Eigen::Ref<const Eigen::VectorXd> f_grid,
820  const Numeric& fmin,
821  const Numeric& fmax) {
822  auto nf = f_grid.size();
823 
824  const bool need_cutoff = (fmax > fmin);
825  if (need_cutoff) {
826  // Find range of simulations
827  start_cutoff = 0;
828  Index i_f_max = nf - 1;
829 
830  // Loop over positions to compute the line shape cutoff point
831  while (start_cutoff < nf and fmin > f_grid[start_cutoff])
832  ++start_cutoff;
833  while (i_f_max >= start_cutoff and fmax < f_grid[i_f_max])
834  --i_f_max;
835 
836  // The extent is one more than the difference between the indices of interest
837  nelem_cutoff = i_f_max - start_cutoff + 1; // min is 0, max is nf
838  } else {
839  start_cutoff = 0;
840  nelem_cutoff = nf;
841  }
842 }
843 
845  Eigen::Ref<Eigen::VectorXcd> F,
846  Eigen::Ref<Eigen::MatrixXcd> dF,
847  Eigen::Ref<Eigen::VectorXcd> N,
848  Eigen::Ref<Eigen::MatrixXcd> dN,
849  const Numeric& r1,
850  const Numeric& r2,
851  const Numeric& g1,
852  const Numeric& g2,
853  const Numeric& A21,
854  const Numeric& F0,
855  const Numeric& T,
856  const AbsorptionLines& band,
857  const Index& line_ind,
858  const ArrayOfRetrievalQuantity& derivatives_data,
859  const ArrayOfIndex& derivatives_data_position) {
860  // Size of the problem
861  auto nppd = derivatives_data_position.nelem();
862 
863  // Physical constants
864  constexpr Numeric c0 = 2.0 * Constant::h / Constant::pow2(Constant::c);
865  constexpr Numeric c1 = Constant::h / (4 * Constant::pi);
866 
867  // Constants based on input
868  const Numeric c2 = c0 * F0 * F0 * F0;
869  const Numeric c3 = c1 * F0;
870  const Numeric x = g2 / g1;
871 
872  /*
873  Einstein 'active' coefficients are as follows:
874  const Numeric B21 = A21 / c2;
875  const Numeric B12 = x * B21;
876  */
877 
878  // Absorption strength
879  const Numeric k = c3 * (r1 * x - r2) * (A21 / c2);
880 
881  // Emission strength
882  const Numeric e = c3 * r2 * A21;
883 
884  // Planck function of this line
885  const Numeric exp_T = std::exp(Constant::h / Constant::k * F0 / T),
886  b = c2 / (exp_T - 1);
887 
888  // Ratio between emission and absorption constant
889  const Numeric ratio = e / b - k;
890 
891  // Constants ALMOST everywhere inside these loops
892  dN.noalias() = dF * ratio;
893  dF *= k;
894 
895  // Partial derivatives
896  for (auto iq = 0; iq < nppd; iq++) {
897  const auto& deriv = derivatives_data[derivatives_data_position[iq]];
898 
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
903  Absorption::id_in_line(band, deriv.QuantumIdentity(), line_ind)) {
904  Numeric done_over_b_df0 =
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;
908 
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) {
912  if (Absorption::id_in_line_lower(band, deriv.QuantumIdentity(), line_ind)) {
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;
917  } else if (Absorption::id_in_line_upper(band, deriv.QuantumIdentity(), line_ind)) {
918  const Numeric dk_dr1 = c3 * x * A21 / c2;
919  dF.col(iq).noalias() = F * dk_dr1;
920  }
921  }
922  }
923 
924  // Set source function to be the relative amount emitted by the line divided by the Planck function
925  N.noalias() = F * ratio;
926 
927  // Set absorption
928  F *= k;
929 }
930 
931 void Linefunctions::set_htp(Eigen::Ref<Eigen::VectorXcd> F,
932  Eigen::Ref<Eigen::MatrixXcd> dF,
933  const Eigen::Ref<const Eigen::VectorXd> f_grid,
934  const Numeric& zeeman_df_si,
935  const Numeric& magnetic_magnitude_si,
936  const Numeric& F0_noshift_si,
937  const Numeric& GD_div_F0_si,
938  const LineShape::Output& lso_si,
939  const AbsorptionLines& band,
940  const Index& line_ind,
941  const ArrayOfRetrievalQuantity& derivatives_data,
942  const ArrayOfIndex& derivatives_data_position,
943  const Numeric& dGD_div_F0_dT_si,
944  const LineShape::Output& dT_si,
945  const LineShape::Output& dVMR_si) {
946  using Constant::inv_sqrt_pi;
947  using Constant::pi;
948  using Constant::pow2;
949  using Constant::pow3;
950  using Constant::pow4;
951  using Constant::sqrt_ln_2;
952  using Constant::sqrt_pi;
954  using std::abs;
955  using std::imag;
956  using std::real;
957  using std::sqrt;
958 
959  // Convert to CGS for using original function
960  const Numeric sg0 =
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);
966 
967  // General normalization
968  const Numeric cte = sqrt_ln_2 / GamD;
969  constexpr Complex iz(0, 1);
970 
971  // Calculating the different parameters
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;
976  const Complex Y = pow2(1 / (2 * cte * c2t));
977  const Complex sqrtY = sqrt(Y);
978 
979  // For all frequencies
980  for (auto iv = 0; iv < f_grid.size(); iv++) {
981  const Complex X = (iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) / c2t;
982  const Complex sqrtXY = sqrt(X + Y);
983  const Complex sqrtX = sqrt(X);
984 
985  // Declare terms required for the derivatives as external to the if-statement
986  Complex Z1, Z2, Zb, W1, W2, Wb;
987  Complex Aterm, Bterm;
988  if (abs(c2t) ==
989  0) { // If this method does not require G2, D2, or ETA==1, then FVC matters.
990  Z1 = (iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) * cte;
991  W1 = w(iz * Z1);
992 
993  Aterm = sqrt_pi * cte * W1;
994  if (abs(Z1) <=
995  4e3) // For very large Z1 (i.e., very large broadening or very large distance from line-center
996  Bterm = sqrt_pi * cte * ((1 - pow2(Z1)) * W1 + Z1 * inv_sqrt_pi);
997  else
998  Bterm = cte * (sqrt_pi * W1 + 0.5 / Z1 - 0.75 / pow3(Z1));
999  } else if (
1000  abs(X) <=
1001  3e-8 *
1002  abs(Y)) { // If this method is executed very close to the line center
1003  Z1 = (iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) * cte;
1004  Z2 = sqrtXY + sqrtY;
1005 
1006  W1 = w(iz * Z1);
1007  W2 = w(iz * Z2);
1008 
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) /
1012  c2t;
1013  } else if (
1014  abs(Y) <=
1015  1e-15 *
1016  abs(X)) { // If this method is executed very far from the line center
1017  Z1 = sqrtXY;
1018 
1019  W1 = w(iz * Z1);
1020 
1021  if (abs(sqrtX) <= 4e3) { // If X is small still
1022  Zb = sqrtX;
1023  Wb = w(iz * Zb);
1024 
1025  Aterm = (2 * sqrt_pi / c2t) * (inv_sqrt_pi - sqrtX * Wb);
1026  Bterm =
1027  (1 / c2t) *
1028  (-1 + 2 * sqrt_pi * (1 - X - 2 * Y) * (inv_sqrt_pi - sqrtX * Wb) +
1029  2 * sqrt_pi * sqrtXY * W1);
1030  } else {
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);
1034  }
1035  } else { // General calculations
1036  Z1 = sqrtXY - sqrtY;
1037  Z2 = Z1 + 2 * sqrtY;
1038 
1039  // NOTE: the region of w might matter according to original code! So this might need changing...
1040  W1 = w(iz * Z1);
1041  W2 = w(iz * Z2);
1042 
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) /
1046  c2t;
1047  }
1048 
1049  F(iv) = Aterm / (pi * (((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm +
1050  Bterm * c2 * lso.ETA + 1));
1051 
1052  for (auto iq = 0; iq < derivatives_data_position.nelem(); iq++) {
1053  const RetrievalQuantity& rt =
1054  derivatives_data[derivatives_data_position[iq]];
1055 
1056  Numeric dcte = 0;
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) *
1060  cte; // for T
1061  dc0 = Complex(dT.G0, -dT.D0); // for T
1062  dc2 = Complex(dT.G2, -dT.D2); // for T
1063  dc0t = (-(c0 - 1.5 * c2) * dT.ETA + dT.FVC) +
1064  ((1 - lso.ETA) * (dc0 - 1.5 * dc2)); // for T
1065  dc2t = (-c2 * dT.ETA) + ((1 - lso.ETA) * dc2); // for T
1066  } else if (rt == JacPropMatType::VMR and
1067  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1068  dc0 = Complex(dV.G0, -dV.D0); // for VMR
1069  dc2 = Complex(dV.G2, -dV.D2); // for VMR
1070  dc0t = (-(c0 - 1.5 * c2) * dV.ETA + dV.FVC) +
1071  ((1 - lso.ETA) * (dc0 - 1.5 * dc2)); // for VMR
1072  dc2t = (-c2 * dV.ETA) + ((1 - lso.ETA) * dc2); // for VMR
1073  } else if (is_pressure_broadening_G2(rt) and
1074  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1075  dc2 =
1076  1; // for Gam2(T, VMR) and for Shift2(T, VMR), the derivative is wrt c2 for later computations
1077  dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1078  dc2t = ((1 - lso.ETA) * dc2);
1079  } else if (is_pressure_broadening_D2(rt) and
1080  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1081  dc2 =
1082  -iz; // for Gam2(T, VMR) and for Shift2(T, VMR), the derivative is wrt c2 for later computations
1083  dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1084  dc2t = ((1 - lso.ETA) * dc2);
1085  } else if (is_pressure_broadening_G0(rt) and
1086  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1087  dc0 =
1088  1; // for Gam0(T, VMR) and for Shift0(T, VMR), the derivative is wrt c0 for later computations
1089  dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1090  } else if (is_pressure_broadening_D0(rt) and
1091  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1092  dc0 =
1093  -iz; // for Gam0(T, VMR) and for Shift0(T, VMR), the derivative is wrt c0 for later computations
1094  dc0t = ((1 - lso.ETA) * (dc0 - 1.5 * dc2));
1095  } else if (is_pressure_broadening_FVC(rt) and
1096  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1097  dc0t = (1); // for FVC(T, VMR)
1098  } else if (is_pressure_broadening_ETA(rt) and
1099  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind)) {
1100  dc0t = (-c0 + 1.5 * c2); // for eta(T, VMR)
1101  dc2t = (-c2); // for eta(T, VMR)
1102  }
1103 
1104  const Complex dY = (-2 * dcte / cte - 2 * dc2t / c2t) * Y; // for all
1105 
1106  Complex dX;
1107  if (is_magnetic_parameter(rt))
1108  dX = -iz * freq2kaycm(zeeman_df_si) / c2t; // for H
1109  else if (rt == JacPropMatType::LineCenter and
1110  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
1111  dX = -iz / c2t; // for sg0
1112  else if (is_frequency_parameter(rt))
1113  dX = iz / c2t; // for sg
1114  else
1115  dX =
1116  (-(iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) * dc2t + c2t * dc0t) /
1117  pow2(c2t); // for c0t and c2t
1118 
1119  Complex dAterm, dBterm;
1120  if (abs(c2t) == 0) {
1121  Complex dZ1;
1122  if (is_magnetic_parameter(rt))
1123  dZ1 = -iz * cte * freq2kaycm(zeeman_df_si); // for H
1124  else if (rt == JacPropMatType::LineCenter and
1125  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
1126  dZ1 = -iz * cte; // for sg0
1127  else if (is_frequency_parameter(rt))
1128  dZ1 = iz * cte; // for sg
1129  else
1130  dZ1 = (iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) * dcte +
1131  cte * dc0t; // for c0t
1132 
1133  const Complex dW1 = iz * dZ1 * dw(Z1, W1); // NEED TO CHECK DW!
1134 
1135  dAterm = sqrt_pi * (W1 * dcte + cte * dW1); // for all
1136 
1137  if (abs(Z1) <= 4e3)
1138  dBterm =
1139  -(sqrt_pi * ((pow2(Z1) - 1) * dW1 + 2 * W1 * Z1 * dZ1) - dZ1) *
1140  cte -
1141  (sqrt_pi * (pow2(Z1) - 1) * W1 - Z1) * dcte; // for all
1142  else
1143  dBterm =
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) *
1146  cte) /
1147  pow4(Z1); // for all
1148  } else if (abs(X) <= 3e-8 * abs(Y)) {
1149  Complex dZ1;
1150  if (is_magnetic_parameter(rt))
1151  dZ1 = -iz * cte * freq2kaycm(zeeman_df_si); // for H
1152  else if (rt == JacPropMatType::LineCenter and
1153  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
1154  dZ1 = -iz * cte; // for sg0
1155  else if (is_frequency_parameter(rt))
1156  dZ1 = iz * cte; // for sg
1157  else
1158  dZ1 = (iz * (freq2kaycm(f_grid[iv]) - sg0) + c0t) * dcte +
1159  cte * dc0t; // for c0t
1160 
1161  const Complex dZ2 = dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1162  dY / (2 * sqrtXY); // for all
1163 
1164  const Complex dW1 = iz * dZ1 * dw(Z1, W1); // NEED TO CHECK DW!
1165  const Complex dW2 = iz * dZ2 * dw(Z2, W2); // NEED TO CHECK DW!
1166 
1167  dAterm = sqrt_pi * ((W1 - W2) * dcte + (dW1 - dW2) * cte); // for all
1168 
1169  dBterm = (sqrt_pi *
1170  (((pow2(Z1) - 1) * W1 - (pow2(Z2) - 1) * W2) * dY +
1171  2 *
1172  (-(pow2(Z1) - 1) * dW1 + (pow2(Z2) - 1) * dW2 -
1173  2 * W1 * Z1 * dZ1 + 2 * W2 * Z2 * dZ2) *
1174  Y) *
1175  c2t +
1176  2 *
1177  (sqrt_pi * (pow2(Z1) - 1) * W1 -
1178  sqrt_pi * (pow2(Z2) - 1) * W2 + 2 * sqrtY) *
1179  Y * dc2t) /
1180  (4 * Y * sqrtY * pow2(c2t)); // for all
1181  } else if (abs(Y) <= 1e-15 * abs(X)) {
1182  const Complex dZ1 = (dX + dY) / (2 * sqrtXY); // for all
1183  const Complex dW1 = iz * dZ1 * dw(Z1, W1); // NEED TO CHECK DW!
1184  if (abs(sqrtX) <= 4e3) {
1185  const Complex dZb = dX / (2 * sqrtX); // for all
1186  const Complex dWb = iz * dZb * dw(Zb, Wb); // NEED TO CHECK DW!
1187 
1188  dAterm = (-sqrt_pi * (Wb * dX + 2 * X * dWb) * c2t +
1189  2 * (sqrt_pi * Wb * sqrtX - 1) * sqrtX * dc2t) /
1190  (sqrtX * pow2(c2t)); // for all
1191 
1192  dBterm =
1193  (-sqrtXY *
1194  (2 * (sqrt_pi * Wb * sqrtX - 1) * (X + 2 * Y - 1) +
1195  2 * sqrt_pi * sqrtXY * W1 - 1) *
1196  sqrtX * dc2t +
1197  (2 *
1198  ((sqrt_pi * Wb * sqrtX - 1) * (dX + 2 * dY) +
1199  sqrt_pi * sqrtXY * dW1) *
1200  sqrtXY * sqrtX +
1201  sqrt_pi * (Wb * dX + 2 * X * dWb) * sqrtXY * (X + 2 * Y - 1) +
1202  sqrt_pi * (dX + dY) * W1 * sqrtX) *
1203  c2t) /
1204  (sqrtXY * sqrtX * pow2(c2t)); // for all
1205  } else {
1206  dAterm = ((-X + 3.0) * c2t * dX - (X - 1.5) * X * dc2t) /
1207  (pow3(X) * pow2(c2t)); // for all
1208 
1209  dBterm = (((-2 * sqrt_pi * sqrtXY * W1 + 1) * pow2(X) +
1210  (X - 1.5) * (X + 2 * Y - 1)) *
1211  sqrtXY * X * dc2t +
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)) *
1216  c2t) /
1217  (sqrtXY * pow3(X) * pow2(c2t)); // for all
1218  }
1219  } else {
1220  const Complex dZ1 = -dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1221  dY / (2 * sqrtXY); // for all
1222  const Complex dZ2 = dY / (2 * sqrtY) + dX / (2 * sqrtXY) +
1223  dY / (2 * sqrtXY); // for all
1224 
1225  const Complex dW1 = iz * dZ1 * dw(Z1, W1); // NEED TO CHECK DW!
1226  const Complex dW2 = iz * dZ2 * dw(Z2, W2); // NEED TO CHECK DW!
1227 
1228  dAterm = sqrt_pi * ((W1 - W2) * dcte + (dW1 - dW2) * cte); // for all
1229 
1230  dBterm = (sqrt_pi *
1231  (((pow2(Z1) - 1) * W1 - (pow2(Z2) - 1) * W2) * dY +
1232  2 *
1233  (-(pow2(Z1) - 1) * dW1 + (pow2(Z2) - 1) * dW2 -
1234  2 * W1 * Z1 * dZ1 + 2 * W2 * Z2 * dZ2) *
1235  Y) *
1236  c2t +
1237  2 *
1238  (sqrt_pi * (pow2(Z1) - 1) * W1 -
1239  sqrt_pi * (pow2(Z2) - 1) * W2 + 2 * sqrtY) *
1240  Y * dc2t) /
1241  (4 * Y * sqrtY * pow2(c2t)); // for all
1242  }
1243 
1244  Numeric dFVC = 0, dETA = 0;
1245  if (rt == JacPropMatType::Temperature) {
1246  dETA = dT.ETA;
1247  dFVC = dT.FVC;
1248  } else if (rt == JacPropMatType::VMR) {
1249  dETA = dV.ETA;
1250  dFVC = dV.FVC;
1251  } else if (is_pressure_broadening_ETA(rt) and
1252  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
1253  dETA = 1;
1254  else if (is_pressure_broadening_FVC(rt) and
1255  Absorption::id_in_line(band, rt.QuantumIdentity(), line_ind))
1256  dFVC = 1;
1257 
1258  dF.col(iq)(iv) =
1259  ((((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm + Bterm * c2 * lso.ETA +
1260  1) *
1261  dAterm -
1262  (((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * dAterm +
1263  ((c0 - 1.5 * c2) * dETA + (dc0 - 1.5 * dc2) * lso.ETA - dFVC) *
1264  Aterm +
1265  Bterm * c2 * dETA + Bterm * lso.ETA * dc2 + c2 * lso.ETA * dBterm) *
1266  Aterm) /
1267  (pi * pow2(((c0 - 1.5 * c2) * lso.ETA - lso.FVC) * Aterm +
1268  Bterm * c2 * lso.ETA + 1)); // for all
1269  }
1270  }
1271 
1272  // Convert back to ARTS
1273  F = F.unaryExpr(&pCqSDHC_to_arts);
1274  for (auto iq = 0; iq < derivatives_data_position.nelem(); iq++) {
1275  const RetrievalQuantity& rt =
1276  derivatives_data[derivatives_data_position[iq]];
1277 
1278  if (rt == JacPropMatType::LineCenter or is_frequency_parameter(rt) or
1281  dF.col(iq) = dF.col(iq).unaryExpr(&pCqSDHC_to_arts_freq_deriv);
1282  else if (is_pressure_broadening_G2(rt))
1283  dF.col(iq) = dF.col(iq).unaryExpr(&pCqSDHC_to_arts_G2_deriv);
1284  else if (is_pressure_broadening_D2(rt))
1285  dF.col(iq) = dF.col(iq).unaryExpr(&pCqSDHC_to_arts_D2_deriv);
1286  else
1287  dF.col(iq) = dF.col(iq).unaryExpr(&pCqSDHC_to_arts);
1288  }
1289 }
1290 
1292  InternalData& scratch,
1293  InternalData& sum,
1294  const ConstVectorView f_grid,
1295  const AbsorptionLines& band,
1296  const ArrayOfRetrievalQuantity& derivatives_data,
1297  const ArrayOfIndex& derivatives_data_active,
1298  const Vector& vmrs,
1299  const EnergyLevelMap& nlte,
1300  const Numeric& P,
1301  const Numeric& T,
1302  const Numeric& isot_ratio,
1303  const Numeric& H,
1304  const Numeric& DC,
1305  const Numeric& dDCdT,
1306  const Numeric& QT,
1307  const Numeric& dQTdT,
1308  const Numeric& QT0,
1309  const bool no_negatives,
1310  const bool zeeman,
1311  const Zeeman::Polarization zeeman_polarization)
1312 {
1313  const Index nj = derivatives_data_active.nelem();
1314  const bool do_temperature = do_temperature_jacobian(derivatives_data);
1315 
1316  // Sum up variable reset
1317  sum.SetZero();
1318 
1319  if (band.NumLines() == 0 or (Absorption::relaxationtype_relmat(band.Population()) and band.LinemixingLimit() > P)) {
1320  return; // No line-by-line computations required/wanted
1321  }
1322 
1323  // Cutoff for Eigen-library types
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;
1330 
1331  // Frequency grid as Eigen type
1332  const auto f_full = MapToEigen(f_grid);
1333 
1334  // Cut off range
1335  const Numeric fmean = (band.Cutoff() == Absorption::CutoffType::BandFixedFrequency) ? band.F_mean() : 0;
1336  Numeric fcut_upp, fcut_low;
1337  Index start, nelem;
1338  fcut_upp = band.CutoffFreq(0);
1339  fcut_low = band.CutoffFreqMinus(0, fmean);
1340  find_cutoff_ranges(start, nelem, f_full, fcut_low, fcut_upp);
1341  fc[0] = fcut_upp;
1342 
1343  // VMR Jacobian check
1344  auto do_vmr = do_vmr_jacobian(derivatives_data, band.QuantumIdentity());
1345 
1346  // Placeholder nothingness
1347  constexpr LineShape::Output empty_output = {0, 0, 0, 0, 0, 0, 0, 0, 0};
1348 
1349  for (Index i=0; i<band.NumLines(); i++) {
1350 
1351  // Select the range of cutoff if different for each line
1352  if (band.Cutoff() == Absorption::CutoffType::LineByLineOffset and i>0) {
1353  fcut_upp = band.CutoffFreq(i);
1354  fcut_low = band.CutoffFreqMinus(i, fmean);
1355  find_cutoff_ranges(start, nelem, f_full, fcut_low, fcut_upp);
1356  fc[0] = fcut_upp;
1357  }
1358 
1359  // Relevant range FIXME: By Band and no-cutoff does not need this...
1360  auto F = scratch.F.segment(start, nelem);
1361  auto N = scratch.N.segment(start, nelem);
1362  auto dF = scratch.dF.middleRows(start, nelem);
1363  auto dN = scratch.dN.middleRows(start, nelem);
1364  auto data = scratch.data.middleRows(start, nelem);
1365  const auto f = f_full.middleRows(start, nelem);
1366 
1367  // Pressure broadening and line mixing terms
1368  const auto X = band.ShapeParameters(i, T, P, vmrs);
1369 
1370  // Partial derivatives for temperature
1371  const auto dXdT = do_temperature ?
1372  band.ShapeParameters_dT(i, T, P, vmrs) : empty_output;
1373 
1374  // Partial derivatives for VMR of self (function works for any species but only do self for now)
1375  const auto dXdVMR = do_vmr.test ?
1376  band.ShapeParameters_dVMR(i, T, P, do_vmr.qid) : empty_output;
1377 
1378  // Zeeman lines if necessary
1379  const Index nz = zeeman ?
1380  band.ZeemanCount(i, zeeman_polarization) : 1;
1381 
1382  for (Index iz=0; iz<nz; iz++) {
1383 
1384  // Zeeman values for this sub-line
1385  const Numeric Sz = zeeman ?
1386  band.ZeemanStrength(i, zeeman_polarization, iz) : 1;
1387  const Numeric dfdH = zeeman ?
1388  band.ZeemanSplitting(i, zeeman_polarization, iz) : 0;
1389 
1390  // Set the line shape and its derivatives
1391  switch (band.LineShapeType()) {
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);
1396  break;
1397  case LineShape::Type::HTP:
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);
1402  break;
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);
1407  break;
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);
1412  break;
1413  }
1414 
1415  // Remove the cutoff values
1416  if (band.Cutoff() not_eq Absorption::CutoffType::None) {
1417  F.array() -= Fc[0];
1418  for (Index ij = 0; ij < nj; ij++) {
1419  dF.col(ij).array() -= dFc[ij];
1420  }
1421  }
1422 
1423  // Set the mirrored line shape
1424  const bool with_mirroring =
1425  band.Mirroring() not_eq Absorption::MirroringType::None and
1427  switch (band.Mirroring()) {
1428  case Absorption::MirroringType::None:
1430  break;
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);
1435  break;
1436  case Absorption::MirroringType::SameAsLineShape:
1437  switch (band.LineShapeType()) {
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);
1442  break;
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);
1447  break;
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);
1452  break;
1453  case LineShape::Type::HTP:
1454  case LineShape::Type::SDVP:
1455  // WARNING: This mirroring is not tested and it might require, e.g., FVC to be treated differently
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);
1459  break;
1460  }
1461  break;
1462  }
1463 
1464  // Remove the mirrored cutoff values
1465  if (band.Cutoff() not_eq Absorption::CutoffType::None and with_mirroring) {
1466  N.array() -= Nc[0];
1467  for (Index ij = 0; ij < nj; ij++) {
1468  dN.col(ij).array() -= dNc[ij];
1469  }
1470  }
1471 
1472  // Mirror and and line mixing is added together (because of conjugate)
1473  if (band.LineShapeType() not_eq LineShape::Type::DP) {
1474  apply_linemixing_scaling_and_mirroring(F, dF, N, dN, X, with_mirroring, band, i, derivatives_data, derivatives_data_active, dXdT, dXdVMR);
1475 
1476  // Apply line mixing and pressure broadening partial derivatives
1477  apply_lineshapemodel_jacobian_scaling(dF, band, i, derivatives_data, derivatives_data_active, T, P, vmrs);
1478  }
1479 
1480  // Normalize the lines
1481  switch (band.Normalization()) {
1482  case Absorption::NormalizationType::None:
1483  break;
1484  case Absorption::NormalizationType::VVH:
1485  apply_VVH_scaling(F, dF, data, f, band.F0(i), T, band, i, derivatives_data, derivatives_data_active);
1486  break;
1487  case Absorption::NormalizationType::VVW:
1488  apply_VVW_scaling(F, dF, f, band.F0(i), band, i, derivatives_data, derivatives_data_active);
1489  break;
1491  apply_rosenkranz_quadratic_scaling(F, dF, f, band.F0(i), T, band, i, derivatives_data, derivatives_data_active);
1492  break;
1493  }
1494 
1495  // Apply line strength by whatever method is necessary
1496  switch (band.Population()) {
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);
1501  break;
1502  case Absorption::PopulationType::ByNLTEVibrationalTemperatures: {
1503  auto nlte_data = nlte.get_vibtemp_params(band, i, T);
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);
1505  } break;
1506  case Absorption::PopulationType::ByNLTEPopulationDistribution: {
1507  auto nlte_data = nlte.get_ratio_params(band, i);
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);
1509  } break;
1510  }
1511 
1512  // Zeeman-adjusted strength
1513  if (zeeman) {
1514  F *= Sz;
1515  N *= Sz;
1516  dF *= Sz;
1517  dN *= Sz;
1518  }
1519 
1520  // Sum up the contributions
1521  sum.F.segment(start, nelem).noalias() += F;
1522  sum.N.segment(start, nelem).noalias() += N;
1523  sum.dF.middleRows(start, nelem).noalias() += dF;
1524  sum.dN.middleRows(start, nelem).noalias() += dN;
1525  }
1526  }
1527 
1528  // Set negative values to zero incase this is requested
1529  if (no_negatives) {
1530  auto reset_zeroes = (sum.F.array().real() < 0);
1531 
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);
1538  }
1539 }
Linefunctions::InternalData::dF
Eigen::MatrixXcd dF
Definition: linefunctions.h:522
do_vmr_jacobian
jacobianVMRcheck do_vmr_jacobian(const ArrayOfRetrievalQuantity &js, const QuantumIdentifier &line_qid) noexcept
Returns the required info for VMR Jacobian.
Definition: jacobian.cc:1283
Absorption::Lines::Cutoff
CutoffType Cutoff() const noexcept
Returns cutoff style.
Definition: absorptionlines.h:1128
Linefunctions::InternalData::datac
Eigen::Matrix< Complex, 1, Linefunctions::ExpectedDataSize()> datac
Definition: linefunctions.h:531
EnergyLevelMap::get_ratio_params
Output2 get_ratio_params(const AbsorptionLines &band, const Index &line_index) const
Get the output required for Population::NLTE.
Definition: energylevelmap.cc:31
fac
Numeric fac(const Index n)
fac
Definition: math_funcs.cc:63
Absorption::Lines::g_low
Numeric g_low(size_t k) const noexcept
Lower level statistical weight.
Definition: absorptionlines.h:1052
is_pressure_broadening_Y
bool is_pressure_broadening_Y(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a Y derivative.
Absorption::SingleLine::F0
Numeric F0() const noexcept
Central frequency.
Definition: absorptionlines.h:342
EnergyLevelMap::get_vibtemp_params
Output4 get_vibtemp_params(const AbsorptionLines &band, const Index &line_index, const Numeric T) const
Get the output required for Population::NLTE-VibrationalTemperatures.
Definition: energylevelmap.cc:58
Linefunctions::apply_rosenkranz_quadratic_scaling
void apply_rosenkranz_quadratic_scaling(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &F0, const Numeric &T, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex())
Applies Rosenkranz quadratic normalization to already set line shape.
Definition: linefunctions.cc:481
Linefunctions::apply_VVH_scaling
void apply_VVH_scaling(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::Matrix< Complex, Eigen::Dynamic, ExpectedDataSize()>> data, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &F0, const Numeric &T, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex())
Applies Van Vleck and Huber normalization to already set line shape.
Definition: linefunctions.cc:531
conj
constexpr Complex conj(Complex c)
the conjugate of c
Definition: complex.h:65
w
Complex w(Complex z) noexcept
The Faddeeva function.
Definition: linefunctions.cc:42
Absorption::Lines::LineShapeType
LineShape::Type LineShapeType() const noexcept
Returns lineshapetype style.
Definition: absorptionlines.h:1180
Absorption::Lines::F_mean
Numeric F_mean() const noexcept
Mean frequency by weight of line strengt.
Definition: absorptionlines.h:983
Absorption::SingleLine::I0
Numeric I0() const noexcept
Reference line strength.
Definition: absorptionlines.h:348
Linefunctions::InternalData::N
Eigen::VectorXcd N
Definition: linefunctions.h:521
Absorption::Lines::NumLines
Index NumLines() const noexcept
Number of lines.
Definition: absorptionlines.h:852
Linefunctions::apply_linemixing_scaling_and_mirroring
void apply_linemixing_scaling_and_mirroring(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, const Eigen::Ref< Eigen::VectorXcd > Fm, const Eigen::Ref< Eigen::MatrixXcd > dFm, const LineShape::Output &lso, const bool with_mirroring, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const LineShape::Output &dT={0, 0, 0, 0, 0, 0, 0, 0, 0}, const LineShape::Output &dVMR={0, 0, 0, 0, 0, 0, 0, 0, 0})
Applies line mixing scaling to already set lineshape and line mirror.
Definition: linefunctions.cc:420
Absorption::id_in_line
bool id_in_line(const Lines &band, const QuantumIdentifier &id, size_t line_index)
Checks if the external quantum identifier match a line's ID.
Definition: absorptionlines.cc:2840
Linefunctions::set_cross_section_of_band
void set_cross_section_of_band(InternalData &scratch, InternalData &sum, const ConstVectorView f_grid, const AbsorptionLines &band, const ArrayOfRetrievalQuantity &derivatives_data, const ArrayOfIndex &derivatives_data_active, const Vector &vmrs, const EnergyLevelMap &nlte, const Numeric &P, const Numeric &T, const Numeric &isot_ratio, const Numeric &H, const Numeric &DC, const Numeric &dDCdT, const Numeric &QT, const Numeric &dQTdT, const Numeric &QT0, const bool no_negatives=false, const bool zeeman=false, const Zeeman::Polarization zeeman_polarization=Zeeman::Polarization::Pi)
Computes the cross-section of an absorption band.
Definition: linefunctions.cc:1291
Linefunctions::InternalData::SetZero
void SetZero()
Definition: linefunctions.h:547
RetrievalQuantity::QuantumIdentity
const QuantumIdentifier & QuantumIdentity() const
Returns the identity of this Jacobian.
Definition: jacobian.h:311
Absorption::Lines::T0
Numeric T0() const noexcept
Returns reference temperature.
Definition: absorptionlines.h:1311
Absorption::Lines::ZeemanSplitting
Numeric ZeemanSplitting(size_t k, Zeeman::Polarization type, Index i) const noexcept
Returns the splitting of a Zeeman split line.
Definition: absorptionlines.h:947
dabsorption_nlte_rate_dT
Numeric dabsorption_nlte_rate_dT(const Numeric &gamma, const Numeric &T, const Numeric &F0, const Numeric &El, const Numeric &Eu, const Numeric &r_upp, const Numeric &r_low)
Computes temperature derivatives of (r_low - r_upp * gamma) / (1 - gamma)
Definition: linescaling.cc:201
linefunctions.h
Stuff related to lineshape functions.
LineShape::Type::DP
@ DP
oem::Matrix
invlib::Matrix< ArtsMatrix > Matrix
invlib wrapper type for ARTS matrices.
Definition: oem.h:34
Linefunctions::ExpectedDataSize
constexpr Index ExpectedDataSize()
Size required for data buffer.
Definition: linefunctions.h:42
abs
#define abs(x)
Definition: legacy_continua.cc:20626
Absorption::NormalizationType
NormalizationType
Describes the type of normalization line effects.
Definition: absorptionlines.h:97
Absorption::Lines::CutoffFreq
Numeric CutoffFreq(size_t k) const noexcept
Returns cutoff frequency or maximum value.
Definition: absorptionlines.h:1281
Conversion::freq2kaycm
constexpr Numeric freq2kaycm(T x)
Definition: constants.h:387
Absorption::Lines::Normalization
NormalizationType Normalization() const noexcept
Returns normalization style.
Definition: absorptionlines.h:1102
Linefunctions::lte_linestrength
Numeric lte_linestrength(Numeric S0, Numeric E0, Numeric F0, Numeric QT0, Numeric T0, Numeric QT, Numeric T)
Gets the local thermodynamic equilibrium line strength.
Definition: linefunctions.cc:618
Linefunctions::InternalData::Fc
Eigen::Matrix< Complex, 1, 1 > Fc
Definition: linefunctions.h:525
pCqSDHC_to_arts_G2_deriv
constexpr Complex pCqSDHC_to_arts_G2_deriv(Complex x) noexcept
Conversion from CGS-style lineshape to ARTS G2 derivative.
Definition: linefunctions.cc:68
dboltzman_ratio_dT
Numeric dboltzman_ratio_dT(const Numeric &boltzmann_ratio, const Numeric &T, const Numeric &E0)
Computes temperature derivatives exp(E0/k (T - T0) / (T * T0))
Definition: linescaling.cc:166
Linefunctions::InternalData::dNc
Eigen::Matrix< Complex, 1, Eigen::Dynamic > dNc
Definition: linefunctions.h:528
is_pressure_broadening_FVC
bool is_pressure_broadening_FVC(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a FVC derivative.
data
G0 G2 FVC Y DV Numeric Numeric Numeric Zeeman LowerQuantumNumbers void * data
Definition: arts_api_classes.cc:232
Linefunctions::InternalData
Definition: linefunctions.h:518
Linefunctions::apply_lineshapemodel_jacobian_scaling
void apply_lineshapemodel_jacobian_scaling(Eigen::Ref< Eigen::MatrixXcd > dF, const AbsorptionLines &band, const Index &line_ind, const ArrayOfRetrievalQuantity &derivatives_data, const ArrayOfIndex &derivatives_data_position, const Numeric &T, const Numeric &P, const Vector &vmrs)
Applies the line-by-line pressure broadening jacobian for the matching lines.
Definition: linefunctions.cc:787
Absorption::id_in_line_lower
bool id_in_line_lower(const Lines &band, const QuantumIdentifier &id, size_t line_index)
Checks if the external quantum identifier match a line's ID.
Definition: absorptionlines.cc:2904
is_pressure_broadening_DV
bool is_pressure_broadening_DV(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a DV derivative.
do_temperature_jacobian
bool do_temperature_jacobian(const ArrayOfRetrievalQuantity &js) noexcept
Returns if the array wants the temperature derivative.
Definition: jacobian.cc:1279
Linefunctions::set_doppler
void set_doppler(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::Matrix< Complex, Eigen::Dynamic, ExpectedDataSize()>> data, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &zeeman_df, const Numeric &magnetic_magnitude, const Numeric &F0_noshift, const Numeric &GD_div_F0, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const Numeric &dGD_div_F0_dT=0.0)
Sets the Doppler line shape.
Definition: linefunctions.cc:375
Absorption::SingleLine
Computations and data for a single absorption line.
Definition: absorptionlines.h:246
Absorption::id_in_line_upper
bool id_in_line_upper(const Lines &band, const QuantumIdentifier &id, size_t line_index)
Checks if the external quantum identifier match a line's ID.
Definition: absorptionlines.cc:2877
Absorption::NormalizationType::None
@ None
Complex
std::complex< Numeric > Complex
Definition: complex.h:33
dabsorption_nlte_rate_dF0
Numeric dabsorption_nlte_rate_dF0(const Numeric &gamma, const Numeric &T, const Numeric &r_upp, const Numeric &r_low)
Computes frequency derivative of (r_low - r_upp * gamma) / (1 - gamma)
Definition: linescaling.cc:228
Linefunctions::apply_VVW_scaling
void apply_VVW_scaling(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &F0, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex())
Applies Van Vleck and Weiskopf normalization to already set line shape.
Definition: linefunctions.cc:580
Absorption::Lines::ShapeParameter_dInternal
Numeric ShapeParameter_dInternal(size_t k, Numeric T, Numeric P, const Vector &vmrs, const RetrievalQuantity &derivative) const noexcept
Line shape parameter internal derivative.
Definition: absorptionlines.cc:131
absorption_nlte_ratio
Numeric absorption_nlte_ratio(const Numeric &gamma, const Numeric &r_upp, const Numeric &r_low)
Computes (r_low - r_upp * gamma) / (1 - gamma)
Definition: linescaling.cc:195
Absorption::Lines::CutoffFreqMinus
Numeric CutoffFreqMinus(size_t k, Numeric fmean) const noexcept
Returns negative cutoff frequency or lowest value.
Definition: absorptionlines.h:1298
sqrt
Numeric sqrt(const Rational r)
Square root.
Definition: rational.h:620
Linefunctions::set_voigt
void set_voigt(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::Matrix< Complex, Eigen::Dynamic, ExpectedDataSize()>> data, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &zeeman_df, const Numeric &magnetic_magnitude, const Numeric &F0_noshift, const Numeric &GD_div_F0, const LineShape::Output &lso, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const Numeric &dGD_div_F0_dT=0.0, const LineShape::Output &dT={0, 0, 0, 0, 0, 0, 0, 0, 0}, const LineShape::Output &dVMR={0, 0, 0, 0, 0, 0, 0, 0, 0})
Sets the Voigt line shape.
Definition: linefunctions.cc:298
is_pressure_broadening_D0
bool is_pressure_broadening_D0(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a D0 derivative.
is_pressure_broadening_G0
bool is_pressure_broadening_G0(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a G0 derivative.
Array
This can be used to make arrays out of anything.
Definition: array.h:108
Absorption::nelem
Index nelem(const Lines &l)
Number of lines.
Definition: absorptionlines.h:1820
Constant::pow4
constexpr T pow4(T x)
power of four
Definition: constants.h:76
MapToEigen
ComplexConstMatrixViewMap MapToEigen(const ConstComplexMatrixView &A)
Definition: complex.cc:1655
dboltzman_ratio_dT_div_boltzmann_ratio
constexpr Numeric dboltzman_ratio_dT_div_boltzmann_ratio(Numeric T, Numeric E0)
Computes temperature derivatives exp(E0/k (T - T0) / (T * T0)) / exp(E0/c (T - T0) / (T * T0))
Definition: linescaling.h:159
Absorption::CutoffType::None
@ None
dw
constexpr Complex dw(Complex z, Complex w) noexcept
The Faddeeva function partial derivative.
Definition: linefunctions.cc:45
is_magnetic_parameter
bool is_magnetic_parameter(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a magnetic parameter.
Definition: jacobian.cc:1128
Conversion::hitran2arts_linestrength
constexpr T hitran2arts_linestrength(T x)
Definition: constants.h:462
Absorption::Lines::Mirroring
MirroringType Mirroring() const noexcept
Returns mirroring style.
Definition: absorptionlines.h:1076
Linefunctions::InternalData::dFc
Eigen::Matrix< Complex, 1, Eigen::Dynamic > dFc
Definition: linefunctions.h:527
Absorption::Lines::ShapeParameters_dVMR
LineShape::Output ShapeParameters_dVMR(size_t k, Numeric T, Numeric P, const QuantumIdentifier &vmr_qid) const noexcept
Line shape parameters vmr derivative.
Definition: absorptionlines.cc:102
Linefunctions::apply_linestrength_scaling_by_lte
void apply_linestrength_scaling_by_lte(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::VectorXcd > N, Eigen::Ref< Eigen::MatrixXcd > dN, const Absorption::SingleLine &line, const Numeric &T, const Numeric &T0, const Numeric &isotopic_ratio, const Numeric &QT, const Numeric &QT0, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const Numeric &dQT_dT=0.0)
Applies linestrength to already set line shape by LTE population type.
Definition: linefunctions.cc:632
Absorption::SingleLine::E0
Numeric E0() const noexcept
Lower level energy.
Definition: absorptionlines.h:345
Absorption::Lines::ShapeParameters
LineShape::Output ShapeParameters(size_t k, Numeric T, Numeric P, const Vector &vmrs) const noexcept
Line shape parameters.
Definition: absorptionlines.cc:68
LineShape::mirroredOutput
constexpr Output mirroredOutput(Output x) noexcept
Output to be used by mirroring calls.
Definition: lineshapemodel.h:895
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:33
Linefunctions::InternalData::Nc
Eigen::Matrix< Complex, 1, 1 > Nc
Definition: linefunctions.h:526
boltzman_ratio
Numeric boltzman_ratio(const Numeric &T, const Numeric &T0, const Numeric &E0)
Computes exp(E0/c (T - T0) / (T * T0))
Definition: linescaling.cc:159
Linefunctions::set_htp
void set_htp(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &zeeman_df, const Numeric &magnetic_magnitude, const Numeric &F0_noshift, const Numeric &GD_div_F0, const LineShape::Output &lso, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const Numeric &dGD_div_F0_dT=0.0, const LineShape::Output &dT={0, 0, 0, 0, 0, 0, 0, 0, 0}, const LineShape::Output &dVMR={0, 0, 0, 0, 0, 0, 0, 0, 0})
Sets the HTP line shape.
Definition: linefunctions.cc:931
Linefunctions::dDopplerConstant_dT
Numeric dDopplerConstant_dT(const Numeric &T, const Numeric &dc)
Returns the temperature derivative of the frequency-independent part of the Doppler broadening.
Definition: linefunctions.cc:811
pCqSDHC_to_arts
constexpr Complex pCqSDHC_to_arts(Complex x) noexcept
Conversion from CGS-style lineshape to ARTS.
Definition: linefunctions.cc:50
LineShape::Output::Y
Numeric Y
Definition: lineshapemodel.h:882
stimulated_relative_emission
Numeric stimulated_relative_emission(const Numeric &gamma, const Numeric &gamma_ref)
Computes (1 - gamma) / (1 - gamma_ref)
Definition: linescaling.cc:128
pCqSDHC_to_arts_freq_deriv
constexpr Complex pCqSDHC_to_arts_freq_deriv(Complex x) noexcept
Conversion from CGS-style lineshape to ARTS for frequncy derivatives.
Definition: linefunctions.cc:59
Absorption::MirroringType
MirroringType
Describes the type of mirroring line effects.
Definition: absorptionlines.h:49
ARTS::Var::f_grid
Vector f_grid(Workspace &ws) noexcept
Definition: autoarts.h:3449
Absorption::Lines
Definition: absorptionlines.h:547
EnergyLevelMap
Definition: energylevelmap.h:60
Absorption::MirroringType::None
@ None
Absorption::Lines::F0
Numeric F0(size_t k) const noexcept
Central frequency.
Definition: absorptionlines.h:970
is_pressure_broadening_G
bool is_pressure_broadening_G(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a G derivative.
LineShape::Type
Type
Type of line shape to compute.
Definition: lineshapemodel.h:788
Linefunctions::InternalData::dN
Eigen::MatrixXcd dN
Definition: linefunctions.h:523
pCqSDHC_to_arts_D2_deriv
constexpr Complex pCqSDHC_to_arts_D2_deriv(Complex x) noexcept
Conversion from CGS-style lineshape to ARTS D2 derivative.
Definition: linefunctions.cc:78
is_pressure_broadening_D2
bool is_pressure_broadening_D2(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a D2 derivative.
Absorption::Lines::Population
PopulationType Population() const noexcept
Returns population style.
Definition: absorptionlines.h:1152
LineShape::Output
Main output of Model.
Definition: lineshapemodel.h:875
Absorption::Lines::ShapeParameters_dT
LineShape::Output ShapeParameters_dT(size_t k, Numeric T, Numeric P, const Vector &vmrs) const noexcept
Line shape parameters temperature derivatives.
Definition: absorptionlines.cc:76
Absorption::Lines::Line
SingleLine & Line(Index) noexcept
Returns a single line.
Definition: absorptionlines.cc:2709
N
#define N
Definition: rng.cc:164
LineShape::Output::G
Numeric G
Definition: lineshapemodel.h:883
Absorption::Lines::g_upp
Numeric g_upp(size_t k) const noexcept
Upper level statistical weight.
Definition: absorptionlines.h:1066
Absorption::Lines::ZeemanStrength
Numeric ZeemanStrength(size_t k, Zeeman::Polarization type, Index i) const noexcept
Returns the strength of a Zeeman split line.
Definition: absorptionlines.h:929
dabsorption_nlte_rate_dTl
Numeric dabsorption_nlte_rate_dTl(const Numeric &gamma, const Numeric &T, const Numeric &Tl, const Numeric &El, const Numeric &r_low)
Computes lower state temperature derivative of (r_low - r_upp * gamma) / (1 - gamma)
Definition: linescaling.cc:239
dstimulated_relative_emission_dT
Numeric dstimulated_relative_emission_dT(const Numeric &gamma, const Numeric &gamma_ref, const Numeric &F0, const Numeric &T)
Computes temperature derivative of (1 - gamma) / (1 - gamma_ref)
Definition: linescaling.cc:133
dstimulated_relative_emission_dF0
Numeric dstimulated_relative_emission_dF0(const Numeric &gamma, const Numeric &gamma_ref, const Numeric &T, const Numeric &T0)
Computes frequency derivative of (1 - gamma) / (1 - gamma_ref)
Definition: linescaling.cc:144
LineShape::Output::DV
Numeric DV
Definition: lineshapemodel.h:884
Constant::pow3
constexpr T pow3(T x)
power of three
Definition: constants.h:70
is_pressure_broadening_G2
bool is_pressure_broadening_G2(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a G0 derivative.
Linefunctions::DopplerConstant
Numeric DopplerConstant(Numeric T, Numeric mass)
Returns the frequency-independent part of the Doppler broadening.
Definition: linefunctions.cc:807
constants.h
Constants of physical expressions as constexpr.
Linefunctions::find_cutoff_ranges
void find_cutoff_ranges(Index &start_cutoff, Index &nelem_cutoff, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &fmin, const Numeric &fmax)
Sets cutoff frequency indices.
Definition: linefunctions.cc:816
Zeeman::start
constexpr Rational start(Rational Ju, Rational Jl, Polarization type) noexcept
Gives the lowest M for a polarization type of this transition.
Definition: zeemandata.h:77
Linefunctions::apply_linestrength_scaling_by_vibrational_nlte
void apply_linestrength_scaling_by_vibrational_nlte(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::VectorXcd > N, Eigen::Ref< Eigen::MatrixXcd > dN, const Absorption::SingleLine &line, const Numeric &T, const Numeric &T0, const Numeric &Tu, const Numeric &Tl, const Numeric &Evu, const Numeric &Evl, const Numeric &isotopic_ratio, const Numeric &QT, const Numeric &QT0, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const Numeric &dQT_dT=0.0)
Applies linestrength to already set line shape by vibrational level temperatures.
Definition: linefunctions.cc:684
linescaling.h
Constains various line scaling functions.
Absorption::Lines::QuantumIdentity
const QuantumIdentifier & QuantumIdentity() const noexcept
Returns identity status.
Definition: absorptionlines.h:1381
dabsorption_nlte_rate_dTu
Numeric dabsorption_nlte_rate_dTu(const Numeric &gamma, const Numeric &T, const Numeric &Tu, const Numeric &Eu, const Numeric &r_upp)
Computes upper state temperature derivative of (r_low - r_upp * gamma) / (1 - gamma)
Definition: linescaling.cc:252
Absorption::relaxationtype_relmat
bool relaxationtype_relmat(PopulationType in)
Definition: absorptionlines.h:201
Absorption::Lines::ZeemanCount
Index ZeemanCount(size_t k, Zeeman::Polarization type) const noexcept
Returns the number of Zeeman split lines.
Definition: absorptionlines.h:911
Linefunctions::set_lineshape
void set_lineshape(Eigen::Ref< Eigen::VectorXcd > F, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Absorption::SingleLine &line, const Numeric &temperature, const Numeric &zeeman_df, const Numeric &magnetic_magnitude, const Numeric &doppler_constant, const LineShape::Output &lso, const LineShape::Type lineshape_type, const Absorption::MirroringType mirroring_type, const Absorption::NormalizationType norm_type)
Sets the lineshape normalized to unity.
Definition: linefunctions.cc:87
Constant::pow2
constexpr T pow2(T x)
power of two
Definition: constants.h:64
LineShape::Output::G0
Numeric G0
Definition: lineshapemodel.h:876
ARTS::Var::x
Vector x(Workspace &ws) noexcept
Definition: autoarts.h:7346
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:39
Absorption::Lines::LinemixingLimit
Numeric LinemixingLimit() const noexcept
Returns line mixing limit.
Definition: absorptionlines.h:1331
F0
G0 G2 FVC Y DV F0
Definition: arts_api_classes.cc:218
RetrievalQuantity
Deals with internal derivatives, Jacobian definition, and OEM calculations.
Definition: jacobian.h:120
Linefunctions::apply_linestrength_from_nlte_level_distributions
void apply_linestrength_from_nlte_level_distributions(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::VectorXcd > N, Eigen::Ref< Eigen::MatrixXcd > dN, const Numeric &r1, const Numeric &r2, const Numeric &g1, const Numeric &g2, const Numeric &A21, const Numeric &F0, const Numeric &T, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex())
Applies non-lte linestrength to already set line shape.
Definition: linefunctions.cc:844
LineShape::Output::D0
Numeric D0
Definition: lineshapemodel.h:877
Linefunctions::set_lorentz
void set_lorentz(Eigen::Ref< Eigen::VectorXcd > F, Eigen::Ref< Eigen::MatrixXcd > dF, Eigen::Ref< Eigen::Matrix< Complex, Eigen::Dynamic, ExpectedDataSize()>> data, const Eigen::Ref< const Eigen::VectorXd > f_grid, const Numeric &zeeman_df, const Numeric &magnetic_magnitude, const Numeric &F0_noshift, const LineShape::Output &lso, const AbsorptionLines &band=AbsorptionLines(), const Index &line_ind=0, const ArrayOfRetrievalQuantity &derivatives_data=ArrayOfRetrievalQuantity(), const ArrayOfIndex &derivatives_data_position=ArrayOfIndex(), const LineShape::Output &dT={0, 0, 0, 0, 0, 0, 0, 0, 0}, const LineShape::Output &dVMR={0, 0, 0, 0, 0, 0, 0, 0, 0})
Sets the Lorentz line shape.
Definition: linefunctions.cc:234
oem::LM
invlib::LevenbergMarquardt< Numeric, CovarianceMatrix, Std > LM
Levenberg-Marquardt (LM) optimization using normed ARTS QR solver.
Definition: oem.h:173
Vector
The Vector class.
Definition: matpackI.h:860
is_frequency_parameter
bool is_frequency_parameter(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a frequency parameter in propagation matrix calculations.
Definition: jacobian.cc:1120
Absorption::Lines::A
Numeric A(size_t k) const noexcept
Einstein spontaneous emission.
Definition: absorptionlines.h:1038
stimulated_emission
Numeric stimulated_emission(Numeric T, Numeric F0)
Computes exp(-hf/kT)
Definition: linescaling.cc:110
ConstVectorView
A constant view of a Vector.
Definition: matpackI.h:476
is_pressure_broadening_ETA
bool is_pressure_broadening_ETA(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a ETA derivative.
Array::nelem
Index nelem() const
Number of elements.
Definition: array.h:195
is_lineshape_parameter
bool is_lineshape_parameter(const RetrievalQuantity &t) noexcept
Returns if the Retrieval quantity is a G0, D0, G2, D2, ETA, FVC, Y, G, DV derivative.
Definition: jacobian.cc:1187
Absorption::PopulationType::ByLTE
@ ByLTE
Zeeman::Polarization
Polarization
Zeeman polarization selection.
Definition: zeemandata.h:43
Linefunctions::InternalData::F
Eigen::VectorXcd F
Definition: linefunctions.h:520
Linefunctions::InternalData::data
Eigen::Matrix< Complex, Eigen::Dynamic, Linefunctions::ExpectedDataSize()> data
Definition: linefunctions.h:530
E0
G0 G2 FVC Y DV Numeric E0
Definition: arts_api_classes.cc:220
LineShape::si2cgs
constexpr Output si2cgs(Output x) noexcept
Output turned from SI to CGS units.
Definition: lineshapemodel.h:905
LineShape::vmrs
Vector vmrs(const ConstVectorView &atmospheric_vmrs, const ArrayOfArrayOfSpeciesTag &atmospheric_species, const QuantumIdentifier &self, const ArrayOfSpeciesTag &lineshape_species, bool self_in_list, bool bath_in_list, Type type)
Returns a VMR vector for this model's main calculations.
Definition: lineshapemodel.cc:474