ARTS  2.3.1285(git:92a29ea9-dirty)
Linefunctions Namespace Reference

Line functions related to line shapes and line strength. More...

Classes

class  InternalData
 

Functions

constexpr Index ExpectedDataSize ()
 Size required for data buffer. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
Numeric lte_linestrength (Numeric S0, Numeric E0, Numeric F0, Numeric QT0, Numeric T0, Numeric QT, Numeric T)
 Gets the local thermodynamic equilibrium line strength. More...
 
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. More...
 
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. More...
 
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. More...
 
Numeric DopplerConstant (Numeric T, Numeric mass)
 Returns the frequency-independent part of the Doppler broadening. More...
 
Numeric dDopplerConstant_dT (const Numeric &T, const Numeric &dc)
 Returns the temperature derivative of the frequency-independent part of the Doppler broadening. More...
 
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. More...
 
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. More...
 
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. More...
 

Detailed Description

Line functions related to line shapes and line strength.

Function Documentation

◆ apply_linemixing_scaling_and_mirroring()

void Linefunctions::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.

Equation: with_mirroring: F := (1+G-iY) * F + (1+G+iY) * Fm else: F := (1+G-iY) * F

and appropriate derivatives

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in]FmMirrored lineshape. Must be right size
[in]dFmMirrored lineshape derivative. Must be right size
[in]lsoLine shape parameters
[in]with_mirroringMirror lineshape check
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dTTemperature derivatives of line shape parameters
[in]dVMRVMR derivatives of line shape parameters

Definition at line 420 of file linefunctions.cc.

References conj(), LineShape::Output::G, Array< base >::nelem(), Temperature, and LineShape::Output::Y.

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ apply_lineshapemodel_jacobian_scaling()

void Linefunctions::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.

Parameters
[in,out]dFLineshape derivative. Must be right size
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]TAtmospheric temperature
[in]PAtmospheric pressure
[in]vmrsVMRs for line shape broadeners

Definition at line 787 of file linefunctions.cc.

References Absorption::id_in_line(), is_lineshape_parameter(), Array< base >::nelem(), RetrievalQuantity::QuantumIdentity(), and Absorption::Lines::ShapeParameter_dInternal().

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ apply_linestrength_from_nlte_level_distributions()

void Linefunctions::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.

Works on ratio-inputs, meaning that the total distribution does not have to be known

Cannot support partial derivatives at this point due to ARTS not possessing its own NLTE ratio calculation agenda

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]NSource lineshape
[in,out]dNSource lineshape derivative
[in]r1Ratio of molecules at energy level 1
[in]r2Ratio of molecules at energy level 2
[in]g1Statistical weight of energy level 1
[in]g2Statistical weight of energy level 2
[in]A21Einstein coefficient for the transition from energy level 2 to energy level 1
[in]F0Central frequency without any shifts
[in]TAtmospheric temperature
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF

Definition at line 844 of file linefunctions.cc.

References F0, Array< base >::nelem(), and Constant::pow2().

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ apply_linestrength_scaling_by_lte()

void Linefunctions::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.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]NSource lineshape
[in,out]dNSource lineshape derivative
[in]lineThe absorption line
[in]TThe atmospheric temperature
[in]T0The reference temperature
[in]isotopic_ratioThe ratio of the isotopologue in the atmosphere
[in]QTPartition function at atmospheric temperature of level
[in]QT0Partition function at reference temperature
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dQT_dTTemperature derivative of QT

Definition at line 632 of file linefunctions.cc.

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ apply_linestrength_scaling_by_vibrational_nlte()

void Linefunctions::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.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]NSource lineshape
[in,out]dNSource lineshape derivative
[in]lineThe absorption line
[in]TThe atmospheric temperature
[in]T0The reference temperature
[in]TuThe upper state vibrational temperature; must be T if level is LTE
[in]TlThe lower state vibrational temperature; must be T if level is LTE
[in]EvuThe upper state vibrational energy; if set funny, yields funny results
[in]EvlThe lower state vibrational energy; if set funny, yields funny results
[in]isotopic_ratioThe ratio of the isotopologue in the atmosphere
[in]QTPartition function at atmospheric temperature of level
[in]QT0Partition function at reference temperature
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dQT_dTTemperature derivative of QT

Definition at line 684 of file linefunctions.cc.

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ apply_rosenkranz_quadratic_scaling()

void Linefunctions::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.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in]f_gridFrequency grid of computations
[in]F0Central frequency without any shifts
[in]TAtmospheric temperature at level
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF

Definition at line 481 of file linefunctions.cc.

References do_temperature_jacobian(), F0, Absorption::id_in_line(), is_frequency_parameter(), LineCenter, Array< base >::nelem(), Absorption::Lines::QuantumIdentity(), and Temperature.

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ apply_VVH_scaling()

void Linefunctions::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.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]dataBlock of allocated memory. Output nonsensical
[in]f_gridFrequency grid of computations
[in]F0Central frequency without any shifts
[in]TAtmospheric temperature at level
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF

Definition at line 531 of file linefunctions.cc.

References data, Absorption::id_in_line(), is_frequency_parameter(), LineCenter, Array< base >::nelem(), Absorption::Lines::QuantumIdentity(), and Temperature.

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ apply_VVW_scaling()

void Linefunctions::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.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in]f_gridFrequency grid of computations
[in]F0Central frequency without any shifts
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF

Definition at line 580 of file linefunctions.cc.

References F0, fac(), Absorption::id_in_line(), is_frequency_parameter(), LineCenter, Array< base >::nelem(), and Absorption::Lines::QuantumIdentity().

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ dDopplerConstant_dT()

Numeric Linefunctions::dDopplerConstant_dT ( const Numeric T,
const Numeric dc 
)

Returns the temperature derivative of the frequency-independent part of the Doppler broadening.

Parameters
[in]TAtmospheric temperature at level
[in]dcOutput of Linefunctions::DopplerConstant(T, mass)
Returns
Doppler broadening constant temperature derivative

Definition at line 811 of file linefunctions.cc.

Referenced by ExpectedDataSize(), xsec_species(), and zeeman_on_the_fly().

◆ DopplerConstant()

Numeric Linefunctions::DopplerConstant ( Numeric  T,
Numeric  mass 
)

Returns the frequency-independent part of the Doppler broadening.

Parameters
[in]TAtmospheric temperature at level
[in]massMass of molecule under consideration
Returns
Doppler broadening constant

Definition at line 807 of file linefunctions.cc.

References sqrt().

Referenced by ExpectedDataSize(), Absorption::PredefinedModel::makarov2020_o2_lines_mpm(), xsec_species(), and zeeman_on_the_fly().

◆ ExpectedDataSize()

◆ find_cutoff_ranges()

void Linefunctions::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.

Parameters
[out]start_cutoffStart pos of cutoff frequency
[out]end_cutoffEnd pos of cutoff frequency
[in]f_gridFrequency grid of computations
[in]fminMinimum frequency
[in]fmaxMaximum frequency

Definition at line 816 of file linefunctions.cc.

Referenced by ExpectedDataSize(), and set_cross_section_of_band().

◆ lte_linestrength()

Numeric Linefunctions::lte_linestrength ( Numeric  S0,
Numeric  E0,
Numeric  F0,
Numeric  QT0,
Numeric  T0,
Numeric  QT,
Numeric  T 
)

Gets the local thermodynamic equilibrium line strength.

Parameters
[in]S0Reference linestrength
[in]E0Reference line lower energy
[in]F0Central frequency without any shifts
[in]QT0Partition function at reference temperature
[in]T0Reference temperature
[in]QTPartition function at atmospheric temperature
[in]TAtmospheric temperature at level
Returns
The local thermodynamic equilibrium line strength

Definition at line 618 of file linefunctions.cc.

Referenced by ExpectedDataSize().

◆ set_cross_section_of_band()

void Linefunctions::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.

Parameters
[in,out]scratchData that is overwritten by every line
[in,out]sunData that is set to zero then added onto by every line
[in]f_gridAs WSV
[in]bandThe absorption band
[in]derivatives_dataDerivatives
[in]derivatives_data_activeDerivatives that are active
[in]vmrsThe VMRs of this band's broadening species
[in]nlteA map of NLTE energy levels
[in]PThe pressure
[in]TThe temperature
[in]isot_ratioThe band isotopic ratio
[in]HThe strength of the magnetic field
[in]DCAs per DopplerConstant
[in]dDCdTTemperature derivative of DC
[in]QTThe partition function at the temperature
[in]dQTdTTemperature derivative of QT
[in]QT0The partition function at the band reference temperature
[in]no_negativesCheck sum.F before output of any real negative values, and removes them if present
[in]zeemanAttempts adding up the fine Zeeman lines
[in]zeeman_polarizationThe polarization of Zeeman model (to know how many Zeeman lines there will be)

Definition at line 1291 of file linefunctions.cc.

References apply_linemixing_scaling_and_mirroring(), apply_lineshapemodel_jacobian_scaling(), apply_linestrength_from_nlte_level_distributions(), apply_linestrength_scaling_by_lte(), apply_linestrength_scaling_by_vibrational_nlte(), apply_rosenkranz_quadratic_scaling(), apply_VVH_scaling(), apply_VVW_scaling(), Absorption::BandFixedFrequency, Absorption::ByHITRANFullRelmat, Absorption::ByHITRANRosenkranzRelmat, Absorption::ByLTE, Absorption::ByNLTEPopulationDistribution, Absorption::ByNLTEVibrationalTemperatures, Absorption::Lines::Cutoff(), Absorption::Lines::CutoffFreq(), Absorption::Lines::CutoffFreqMinus(), data, Linefunctions::InternalData::data, Linefunctions::InternalData::datac, Linefunctions::InternalData::dF, Linefunctions::InternalData::dFc, dN, Linefunctions::InternalData::dN, Linefunctions::InternalData::dNc, do_temperature_jacobian(), do_vmr_jacobian(), LineShape::DP, Linefunctions::InternalData::F, Absorption::Lines::F0(), Absorption::Lines::F_mean(), Linefunctions::InternalData::Fc, find_cutoff_ranges(), EnergyLevelMap::get_ratio_params(), EnergyLevelMap::get_vibtemp_params(), LineShape::HTP, i, Absorption::Lines::Line(), Absorption::LineByLineOffset, Absorption::Lines::LinemixingLimit(), Absorption::Lines::LineShapeType(), Absorption::Lorentz, LineShape::LP, Absorption::Manual, MapToEigen(), LineShape::mirroredOutput(), Absorption::Lines::Mirroring(), N, Linefunctions::InternalData::N, Linefunctions::InternalData::Nc, Array< base >::nelem(), Absorption::nelem(), Absorption::None, Absorption::Lines::Normalization(), Absorption::Lines::NumLines(), Absorption::Lines::Population(), Absorption::Lines::QuantumIdentity(), Absorption::relaxationtype_relmat(), Absorption::RosenkranzQuadratic, Absorption::SameAsLineShape, LineShape::SDVP, set_doppler(), set_htp(), set_lorentz(), set_voigt(), Linefunctions::InternalData::SetZero(), Absorption::Lines::ShapeParameters(), Absorption::Lines::ShapeParameters_dT(), Absorption::Lines::ShapeParameters_dVMR(), Zeeman::start(), Absorption::Lines::T0(), LineShape::VP, Absorption::VVH, Absorption::VVW, Absorption::Lines::ZeemanCount(), Absorption::Lines::ZeemanSplitting(), and Absorption::Lines::ZeemanStrength().

Referenced by Linefunctions::InternalData::SetZero(), xsec_species(), and zeeman_on_the_fly().

◆ set_doppler()

void Linefunctions::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.

Normalization is unity.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]dataBlock of allocated memory. Output nonsensical
[in]f_gridFrequency grid of computations
[in]zeeman_dfZeeman shift parameter for the line
[in]magnetic_magnitudeAbsolute strength of the magnetic field
[in]F0_noshiftCentral frequency without any shifts
[in]GD_div_F0Frequency-independent part of the Doppler broadening
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dGD_div_F0_dTTemperature derivative of GD_div_F0

Definition at line 375 of file linefunctions.cc.

References data, F0, Absorption::id_in_line(), is_frequency_parameter(), LineCenter, Array< base >::nelem(), Absorption::Lines::QuantumIdentity(), Temperature, and VMR.

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ set_htp()

void Linefunctions::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.

Normalization is unity.

Note: We are not experienced with this line shape and cannot tell what parameters depends on what input. It is therefore likely that this function will have to be adapted in the future

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in]f_gridFrequency grid of computations
[in]zeeman_dfZeeman shift parameter for the line
[in]magnetic_magnitudeAbsolute strength of the magnetic field
[in]F0_noshiftCentral frequency without any shifts
[in]GD_div_F0Frequency-independent part of the Doppler broadening
[in]lsoLine shape parameters
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dGD_div_F0_dTTemperature derivative of GD_div_F0
[in]dTTemperature derivatives of line shape parameters
[in]dVMRVMR derivatives of line shape parameters

Definition at line 931 of file linefunctions.cc.

References abs, Conversion::freq2kaycm(), Constant::pow2(), Constant::pow3(), Constant::pow4(), and sqrt().

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ set_lineshape()

void Linefunctions::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.

No line mixing or linestrength is computed.

Parameters
[in,out]FLineshape. Must be right size
[in]f_gridFrequency grid of computations
[in]lineThe absortion line
[in]temperatureAtmospheric temperature
[in]zeeman_dfZeeman splitting coefficient
[in]magnetic_magnitudeStrength of local magnetic field
[in]lsoLine shape parameters
[in]lineshape_typeLine shape scheme
[in]mirroring_typeMirroring scheme
[in]norm_typeNormalization scheme

Definition at line 87 of file linefunctions.cc.

References apply_rosenkranz_quadratic_scaling(), apply_VVH_scaling(), apply_VVW_scaling(), data, LineShape::DP, ExpectedDataSize(), Absorption::SingleLine::F0(), LineShape::HTP, Absorption::Lorentz, LineShape::LP, Absorption::Manual, LineShape::mirroredOutput(), Absorption::None, Absorption::RosenkranzQuadratic, Absorption::SameAsLineShape, LineShape::SDVP, set_doppler(), set_htp(), set_lorentz(), set_voigt(), LineShape::VP, Absorption::VVH, and Absorption::VVW.

Referenced by ExpectedDataSize().

◆ set_lorentz()

void Linefunctions::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.

Normalization is unity.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]dataBlock of allocated memory. Output nonsensical
[in]f_gridFrequency grid of computations
[in]zeeman_dfZeeman shift parameter for the line
[in]magnetic_magnitudeAbsolute strength of the magnetic field
[in]F0_noshiftCentral frequency without any shifts
[in]lsoLine shape parameters
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dTTemperature derivatives of line shape parameters
[in]dVMRVMR derivatives of line shape parameters

Definition at line 234 of file linefunctions.cc.

References LineShape::Output::D0, data, LineShape::Output::DV, dw(), F0, LineShape::Output::G0, Absorption::id_in_line(), is_frequency_parameter(), is_magnetic_parameter(), is_pressure_broadening_D0(), is_pressure_broadening_DV(), is_pressure_broadening_G0(), LineCenter, Array< base >::nelem(), Absorption::Lines::QuantumIdentity(), Temperature, and VMR.

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().

◆ set_voigt()

void Linefunctions::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.

Normalization is unity.

Parameters
[in,out]FLineshape. Must be right size
[in,out]dFLineshape derivative. Must be right size
[in,out]dataBlock of allocated memory. Output nonsensical
[in]f_gridFrequency grid of computations
[in]zeeman_dfZeeman shift parameter for the line
[in]magnetic_magnitudeAbsolute strength of the magnetic field
[in]F0_noshiftCentral frequency without any shifts
[in]GD_div_F0Frequency-independent part of the Doppler broadening
[in]lsoLine shape parameters
[in]bandThe absorption lines
[in]line_indThe current line's ID
[in]derivatives_dataThe derivatives in dF
[in]derivatives_data_positionThe derivatives positions in dF
[in]dGD_div_F0_dTTemperature derivative of GD_div_F0
[in]dTTemperature derivatives of line shape parameters
[in]dVMRVMR derivatives of line shape parameters

Definition at line 298 of file linefunctions.cc.

References LineShape::Output::D0, data, LineShape::Output::DV, dw(), F0, fac(), LineShape::Output::G0, Absorption::id_in_line(), is_frequency_parameter(), is_magnetic_parameter(), is_pressure_broadening_D0(), is_pressure_broadening_DV(), is_pressure_broadening_G0(), LineCenter, Array< base >::nelem(), Absorption::Lines::QuantumIdentity(), Temperature, VMR, and w().

Referenced by ExpectedDataSize(), set_cross_section_of_band(), and set_lineshape().