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 "auto_md.h"
#include "messages.h"

Go to the source code of this file.

Functions

void define_species_map ()
 Define the species data map.
ostream & operator<< (ostream &os, const LineRecord &lr)
 Output operator for LineRecord.
template<class T>
void extract (T &x, String &line, Index n)
 Extract something from a catalogue line.
ostream & operator<< (ostream &os, const OneTag &ot)
 Output operator for OneTag.
void get_tagindex_for_Strings (ArrayOfIndex &tags1_index, const TagGroups &tags1, const ArrayOfString &tags2_Strings)
 Returns the index among some tag groups for an array of tag Strings.
void get_tag_group_index_for_tag_group (Index &tgs1_index, const TagGroups &tgs1, const Array< OneTag > &tg2)
 Returns the index of the tag group tg2 within the array of tag groups tgs1.
String get_tag_group_name (const Array< OneTag > &tg)
 Print the name of a tag group.
void write_lines_to_stream (ostream &os, const ArrayOfLineRecord &lines)
 A helper function that writes lines in a line list to a stream.
bool is_sorted (ConstVectorView x)
 Checks if a vector is sorted in ascending order.
void xsec_species (MatrixView xsec, ConstVectorView f_mono, ConstVectorView p_abs, ConstVectorView t_abs, ConstVectorView h2o_abs, ConstVectorView vmr, const ArrayOfLineRecord &lines, const Index ind_ls, const Index ind_lsn, const Numeric cutoff)
 Calculate line absorption cross sections for one tag group.
void refr_index_BoudourisDryAir (Vector &refr_index, ConstVectorView p_abs, ConstVectorView t_abs)
 Calculates the refractive index for dry air at microwave frequncies following Boudouris 1963.
void refr_index_Boudouris (Vector &refr_index, ConstVectorView p_abs, ConstVectorView t_abs, ConstVectorView h2o_abs)
 Calculates the refractive index at microwave frequncies following Boudouris 1963.
Numeric wavenumber_to_joule (Numeric e)
 A little helper function to convert energy from units of wavenumber (cm^-1) to Joule (J).
void convHitranIERF (Numeric &mdf, const Index &df)
void convHitranIERSH (Numeric &mdh, const Index &dh)
void convMytranIER (Numeric &mdh, const Index &dh)

Variables

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


Detailed Description

Physical absorption routines.

The absorption workspace methods are in file m_abs.cc

Author:
Stefan Buehler

Definition in file absorption.cc.


Function Documentation

void convHitranIERF ( Numeric mdf,
const Index df 
)

Definition at line 2901 of file absorption.cc.

void convHitranIERSH ( Numeric mdh,
const Index dh 
)

Definition at line 2948 of file absorption.cc.

void convMytranIER ( Numeric mdh,
const Index dh 
)

Definition at line 3007 of file absorption.cc.

void define_species_map (  ) 

Define the species data map.

Author:
Stefan Buehler

Definition at line 64 of file absorption.cc.

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

Extract something from a catalogue line.

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

Return values:
x What was extracted from the beginning of the line.
line What was extracted is also cut away from line.
Parameters:
n The width of the stuff to extract.
Author:
Stefan Buehler

Definition at line 127 of file absorption.cc.

void get_tag_group_index_for_tag_group ( Index tgs1_index,
const TagGroups tgs1,
const Array< OneTag > &  tg2 
)

Returns the index of the tag group tg2 within the array of tag groups tgs1.

Slightly modified copy of get_tagindex_for_Strings.

Exceptions:
runtime_error Could not find tg2 in tgs1.
Return values:
tgs1_index Index in tgs1 for tg2
Parameters:
tgs1 The tags groups to search in.
tg2 The tag group for which the index shall be found.
Author:
Patrick Eriksson, Axel von Engeln, and Stefan Buehler
Date:
2001-01-31

Definition at line 2233 of file absorption.cc.

String get_tag_group_name ( const Array< OneTag > &  tg  ) 

Print the name of a tag group.

A tag group consists of several elementary OneTags. This functions returns a String with the name of the entire tag group. This is nice for informational output messages, for example in the absorption routines.

Parameters:
tg The tag group in question.
Returns:
The full name of the tag group, as it could occur in the controlfile.
Author:
Stefan Buehler
Date:
2001-03-13

Definition at line 2286 of file absorption.cc.

void get_tagindex_for_Strings ( ArrayOfIndex tags1_index,
const TagGroups tags1,
const ArrayOfString tags2_Strings 
)

Returns the index among some tag groups for an array of tag Strings.

{verbatim} For example, if tags1 correspond to the definition ["O3","H2O-161,H2O-162"] and the tag Strings are ["H2O-161,H2O-162","O3"] the tags1_index becomes [2,1] {verbatim}

Exceptions:
runtime_error Some String is not a valid tag item.
runtime_error Not all Strings are not found among the tags.
Return values:
tags1_index Index in tags1 for tags2_Strings
Parameters:
tags1 The tags to search in.
tags2_Strings The tag Strings for which indeces shall be found.
Author:
Patrick Eriksson
Date:
2000-12-06

Definition at line 2176 of file absorption.cc.

bool is_sorted ( ConstVectorView  x  ) 

Checks if a vector is sorted in ascending order.

Duplicated values are allowed.

We need this to make sure that the frequency grid is sorted in the case of lineshape with cutoff. Duplicate values are allowed.

THIS FUNCTION IS TAKEN FROM ARTS-1-1.

Parameters:
x A vector.
Returns:
True if sorted.

Definition at line 2340 of file absorption.cc.

ostream& operator<< ( ostream &  os,
const OneTag ot 
)

Output operator for OneTag.

Author:
Stefan Buehler

Definition at line 2147 of file absorption.cc.

ostream& operator<< ( ostream &  os,
const LineRecord lr 
)

Output operator for LineRecord.

The result should look like a catalogue line.

Author:
Stefan Buehler

Definition at line 76 of file absorption.cc.

void refr_index_Boudouris ( Vector refr_index,
ConstVectorView  p_abs,
ConstVectorView  t_abs,
ConstVectorView  h2o_abs 
)

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_index refractive index
Parameters:
p_abs absorption pressure grid
t_abs temperatures at p_abs
h2o_abs H2O vmr at p_abs
Author:
Patrick Eriksson
Date:
2001-02-16

Definition at line 2842 of file absorption.cc.

void refr_index_BoudourisDryAir ( Vector refr_index,
ConstVectorView  p_abs,
ConstVectorView  t_abs 
)

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_index refractive index
Parameters:
p_abs absorption pressure grid
t_abs temperatures at p_abs
Author:
Patrick Eriksson
Date:
2001-02-16

Definition at line 2810 of file absorption.cc.

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:
e Energy in cm^-1.
Author:
Stefan Buehler
Date:
2001-06-26

Definition at line 2879 of file absorption.cc.

void write_lines_to_stream ( ostream &  os,
const ArrayOfLineRecord lines 
)

A helper function that writes lines in a line list to a stream.

Return values:
os The stream to write to.
Parameters:
lines The line list to write.
Date:
2000-06-12
Author:
Stefan Buehler

Definition at line 2309 of file absorption.cc.

void xsec_species ( MatrixView  xsec,
ConstVectorView  f_mono,
ConstVectorView  p_abs,
ConstVectorView  t_abs,
ConstVectorView  h2o_abs,
ConstVectorView  vmr,
const ArrayOfLineRecord lines,
const Index  ind_ls,
const Index  ind_lsn,
const Numeric  cutoff 
)

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 lines_per_tgCreateFromLines, so it is only verified with assert. Also, the input vectors p_abs, and t_abs 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:
xsec Cross section of one tag group.
Parameters:
f_mono Frequency grid.
p_abs Pressure grid.
t_abs Temperatures associated with p_abs.
h2o_abs Total volume mixing ratio of water vapor.
vmr Volume mixing ratio of the calculated species.
lines The spectroscopic line list.
ind_ls Lineshape specifications.
Author:
Stefan Buehler and Axel von Engeln
Date:
2001-01-11

Definition at line 2379 of file absorption.cc.


Variable Documentation

std::map<String, Index> SpeciesMap

The map associated with species_data.

Definition at line 42 of file absorption.cc.


Generated on Wed Feb 4 08:17:21 2009 for ARTS by  doxygen 1.5.6