ARTS
2.3.1285(git:92a29ea9-dirty)
|
Line mixing calculation implementation. More...
Go to the source code of this file.
Classes | |
class | AdiabaticFactor |
Adiabatic factor computations. More... | |
class | BasisRate |
Basis rate of transitions. More... | |
struct | OffDiagonalElementOutput |
Struct to help keep position of the two matched outputs clear. More... | |
struct | SecondOrderLineMixingCoeffs |
To keep track of (y0 + y1 * (T0 / T - 1.)) * pow(T0 / T, n) More... | |
Namespaces | |
Molecule | |
Molecular constant. | |
Molecule::O2_66 | |
O2-66 constants. | |
Molecule::CO2_626 | |
CO-626 constants. | |
OffDiagonalElement | |
Methods to compute off diagonal elements. | |
Enumerations | |
enum | Molecule::Name { Molecule::Name::O2_66, Molecule::Name::CO2_626 } |
Name of molecules. More... | |
enum | OffDiagonalElement::Type { OffDiagonalElement::Type::CO2_IR, OffDiagonalElement::Type::O2_66_MW } |
Type of off diagonal element computations. More... | |
enum | RedPoleType { RedPoleType::ElectricRoVibDipole, RedPoleType::MagneticQuadrapole } |
Type of reduced dipole. More... | |
Functions | |
template<class T > | |
Numeric | rotational_energy_hund_b_molecule (T J, T N, T J2, Numeric B, Numeric D, Numeric H, Numeric gamma, Numeric gamma_D, Numeric gamma_H, Numeric lambda, Numeric lambda_D, Numeric lambda_H) |
Compute the rotational energy of a Hund b case molecule. More... | |
template<class T > | |
Numeric | Molecule::O2_66::hamiltonian_freq (T J, int dcol=0, int drow=0) |
Hamiltonian frequency. More... | |
template<class T > | |
constexpr Numeric | Molecule::CO2_626::hamiltonian_freq (T J) |
Hamiltonian frequency. More... | |
OffDiagonalElementOutput | OffDiagonalElement::CO2_IR (const Rational &Jku, const Rational &Jju, const Rational &Jkl, const Rational &Jjl, const Rational &l2ku, const Rational &l2ju, const Rational &l2kl, const Rational &l2jl, const Numeric &j_rho, const Numeric &k_rho, const BasisRate &br, const AdiabaticFactor &af, const Numeric &T, const Numeric &B0, const Numeric &main_mass, const Numeric &collider_mass) |
CO2 IR off diagonal element computer. More... | |
OffDiagonalElementOutput | OffDiagonalElement::O2_66_MW (const Rational &J1u, const Rational &N1u, const Rational &J1l, const Rational &N1l, const Rational &J2u, const Rational &N2u, const Rational &J2l, const Rational &N2l, const Numeric &rho1, const Numeric &rho2, const Numeric &T, const Numeric &collider_mass) |
O2-66 MW off diagonal element computer. More... | |
Matrix | hartmann_ecs_interface (const AbsorptionLines &abs_lines, const ArrayOfSpeciesTag &collider_species, const Vector &collider_species_vmr, const SpeciesAuxData &partition_functions, const Numeric &T, const Index &size) |
Energy corrected sudden relaxation matrix using Hartmann's method. More... | |
Vector | population_density_vector (const AbsorptionLines &abs_lines, const SpeciesAuxData &partition_functions, const Numeric &T) |
Compute the population density. More... | |
Vector | dipole_vector (const AbsorptionLines &abs_lines, const SpeciesAuxData &partition_functions) |
Dipole vector. More... | |
Vector | reduced_dipole_vector (const AbsorptionLines &abs_lines, const RedPoleType type) |
Reduced dipole vector. More... | |
Vector | rosenkranz_scaling_second_order (const AbsorptionLines &abs_lines, const Matrix &W, const Vector &d0) |
Computes G for Rosenkranz's line mixing coefficients. More... | |
Vector | rosenkranz_shifting_second_order (const AbsorptionLines &abs_lines, const Matrix &W) |
Computes DV for Rosenkranz's line mixing coefficients. More... | |
Vector | rosenkranz_first_order (const AbsorptionLines &abs_lines, const Matrix &W, const Vector &d0) |
Computes Y for Rosenkranz's line mixing coefficients. More... | |
SecondOrderLineMixingCoeffs | compute_2nd_order_lm_coeff (ConstVectorView y, ConstVectorView x, const Numeric exp, const Numeric x0) |
Fit to a second order line mixing formula the input. More... | |
ComplexVector | equivalent_linestrengths (const Vector &population, const Vector &dipole, const Eigen::ComplexEigenSolver< Eigen::MatrixXcd > &M) |
Equivalent line strengths. More... | |
Numeric | total_linestrengths (const Vector &population, const Vector &dipole) |
Sum of line strengths. More... | |
Matrix | CO2_ir_training (const ArrayOfRational &Ji, const ArrayOfRational &Jf, const ArrayOfRational &l2i, const ArrayOfRational &l2f, const Vector &F0, const Vector &d, const Vector &rho, const Vector &gamma, const Numeric &T) |
CO2 IR training algorithm for linearization. More... | |
void | relmatInAir (Matrix &relmat, const AbsorptionLines &abs_lines, const SpeciesAuxData &partition_functions, const Index &wigner_initialized, const Numeric &temperature) |
Compute the relaxation matrix in air mixture. More... | |
Variables | |
constexpr Numeric | Molecule::O2_66::g_S = 2.002084 |
constexpr Numeric | Molecule::O2_66::ge_l = 2.77e-3 |
constexpr Numeric | Molecule::O2_66::g_r = -1.16e-4 |
constexpr Numeric | Molecule::O2_66::B = 43100.44276e6 |
constexpr Numeric | Molecule::O2_66::D = 145.1271e3 |
constexpr Numeric | Molecule::O2_66::H = 49e-3 |
constexpr Numeric | Molecule::O2_66::lambda = 59501.3438e6 |
constexpr Numeric | Molecule::O2_66::lambda_D = 58.3680e3 |
constexpr Numeric | Molecule::O2_66::lambda_H = 290.8e-3 |
constexpr Numeric | Molecule::O2_66::gamma = -252.58634e6 |
constexpr Numeric | Molecule::O2_66::gamma_D = -243.42 |
constexpr Numeric | Molecule::O2_66::gamma_H = -1.46e-3 |
constexpr Numeric | Molecule::O2_66::mass = 31.989830 |
constexpr Numeric | Molecule::CO2_626::B = Conversion::kaycm2freq(0.39021) |
Line mixing calculation implementation.
This file contains only experimental code to test the F90-routines. It might evolve to replace it at some point but that is beyond present ability
Definition in file linemixing.h.
|
strong |
Type of reduced dipole.
Enumerator | |
---|---|
ElectricRoVibDipole | |
MagneticQuadrapole |
Definition at line 419 of file linemixing.h.
Matrix CO2_ir_training | ( | const ArrayOfRational & | Ji, |
const ArrayOfRational & | Jf, | ||
const ArrayOfRational & | l2i, | ||
const ArrayOfRational & | l2f, | ||
const Vector & | F0, | ||
const Vector & | d, | ||
const Vector & | rho, | ||
const Vector & | gamma, | ||
const Numeric & | T | ||
) |
CO2 IR training algorithm for linearization.
[in] | Ji | J init for all lines |
[in] | Jf | J final for all lines |
[in] | l2i | l2 init for all lines |
[in] | l2f | l2 final for all lines |
[in] | F0 | Central frequency |
[in] | d | Dipole for each line |
[in] | rho | The population density for each line |
[in] | gamma | Pressure broadening for each line |
[in] | T | Temperature |
Definition at line 665 of file linemixing.cc.
References i, joker, linearized_relaxation_matrix(), lsf(), M, n, Array< base >::nelem(), and params.
SecondOrderLineMixingCoeffs compute_2nd_order_lm_coeff | ( | ConstVectorView | y, |
ConstVectorView | x, | ||
const Numeric | exp, | ||
const Numeric | x0 | ||
) |
Fit to a second order line mixing formula the input.
Finds best fit [c0, c1] of y(x) = (c0 + c1 * (x0 / x - 1.)) * pow(x0 / x, exp)
[in] | y | Y-axis |
[in] | x | X-axis |
[in] | exp | Exponent |
[in] | x0 | Zero-value of x0 |
Definition at line 835 of file linemixing.cc.
References i, n, and ConstVectorView::nelem().
Vector dipole_vector | ( | const AbsorptionLines & | abs_lines, |
const SpeciesAuxData & | partition_functions | ||
) |
Dipole vector.
[in] | abs_lines | One band of lines |
[in] | partition_functions | Method to compute the partition function |
Definition at line 371 of file linemixing.cc.
References n, Absorption::Lines::NumLines(), population_density_vector(), and Absorption::Lines::T0().
Referenced by abs_lines_per_speciesSetLineMixingFromRelmat(), abs_xsec_per_speciesAddLineMixedLines(), abs_xsec_per_speciesAddLineMixedLinesInAir(), and hartmann_ecs_interface().
ComplexVector equivalent_linestrengths | ( | const Vector & | population, |
const Vector & | dipole, | ||
const Eigen::ComplexEigenSolver< Eigen::MatrixXcd > & | M | ||
) |
Equivalent line strengths.
[in] | population | The population density for each line |
[in] | dipole | Dipole for each line |
[in] | M | Solver |
Definition at line 884 of file linemixing.cc.
References Molecule::O2_66::B, i, n, and ConstVectorView::nelem().
Referenced by abs_xsec_per_speciesAddLineMixedLines(), and abs_xsec_per_speciesAddLineMixedLinesInAir().
Matrix hartmann_ecs_interface | ( | const AbsorptionLines & | abs_lines, |
const ArrayOfSpeciesTag & | collider_species, | ||
const Vector & | collider_species_vmr, | ||
const SpeciesAuxData & | partition_functions, | ||
const Numeric & | T, | ||
const Index & | size | ||
) |
Energy corrected sudden relaxation matrix using Hartmann's method.
[in] | abs_lines | One band of lines |
[in] | collider_species | Species tag of collider |
[in] | collider_species_vmr | VMR of collider |
[in] | partition_functions | Method to compute the partition function |
[in] | T | Temperature |
[in] | size | Number of elements |
Definition at line 501 of file linemixing.cc.
References dipole_vector(), hartmann_relaxation_matrix(), and population_density_vector().
Referenced by relmatInAir().
Vector population_density_vector | ( | const AbsorptionLines & | abs_lines, |
const SpeciesAuxData & | partition_functions, | ||
const Numeric & | T | ||
) |
Compute the population density.
[in] | abs_lines | One band of lines |
[in] | partition_functions | Method to compute the partition function |
[in] | T | Temperature |
Definition at line 342 of file linemixing.cc.
References SpeciesAuxData::getParam(), SpeciesAuxData::getParamType(), n, Absorption::Lines::NumLines(), Absorption::Lines::QuantumIdentity(), and single_partition_function().
Referenced by abs_xsec_per_speciesAddLineMixedLines(), abs_xsec_per_speciesAddLineMixedLinesInAir(), dipole_vector(), and hartmann_ecs_interface().
Vector reduced_dipole_vector | ( | const AbsorptionLines & | band, |
const RedPoleType | type | ||
) |
Reduced dipole vector.
[in] | abs_lines | One band of lines |
[in] | type | Type of reduced dipole |
Computes the dipole moment
abs_lines | all absorption lines of the band |
Definition at line 401 of file linemixing.cc.
References ElectricRoVibDipole, i, J, l2, Absorption::Lines::LowerQuantumNumber(), MagneticQuadrapole, N, n, Absorption::Lines::NumLines(), Absorption::reduced_magnetic_quadrapole(), Absorption::reduced_rovibrational_dipole(), and Absorption::Lines::UpperQuantumNumber().
Referenced by normalize_relaxation_matrix().
void relmatInAir | ( | Matrix & | relmat, |
const AbsorptionLines & | abs_lines, | ||
const SpeciesAuxData & | partition_functions, | ||
const Index & | wigner_initialized, | ||
const Numeric & | temperature | ||
) |
Compute the relaxation matrix in air mixture.
[out] | relmat | Relaxation matrix |
[in] | abs_lines | Absorption band |
[in] | partition_functions | Partition functions |
[in] | wigner_initialized | Indication of the Wigner state |
[in] | temperature | Atmospheric temperature |
Definition at line 915 of file linemixing.cc.
References Absorption::Lines::BroadeningSpecies(), Absorption::ByRelmatHartmannLTE, hartmann_ecs_interface(), SpeciesTag::Isotopologue(), Absorption::Lines::Isotopologue(), Array< base >::nelem(), Absorption::Lines::Population(), SpeciesTag::Species(), and Absorption::Lines::Species().
Referenced by abs_xsec_per_speciesAddLineMixedLinesInAir().
Vector rosenkranz_first_order | ( | const AbsorptionLines & | abs_lines, |
const Matrix & | W, | ||
const Vector & | d0 | ||
) |
Computes Y for Rosenkranz's line mixing coefficients.
[in] | abs_lines | One band of lines |
[in] | W | Relaxation Matrix |
[in] | d0 | Dipole vector |
Definition at line 431 of file linemixing.cc.
References Absorption::Lines::F0(), i, n, and Absorption::Lines::NumLines().
Referenced by abs_lines_per_speciesSetLineMixingFromRelmat().
Vector rosenkranz_scaling_second_order | ( | const AbsorptionLines & | abs_lines, |
const Matrix & | W, | ||
const Vector & | d0 | ||
) |
Computes G for Rosenkranz's line mixing coefficients.
[in] | abs_lines | One band of lines |
[in] | W | Relaxation Matrix |
[in] | d0 | Dipole vector |
Definition at line 466 of file linemixing.cc.
References Absorption::Lines::F0(), i, n, Absorption::Lines::NumLines(), Constant::pow2(), and r.
Referenced by abs_lines_per_speciesSetLineMixingFromRelmat().
Vector rosenkranz_shifting_second_order | ( | const AbsorptionLines & | abs_lines, |
const Matrix & | W | ||
) |
Computes DV for Rosenkranz's line mixing coefficients.
[in] | abs_lines | One band of lines |
[in] | W | Relaxation Matrix |
Definition at line 449 of file linemixing.cc.
References Absorption::Lines::F0(), i, n, and Absorption::Lines::NumLines().
Referenced by abs_lines_per_speciesSetLineMixingFromRelmat().
|
inline |
Compute the rotational energy of a Hund b case molecule.
[in] | J | Rotational angular momentum w/ spin |
[in] | N | Rotational angular momentum w/0 spin |
[in] | J2 | Rotational angular momentum w/ spin for other level |
[in] | B | Main rotation constant |
[in] | D | Main rotation constant squared |
[in] | H | Main rotation constant cubic |
[in] | gamma | Second rotation constant offset |
[in] | gamma_D | Second rotation constant |
[in] | gamma_H | Second rotation constant squared |
[in] | lambda | Third energy constant offset |
[in] | lambda_D | Third energy constant |
[in] | lambda_H | Third energy constant squared |
Definition at line 54 of file linemixing.h.
References sqrt().
Referenced by Molecule::O2_66::hamiltonian_freq().
Sum of line strengths.
[in] | population | The population density for each line |
[in] | dipole | Dipole for each line |
Definition at line 906 of file linemixing.cc.
References i, n, and ConstVectorView::nelem().