ARTS  2.2.66
absorption.cc File Reference

Physical absorption routines. More...

#include "arts.h"
#include <map>
#include <cfloat>
#include <algorithm>
#include <cmath>
#include <cstdlib>
#include "file.h"
#include "absorption.h"
#include "math_funcs.h"
#include "messages.h"
#include "logic.h"
#include "interpolation_poly.h"
#include "global_data.h"

Go to the source code of this file.

Functions

void checkIsotopologueRatios (const ArrayOfArrayOfSpeciesTag &abs_species, const SpeciesAuxData &isoratios)
 Check that isotopologue ratios for the given species are correctly defined. More...
 
void fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData (SpeciesAuxData &sad)
 Fill SpeciesAuxData with default isotopologue ratios from species data. More...
 
void define_species_map ()
 Define the species data map. More...
 
ostream & operator<< (ostream &os, const SpeciesRecord &sr)
 Output operator for SpeciesRecord. More...
 
ostream & operator<< (ostream &os, const SpeciesAuxData &sad)
 Output operator for SpeciesAuxData. More...
 
void find_broad_spec_locations (ArrayOfIndex &broad_spec_locations, const ArrayOfArrayOfSpeciesTag &abs_species, const Index this_species)
 Find the location of all broadening species in abs_species. More...
 
void calc_gamma_and_deltaf_artscat4 (Numeric &gamma, Numeric &deltaf, const Numeric p, const Numeric t, ConstVectorView vmrs, const Index this_species, const ArrayOfIndex &broad_spec_locations, const LineRecord &l_l, const Verbosity &verbosity)
 Calculate line width and pressure shift for artscat4. More...
 
void xsec_species (MatrixView xsec_attenuation, MatrixView xsec_phase, ConstVectorView f_grid, ConstVectorView abs_p, ConstVectorView abs_t, ConstMatrixView all_vmrs, const ArrayOfArrayOfSpeciesTag &abs_species, const Index this_species, const ArrayOfLineRecord &abs_lines, const Index ind_ls, const Index ind_lsn, const Numeric cutoff, const SpeciesAuxData &isotopologue_ratios, const Verbosity &verbosity)
 Calculate line absorption cross sections for one tag group. More...
 
Numeric wavenumber_to_joule (Numeric e)
 A little helper function to convert energy from units of wavenumber (cm^-1) to Joule (J). More...
 
Index species_index_from_species_name (String name)
 Return species index for given species name. More...
 
String species_name_from_species_index (const Index spec_ind)
 Return species name for given species index. More...
 
void convHitranIERF (Numeric &mdf, const Index &df)
 
void convHitranIERSH (Numeric &mdh, const Index &dh)
 
void convMytranIER (Numeric &mdh, const Index &dh)
 
ostream & operator<< (ostream &os, const LineshapeSpec &lsspec)
 
void xsec_species_line_mixing_wrapper (MatrixView xsec_attenuation, MatrixView xsec_phase, const ArrayOfArrayOfLineMixingRecord &line_mixing_data, const ArrayOfArrayOfIndex &line_mixing_data_lut, ConstVectorView f_grid, ConstVectorView abs_p, ConstVectorView abs_t, ConstMatrixView all_vmrs, const ArrayOfArrayOfSpeciesTag &abs_species, const Index this_species, const ArrayOfLineRecord &abs_lines, const Index ind_ls, const Index ind_lsn, const Numeric cutoff, const SpeciesAuxData &isotopologue_ratios, const Verbosity &verbosity)
 This will work as a wrapper for linemixing when abs_species contain relevant data. More...
 
void xsec_species_line_mixing_2nd_order (MatrixView xsec_attenuation, MatrixView xsec_phase, const ArrayOfArrayOfLineMixingRecord &line_mixing_data, const ArrayOfArrayOfIndex &line_mixing_data_lut, ConstVectorView f_grid, ConstVectorView abs_p, ConstVectorView abs_t, ConstMatrixView all_vmrs, const ArrayOfArrayOfSpeciesTag &abs_species, const Index this_species, const ArrayOfLineRecord &abs_lines, const Index ind_ls, const Index ind_lsn, const Numeric cutoff, const SpeciesAuxData &isotopologue_ratios, const Verbosity &verbosity)
 This is the second order line mixing correction as presented by: Makarov, D.S, M.Yu. More...
 

Variables

std::map< String, IndexSpeciesMap
 The map associated with species_data. More...
 

Detailed Description

Physical absorption routines.

The absorption workspace methods are in file m_abs.cc

This is the file from arts-1-0, back-ported to arts-1-1.

Author
Stefan Buehler and Axel von Engeln

Definition in file absorption.cc.

Function Documentation

◆ calc_gamma_and_deltaf_artscat4()

void calc_gamma_and_deltaf_artscat4 ( Numeric gamma,
Numeric deltaf,
const Numeric  p,
const Numeric  t,
ConstVectorView  vmrs,
const Index  this_species,
const ArrayOfIndex broad_spec_locations,
const LineRecord l_l,
const Verbosity verbosity 
)

Calculate line width and pressure shift for artscat4.

Return values
gammaLine width [Hz].
deltafPressure shift [Hz].
Parameters
pPressure [Pa].
tTemperature [K].
vmrsVector of VMRs for different species [dimensionless].
this_speciesIndex of current species in vmrs.
broad_spec_locationsHas length of number of allowed broadening species (6 in artscat-4). Gives for each species the position in vmrs, or negative if it should be ignored. See function find_broad_spec_locations for details.
l_lSpectral line data record (a single line).
verbosityVerbosity flag.
Author
Stefan Buehler
Date
2012-09-05

Definition at line 432 of file absorption.cc.

References abs, LineRecord::BroadSpecName(), CREATE_OUT2, LineRecord::Delta_foreign(), LineRecord::Gamma_foreign(), LineRecord::N_foreign(), LineRecord::NBroadSpec(), Array< base >::nelem(), LineRecord::Nself(), LineRecord::Sgam(), and LineRecord::Ti0().

Referenced by xsec_species().

◆ checkIsotopologueRatios()

void checkIsotopologueRatios ( const ArrayOfArrayOfSpeciesTag abs_species,
const SpeciesAuxData isoratios 
)

Check that isotopologue ratios for the given species are correctly defined.

Definition at line 245 of file absorption.cc.

References SpeciesAuxData::getParam(), SpeciesAuxData::getParams(), iso(), SpeciesRecord::Isotopologue(), SpeciesRecord::Name(), Array< base >::nelem(), and global_data::species_data.

Referenced by abs_xsec_per_speciesAddLines(), and propmat_clearskyAddZeeman().

◆ convHitranIERF()

void convHitranIERF ( Numeric mdf,
const Index df 
)

◆ convHitranIERSH()

void convHitranIERSH ( Numeric mdh,
const Index dh 
)

◆ convMytranIER()

void convMytranIER ( Numeric mdh,
const Index dh 
)

Definition at line 1433 of file absorption.cc.

Referenced by LineRecord::ReadFromMytran2Stream().

◆ define_species_map()

void define_species_map ( )

Define the species data map.

Define the species data map.

Author
Stefan Buehler

Definition at line 322 of file absorption.cc.

References global_data::species_data, and SpeciesMap.

Referenced by main().

◆ fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData()

void fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData ( SpeciesAuxData sad)

Fill SpeciesAuxData with default isotopologue ratios from species data.

Definition at line 305 of file absorption.cc.

References SpeciesAuxData::initParams(), SpeciesAuxData::setParam(), and global_data::species_data.

Referenced by isotopologue_ratiosInitFromBuiltin().

◆ find_broad_spec_locations()

void find_broad_spec_locations ( ArrayOfIndex broad_spec_locations,
const ArrayOfArrayOfSpeciesTag abs_species,
const Index  this_species 
)

Find the location of all broadening species in abs_species.

Set to -1 if not found. The length of array broad_spec_locations is the number of allowed broadening species (in ARTSCAT-4 N2, O2, H2O, CO2, H2, H2). The value means:

-1 = not in abs_species

-2 = in abs_species, but should be ignored because it is identical to Self

N = species is number N in abs_species

The catalogue contains also broadening parameters for "Self", but we ignore this here, since we know the position of species "Self" anyway.

Return values
broad_spec_locationsSee above.
Parameters
abs_speciesList of absorption species
this_speciesIndex of the current species in abs_species
Author
Stefan Buehler
Date
2012-09-04

Definition at line 372 of file absorption.cc.

References LineRecord::BroadSpecSpecIndex(), find_first_species_tg(), LineRecord::NBroadSpec(), and Array< base >::nelem().

Referenced by xsec_species().

◆ operator<<() [1/3]

ostream& operator<< ( ostream &  os,
const LineshapeSpec lsspec 
)

◆ operator<<() [2/3]

ostream& operator<< ( ostream &  os,
const SpeciesAuxData sad 
)

Output operator for SpeciesAuxData.

Author
Oliver Lemke

Definition at line 343 of file absorption.cc.

References SpeciesAuxData::getParams().

◆ operator<<() [3/3]

ostream& operator<< ( ostream &  os,
const SpeciesRecord sr 
)

Output operator for SpeciesRecord.

Incomplete version: only writes SpeciesName.

Author
Jana Mendrok

Definition at line 333 of file absorption.cc.

References SpeciesRecord::Isotopologue(), and SpeciesRecord::Name().

◆ species_index_from_species_name()

Index species_index_from_species_name ( String  name)

Return species index for given species name.

This is useful in connection with other functions that need a species index.

See also
find_first_species_tg.
Parameters
nameSpecies name.
Returns
Species index, -1 means not found.
Author
Stefan Buehler
Date
2003-01-13

Definition at line 1260 of file absorption.cc.

References SpeciesMap, and my_basic_string< charT >::trim().

Referenced by abs_h2oSet(), abs_lookupSetup(), abs_lookupSetupBatch(), abs_lookupSetupWide(), abs_lookupTestAccMC(), abs_lookupTestAccuracy(), abs_n2Set(), abs_o2Set(), LineRecord::BroadSpecSpecIndex(), choose_abs_nls(), GasAbsLookup::Extract(), find_nonlinear_continua(), refr_index_airMWgeneral(), refr_index_airThayer(), CIARecord::SetMoleculeName(), SpeciesTag::SpeciesTag(), xml_read_from_stream(), and z_fieldFromHSE().

◆ species_name_from_species_index()

String species_name_from_species_index ( const Index  spec_ind)

Return species name for given species index.

This is useful in connection with other functions that use a species index.

Does an assertion that the index really corresponds to a species.

Parameters
spec_indSpecies index.
Returns
Species name
Author
Stefan Buehler
Date
2013-01-04

Definition at line 1301 of file absorption.cc.

References SpeciesRecord::Name(), and global_data::species_data.

Referenced by abs_cia_dataReadFromCIA(), abs_cia_dataReadFromXML(), get_species_name(), CIARecord::MoleculeName(), SpeciesTag::Name(), and xml_write_to_stream().

◆ wavenumber_to_joule()

Numeric wavenumber_to_joule ( Numeric  e)

A little helper function to convert energy from units of wavenumber (cm^-1) to Joule (J).

This is used when reading HITRAN or JPL catalogue files, which have the lower state energy in cm^-1.

Returns
Energy in J.
Parameters
eEnergy in cm^-1.
Author
Stefan Buehler
Date
2001-06-26

Definition at line 1227 of file absorption.cc.

References PLANCK_CONST, and SPEED_OF_LIGHT.

Referenced by LineRecord::ReadFromHitran2001Stream(), LineRecord::ReadFromHitran2004Stream(), LineRecord::ReadFromJplStream(), and LineRecord::ReadFromMytran2Stream().

◆ xsec_species()

void xsec_species ( MatrixView  xsec_attenuation,
MatrixView  xsec_phase,
ConstVectorView  f_grid,
ConstVectorView  abs_p,
ConstVectorView  abs_t,
ConstMatrixView  all_vmrs,
const ArrayOfArrayOfSpeciesTag abs_species,
const Index  this_species,
const ArrayOfLineRecord abs_lines,
const Index  ind_ls,
const Index  ind_lsn,
const Numeric  cutoff,
const SpeciesAuxData isotopologue_ratios,
const Verbosity verbosity 
)

Calculate line absorption cross sections for one tag group.

All lines in the line list must belong to the same species. This must be ensured by abs_lines_per_speciesCreateFromLines, so it is only verified with assert. Also, the input vectors abs_p, and abs_t must all have the same dimension.

This is mainly a copy of abs_species which is removed now, with the difference that the vmrs are removed from the absorption coefficient calculation. (the vmr is still used for the self broadening)

Continua are not handled by this function, you have to call xsec_continuum_tag for those.

Return values
xsecCross section of one tag group. This is now the true absorption cross section in units of m^2.
Parameters
f_gridFrequency grid.
abs_pPressure grid.
abs_tTemperatures associated with abs_p.
all_vmrsGas volume mixing ratios [nspecies, np].
abs_speciesSpecies tags for all species.
this_speciesIndex of the current species in abs_species.
abs_linesThe spectroscopic line list.
ind_lsIndex to lineshape function.
ind_lsnIndex to lineshape norm.
cutoffLineshape cutoff.
isotopologue_ratiosIsotopologue ratios.
Author
Stefan Buehler and Axel von Engeln
Date
2001-01-11

Changed from pseudo cross sections to true cross sections

Author
Stefan Buehler
Date
2007-08-08

Adapted to new Perrin line parameters, treating broadening by different gases explicitly

Author
Stefan Buehler
Date
2012-09-03

Definition at line 594 of file absorption.cc.

References LineRecord::Agam(), arts_omp_get_max_threads(), arts_omp_get_thread_num(), arts_omp_in_parallel(), AVOGADROS_NUMB, BOLTZMAN_CONST, calc_gamma_and_deltaf_artscat4(), IsotopologueRecord::CalculatePartitionFctRatio(), LineRecord::Elow(), LineRecord::F(), fac(), find_broad_spec_locations(), SpeciesAuxData::getParam(), LineRecord::I0(), is_sorted(), LineRecord::Isotopologue(), LineRecord::IsotopologueData(), joker, global_data::lineshape_data, global_data::lineshape_norm_data, IsotopologueRecord::Mass(), LineRecord::Nair(), ConstMatrixView::ncols(), Array< base >::nelem(), ConstVectorView::nelem(), ConstMatrixView::nrows(), LineRecord::Nself(), PLANCK_CONST, LineRecord::Psf(), LineRecord::Sgam(), LineRecord::Species(), SPEED_OF_LIGHT, LineRecord::Tgam(), LineRecord::Ti0(), and LineRecord::Version().

Referenced by abs_xsec_per_speciesAddLines(), xsec_species_line_mixing_2nd_order(), and xsec_species_line_mixing_wrapper().

◆ xsec_species_line_mixing_2nd_order()

void xsec_species_line_mixing_2nd_order ( MatrixView  xsec_attenuation,
MatrixView  xsec_phase,
const ArrayOfArrayOfLineMixingRecord line_mixing_data,
const ArrayOfArrayOfIndex line_mixing_data_lut,
ConstVectorView  f_grid,
ConstVectorView  abs_p,
ConstVectorView  abs_t,
ConstMatrixView  all_vmrs,
const ArrayOfArrayOfSpeciesTag abs_species,
const Index  this_species,
const ArrayOfLineRecord abs_lines,
const Index  ind_ls,
const Index  ind_lsn,
const Numeric  cutoff,
const SpeciesAuxData isotopologue_ratios,
const Verbosity verbosity 
)

This is the second order line mixing correction as presented by: Makarov, D.S, M.Yu.

Tretyakov and P.W. Rosenkranz, 2011, 60-GHz oxygen band: Precise experimental profiles and extended absorption modeling in a wide temperature range, JQSRT 112, 1420-1428.

Return values
xsec_attenuationCross section of one tag group. This is now the true attenuation cross section in units of m^2.
xsec_phaseCross section of one tag group. This is now the true phase cross section in units of m^2.
Parameters
line_mixing_dataLine mixing data for specific type of line mixing.
line_mixing_data_lutLine mixing data lookup table.
f_gridFrequency grid.
abs_pPressure grid.
abs_tTemperatures associated with abs_p.
all_vmrsGas volume mixing ratios [nspecies, np].
abs_speciesSpecies tags for all species.
this_speciesIndex of the current species in abs_species.
abs_linesThe spectroscopic line list.
ind_lsIndex to lineshape function.
ind_lsnIndex to lineshape norm.
cutoffLineshape cutoff.
isotopologue_ratiosIsotopologue ratios.
Author
Richard Larsson
Date
2013-04-24

Definition at line 1598 of file absorption.cc.

References joker, global_data::lineshape_data, ll, Array< base >::nelem(), ConstVectorView::nelem(), ConstMatrixView::nrows(), and xsec_species().

Referenced by xsec_species_line_mixing_wrapper().

◆ xsec_species_line_mixing_wrapper()

void xsec_species_line_mixing_wrapper ( MatrixView  xsec_attenuation,
MatrixView  xsec_phase,
const ArrayOfArrayOfLineMixingRecord line_mixing_data,
const ArrayOfArrayOfIndex line_mixing_data_lut,
ConstVectorView  f_grid,
ConstVectorView  abs_p,
ConstVectorView  abs_t,
ConstMatrixView  all_vmrs,
const ArrayOfArrayOfSpeciesTag abs_species,
const Index  this_species,
const ArrayOfLineRecord abs_lines,
const Index  ind_ls,
const Index  ind_lsn,
const Numeric  cutoff,
const SpeciesAuxData isotopologue_ratios,
const Verbosity verbosity 
)

This will work as a wrapper for linemixing when abs_species contain relevant data.

The funciton will only pass on arguments to xsec_species if there is no linemixing.

Return values
xsec_attenuationCross section of one tag group. This is now the true attenuation cross section in units of m^2.
xsec_phaseCross section of one tag group. This is now the true phase cross section in units of m^2.
Parameters
line_mixing_dataLine mixing data for specific type of line mixing.
line_mixing_data_lutLine mixing data lookup table.
f_gridFrequency grid.
abs_pPressure grid.
abs_tTemperatures associated with abs_p.
all_vmrsGas volume mixing ratios [nspecies, np].
abs_speciesSpecies tags for all species.
this_speciesIndex of the current species in abs_species.
abs_linesThe spectroscopic line list.
ind_lsIndex to lineshape function.
ind_lsnIndex to lineshape norm.
cutoffLineshape cutoff.
isotopologue_ratiosIsotopologue ratios.
Author
Richard Larsson
Date
2013-04-24

Definition at line 1530 of file absorption.cc.

References SpeciesTag::LINE_MIXING_TYPE_2NDORDER, SpeciesTag::LINE_MIXING_TYPE_NONE, xsec_species(), and xsec_species_line_mixing_2nd_order().

Referenced by abs_xsec_per_speciesAddLines(), and xsec_species_line_mixing_wrapper_with_zeeman().

Variable Documentation

◆ SpeciesMap

std::map<String, Index> SpeciesMap

The map associated with species_data.

Definition at line 50 of file absorption.cc.

Referenced by define_species_map(), and species_index_from_species_name().