ARTS  2.0.49
absorption.cc File Reference

Physical absorption routines. More...

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

Go to the source code of this file.

Functions

void define_species_map ()
 Define the species data map. More...
 
ostream & operator<< (ostream &os, const LineRecord &lr)
 Output operator for LineRecord. More...
 
template<class T >
void extract (T &x, String &line, Index n)
 Extract something from a catalogue line. More...
 
void xsec_species (MatrixView xsec, ConstVectorView f_grid, ConstVectorView abs_p, ConstVectorView abs_t, ConstVectorView abs_h2o_orig, ConstVectorView vmr, const ArrayOfLineRecord &abs_lines, const Index ind_ls, const Index ind_lsn, const Numeric cutoff, const Verbosity &verbosity)
 Calculate line absorption cross sections for one tag group. More...
 
void refr_index_BoudourisDryAir (Vector &refr_index, ConstVectorView abs_p, ConstVectorView abs_t)
 Calculates the refractive index for dry air at microwave frequncies following Boudouris 1963. More...
 
void refr_index_Boudouris (Vector &refr_index, ConstVectorView abs_p, ConstVectorView abs_t, ConstVectorView abs_h2o)
 Calculates the refractive index at microwave frequncies following Boudouris 1963. More...
 
Numeric wavenumber_to_joule (Numeric e)
 A little helper function to convert energy from units of wavenumber (cm^-1) to Joule (J). 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 &)
 

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

◆ 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 2973 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 70 of file absorption.cc.

References species_data, and SpeciesMap.

Referenced by main().

◆ extract()

template<class T >
void extract ( T &  x,
String line,
Index  n 
)

Extract something from a catalogue line.

This is just a small helper function to safe some typing.

Return values
xWhat was extracted from the beginning of the line.
lineWhat was extracted is also cut away from line.
Parameters
nThe width of the stuff to extract.
Author
Stefan Buehler

Definition at line 179 of file absorption.cc.

Referenced by LineRecord::ReadFromArtscat3Stream(), LineRecord::ReadFromArtscat4Stream(), LineRecord::ReadFromHitran2004Stream(), LineRecord::ReadFromHitranStream(), LineRecord::ReadFromJplStream(), and LineRecord::ReadFromMytran2Stream().

◆ operator<<() [1/2]

◆ operator<<() [2/2]

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

Definition at line 3034 of file absorption.cc.

◆ refr_index_Boudouris()

void refr_index_Boudouris ( Vector refr_index,
ConstVectorView  abs_p,
ConstVectorView  abs_t,
ConstVectorView  abs_h2o 
)

Calculates the refractive index at microwave frequncies following Boudouris 1963.

The expression is also found in Chapter 5 of the Janssen book.

Return values
refr_indexrefractive index
Parameters
abs_pabsorption pressure grid
abs_ttemperatures at abs_p
abs_h2oH2O vmr at abs_p
Author
Patrick Eriksson
Date
2001-02-16

Definition at line 2808 of file absorption.cc.

References ConstVectorView::nelem(), and Vector::resize().

◆ refr_index_BoudourisDryAir()

void refr_index_BoudourisDryAir ( Vector refr_index,
ConstVectorView  abs_p,
ConstVectorView  abs_t 
)

Calculates the refractive index for dry air at microwave frequncies following Boudouris 1963.

The expression is also found in Chapter 5 of the Janssen book.

The atmosphere is assumed to have no water vapour.

Return values
refr_indexrefractive index
Parameters
abs_pabsorption pressure grid
abs_ttemperatures at abs_p
Author
Patrick Eriksson
Date
2001-02-16

Definition at line 2776 of file absorption.cc.

References ConstVectorView::nelem(), and Vector::resize().

◆ 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 2845 of file absorption.cc.

References PLANCK_CONST, and SPEED_OF_LIGHT.

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

◆ xsec_species()

void xsec_species ( MatrixView  xsec,
ConstVectorView  f_grid,
ConstVectorView  abs_p,
ConstVectorView  abs_t,
ConstVectorView  abs_h2o_orig,
ConstVectorView  vmr,
const ArrayOfLineRecord abs_lines,
const Index  ind_ls,
const Index  ind_lsn,
const Numeric  cutoff,
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.
abs_h2o_origTotal volume mixing ratio of water vapor.
vmrVolume mixing ratio of the calculated species.
abs_linesThe spectroscopic line list.
ind_lsIndex to lineshape function.
ind_lsnIndex to lineshape norm.
cutoffLineshape cutoff.
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

Definition at line 2240 of file absorption.cc.

References IsotopeRecord::Abundance(), LineRecord::Agam(), arts_omp_get_max_threads(), arts_omp_get_thread_num(), arts_omp_in_parallel(), LineRecord::Aux(), AVOGADROS_NUMB, BOLTZMAN_CONST, IsotopeRecord::CalculatePartitionFctRatio(), CREATE_OUT0, LineRecord::Elow(), exit_or_rethrow(), LineRecord::F(), fac(), LineRecord::I0(), is_sorted(), LineRecord::IsotopeData(), joker, lineshape_data, lineshape_norm_data, IsotopeRecord::Mass(), LineRecord::Nair(), LineRecord::Naux(), ConstMatrixView::ncols(), Array< base >::nelem(), ConstVectorView::nelem(), ConstMatrixView::nrows(), LineRecord::Nself(), PLANCK_CONST, LineRecord::Psf(), LineRecord::Sgam(), SPEED_OF_LIGHT, LineRecord::Tgam(), and LineRecord::Ti0().

Referenced by abs_xsec_per_speciesAddLines().

Variable Documentation

◆ SpeciesMap

std::map<String, Index> SpeciesMap

The map associated with species_data.

Definition at line 45 of file absorption.cc.

Referenced by define_species_map(), and species_index_from_species_name().