ARTS  2.2.66
absorption.h
Go to the documentation of this file.
1 /* Copyright (C) 2000-2012
2  Stefan Buehler <sbuehler@ltu.se>
3  Axel von Engeln <engeln@uni-bremen.de>
4 
5  This program is free software; you can redistribute it and/or modify it
6  under the terms of the GNU General Public License as published by the
7  Free Software Foundation; either version 2, or (at your option) any
8  later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License
16  along with this program; if not, write to the Free Software
17  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
18  USA. */
19 
28 #ifndef absorption_h
29 #define absorption_h
30 
31 #include <stdexcept>
32 #include <cmath>
33 #include "matpackI.h"
34 #include "array.h"
35 #include "mystring.h"
36 #include "make_array.h"
37 #include "messages.h"
38 #include "abs_species_tags.h"
39 #include "linerecord.h"
40 #include "linemixingrecord.h"
41 
42 
45 typedef void (*lsf_type)(Vector&,
46  Vector&,
47  Vector&,
48  const Numeric,
49  const Numeric,
50  const Numeric,
52 
59 public:
60 
63  mdescription(),
64  mphase(),
65  mfunction()
66  { /* Nothing to do here. */ }
67 
69  LineshapeRecord(const String& name,
70  const String& description,
71  lsf_type function,
72  const bool phase)
73  : mname(name),
74  mdescription(description),
75  mphase(phase),
76  mfunction(function)
77  { /* Nothing to do here. */ }
79  const String& Name() const { return mname; }
81  const String& Description() const { return mdescription; }
83  lsf_type Function() const { return mfunction; }
85  bool Phase() const { return mphase; }
86 private:
89  bool mphase;
91 
92 };
93 
96 typedef void (*lsnf_type)(Vector&,
97  const Numeric,
99  const Numeric);
100 
108 public:
109 
112  mdescription(),
113  mfunction()
114  { /* Nothing to do here. */ }
115 
118  const String& description,
119  lsnf_type function)
120  : mname(name),
121  mdescription(description),
122  mfunction(function)
123  { /* Nothing to do here. */ }
125  const String& Name() const { return mname; }
127  const String& Description() const { return mdescription; }
129  lsnf_type Function() const { return mfunction; }
130 private:
134 };
135 
142 public:
143 
146  mind_lsn(-1),
147  mcutoff(0.)
148  { /* Nothing to do here. */ }
149 
151  LineshapeSpec(const Index& ind_ls,
152  const Index& ind_lsn,
153  const Numeric& cutoff)
154  : mind_ls(ind_ls),
155  mind_lsn(ind_lsn),
156  mcutoff(cutoff)
157  { /* Nothing to do here. */ }
158 
160  const Index& Ind_ls() const { return mind_ls; }
162  void SetInd_ls( Index ind_ls ) { mind_ls = ind_ls; }
163 
165  const Index& Ind_lsn() const { return mind_lsn; }
167  void SetInd_lsn( Index ind_lsn ) { mind_lsn = ind_lsn; }
168 
172  const Numeric& Cutoff() const { return mcutoff; }
174  void SetCutoff( Numeric cutoff ) { mcutoff = cutoff; }
175 private:
179 };
180 
181 ostream& operator<< (ostream& os, const LineshapeSpec& lsspec);
182 
186 
187 
188 
192 public:
193 
196  mabundance(0.),
197  mmass(0.),
198  mmytrantag(-1),
199  mhitrantag(-1),
200  mjpltags(),
201  mqcoeff(),
202  mqcoefftype(-1),
203  mqcoeffgrid(),
205  { /* Nothing left to do here. */ }
206 
210  mname(x.mname),
212  mmass(x.mmass),
215  mjpltags(x.mjpltags),
216  mqcoeff(),
217  mqcoefftype(),
218  mqcoeffgrid(),
220  { /* Nothing left to do here. */ }
221 
224  const Numeric& abundance,
225  const Numeric& mass,
226  const Index& mytrantag,
227  const Index& hitrantag,
228  const MakeArray<Index>& jpltags) :
229  mname(name),
230  mabundance(abundance),
231  mmass(mass),
232  mmytrantag(mytrantag),
233  mhitrantag(hitrantag),
234  mjpltags(jpltags),
235  mqcoeff(),
237  mqcoeffgrid(),
239  {
240  // With Matpack, initialization of mjpltags from jpltags should now work correctly.
241 
242  // Some consistency checks whether the given data makes sense.
243 #ifndef NDEBUG
244  {
245  /* 1. All the tags must be positive or -1 */
246  assert( (0<mmytrantag) || (-1==mmytrantag) );
247  assert( (0<mhitrantag) || (-1==mhitrantag) );
248  for ( Index i=0; i<mjpltags.nelem(); ++i )
249  assert( (0<mjpltags[i]) || (-1==mjpltags[i]) );
250  }
251 #endif // ifndef NDEBUG
252  }
253 
255  const String& Name() const { return mname; }
257  const Numeric& Abundance() const { return mabundance; }
260  const Numeric& Mass() const { return mmass; }
262  const Index& MytranTag() const { return mmytrantag; }
264  const Index& HitranTag() const { return mhitrantag; }
268  const ArrayOfIndex& JplTags() const { return mjpltags; }
269 
271 
274  bool isContinuum() const { return mname.length() && !isdigit(mname[0]); }
275 
276  void SetPartitionFctCoeff( const ArrayOfNumeric& qcoeff, const Index& qcoefftype )
277  {
278  mqcoeff = qcoeff;
279  mqcoefftype = qcoefftype;
280  }
281 
283 
296  Numeric actual_temperature ) const
297  {
298 
299  Numeric qcoeff_at_t_ref, qtemp;
300 
301  switch(mqcoefftype)
302  {
303  case PF_FROMCOEFF:
304  qcoeff_at_t_ref =
305  CalculatePartitionFctAtTempFromCoeff( reference_temperature );
306  qtemp =
307  CalculatePartitionFctAtTempFromCoeff( actual_temperature );
308  break;
309  case PF_FROMTEMP:
310  qcoeff_at_t_ref =
311  CalculatePartitionFctAtTempFromData( actual_temperature );
312  qtemp =
313  CalculatePartitionFctAtTempFromData( actual_temperature );
314  break;
315  default:
316  throw runtime_error("The partition functions are incorrect.\n");
317  break;
318  }
319 /* cout << "ref_t: " << reference_temperature << ", act_t:" <<
320  actual_temperature << "\n";
321  cout << "ref_q: " << qcoeff_at_t_ref << ", act_q:" <<
322  qtemp << "\n";
323 */
324  if ( qtemp > 0. )
325  return qcoeff_at_t_ref / qtemp;
326  else
327  {
328  ostringstream os;
329  os << "Partition function of "
330  << "Isotopologue " << mname
331 // << " is unknown.";
332  << " at T=" << actual_temperature << "K is zero or negative.";
333  throw runtime_error(os.str());
334  }
335  }
336 
337 
338  enum {
339  PF_FROMCOEFF, // Partition function will be from coefficients
340  PF_FROMTEMP, // Partition function will be from temperature field
341  PF_NOTHING // This will be the designated starter value
342  };
343 
344 private:
345 
346  // calculate the partition fct at a certain temperature
347  // this is only the prototyping
350 
361 };
362 
363 
368 public:
369 
372  mdegfr(-1),
373  misotopologue() { /* Nothing to do here */ }
374 
376  SpeciesRecord(const char name[],
377  const Index degfr,
378  const MakeArray<IsotopologueRecord>& isotopologue)
379  : mname(name),
380  mdegfr(degfr),
381  misotopologue(isotopologue)
382  {
383 
384  // Thanks to Matpack, initialization of misotopologue with isotopologue
385  // should now work correctly.
386 
387 #ifndef NDEBUG
388  {
389  /* Check that the isotopologues are correctly sorted. */
390  for ( Index i=0; i<misotopologue.nelem()-1; ++i )
391  {
392  assert(isnan(misotopologue[i].Abundance()) || isnan(misotopologue[i+1].Abundance())
393  || misotopologue[i].Abundance() >= misotopologue[i+1].Abundance());
394  }
395 
396  /* Check that the Mytran tags are correctly sorted. */
397  for ( Index i=0; i<misotopologue.nelem()-1; ++i )
398  {
399  if ( (0<misotopologue[i].MytranTag()) && (0<misotopologue[i+1].MytranTag()) )
400  {
401  assert( misotopologue[i].MytranTag() < misotopologue[i+1].MytranTag() );
402 
403  // Also check that the tags have the same base number:
404  assert( misotopologue[i].MytranTag()/10 == misotopologue[i].MytranTag()/10 );
405  }
406  }
407 
408  /* Check that the Hitran tags are correctly sorted. */
409  for ( Index i=0; i<misotopologue.nelem()-1; ++i )
410  {
411  if ( (0<misotopologue[i].HitranTag()) && (0<misotopologue[i+1].HitranTag()) )
412  {
413 // assert( misotopologue[i].HitranTag() < misotopologue[i+1].HitranTag() );
414 
415  // Also check that the tags have the same base number:
416  assert( misotopologue[i].HitranTag()/10 == misotopologue[i+1].HitranTag()/10 );
417  }
418  }
419  }
420 #endif // #ifndef NDEBUG
421  }
422 
423  const String& Name() const { return mname; }
424  Index Degfr() const { return mdegfr; }
427 
428 private:
435 };
436 
437 
440 {
441 public:
444 
446  void initParams(Index nparams);
447 
449  Numeric getParam(Index species, Index isotopologue, Index col) const
450  {
451  return mparams[species](isotopologue, col);
452  }
453 
455  void setParam(Index species, Index isotopologue, Index col, Numeric v)
456  {
457  mparams[species](isotopologue, col) = v;
458  }
459 
461  const ArrayOfMatrix& getParams() const { return mparams; }
462 
464  bool ReadFromStream(String& artsid, istream& is, Index nparams,
465  const Verbosity& verbosity);
466 
467 private:
469 };
470 
471 
473 void checkIsotopologueRatios(const ArrayOfArrayOfSpeciesTag& abs_species,
474  const SpeciesAuxData& sad);
475 
478 
479 
480 // is needed to map jpl tags/arts identifier to the species/isotopologue data within arts
482 public:
484  SpecIsoMap(const Index& speciesindex,
485  const Index& isotopologueindex)
486  : mspeciesindex(speciesindex),
487  misotopologueindex(isotopologueindex)
488  {}
489 
490  // Return the index to the species
491  const Index& Speciesindex() const { return mspeciesindex; }
492  // Return the index to the isotopologue
493  const Index& Isotopologueindex() const { return misotopologueindex; }
494 
495 private:
498 };
499 
500 
501 
506 ostream& operator<< (ostream& os, const SpeciesRecord& sr);
507 
510 ostream& operator<< (ostream& os, const SpeciesAuxData& sad);
511 
512 
513 
517 void define_species_map();
518 
519 
520 void xsec_species(MatrixView xsec_attenuation,
521  MatrixView xsec_phase,
522  ConstVectorView f_grid,
523  ConstVectorView abs_p,
524  ConstVectorView abs_t,
525  ConstMatrixView all_vmrs,
526  const ArrayOfArrayOfSpeciesTag& abs_species,
527  const Index this_species,
528  const ArrayOfLineRecord& abs_lines,
529  const Index ind_ls,
530  const Index ind_lsn,
531  const Numeric cutoff,
532  const SpeciesAuxData& isotopologue_ratios,
533  const Verbosity& verbosity );
534 
535 
536 void xsec_species_line_mixing_wrapper( MatrixView xsec_attenuation,
537  MatrixView xsec_phase,
538  const ArrayOfArrayOfLineMixingRecord& line_mixing_data,
539  const ArrayOfArrayOfIndex& line_mixing_data_lut,
540  ConstVectorView f_grid,
541  ConstVectorView abs_p,
542  ConstVectorView abs_t,
543  ConstMatrixView all_vmrs,
544  const ArrayOfArrayOfSpeciesTag& abs_species,
545  const Index this_species,
546  const ArrayOfLineRecord& abs_lines,
547  const Index ind_ls,
548  const Index ind_lsn,
549  const Numeric cutoff,
550  const SpeciesAuxData& isotopologue_ratios,
551  const Verbosity& verbosity );
552 
553 
554 void xsec_species_line_mixing_2nd_order( MatrixView xsec_attenuation,
555  MatrixView xsec_phase,
556  const ArrayOfArrayOfLineMixingRecord& line_mixing_data,
557  const ArrayOfArrayOfIndex& line_mixing_data_lut,
558  ConstVectorView f_grid,
559  ConstVectorView abs_p,
560  ConstVectorView abs_t,
561  ConstMatrixView all_vmrs,
562  const ArrayOfArrayOfSpeciesTag& abs_species,
563  const Index this_species,
564  const ArrayOfLineRecord& abs_lines,
565  const Index ind_ls,
566  const Index ind_lsn,
567  const Numeric cutoff,
568  const SpeciesAuxData& isotopologue_ratios,
569  const Verbosity& verbosity );
570 
571 
572 // A helper function for energy conversion:
574 
575 
576 //======================================================================
577 // Functions related to species
578 //======================================================================
579 
581 
583 
584 
585 //======================================================================
586 // Functions to convert the accuracy index
587 //======================================================================
588 
589 // ********* for HITRAN database *************
590 // convert index for the frequency accuracy.
591 void convHitranIERF(
592  Numeric& mdf,
593  const Index& df
594  );
595 
596 // convert to percents index for intensity and halfwidth accuracy.
597 
598 void convHitranIERSH(
599  Numeric& mdh,
600  const Index& dh
601  );
602 
603 // ********* for MYTRAN database *************
604 // convert index for the halfwidth accuracy.
605 void convMytranIER(
606  Numeric& mdh,
607  const Index & dh
608  );
609 
610 
611 // Functions to set abs_n2 and abs_h2o:
612 
613 void abs_n2Set(Vector& abs_n2,
614  const ArrayOfArrayOfSpeciesTag& abs_species,
615  const Matrix& abs_vmrs,
616  const Verbosity&);
617 
618 void abs_h2oSet(Vector& abs_h2o,
619  const ArrayOfArrayOfSpeciesTag& abs_species,
620  const Matrix& abs_vmrs,
621  const Verbosity&);
622 
623 #endif // absorption_h
IsotopologueRecord::IsotopologueRecord
IsotopologueRecord(const String &name, const Numeric &abundance, const Numeric &mass, const Index &mytrantag, const Index &hitrantag, const MakeArray< Index > &jpltags)
Constructor that sets the values.
Definition: absorption.h:223
fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData
void fillSpeciesAuxDataWithIsotopologueRatiosFromSpeciesData(SpeciesAuxData &sad)
Fill SpeciesAuxData with default isotopologue ratios from species data.
Definition: absorption.cc:305
Matrix
The Matrix class.
Definition: matpackI.h:788
IsotopologueRecord::Mass
const Numeric & Mass() const
Mass of the isotopologue.
Definition: absorption.h:260
LineshapeNormRecord::Description
const String & Description() const
Return the description text.
Definition: absorption.h:127
IsotopologueRecord::IsotopologueRecord
IsotopologueRecord(const IsotopologueRecord &x)
Copy constructor.
Definition: absorption.h:209
SpeciesRecord::Degfr
Index Degfr() const
Definition: absorption.h:424
MatrixView
The MatrixView class.
Definition: matpackI.h:679
ArrayOfLineshapeSpec
Array< LineshapeSpec > ArrayOfLineshapeSpec
Holds a list of lineshape specifications: function, normalization, cutoff.
Definition: absorption.h:185
IsotopologueRecord::SetPartitionFctCoeff
void SetPartitionFctCoeff(const ArrayOfNumeric &qcoeff, const Index &qcoefftype)
Definition: absorption.h:276
SpecIsoMap
Definition: absorption.h:481
LineshapeRecord::Phase
bool Phase() const
Returns true if lineshape function calculates phase information.
Definition: absorption.h:85
LineshapeNormRecord::mfunction
lsnf_type mfunction
Pointer to lineshape normalization function.
Definition: absorption.h:133
SpeciesAuxData::SpeciesAuxData
SpeciesAuxData()
Default constructor.
Definition: absorption.h:443
SpecIsoMap::misotopologueindex
Index misotopologueindex
Definition: absorption.h:497
LineshapeRecord::mfunction
lsf_type mfunction
Pointer to lineshape function.
Definition: absorption.h:90
LineshapeRecord::Description
const String & Description() const
Return the description text.
Definition: absorption.h:81
LineshapeNormRecord::Function
lsnf_type Function() const
Return pointer to lineshape normalization function.
Definition: absorption.h:129
LineshapeSpec::Ind_ls
const Index & Ind_ls() const
Return the index of this lineshape.
Definition: absorption.h:160
LineshapeNormRecord::Name
const String & Name() const
Return the name of this lineshape.
Definition: absorption.h:125
IsotopologueRecord::mhitrantag
Index mhitrantag
Definition: absorption.h:355
SpeciesRecord::misotopologue
Array< IsotopologueRecord > misotopologue
Isotopologue data.
Definition: absorption.h:434
LineshapeRecord::mphase
bool mphase
Does this lineshape calculate phase information?
Definition: absorption.h:89
IsotopologueRecord::PF_FROMCOEFF
@ PF_FROMCOEFF
Definition: absorption.h:339
LineshapeNormRecord::LineshapeNormRecord
LineshapeNormRecord(const String &name, const String &description, lsnf_type function)
Initializing constructor, used to build the lookup table.
Definition: absorption.h:117
SpecIsoMap::Isotopologueindex
const Index & Isotopologueindex() const
Definition: absorption.h:493
LineshapeRecord
Lineshape related information.
Definition: absorption.h:58
SpeciesAuxData::ReadFromStream
bool ReadFromStream(String &artsid, istream &is, Index nparams, const Verbosity &verbosity)
Read parameters from input stream.
Definition: absorption.cc:104
species_index_from_species_name
Index species_index_from_species_name(String name)
Return species index for given species name.
Definition: absorption.cc:1260
IsotopologueRecord::mqcoeffgrid
Vector mqcoeffgrid
Definition: absorption.h:359
SpeciesAuxData::mparams
ArrayOfMatrix mparams
Definition: absorption.h:468
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.
Definition: absorption.cc:1530
LineshapeRecord::Function
lsf_type Function() const
Return pointer to lineshape function.
Definition: absorption.h:83
array.h
This file contains the definition of Array.
LineshapeRecord::LineshapeRecord
LineshapeRecord(const String &name, const String &description, lsf_type function, const bool phase)
Initializing constructor, used to build the lookup table.
Definition: absorption.h:69
SpeciesRecord::Isotopologue
Array< IsotopologueRecord > & Isotopologue()
Definition: absorption.h:426
LineshapeSpec::Cutoff
const Numeric & Cutoff() const
Return the cutoff frequency (in Hz).
Definition: absorption.h:172
SpeciesAuxData
Auxiliary data for isotopologues.
Definition: absorption.h:440
SpecIsoMap::SpecIsoMap
SpecIsoMap(const Index &speciesindex, const Index &isotopologueindex)
Definition: absorption.h:484
IsotopologueRecord::mmass
Numeric mmass
Definition: absorption.h:353
matpackI.h
LineshapeSpec::Ind_lsn
const Index & Ind_lsn() const
Return the index of the normalization factor.
Definition: absorption.h:165
IsotopologueRecord::HitranTag
const Index & HitranTag() const
HITRAN-96 tag numbers for all isotopologues.
Definition: absorption.h:264
convHitranIERF
void convHitranIERF(Numeric &mdf, const Index &df)
Definition: absorption.cc:1327
Array
This can be used to make arrays out of anything.
Definition: array.h:107
convMytranIER
void convMytranIER(Numeric &mdh, const Index &dh)
Definition: absorption.cc:1433
LineshapeRecord::mdescription
String mdescription
Short description.
Definition: absorption.h:88
abs_h2oSet
void abs_h2oSet(Vector &abs_h2o, const ArrayOfArrayOfSpeciesTag &abs_species, const Matrix &abs_vmrs, const Verbosity &)
abs_h2oSet.
Definition: m_abs.cc:1287
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).
Definition: absorption.cc:1227
operator<<
ostream & operator<<(ostream &os, const LineshapeSpec &lsspec)
Definition: absorption.cc:1494
SpeciesAuxData::getParams
const ArrayOfMatrix & getParams() const
Return a constant reference to the parameters.
Definition: absorption.h:461
SpeciesAuxData::getParam
Numeric getParam(Index species, Index isotopologue, Index col) const
Get single parameter value.
Definition: absorption.h:449
convHitranIERSH
void convHitranIERSH(Numeric &mdh, const Index &dh)
Definition: absorption.cc:1374
messages.h
Declarations having to do with the four output streams.
SpeciesRecord::SpeciesRecord
SpeciesRecord()
Default constructor.
Definition: absorption.h:371
IsotopologueRecord::mqcoefftype
Index mqcoefftype
Definition: absorption.h:358
my_basic_string< char >
IsotopologueRecord::mqcoeff
Vector mqcoeff
Definition: absorption.h:357
SpeciesAuxData::initParams
void initParams(Index nparams)
Resize according to builtin isotopologues in species data.
Definition: absorption.cc:90
IsotopologueRecord::CalculatePartitionFctRatio
Numeric CalculatePartitionFctRatio(Numeric reference_temperature, Numeric actual_temperature) const
Calculate partition function ratio.
Definition: absorption.h:295
IsotopologueRecord::Abundance
const Numeric & Abundance() const
Normal abundance ( = isotopologue ratio).
Definition: absorption.h:257
define_species_map
void define_species_map()
Define the species data map.
Definition: absorption.cc:322
LineshapeSpec::SetInd_ls
void SetInd_ls(Index ind_ls)
Set it.
Definition: absorption.h:162
linemixingrecord.h
LineMixingRecord class for storing line mixing data.
IsotopologueRecord::JplTags
const ArrayOfIndex & JplTags() const
JPL tag numbers for all isotopologues.
Definition: absorption.h:268
lsf_type
void(* lsf_type)(Vector &, Vector &, Vector &, const Numeric, const Numeric, const Numeric, ConstVectorView)
The type that is used to store pointers to lineshape functions.
Definition: absorption.h:45
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:29
Verbosity
Definition: messages.h:50
SpecIsoMap::Speciesindex
const Index & Speciesindex() const
Definition: absorption.h:491
IsotopologueRecord::mabundance
Numeric mabundance
Definition: absorption.h:352
make_array.h
Implements the class MakeArray, which is a derived class of Array, allowing explicit initialization.
IsotopologueRecord::PF_NOTHING
@ PF_NOTHING
Definition: absorption.h:341
LineshapeSpec::SetCutoff
void SetCutoff(Numeric cutoff)
Set it.
Definition: absorption.h:174
abs_species_tags.h
Header file for stuff related to absorption species tags.
SpecIsoMap::mspeciesindex
Index mspeciesindex
Definition: absorption.h:496
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.
Definition: absorption.cc:1598
SpeciesRecord
Contains the lookup data for one species.
Definition: absorption.h:367
LineshapeSpec
Lineshape related specification like which lineshape to use, the normalizationfactor,...
Definition: absorption.h:141
IsotopologueRecord::CalculatePartitionFctAtTempFromCoeff
Numeric CalculatePartitionFctAtTempFromCoeff(Numeric temperature) const
Definition: absorption.cc:56
ConstMatrixView
A constant view of a Matrix.
Definition: matpackI.h:596
IsotopologueRecord::PF_FROMTEMP
@ PF_FROMTEMP
Definition: absorption.h:340
SpeciesRecord::mdegfr
Index mdegfr
Degrees of freedom.
Definition: absorption.h:432
LineshapeRecord::mname
String mname
Name of the function (e.g., Lorentz).
Definition: absorption.h:87
lsnf_type
void(* lsnf_type)(Vector &, const Numeric, ConstVectorView, const Numeric)
The type that is used to store pointers to lineshape normalization functions.
Definition: absorption.h:96
IsotopologueRecord::Name
const String & Name() const
Isotopologue name.
Definition: absorption.h:255
SpeciesRecord::Isotopologue
const Array< IsotopologueRecord > & Isotopologue() const
Definition: absorption.h:425
SpeciesRecord::mname
String mname
Species name.
Definition: absorption.h:430
LineshapeRecord::Name
const String & Name() const
Return the name of this lineshape.
Definition: absorption.h:79
IsotopologueRecord::mmytrantag
Index mmytrantag
Definition: absorption.h:354
LineshapeSpec::SetInd_lsn
void SetInd_lsn(Index ind_lsn)
Set it.
Definition: absorption.h:167
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.
Definition: absorption.cc:594
IsotopologueRecord::mname
String mname
Definition: absorption.h:351
IsotopologueRecord
Contains the lookup data for one isotopologue.
Definition: absorption.h:191
IsotopologueRecord::MytranTag
const Index & MytranTag() const
MYTRAN2 tag numbers for all isotopologues.
Definition: absorption.h:262
SpeciesRecord::Name
const String & Name() const
Definition: absorption.h:423
LineshapeNormRecord::mdescription
String mdescription
Short description.
Definition: absorption.h:132
checkIsotopologueRatios
void checkIsotopologueRatios(const ArrayOfArrayOfSpeciesTag &abs_species, const SpeciesAuxData &sad)
Check that isotopologue ratios for the given species are correctly defined.
Definition: absorption.cc:245
SpecIsoMap::SpecIsoMap
SpecIsoMap()
Definition: absorption.h:483
IsotopologueRecord::mjpltags
ArrayOfIndex mjpltags
Definition: absorption.h:356
LineshapeNormRecord::LineshapeNormRecord
LineshapeNormRecord()
Default constructor.
Definition: absorption.h:111
SpeciesAuxData::setParam
void setParam(Index species, Index isotopologue, Index col, Numeric v)
Set parameter.
Definition: absorption.h:455
LineshapeSpec::LineshapeSpec
LineshapeSpec(const Index &ind_ls, const Index &ind_lsn, const Numeric &cutoff)
Initializing constructor.
Definition: absorption.h:151
IsotopologueRecord::IsotopologueRecord
IsotopologueRecord()
Default constructor.
Definition: absorption.h:195
LineshapeSpec::LineshapeSpec
LineshapeSpec()
Default constructor.
Definition: absorption.h:145
species_name_from_species_index
String species_name_from_species_index(const Index spec_ind)
Return species name for given species index.
Definition: absorption.cc:1301
linerecord.h
LineRecord class for managing line catalog data.
LineshapeSpec::mind_ls
Index mind_ls
Definition: absorption.h:176
MakeArray
Explicit construction of Arrays.
Definition: make_array.h:52
abs_n2Set
void abs_n2Set(Vector &abs_n2, const ArrayOfArrayOfSpeciesTag &abs_species, const Matrix &abs_vmrs, const Verbosity &)
abs_n2Set.
Definition: m_abs.cc:1315
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:35
IsotopologueRecord::CalculatePartitionFctAtTempFromData
Numeric CalculatePartitionFctAtTempFromData(Numeric temperature) const
Definition: absorption.cc:79
LineshapeNormRecord::mname
String mname
Name of the function (e.g., linear).
Definition: absorption.h:131
Vector
The Vector class.
Definition: matpackI.h:556
LineshapeRecord::LineshapeRecord
LineshapeRecord()
Default constructor.
Definition: absorption.h:62
SpeciesRecord::SpeciesRecord
SpeciesRecord(const char name[], const Index degfr, const MakeArray< IsotopologueRecord > &isotopologue)
The constructor used in define_species_data.
Definition: absorption.h:376
LineshapeSpec::mcutoff
Numeric mcutoff
Definition: absorption.h:178
ConstVectorView
A constant view of a Vector.
Definition: matpackI.h:292
LineshapeNormRecord
Lineshape related normalization function information.
Definition: absorption.h:107
Array::nelem
Index nelem() const
Number of elements.
Definition: array.h:176
LineshapeSpec::mind_lsn
Index mind_lsn
Definition: absorption.h:177
mystring.h
This file contains the definition of String, the ARTS string class.
IsotopologueRecord::mqcoeffinterporder
Index mqcoeffinterporder
Definition: absorption.h:360
IsotopologueRecord::isContinuum
bool isContinuum() const
Check if isotopologue is actually a continuum.
Definition: absorption.h:274