Go to the documentation of this file.
67 result += *it * exponent;
68 exponent *= temperature;
114 static map<String, SpecIsoMap> ArtsMap;
117 static bool hinit =
false;
122 out3 <<
" ARTS index table:\n";
135 ArtsMap[buf] = indicies;
142 const Index& i2 = ArtsMap[buf].Isotopologueindex();
144 out3 <<
" Arts Identifier = " << buf <<
" Species = "
145 << setw(10) << setiosflags(ios::left)
168 if (is.eof())
return true;
171 if (!is)
throw runtime_error (
"Stream bad.");
180 if (line.
nelem() == 0 && is.eof())
return true;
195 istringstream icecream(line);
199 if (artsid.length() != 0)
206 const map<String, SpecIsoMap>::const_iterator i = ArtsMap.find(artsid);
207 if ( i == ArtsMap.end() )
210 os <<
"ARTS Tag: " << artsid <<
" is unknown.";
211 throw runtime_error(os.str());
221 misotopologue =
id.Isotopologueindex();
228 for (
Index ip = 0; ip < nparams; ip++)
231 params(misotopologue, ip) = p;
234 catch (runtime_error)
236 throw runtime_error(
"Error reading SpeciesAuxData.");
254 os <<
"Number of species in SpeciesAuxData (" << isoratios.
getParams().
nelem()
255 <<
"does not fit builtin species data (" <<
species_data.nelem() <<
").";
256 throw runtime_error(os.str());
267 const Index sp = abs_species[i][0].Species();
276 os <<
"Incorrect number of isotopologues in isotopologue data.\n"
277 <<
"Species: " << this_sd.
Name() <<
".\n"
278 <<
"Number of isotopes in SpeciesAuxData ("
279 << isoratios.
getParams()[i].nrows() <<
") "
280 <<
"does not fit builtin species data (" << this_sd.
Isotopologue().
nelem() <<
").";
281 throw runtime_error(os.str());
293 os <<
"Invalid isotopologue ratio.\n"
294 <<
"Species: " << this_sd.
Name() <<
"-"
297 throw runtime_error(os.str());
374 const Index this_species)
377 assert(this_species>=0);
378 assert(this_species<abs_species.
nelem());
384 broad_spec_locations.resize(nbs);
397 for (
Index i=0; i<nbs; ++i) {
402 if ( isi == abs_species[this_species][0].Species() )
403 broad_spec_locations[i] = -2;
437 const Index this_species,
446 assert(nbs==broad_spec_locations.
nelem());
453 const Numeric p_self = vmrs[this_species] * p;
454 const Numeric p_foreign = p-p_self;
458 Numeric broad_spec_vmr_sum = 0;
461 gamma = l_l.
Sgam() * pow(theta, l_l.
Nself()) * p_self;
471 for (
Index i=0; i<nbs; ++i) {
472 if ( broad_spec_locations[i] < -1 ) {
479 os <<
"Inconsistency in LineRecord, self broadening and line "
481 <<
"should be identical.\n"
484 throw runtime_error(os.str());
486 }
else if ( broad_spec_locations[i] >= 0 ) {
489 broad_spec_vmr_sum += vmrs[broad_spec_locations[i]];
493 * vmrs[broad_spec_locations[i]];
501 * vmrs[broad_spec_locations[i]];
506 if (
abs(vmrs[this_species]+broad_spec_vmr_sum-1) > 0.1
507 && out2.sufficient_priority() )
510 os <<
"Warning: The total VMR of all your defined broadening\n"
511 <<
"species (including \"self\") is "
512 << vmrs[this_species]+broad_spec_vmr_sum
513 <<
", more than 10% " <<
"different from 1.\n";
519 if (broad_spec_vmr_sum != 0.)
521 gamma_foreign /= broad_spec_vmr_sum;
522 deltaf /= broad_spec_vmr_sum;
524 else if (p_self > 0.)
531 gamma_foreign = gamma/p_self;
541 gamma_foreign *= p_foreign;
546 gamma += gamma_foreign;
601 const Index this_species,
643 Vector ls_attenuation(nf+1);
647 const bool cut = (cutoff != -1) ?
true :
false;
658 os <<
"If you use a lineshape function with cutoff, your\n"
659 <<
"frequency grid *f_grid* must be sorted.\n"
660 <<
"(Duplicate values are allowed.)";
661 throw runtime_error(os.str());
666 bool negative =
false;
668 for (
Index i = 0; !negative && i < abs_t.
nelem (); i++)
677 os <<
"abs_t contains at least one negative temperature value.\n"
678 <<
"This is not allowed.";
679 throw runtime_error(os.str());
696 Index ii = (nf+1 < 10) ? 10 : nf+1;
703 if ( abs_t.
nelem() != np )
706 os <<
"Variable abs_t must have the same dimension as abs_p.\n"
707 <<
"abs_t.nelem() = " << abs_t.
nelem() <<
'\n'
708 <<
"abs_p.nelem() = " << np;
709 throw runtime_error(os.str());
714 if ( all_vmrs.
ncols() != np )
717 os <<
"Number of columns of all_vmrs must match abs_p.\n"
718 <<
"all_vmrs.ncols() = " << all_vmrs.
ncols() <<
'\n'
719 <<
"abs_p.nelem() = " << np;
720 throw runtime_error(os.str());
725 if ( all_vmrs.
nrows() != nspecies)
728 os <<
"Number of rows of all_vmrs must match abs_species.\n"
729 <<
"all_vmrs.nrows() = " << all_vmrs.
nrows() <<
'\n'
730 <<
"abs_species.nelem() = " << nspecies;
731 throw runtime_error(os.str());
768 if ( xsec_attenuation.
nrows() != nf || xsec_attenuation.
ncols() != np )
771 os <<
"Variable xsec must have dimensions [f_grid.nelem(),abs_p.nelem()].\n"
772 <<
"[xsec_attenuation.nrows(),xsec_attenuation.ncols()] = [" << xsec_attenuation.
nrows()
773 <<
", " << xsec_attenuation.
ncols() <<
"]\n"
774 <<
"f_grid.nelem() = " << nf <<
'\n'
775 <<
"abs_p.nelem() = " << np;
776 throw runtime_error(os.str());
778 if ( xsec_phase.
nrows() != nf || xsec_phase.
ncols() != np )
781 os <<
"Variable xsec must have dimensions [f_grid.nelem(),abs_p.nelem()].\n"
782 <<
"[xsec_phase.nrows(),xsec_phase.ncols()] = [" << xsec_phase.
nrows()
783 <<
", " << xsec_phase.
ncols() <<
"]\n"
784 <<
"f_grid.nelem() = " << nf <<
'\n'
785 <<
"abs_p.nelem() = " << np;
786 throw runtime_error(os.str());
806 #pragma omp parallel for \
807 if (!arts_omp_in_parallel() \
808 && np >= arts_omp_get_max_threads()) \
809 firstprivate(ls_attenuation, ls_phase, fac, f_local, aux)
810 for (
Index i=0; i<np; ++i )
812 if (failed)
continue;
817 const Numeric vmr_i = all_vmrs(this_species,i);
828 const Numeric p_partial = p_i * vmr_i;
856 Matrix xsec_accum_attenuation(n_lbl_threads, xsec_i_attenuation.
nelem(), 0);
857 Matrix xsec_accum_phase(n_lbl_threads, xsec_i_phase.
nelem(), 0);
862 #pragma omp parallel for \
863 if (!arts_omp_in_parallel() \
864 && nl >= arts_omp_get_max_threads()) \
865 firstprivate(ls_attenuation, ls_phase, fac, f_local, aux)
866 for (
Index l=0; l< nl; ++l )
869 if (failed)
continue;
885 f_local[
Range(0,nf)] = f_grid;
904 os <<
"Unknown spectral line catalogue version (artscat-"
906 <<
"Allowed are artscat-3 and artscat-4.";
907 throw runtime_error(os.str());
946 intensity *= part_fct_ratio * nom / denom;
958 intensity = intensity * mafac / sinh(mafac);
976 broad_spec_locations,
1019 const Numeric theta = Tgam / t_i;
1020 const Numeric theta_Nair = pow(theta, Nair);
1022 gamma = Agam * theta_Nair * (p_i - p_partial)
1023 + Sgam * pow(theta, Nself) * p_partial;
1030 deltaf = Psf * p_i *
1046 const Numeric sigma = F0 * doppler_const *
1057 Index i_f_max = nf-1;
1067 while ( i_f_min < nf && (F0 - cutoff) > f_grid[i_f_min] )
1079 while ( i_f_max >= 0 && (F0 + cutoff) < f_grid[i_f_max] )
1087 f_local[i_f_max] = F0 + cutoff;
1090 nfls = i_f_max - i_f_min + 1;
1118 f_local[
Range(i_f_min,nfls)]);
1122 f_local[
Range(i_f_min,nfls)],
1143 this_ls_attenuation -= ls_attenuation[nfls-1];
1144 if (calc_phase) this_ls_phase -= ls_phase[nfls-1];
1161 const Numeric factors = intensity
1170 this_ls_attenuation *= this_fac;
1171 this_ls_attenuation *= factors;
1172 this_xsec_attenuation += this_ls_attenuation;
1176 this_ls_phase *= this_fac;
1177 this_ls_phase *= factors;
1178 this_xsec_phase += this_ls_phase;
1184 catch (runtime_error e)
1186 #pragma omp critical (xsec_species_fail)
1187 { fail_msg = e.what(); failed =
true; }
1193 if (failed)
continue;
1196 for (
Index j=0; j<xsec_accum_attenuation.
nrows(); ++j)
1198 xsec_i_attenuation += xsec_accum_attenuation(j,
Range(
joker));
1202 for (
Index j=0; j<xsec_accum_phase.
nrows(); ++j)
1204 xsec_i_phase += xsec_accum_phase(j,
Range(
joker));
1208 if (failed)
throw runtime_error(
"Run-time error in function: xsec_species\n" + fail_msg);
1238 return e*lower_energy_const;
1271 map<String, Index>::const_iterator mi =
SpeciesMap.find(name);
1275 mspecies = mi->second;
1309 assert( spec_ind>=0 );
1496 os <<
"LineshapeSpec Index: " << lsspec.
Ind_ls() <<
", NormIndex: " << lsspec.
Ind_lsn()
1497 <<
", Cutoff: " << lsspec.
Cutoff() << endl;
1539 const Index this_species,
1542 const Index ind_lsn,
1548 xsec_species(xsec_attenuation, xsec_phase, f_grid, abs_p, abs_t, all_vmrs,
1549 abs_species, this_species, abs_lines, ind_ls, ind_lsn, cutoff,
1550 isotopologue_ratios, verbosity);
1553 switch(abs_species[this_species][0].LineMixingType())
1557 line_mixing_data_lut, f_grid, abs_p, abs_t, all_vmrs,
1558 abs_species, this_species, abs_lines, ind_ls, ind_lsn,
1559 cutoff, isotopologue_ratios, verbosity);
1607 const Index this_species,
1610 const Index ind_lsn,
1620 os <<
"This is an error message. You are using " <<
lineshape_data[ind_ls].Name() <<
1621 ".\n"<<
"This line shape does not include phase in its calculations and\nis therefore invalid for " <<
1622 "second order line mixing.\n\n";
1623 throw runtime_error(os.str());
1626 const ArrayOfIndex& lut = line_mixing_data_lut[this_species];
1639 const Numeric p = abs_p[jj], t = abs_t[jj];
1644 ll[0] = abs_lines[ii];
1651 const Vector& dat = data[lut[ii]].Data();
1652 if( dat.
nelem() != 10 )
1655 os <<
"This is an error message. You have not defined the linemixing data vector " <<
1656 "properly. It should be of length 10, yet your version is of length: " << dat.
nelem() <<
".\n" <<
1657 "The structure of the line mixing data Vector should be as follows:\n" <<
1658 "data[0] is the first order zeroth phase correction\n" <<
1659 "data[1] is the first order first phase correction\n" <<
1660 "data[2] is the second order zeroth phase correction\n" <<
1661 "data[3] is the second order first phase correction \n" <<
1662 "data[4] is the second order zeroth frequency correction\n" <<
1663 "data[5] is the second order first frequency correction \n" <<
1664 "data[6] is the reference temperature\n" <<
1665 "data[7] is the first order phase correction temperature dependence exponent\n" <<
1666 "data[8] is the second order attenuation correction temperature dependence exponent\n" <<
1667 "data[9] is the second order frequency correction temperature dependence exponent";
1668 throw std::runtime_error(os.str());
1673 * ( ( dat[0] + dat[1] * ( dat[6]/t-1. ) )
1674 * pow( dat[6]/t, dat[7] ) );
1678 * ( ( dat[2] + dat[3] * ( dat[6]/t-1. ) )
1679 * pow( dat[6]/t, dat[8] ) );
1683 * ( ( dat[4] + dat[5] * ( dat[6]/t-1. ) )
1684 * pow( dat[6]/t, dat[9] ) ) +
ll[0].F());
1694 abs_species, this_species,
ll, ind_ls, ind_lsn, cutoff,
1695 isotopologue_ratios, verbosity);
1700 xsec_attenuation(
joker,jj) += phase(
joker,0);
1701 xsec_attenuation(
joker,jj) += attenuation(
joker,0);
1702 attenuation *= a[1];
1703 xsec_attenuation(
joker,jj) += attenuation(
joker,0);
Numeric Ti0() const
Reference temperature for I0 in K:
int arts_omp_get_thread_num()
Wrapper for omp_get_thread_num.
const Numeric & Mass() const
Mass of the isotopologue.
Numeric fac(const Index n)
fac
Numeric Sgam() const
Self broadened width in Hz/Pa:
const Numeric AVOGADROS_NUMB
Global constant, the Avogadro's number [molec/kg].
Declarations required for the calculation of absorption coefficients.
Numeric Gamma_foreign(const Index i) const
ARTSCAT-4 foreign broadening parameters in Hz/Pa :
The constant iterator class for sub vectors.
@ LINE_MIXING_TYPE_2NDORDER
const Index & Ind_ls() const
Return the index of this lineshape.
Numeric Agam() const
Air broadened width in Hz/Pa:
Numeric Nself() const
SGAM temperature exponent (dimensionless):
const Numeric PLANCK_CONST
Global constant, the Planck constant [Js].
static Index BroadSpecSpecIndex(const Index i)
Return the internal species index (index in species_data) of an artscat-4 broadening species,...
Numeric interp(ConstVectorView itw, ConstVectorView a, const GridPos &tc)
Red 1D Interpolate.
void convMytranIER(Numeric &mdh, const Index &dh)
void convHitranIERSH(Numeric &mdh, const Index &dh)
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.
int arts_omp_get_max_threads()
Wrapper for omp_get_max_threads.
std::map< String, Index > SpeciesMap
The map associated with species_data.
bool ReadFromStream(String &artsid, istream &is, Index nparams, const Verbosity &verbosity)
Read parameters from input stream.
Index nrows() const
Returns the number of rows.
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.
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.
const Array< SpeciesRecord > species_data
Species Data.
Numeric N_foreign(const Index i) const
ARTSCAT-4 foreign temperature exponents (dimensionless):
const Numeric & Cutoff() const
Return the cutoff frequency (in Hz).
Auxiliary data for isotopologues.
void extract(T &x, String &line, Index n)
Extract something from the beginning of a string.
const Index & Ind_lsn() const
Return the index of the normalization factor.
This can be used to make arrays out of anything.
Numeric Elow() const
Lower state energy in cm^-1:
Index species_index_from_species_name(String name)
Return species index for given species name.
bool arts_omp_in_parallel()
Wrapper for omp_in_parallel.
void checkIsotopologueRatios(const ArrayOfArrayOfSpeciesTag &abs_species, const SpeciesAuxData &isoratios)
Check that isotopologue ratios for the given species are correctly defined.
Numeric Psf() const
The pressure shift parameter in Hz/Pa.
const ArrayOfMatrix & getParams() const
Return a constant reference to the parameters.
String species_name_from_species_index(const Index spec_ind)
Return species name for given species index.
void iso(Array< IsotopologueRecord >::iterator &ii, String name, const ArrayOfNumeric &coeff, const Index &coefftype)
Numeric getParam(Index species, Index isotopologue, Index col) const
Get single parameter value.
Declarations having to do with the four output streams.
ostream & operator<<(ostream &os, const SpeciesRecord &sr)
Output operator for SpeciesRecord.
Index ncols() const
Returns the number of columns.
Structure to store a grid position for higher order interpolation.
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.
void initParams(Index nparams)
Resize according to builtin isotopologues in species data.
void convHitranIERF(Numeric &mdf, const Index &df)
Numeric CalculatePartitionFctRatio(Numeric reference_temperature, Numeric actual_temperature) const
Calculate partition function ratio.
Index Isotopologue() const
The index of the isotopologue species that this line belongs to.
Index nelem() const
Returns the number of elements.
Input manipulator class for doubles to enable nan and inf parsing.
const Array< LineshapeRecord > lineshape_data
The lookup data for the different lineshapes.
NUMERIC Numeric
The type to use for all floating point numbers.
const IsotopologueRecord & IsotopologueData() const
The matching IsotopologueRecord from species_data.
const Index & Speciesindex() const
void define_species_map()
Define the species data map.
Spectral line catalog data.
void trim()
Trim leading and trailing whitespace.
Index Species() const
The index of the molecular species that this line belongs to.
ConstIterator1D begin() const
Return const iterator to first element.
Contains the lookup data for one species.
Index nelem() const
Number of elements.
Header file for interpolation_poly.cc.
Lineshape related specification like which lineshape to use, the normalizationfactor,...
Numeric CalculatePartitionFctAtTempFromCoeff(Numeric temperature) const
A constant view of a Matrix.
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.
const Array< LineshapeNormRecord > lineshape_norm_data
The lookup data for the different normalization factors to the lineshapes.
const Array< IsotopologueRecord > & Isotopologue() const
void fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData(SpeciesAuxData &sad)
Fill SpeciesAuxData with default isotopologue ratios from species data.
Numeric wavenumber_to_joule(Numeric e)
A little helper function to convert energy from units of wavenumber (cm^-1) to Joule (J).
Header file for logic.cc.
const String & Name() const
Numeric Delta_foreign(const Index i) const
ARTSCAT-4 pressure shift parameters in Hz/Pa :
Index Version() const
Return the version number.
void setParam(Index species, Index isotopologue, Index col, Numeric v)
Set parameter.
ConstIterator1D end() const
Return const iterator behind last element.
const Numeric SPEED_OF_LIGHT
This file contains basic functions to handle ASCII files.
static String BroadSpecName(const Index i)
Return the name of an artscat-4 broadening species, as function of its broadening species index.
INDEX Index
The type to use for all integer numbers and indices.
Numeric CalculatePartitionFctAtTempFromData(Numeric temperature) const
void gridpos_poly(ArrayOfGridPosPoly &gp, ConstVectorView old_grid, ConstVectorView new_grid, const Index order, const Numeric &extpolfac)
Set up grid positions for higher order interpolation.
static Index NBroadSpec()
Return the number of artscat-4 foreign broadening species (6).
Numeric Tgam() const
Reference temperature for AGAM and SGAM in K:
Numeric F() const
The line center frequency in Hz.
bool is_sorted(ConstVectorView x)
Checks if a vector is sorted in ascending order.
A constant view of a Vector.
Index nelem() const
Number of elements.
Numeric I0() const
The line intensity in m^2*Hz at the reference temperature Ti0.
void interpweights(VectorView itw, const GridPos &tc)
Red 1D interpolation weights.
const Numeric BOLTZMAN_CONST
Global constant, the Boltzmann constant [J/K].
The global header file for ARTS.
Numeric Nair() const
AGAM temperature exponent (dimensionless):