ARTS  1.0.222
absorption.cc
Go to the documentation of this file.
1 /* Copyright (C) 2000, 2001 Stefan Buehler <sbuehler@uni-bremen.de>
2  Axel von Engeln <engeln@uni-bremen.de>
3 
4  This program is free software; you can redistribute it and/or modify it
5  under the terms of the GNU General Public License as published by the
6  Free Software Foundation; either version 2, or (at your option) any
7  later version.
8 
9  This program is distributed in the hope that it will be useful,
10  but WITHOUT ANY WARRANTY; without even the implied warranty of
11  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  GNU General Public License for more details.
13 
14  You should have received a copy of the GNU General Public License
15  along with this program; if not, write to the Free Software
16  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
17  USA. */
18 
30 #include "arts.h"
31 #include <map>
32 #include <cfloat>
33 #include <algorithm>
34 #include <cmath>
35 #include "absorption.h"
36 #include "math_funcs.h"
37 #include "auto_md.h"
38 #include "messages.h"
39 
40 
42 std::map<String, Index> SpeciesMap;
43 
44 
45 // member fct of isotoperecord, calculates the partition fct at the
46 // given temperature from the partition fct coefficients (3rd order
47 // polynomial in temperature)
49  temperature ) const
50 {
51  Numeric result = 0.;
52  Numeric exponent = 1.;
53 
54  ArrayOfNumeric::const_iterator it;
55 
56  for (it=mqcoeff.begin(); it != mqcoeff.end(); it++)
57  {
58  result += *it * exponent;
59  exponent *= temperature;
60  }
61  return result;
62 }
63 
65 {
67  extern std::map<String, Index> SpeciesMap;
68 
69  for ( Index i=0 ; i<species_data.nelem() ; ++i)
70  {
71  SpeciesMap[species_data[i].Name()] = i;
72  }
73 }
74 
75 
76 ostream& operator << (ostream& os, const LineRecord& lr)
77 {
78  // Determine the precision, depending on whether Numeric is double
79  // or float:
80  Index precision;
81  switch (sizeof(Numeric)) {
82  case sizeof(float) : precision = FLT_DIG; break;
83  case sizeof(double) : precision = DBL_DIG; break;
84  default: out0 << "Numeric must be double or float\n"; exit(1);
85  }
86 
87  os << "@"
88  << " " << lr.Name ()
89  << " "
90  << setprecision((int)precision)
91  << lr.F ()
92  << " " << lr.Psf ()
93  << " " << lr.I0 ()
94  << " " << lr.Ti0 ()
95  << " " << lr.Elow ()
96  << " " << lr.Agam ()
97  << " " << lr.Sgam ()
98  << " " << lr.Nair ()
99  << " " << lr.Nself ()
100  << " " << lr.Tgam ()
101  << " " << lr.Naux ()
102  << " " << lr.dF ()
103  << " " << lr.dI0 ()
104  << " " << lr.dAgam ()
105  << " " << lr.dSgam ()
106  << " " << lr.dNair ()
107  << " " << lr.dNself()
108  << " " << lr.dPsf ();
109 
110  // Added new lines for the spectroscopic parameters accuracies.
111  for ( Index i=0; i<lr.Naux(); ++i )
112  os << " " << lr.Aux()[i];
113 
114  return os;
115 }
116 
117 
126 template<class T>
127 void extract(T& x,
128  String& line,
129  Index n)
130 {
131  // Initialize output to zero! This is important, because otherwise
132  // the output variable could `remember' old values.
133  x = T(0);
134 
135  // This will contain the short subString with the item to extract.
136  // Make it a String stream, for easy parsing,
137  // extracting subString of width n from line:
138  istringstream item( line.substr(0,n) );
139 
140 // cout << "line = '" << line << "'\n";
141 // cout << "line.substr(0,n) = " << line.substr(0,n) << endl;
142 // cout << "item = " << item.str() << endl;
143 
144  // Shorten line by n:
145  line.erase(0,n);
146 // cout << "line = " << line << endl;
147 
148  // Convert with the aid of String stream item:
149  item >> x;
150 }
151 
152 
153 // This function reads line data in the Hitran 1986-2001 format. For the Hitran
154 // 2004 data format use ReadFromHitran2004Stream.
155 //
156 // 2005/03/29: A check for the right data record length (100 characters) has
157 // been added by Hermann Berg (h.berg@utoronto.ca)
158 //
160 {
161  // Global species lookup data:
163 
164  // This value is used to flag missing data both in species and
165  // isotope lists. Could be any number, it just has to be made sure
166  // that it is neither the index of a species nor of an isotope.
167  const Index missing = species_data.nelem() + 100;
168 
169  // We need a species index sorted by HITRAN tag. Keep this in a
170  // static variable, so that we have to do this only once. The ARTS
171  // species index is hind[<HITRAN tag>].
172  //
173  // Allow for up to 100 species in HITRAN in the future.
174  static Array< Index > hspec(100);
175 
176  // This is an array of arrays for each hitran tag. It contains the
177  // ARTS indices of the HITRAN isotopes.
178  static Array< ArrayOfIndex > hiso(100);
179 
180  // Remeber if this stuff has already been initialized:
181  static bool hinit = false;
182 
183  // Remember, about which missing species we have already issued a
184  // warning:
185  static ArrayOfIndex warned_missing;
186 
187  if ( !hinit )
188  {
189  // Initialize hspec.
190  // The value of missing means that we don't have this species.
191  hspec = missing; // Matpack can set all elements like this.
192 
193  for ( Index i=0; i<species_data.nelem(); ++i )
194  {
195  const SpeciesRecord& sr = species_data[i];
196 
197  // We have to be careful and check for the case that all
198  // HITRAN isotope tags are -1 (this species is missing in HITRAN).
199 
200  if ( 0 < sr.Isotope()[0].HitranTag() )
201  {
202  // The HITRAN tags are stored as species plus isotope tags
203  // (MO and ISO)
204  // in the Isotope() part of the species record.
205  // We can extract the MO part from any of the isotope tags,
206  // so we use the first one. We do this by taking an integer
207  // division by 10.
208 
209  Index mo = sr.Isotope()[0].HitranTag() / 10;
210  // cout << "mo = " << mo << endl;
211  hspec[mo] = i;
212 
213  // Get a nicer to handle array of HITRAN iso tags:
214  Index n_iso = sr.Isotope().nelem();
215  ArrayOfIndex iso_tags;
216  iso_tags.resize(n_iso);
217  for ( Index j=0; j<n_iso; ++j )
218  {
219  iso_tags[j] = sr.Isotope()[j].HitranTag();
220  }
221 
222  // Reserve elements for the isotope tags. How much do we
223  // need? This depends on the largest HITRAN tag that we know
224  // about!
225  // Also initialize the tags to missing.
226  // cout << "iso_tags = " << iso_tags << endl;
227  // cout << "static_cast<Index>(max(iso_tags))%10 + 1 = "
228  // << static_cast<Index>(max(iso_tags))%10 + 1 << endl;
229  hiso[mo].resize( max(iso_tags)%10 + 1 );
230  hiso[mo] = missing; // Matpack can set all elements like this.
231 
232 
233  // Set the isotope tags:
234  for ( Index j=0; j<n_iso; ++j )
235  {
236  if ( 0 < iso_tags[j] ) // ignore -1 elements
237  {
238  // To get the iso tags from HitranTag() we also have to take
239  // modulo 10 to get rid of mo.
240  hiso[mo][iso_tags[j] % 10] = j;
241  }
242  }
243  }
244  }
245 
246 
247  // Print the generated data structures (for debugging):
248  out3 << " HITRAN index table:\n";
249  for ( Index i=0; i<hspec.nelem(); ++i )
250  {
251  if ( missing != hspec[i] )
252  {
253  // The explicit conversion of Name to a c-String is
254  // necessary, because setw does not work correctly for
255  // stl Strings.
256  out3 << " mo = " << i << " Species = "
257  << setw(10) << setiosflags(ios::left)
258  << species_data[hspec[i]].Name().c_str()
259  << "iso = ";
260  for ( Index j=1; j<hiso[i].nelem(); ++j )
261  {
262  if ( missing==hiso[i][j] )
263  out3 << " " << "m";
264  else
265  out3 << " " << species_data[hspec[i]].Isotope()[hiso[i][j]].Name();
266  }
267  out3 << "\n";
268  }
269  }
270 
271  hinit = true;
272  }
273 
274 
275  // This contains the rest of the line to parse. At the beginning the
276  // entire line. Line gets shorter and shorter as we continue to
277  // extract stuff from the beginning.
278  String line;
279 
280  // The first item is the molecule number:
281  Index mo;
282 
283  // Look for more comments?
284  bool comment = true;
285 
286  while (comment)
287  {
288  // Return true if eof is reached:
289  if (is.eof()) return true;
290 
291  // Throw runtime_error if stream is bad:
292  if (!is) throw runtime_error ("Stream bad.");
293 
294  // Read line from file into linebuffer:
295  getline(is,line);
296 
297  // It is possible that we were exactly at the end of the file before
298  // calling getline. In that case the previous eof() was still false
299  // because eof() evaluates only to true if one tries to read after the
300  // end of the file. The following check catches this.
301  if (line.nelem() == 0 && is.eof()) return true;
302 
303  // If the catalogue is in dos encoding, throw away the
304  // additional carriage return
305  if (line[line.nelem () - 1] == 13)
306  {
307  line.erase (line.nelem () - 1, 1);
308  }
309 
310  // Because of the fixed FORTRAN format, we need to break up the line
311  // explicitly in apropriate pieces. Not elegant, but works!
312 
313  // Extract molecule number:
314  mo = 0;
315  // Initialization of mo is important, because mo stays the same
316  // if line is empty.
317  extract(mo,line,2);
318  // cout << "mo = " << mo << endl;
319 
320  // If mo == 0 this is just a comment line:
321  if ( 0 != mo )
322  {
323  // See if we know this species. Exit with an error if the species is unknown.
324  if ( missing != hspec[mo] )
325  {
326  comment = false;
327 
328  // Check if data record has the right number of characters for the
329  // in Hitran 1986-2001 format
330  Index nChar = line.nelem() + 2; // number of characters in data record;
331  if ( nChar != 100 )
332  {
333  ostringstream os;
334  os << "Invalid HITRAN 1986-2001 line data record with " << nChar <<
335  " characters (expected: 100)." << endl << line << " n: " << line.nelem ();
336  throw runtime_error(os.str());
337  }
338 
339  }
340  else
341  {
342  // See if this is already in warned_missing, use
343  // std::count for that:
344  if ( 0 == std::count(warned_missing.begin(),
345  warned_missing.end(),
346  mo) )
347  {
348  out0 << "Error: HITRAN mo = " << mo << " is not "
349  << "known to ARTS.\n";
350  warned_missing.push_back(mo);
351  exit(1);
352  // SAB 08.08.2000 If you want to make the program
353  // continue anyway, just comment out the exit
354  // line.
355  }
356  }
357  }
358  }
359 
360  // Ok, we seem to have a valid species here.
361 
362  // Set mspecies from my cool index table:
363  mspecies = hspec[mo];
364 
365  // Extract isotope:
366  Index iso;
367  extract(iso,line,1);
368  // cout << "iso = " << iso << endl;
369 
370 
371  // Set misotope from the other cool index table.
372  // We have to be careful to issue an error for unknown iso tags. Iso
373  // could be either larger than the size of hiso[mo], or set
374  // explicitly to missing. Unfortunately we have to test both cases.
375  misotope = missing;
376  if ( iso < hiso[mo].nelem() )
377  if ( missing != hiso[mo][iso] )
378  misotope = hiso[mo][iso];
379 
380  // Issue error message if misotope is still missing:
381  if (missing == misotope)
382  {
383  ostringstream os;
384  os << "Species: " << species_data[mspecies].Name()
385  << ", isotope iso = " << iso
386  << " is unknown.";
387  throw runtime_error(os.str());
388  }
389 
390 
391  // Position.
392  {
393  // HITRAN position in wavenumbers (cm^-1):
394  Numeric v;
395  // External constant from constants.cc:
396  extern const Numeric SPEED_OF_LIGHT;
397  // Conversion from wavenumber to Hz. If you multiply a line
398  // position in wavenumber (cm^-1) by this constant, you get the
399  // frequency in Hz.
400  const Numeric w2Hz = SPEED_OF_LIGHT * 100.;
401 
402  // Extract HITRAN postion:
403  extract(v,line,12);
404 
405  // ARTS position in Hz:
406  mf = v * w2Hz;
407 // cout << "mf = " << mf << endl;
408  }
409 
410  // Intensity.
411  {
412  extern const Numeric SPEED_OF_LIGHT; // in [m/s]
413 
414  // HITRAN intensity is in cm-1/(molec * cm-2) at 296 Kelvin.
415  // It already includes the isotpic ratio.
416  // The first cm-1 is the frequency unit (it cancels with the
417  // 1/frequency unit of the line shape function).
418  //
419  // We need to do the following:
420  // 1. Convert frequency from wavenumber to Hz (factor 1e2 * c).
421  // 2. Convert [molec * cm-2] to [molec * m-2] (factor 1e-4).
422  // 3. Take out the isotopic ratio.
423 
424  const Numeric hi2arts = 1e-2 * SPEED_OF_LIGHT;
425 
426  Numeric s;
427 
428  // Extract HITRAN intensity:
429  extract(s,line,10);
430  // Convert to ARTS units (Hz / (molec * m-2) ), or shorter: Hz*m^2
431  mi0 = s * hi2arts;
432  // Take out isotopic ratio:
433  mi0 /= species_data[mspecies].Isotope()[misotope].Abundance();
434  }
435 
436  // Skip transition probability:
437  {
438  Numeric r;
439  extract(r,line,10);
440  }
441 
442 
443  // Air broadening parameters.
444  {
445  // HITRAN parameter is in cm-1/atm at 296 Kelvin
446  // All parameters are HWHM (I hope this is true!)
447  Numeric gam;
448  // External constant from constants.cc: Converts atm to
449  // Pa. Multiply value in atm by this number to get value in Pa.
450  extern const Numeric ATM2PA;
451  // External constant from constants.cc:
452  extern const Numeric SPEED_OF_LIGHT;
453  // Conversion from wavenumber to Hz. If you multiply a value in
454  // wavenumber (cm^-1) by this constant, you get the value in Hz.
455  const Numeric w2Hz = SPEED_OF_LIGHT * 1e2;
456  // Ok, put together the end-to-end conversion that we need:
457  const Numeric hi2arts = w2Hz / ATM2PA;
458 
459  // Extract HITRAN AGAM value:
460  extract(gam,line,5);
461 
462  // ARTS parameter in Hz/Pa:
463  magam = gam * hi2arts;
464 
465  // Extract HITRAN SGAM value:
466  extract(gam,line,5);
467 
468  // ARTS parameter in Hz/Pa:
469  msgam = gam * hi2arts;
470 
471  // If zero, set to agam:
472  if (0==msgam)
473  msgam = magam;
474 
475  // cout << "agam, sgam = " << magam << ", " << msgam << endl;
476  }
477 
478 
479  // Lower state energy.
480  {
481  // HITRAN parameter is in wavenumbers (cm^-1).
482  // We have to convert this to the ARTS unit Joule.
483 
484  // Extract from Catalogue line
485  extract(melow,line,10);
486 
487  // Convert to Joule:
489  }
490 
491 
492  // Temperature coefficient of broadening parameters.
493  {
494  // This is dimensionless, we can also extract directly.
495  extract(mnair,line,4);
496 
497  // Set self broadening temperature coefficient to the same value:
498  mnself = mnair;
499 // cout << "mnair = " << mnair << endl;
500  }
501 
502 
503  // Pressure shift.
504  {
505  // HITRAN value in cm^-1 / atm. So the conversion goes exactly as
506  // for the broadening parameters.
507  Numeric d;
508  // External constant from constants.cc: Converts atm to
509  // Pa. Multiply value in atm by this number to get value in Pa.
510  extern const Numeric ATM2PA;
511  // External constant from constants.cc:
512  extern const Numeric SPEED_OF_LIGHT;
513  // Conversion from wavenumber to Hz. If you multiply a value in
514  // wavenumber (cm^-1) by this constant, you get the value in Hz.
515  const Numeric w2Hz = SPEED_OF_LIGHT * 1e2;
516  // Ok, put together the end-to-end conversion that we need:
517  const Numeric hi2arts = w2Hz / ATM2PA;
518 
519  // Extract HITRAN value:
520  extract(d,line,8);
521 
522  // ARTS value in Hz/Pa
523  mpsf = d * hi2arts;
524  }
525  // Set the accuracies using the definition of HITRAN
526  // indices. If some are missing, they are set to -1.
527 
528  //Skip upper state global quanta index
529  {
530  Index eu;
531  extract(eu,line,3);
532  }
533 
534  //Skip lower state global quanta index
535  {
536  Index el;
537  extract(el,line,3);
538  }
539 
540  //Skip upper state local quanta
541  {
542  Index eul;
543  extract(eul,line,9);
544  }
545 
546  //Skip lower state local quanta
547  {
548  Index ell;
549  extract(ell,line,9);
550  }
551 
552  // Accuracy index for frequency reference
553  {
554  Index df;
555  // Extract HITRAN value:
556  extract(df,line,1);
557  // Convert it to ARTS units (Hz)
558  convHitranIERF(mdf,df);
559  }
560 
561  // Accuracy index for intensity reference
562  {
563  Index di0;
564  // Extract HITRAN value:
565  extract(di0,line,1);
566  convHitranIERSH(mdi0,di0);
567  }
568 
569  // Accuracy index for halfwidth reference
570  {
571  Index dgam;
572  // Extract HITRAN value:
573  extract(dgam,line,1);
574  //Convert to ARTS units (%)
575  convHitranIERSH(mdagam,dgam);
576  // convHitranIERSH(mdsgam,dgam);
577  // convHitranIERSH(mdnair,dgam);
578  // convHitranIERSH(mdnself,dgam);
579  }
580 
581  // Accuracy for pressure shift
582  // This is missing in HITRAN catalogue and it is set to -1.
583  mdpsf =-1;
584 
585  // These were all the parameters that we can extract from
586  // HITRAN. However, we still have to set the reference temperatures
587  // to the appropriate value:
588 
589  // Reference temperature for Intensity in K.
590  // (This is fix for HITRAN)
591  mti0 = 296.0;
592 
593  // Reference temperature for AGAM and SGAM in K.
594  // (This is also fix for HITRAN)
595  mtgam = 296.0;
596 
597  // That's it!
598  return false;
599 }
600 
601 
602 // This function reads line data in the Hitran 2004 format. For the Hitran
603 // 1986-2004 data format use ReadFromHitranStream().
604 //
605 // 2005/03/29: This function was added based on ReadFromHitranStream(). There
606 // is a check for the right data record length (160 characters).
607 // Hermann Berg (h.berg@utoronto.ca)
608 //
610 {
611  // Global species lookup data:
613 
614  // This value is used to flag missing data both in species and
615  // isotope lists. Could be any number, it just has to be made sure
616  // that it is neither the index of a species nor of an isotope.
617  const Index missing = species_data.nelem() + 100;
618 
619  // We need a species index sorted by HITRAN tag. Keep this in a
620  // static variable, so that we have to do this only once. The ARTS
621  // species index is hind[<HITRAN tag>].
622  //
623  // Allow for up to 100 species in HITRAN in the future.
624  static Array< Index > hspec(100);
625 
626  // This is an array of arrays for each hitran tag. It contains the
627  // ARTS indices of the HITRAN isotopes.
628  static Array< ArrayOfIndex > hiso(100);
629 
630  // Remeber if this stuff has already been initialized:
631  static bool hinit = false;
632 
633  // Remember, about which missing species we have already issued a
634  // warning:
635  static ArrayOfIndex warned_missing;
636 
637  if ( !hinit )
638  {
639  // Initialize hspec.
640  // The value of missing means that we don't have this species.
641  hspec = missing; // Matpack can set all elements like this.
642 
643  for ( Index i=0; i<species_data.nelem(); ++i )
644  {
645  const SpeciesRecord& sr = species_data[i];
646 
647  // We have to be careful and check for the case that all
648  // HITRAN isotope tags are -1 (this species is missing in HITRAN).
649 
650  if ( 0 < sr.Isotope()[0].HitranTag() )
651  {
652  // The HITRAN tags are stored as species plus isotope tags
653  // (MO and ISO)
654  // in the Isotope() part of the species record.
655  // We can extract the MO part from any of the isotope tags,
656  // so we use the first one. We do this by taking an integer
657  // division by 10.
658 
659  Index mo = sr.Isotope()[0].HitranTag() / 10;
660  // cout << "mo = " << mo << endl;
661  hspec[mo] = i;
662 
663  // Get a nicer to handle array of HITRAN iso tags:
664  Index n_iso = sr.Isotope().nelem();
665  ArrayOfIndex iso_tags;
666  iso_tags.resize(n_iso);
667  for ( Index j=0; j<n_iso; ++j )
668  {
669  iso_tags[j] = sr.Isotope()[j].HitranTag();
670  }
671 
672  // Reserve elements for the isotope tags. How much do we
673  // need? This depends on the largest HITRAN tag that we know
674  // about!
675  // Also initialize the tags to missing.
676  // cout << "iso_tags = " << iso_tags << endl;
677  // cout << "static_cast<Index>(max(iso_tags))%10 + 1 = "
678  // << static_cast<Index>(max(iso_tags))%10 + 1 << endl;
679  hiso[mo].resize( max(iso_tags)%10 + 1 );
680  hiso[mo] = missing; // Matpack can set all elements like this.
681 
682 
683  // Set the isotope tags:
684  for ( Index j=0; j<n_iso; ++j )
685  {
686  if ( 0 < iso_tags[j] ) // ignore -1 elements
687  {
688  // To get the iso tags from HitranTag() we also have to take
689  // modulo 10 to get rid of mo.
690  hiso[mo][iso_tags[j] % 10] = j;
691  }
692  }
693  }
694  }
695 
696 
697  // Print the generated data structures (for debugging):
698  out3 << " HITRAN index table:\n";
699  for ( Index i=0; i<hspec.nelem(); ++i )
700  {
701  if ( missing != hspec[i] )
702  {
703  // The explicit conversion of Name to a c-String is
704  // necessary, because setw does not work correctly for
705  // stl Strings.
706  out3 << " mo = " << i << " Species = "
707  << setw(10) << setiosflags(ios::left)
708  << species_data[hspec[i]].Name().c_str()
709  << "iso = ";
710  for ( Index j=1; j<hiso[i].nelem(); ++j )
711  {
712  if ( missing==hiso[i][j] )
713  out3 << " " << "m";
714  else
715  out3 << " " << species_data[hspec[i]].Isotope()[hiso[i][j]].Name();
716  }
717  out3 << "\n";
718  }
719  }
720 
721  hinit = true;
722  }
723 
724 
725  // This contains the rest of the line to parse. At the beginning the
726  // entire line. Line gets shorter and shorter as we continue to
727  // extract stuff from the beginning.
728  String line;
729 
730  // The first item is the molecule number:
731  Index mo;
732 
733  // Look for more comments?
734  bool comment = true;
735 
736  while (comment)
737  {
738  // Return true if eof is reached:
739  if (is.eof()) return true;
740 
741  // Throw runtime_error if stream is bad:
742  if (!is) throw runtime_error ("Stream bad.");
743 
744  // Read line from file into linebuffer:
745  getline(is,line);
746 
747  // It is possible that we were exactly at the end of the file before
748  // calling getline. In that case the previous eof() was still false
749  // because eof() evaluates only to true if one tries to read after the
750  // end of the file. The following check catches this.
751  if (line.nelem() == 0 && is.eof()) return true;
752 
753  // If the catalogue is in dos encoding, throw away the
754  // additional carriage return
755  if (line[line.nelem () - 1] == 13)
756  {
757  line.erase (line.nelem () - 1, 1);
758  }
759 
760  // Because of the fixed FORTRAN format, we need to break up the line
761  // explicitly in apropriate pieces. Not elegant, but works!
762 
763  // Extract molecule number:
764  mo = 0;
765  // Initialization of mo is important, because mo stays the same
766  // if line is empty.
767  extract(mo,line,2);
768  // cout << "mo = " << mo << endl;
769 
770  // If mo == 0 this is just a comment line:
771  if ( 0 != mo )
772  {
773  // See if we know this species. Exit with an error if the species is unknown.
774  if ( missing != hspec[mo] )
775  {
776  comment = false;
777 
778  // Check if data record has the right number of characters for the
779  // in Hitran 2004 format
780  Index nChar = line.nelem() + 2; // number of characters in data record;
781  if ( nChar != 160 )
782  {
783  ostringstream os;
784  os << "Invalid HITRAN 2004 line data record with " << nChar <<
785  " characters (expected: 160).";
786  throw runtime_error(os.str());
787  }
788 
789  }
790  else
791  {
792  // See if this is already in warned_missing, use
793  // std::count for that:
794  if ( 0 == std::count(warned_missing.begin(),
795  warned_missing.end(),
796  mo) )
797  {
798  out0 << "Warning: HITRAN molecule number mo = " << mo << " is not "
799  << "known to ARTS.\n";
800  warned_missing.push_back(mo);
801  }
802  }
803  }
804  }
805 
806  // Ok, we seem to have a valid species here.
807 
808  // Set mspecies from my cool index table:
809  mspecies = hspec[mo];
810 
811  // Extract isotope:
812  Index iso;
813  extract(iso,line,1);
814  // cout << "iso = " << iso << endl;
815 
816 
817  // Set misotope from the other cool index table.
818  // We have to be careful to issue an error for unknown iso tags. Iso
819  // could be either larger than the size of hiso[mo], or set
820  // explicitly to missing. Unfortunately we have to test both cases.
821  misotope = missing;
822  if ( iso < hiso[mo].nelem() )
823  if ( missing != hiso[mo][iso] )
824  misotope = hiso[mo][iso];
825 
826  // Issue error message if misotope is still missing:
827  if (missing == misotope)
828  {
829  ostringstream os;
830  os << "Species: " << species_data[mspecies].Name()
831  << ", isotope iso = " << iso
832  << " is unknown.";
833  throw runtime_error(os.str());
834  }
835 
836 
837  // Position.
838  {
839  // HITRAN position in wavenumbers (cm^-1):
840  Numeric v;
841  // External constant from constants.cc:
842  extern const Numeric SPEED_OF_LIGHT;
843  // Conversion from wavenumber to Hz. If you multiply a line
844  // position in wavenumber (cm^-1) by this constant, you get the
845  // frequency in Hz.
846  const Numeric w2Hz = SPEED_OF_LIGHT * 100.;
847 
848  // Extract HITRAN postion:
849  extract(v,line,12);
850 
851  // ARTS position in Hz:
852  mf = v * w2Hz;
853 // cout << "mf = " << mf << endl;
854  }
855 
856  // Intensity.
857  {
858  extern const Numeric SPEED_OF_LIGHT; // in [m/s]
859 
860  // HITRAN intensity is in cm-1/(molec * cm-2) at 296 Kelvin.
861  // It already includes the isotpic ratio.
862  // The first cm-1 is the frequency unit (it cancels with the
863  // 1/frequency unit of the line shape function).
864  //
865  // We need to do the following:
866  // 1. Convert frequency from wavenumber to Hz (factor 1e2 * c).
867  // 2. Convert [molec * cm-2] to [molec * m-2] (factor 1e-4).
868  // 3. Take out the isotopic ratio.
869 
870  const Numeric hi2arts = 1e-2 * SPEED_OF_LIGHT;
871 
872  Numeric s;
873 
874  // Extract HITRAN intensity:
875  extract(s,line,10);
876  // Convert to ARTS units (Hz / (molec * m-2) ), or shorter: Hz*m^2
877  mi0 = s * hi2arts;
878  // Take out isotopic ratio:
879  mi0 /= species_data[mspecies].Isotope()[misotope].Abundance();
880  }
881 
882  // Skip Einstein coefficient
883  {
884  Numeric r;
885  extract(r,line,10);
886  }
887 
888 
889  // Air broadening parameters.
890  {
891  // HITRAN parameter is in cm-1/atm at 296 Kelvin
892  // All parameters are HWHM (I hope this is true!)
893  Numeric gam;
894  // External constant from constants.cc: Converts atm to
895  // Pa. Multiply value in atm by this number to get value in Pa.
896  extern const Numeric ATM2PA;
897  // External constant from constants.cc:
898  extern const Numeric SPEED_OF_LIGHT;
899  // Conversion from wavenumber to Hz. If you multiply a value in
900  // wavenumber (cm^-1) by this constant, you get the value in Hz.
901  const Numeric w2Hz = SPEED_OF_LIGHT * 1e2;
902  // Ok, put together the end-to-end conversion that we need:
903  const Numeric hi2arts = w2Hz / ATM2PA;
904 
905  // Extract HITRAN AGAM value:
906  extract(gam,line,5);
907 
908  // ARTS parameter in Hz/Pa:
909  magam = gam * hi2arts;
910 
911  // Extract HITRAN SGAM value:
912  extract(gam,line,5);
913 
914  // ARTS parameter in Hz/Pa:
915  msgam = gam * hi2arts;
916 
917  // If zero, set to agam:
918  if (0==msgam)
919  msgam = magam;
920 
921  // cout << "agam, sgam = " << magam << ", " << msgam << endl;
922  }
923 
924 
925  // Lower state energy.
926  {
927  // HITRAN parameter is in wavenumbers (cm^-1).
928  // We have to convert this to the ARTS unit Joule.
929 
930  // Extract from Catalogue line
931  extract(melow,line,10);
932 
933  // Convert to Joule:
935  }
936 
937 
938  // Temperature coefficient of broadening parameters.
939  {
940  // This is dimensionless, we can also extract directly.
941  extract(mnair,line,4);
942 
943  // Set self broadening temperature coefficient to the same value:
944  mnself = mnair;
945 // cout << "mnair = " << mnair << endl;
946  }
947 
948 
949  // Pressure shift.
950  {
951  // HITRAN value in cm^-1 / atm. So the conversion goes exactly as
952  // for the broadening parameters.
953  Numeric d;
954  // External constant from constants.cc: Converts atm to
955  // Pa. Multiply value in atm by this number to get value in Pa.
956  extern const Numeric ATM2PA;
957  // External constant from constants.cc:
958  extern const Numeric SPEED_OF_LIGHT;
959  // Conversion from wavenumber to Hz. If you multiply a value in
960  // wavenumber (cm^-1) by this constant, you get the value in Hz.
961  const Numeric w2Hz = SPEED_OF_LIGHT * 1e2;
962  // Ok, put together the end-to-end conversion that we need:
963  const Numeric hi2arts = w2Hz / ATM2PA;
964 
965  // Extract HITRAN value:
966  extract(d,line,8);
967 
968  // ARTS value in Hz/Pa
969  mpsf = d * hi2arts;
970  }
971  // Set the accuracies using the definition of HITRAN
972  // indices. If some are missing, they are set to -1.
973 
974  //Skip upper state global quanta index
975  {
976  Index eu;
977  extract(eu,line,15);
978  }
979 
980  //Skip lower state global quanta index
981  {
982  Index el;
983  extract(el,line,15);
984  }
985 
986  //Skip upper state local quanta
987  {
988  Index eul;
989  extract(eul,line,15);
990  }
991 
992  //Skip lower state local quanta
993  {
994  Index ell;
995  extract(ell,line,15);
996  }
997 
998  // Accuracy index for frequency
999  {
1000  Index df;
1001  // Extract HITRAN value:
1002  extract(df,line,1);
1003  // Convert it to ARTS units (Hz)
1004  convHitranIERF(mdf,df);
1005  }
1006 
1007  // Accuracy index for intensity
1008  {
1009  Index di0;
1010  // Extract HITRAN value:
1011  extract(di0,line,1);
1012  //Convert to ARTS units (%)
1013  convHitranIERSH(mdi0,di0);
1014  }
1015 
1016  // Accuracy index for air-broadened halfwidth
1017  {
1018  Index dagam;
1019  // Extract HITRAN value:
1020  extract(dagam,line,1);
1021  //Convert to ARTS units (%)
1022  convHitranIERSH(mdagam,dagam);
1023  }
1024 
1025  // Accuracy index for self-broadened half-width
1026  {
1027  Index dsgam;
1028  // Extract HITRAN value:
1029  extract(dsgam,line,1);
1030  //Convert to ARTS units (%)
1031  convHitranIERSH(mdsgam,dsgam);
1032  }
1033 
1034  // Accuracy index for temperature-dependence exponent for agam
1035  {
1036  Index dnair;
1037  // Extract HITRAN value:
1038  extract(dnair,line,1);
1039  //Convert to ARTS units (%)
1040  convHitranIERSH(mdnair,dnair);
1041  }
1042 
1043  // Accuracy index for temperature-dependence exponent for sgam
1044  // This is missing in HITRAN catalogue and is set to -1.
1045  mdnself =-1;
1046 
1047  // Accuracy index for pressure shift
1048  {
1049  Index dpsf;
1050  // Extract HITRAN value (given in cm-1):
1051  extract(dpsf,line,1);
1052  // Convert it to ARTS units (Hz)
1053  convHitranIERF(mdpsf,dpsf);
1054  // ARTS wants this error in %
1055  mdpsf = mdpsf / mf;
1056  }
1057 
1058  // These were all the parameters that we can extract from
1059  // HITRAN 2004. However, we still have to set the reference temperatures
1060  // to the appropriate value:
1061 
1062  // Reference temperature for Intensity in K.
1063  // (This is fix for HITRAN 2004)
1064  mti0 = 296.0;
1065 
1066  // Reference temperature for AGAM and SGAM in K.
1067  // (This is also fix for HITRAN 2004)
1068  mtgam = 296.0;
1069 
1070  // That's it!
1071  return false;
1072 }
1073 
1074 
1076 {
1077  // Global species lookup data:
1079 
1080  // This value is used to flag missing data both in species and
1081  // isotope lists. Could be any number, it just has to be made sure
1082  // that it is neither the index of a species nor of an isotope.
1083  const Index missing = species_data.nelem() + 100;
1084 
1085  // We need a species index sorted by MYTRAN tag. Keep this in a
1086  // static variable, so that we have to do this only once. The ARTS
1087  // species index is hind[<MYTRAN tag>]. The value of
1088  // missing means that we don't have this species.
1089  //
1090  // Allow for up to 100 species in MYTRAN in the future.
1091  static Array< Index > hspec(100,missing);
1092 
1093  // This is an array of arrays for each mytran tag. It contains the
1094  // ARTS indices of the MYTRAN isotopes.
1095  static Array< ArrayOfIndex > hiso(100);
1096 
1097  // Remeber if this stuff has already been initialized:
1098  static bool hinit = false;
1099 
1100  // Remember, about which missing species we have already issued a
1101  // warning:
1102  static ArrayOfIndex warned_missing;
1103 
1104  if ( !hinit )
1105  {
1106  for ( Index i=0; i<species_data.nelem(); ++i )
1107  {
1108  const SpeciesRecord& sr = species_data[i];
1109 
1110  // We have to be careful and check for the case that all
1111  // MYTRAN isotope tags are -1 (this species is missing in MYTRAN).
1112 
1113  if ( 0 < sr.Isotope()[0].MytranTag() )
1114  {
1115  // The MYTRAN tags are stored as species plus isotope tags
1116  // (MO and ISO)
1117  // in the Isotope() part of the species record.
1118  // We can extract the MO part from any of the isotope tags,
1119  // so we use the first one. We do this by taking an integer
1120  // division by 10.
1121 
1122  Index mo = sr.Isotope()[0].MytranTag() / 10;
1123  // cout << "mo = " << mo << endl;
1124  hspec[mo] = i;
1125 
1126  // Get a nicer to handle array of MYTRAN iso tags:
1127  Index n_iso = sr.Isotope().nelem();
1128  ArrayOfIndex iso_tags;
1129  iso_tags.resize(n_iso);
1130  for ( Index j=0; j<n_iso; ++j )
1131  {
1132  iso_tags[j] = sr.Isotope()[j].MytranTag();
1133  }
1134 
1135  // Reserve elements for the isotope tags. How much do we
1136  // need? This depends on the largest MYTRAN tag that we know
1137  // about!
1138  // Also initialize the tags to missing.
1139  // cout << "iso_tags = " << iso_tags << endl;
1140  // cout << "static_cast<Index>(max(iso_tags))%10 + 1 = "
1141  // << static_cast<Index>(max(iso_tags))%10 + 1 << endl;
1142  hiso[mo].resize( max(iso_tags)%10 + 1 );
1143  hiso[mo] = missing; // Matpack can set all elements like this.
1144 
1145  // Set the isotope tags:
1146  for ( Index j=0; j<n_iso; ++j )
1147  {
1148  if ( 0 < iso_tags[j] ) // ignore -1 elements
1149  {
1150  // To get the iso tags from MytranTag() we also have to take
1151  // modulo 10 to get rid of mo.
1152  // cout << "iso_tags[j] % 10 = " << iso_tags[j] % 10 << endl;
1153  hiso[mo][iso_tags[j] % 10] = j;
1154  }
1155  }
1156  }
1157  }
1158 
1159 // cout << "hiso = " << hiso << endl << "***********" << endl;
1160 
1161 
1162  // Print the generated data structures (for debugging):
1163  out3 << " MYTRAN index table:\n";
1164  for ( Index i=0; i<hspec.nelem(); ++i )
1165  {
1166  if ( missing != hspec[i] )
1167  {
1168  // The explicit conversion of Name to a c-String is
1169  // necessary, because setw does not work correctly for
1170  // stl Strings.
1171  out3 << " mo = " << i << " Species = "
1172  << setw(10) << setiosflags(ios::left)
1173  << species_data[hspec[i]].Name().c_str()
1174  << "iso = ";
1175  for ( Index j=1; j<hiso[i].nelem(); ++j )
1176  {
1177  if ( missing==hiso[i][j] )
1178  out3 << " " << "m";
1179  else
1180  out3 << " " << species_data[hspec[i]].Isotope()[hiso[i][j]].Name();
1181  }
1182  out3 << "\n";
1183  }
1184  }
1185 
1186  hinit = true;
1187  }
1188 
1189 
1190  // This contains the rest of the line to parse. At the beginning the
1191  // entire line. Line gets shorter and shorter as we continue to
1192  // extract stuff from the beginning.
1193  String line;
1194 
1195  // The first item is the molecule number:
1196  Index mo;
1197 
1198  // Look for more comments?
1199  bool comment = true;
1200 
1201  while (comment)
1202  {
1203  // Return true if eof is reached:
1204  if (is.eof()) return true;
1205 
1206  // Throw runtime_error if stream is bad:
1207  if (!is) throw runtime_error ("Stream bad.");
1208 
1209  // Read line from file into linebuffer:
1210  getline(is,line);
1211 
1212  // It is possible that we were exactly at the end of the file before
1213  // calling getline. In that case the previous eof() was still false
1214  // because eof() evaluates only to true if one tries to read after the
1215  // end of the file. The following check catches this.
1216  if (line.nelem() == 0 && is.eof()) return true;
1217 
1218  // Because of the fixed FORTRAN format, we need to break up the line
1219  // explicitly in apropriate pieces. Not elegant, but works!
1220 
1221  // Extract molecule number:
1222  mo = 0;
1223  // Initialization of mo is important, because mo stays the same
1224  // if line is empty.
1225  extract(mo,line,2);
1226  // cout << "mo = " << mo << endl;
1227 
1228  // If mo == 0 this is just a comment line:
1229  if ( 0 != mo )
1230  {
1231  // See if we know this species. We will give an error if a
1232  // species is not known.
1233  if ( missing != hspec[mo] ) comment = false ;
1234  else
1235  {
1236  // See if this is already in warned_missing, use
1237  // std::count for that:
1238  if ( 0 == std::count(warned_missing.begin(),
1239  warned_missing.end(),
1240  mo) )
1241  {
1242  out0 << "Error: MYTRAN mo = " << mo << " is not "
1243  << "known to ARTS.\n";
1244  warned_missing.push_back(mo);
1245  exit(1);
1246  // SAB 08.08.2000 If you want to make the program
1247  // continue anyway, just comment out the exit
1248  // line.
1249  }
1250  }
1251  }
1252  }
1253 
1254  // Ok, we seem to have a valid species here.
1255 
1256  // Set mspecies from my cool index table:
1257  mspecies = hspec[mo];
1258 
1259  // Extract isotope:
1260  Index iso;
1261  extract(iso,line,1);
1262  // cout << "iso = " << iso << endl;
1263 
1264 
1265  // Set misotope from the other cool index table.
1266  // We have to be careful to issue an error for unknown iso tags. Iso
1267  // could be either larger than the size of hiso[mo], or set
1268  // explicitly to missing. Unfortunately we have to test both cases.
1269  misotope = missing;
1270  if ( iso < hiso[mo].nelem() )
1271  if ( missing != hiso[mo][iso] )
1272  misotope = hiso[mo][iso];
1273 
1274  // Issue error message if misotope is still missing:
1275  if (missing == misotope)
1276  {
1277  ostringstream os;
1278  os << "Species: " << species_data[mspecies].Name()
1279  << ", isotope iso = " << iso
1280  << " is unknown.";
1281  throw runtime_error(os.str());
1282  }
1283 
1284 
1285  // Position.
1286  {
1287  // MYTRAN position in MHz:
1288  Numeric v;
1289 
1290  // Extract MYTRAN postion:
1291  extract(v,line,13);
1292 
1293  // ARTS position in Hz:
1294  mf = v * 1E6;
1295  // cout << "mf = " << mf << endl;
1296  }
1297 
1298  // Accuracy for line position
1299  {
1300  // Extract MYTRAN postion accuracy:
1301  Numeric df;
1302  extract(df,line,8);
1303  // ARTS accuracy of line position in Hz:
1304  mdf = df * 1E6;
1305  }
1306 
1307  // Intensity.
1308  {
1309  extern const Numeric SPEED_OF_LIGHT; // in [m/s]
1310 
1311  // MYTRAN2 intensity is in cm-1/(molec * cm-2) at 296 Kelvin.
1312  // (just like HITRAN, only isotopic ratio is not included)
1313  // The first cm-1 is the frequency unit (it cancels with the
1314  // 1/frequency unit of the line shape function).
1315  //
1316  // We need to do the following:
1317  // 1. Convert frequency from wavenumber to Hz (factor 1e2 * c)
1318  // 2. Convert [molec * cm-2] to [molec * m-2] (factor 1e-4)
1319 
1320  const Numeric hi2arts = 1e-2 * SPEED_OF_LIGHT;
1321 
1322  Numeric s;
1323 
1324  // Extract HITRAN intensity:
1325  extract(s,line,10);
1326 
1327  // Convert to ARTS units (Hz / (molec * m-2) ), or shorter: Hz*m^2
1328  mi0 = s * hi2arts;
1329  }
1330 
1331 
1332  // Air broadening parameters.
1333  {
1334  // MYTRAN parameter is in MHz/Torr at reference temperature
1335  // All parameters are HWHM
1336  Numeric gam;
1337  // External constant from constants.cc: Converts torr to
1338  // Pa. Multiply value in torr by this number to get value in Pa.
1339  extern const Numeric TORR2PA;
1340 
1341  // Extract HITRAN AGAM value:
1342  extract(gam,line,5);
1343 
1344  // ARTS parameter in Hz/Pa:
1345  magam = gam * 1E6 / TORR2PA;
1346 
1347  // Extract MYTRAN SGAM value:
1348  extract(gam,line,5);
1349 
1350  // ARTS parameter in Hz/Pa:
1351  msgam = gam * 1E6 / TORR2PA;
1352 
1353 // cout << "agam, sgam = " << magam << ", " << msgam << endl;
1354  }
1355 
1356 
1357  // Lower state energy.
1358  {
1359  // MYTRAN parameter is in wavenumbers (cm^-1).
1360  // We have to convert this to the ARTS unit Joule.
1361 
1362  // Extract from Catalogue line
1363  extract(melow,line,10);
1364 
1365  // Convert to Joule:
1367  }
1368 
1369 
1370  // Temperature coefficient of broadening parameters.
1371  {
1372  // This is dimensionless, we can also extract directly.
1373  extract(mnair,line,4);
1374 
1375  // Extract the self broadening parameter:
1376  extract(mnself,line,4);
1377 // cout << "mnair = " << mnair << endl;
1378  }
1379 
1380 
1381  // Reference temperature for broadening parameter in K:
1382  {
1383  // correct units, extract directly
1384  extract(mtgam,line,7);
1385  }
1386 
1387 
1388  // Pressure shift.
1389  {
1390  // MYTRAN value in MHz/Torr
1391  Numeric d;
1392  // External constant from constants.cc: Converts torr to
1393  // Pa. Multiply value in torr by this number to get value in Pa.
1394  extern const Numeric TORR2PA;
1395 
1396  // Extract MYTRAN value:
1397  extract(d,line,9);
1398 
1399  // ARTS value in Hz/Pa
1400  mpsf = d * 1E6 / TORR2PA;
1401  }
1402  // Set the accuracies using the definition of MYTRAN accuracy
1403  // indices. If some are missing, they are set to -1.
1404 
1405  //Skip upper state global quanta index
1406  {
1407  Index eu;
1408  extract(eu,line,3);
1409  }
1410 
1411  //Skip lower state global quanta index
1412  {
1413  Index el;
1414  extract(el,line,3);
1415  }
1416 
1417  //Skip upper state local quanta
1418  {
1419  Index eul;
1420  extract(eul,line,9);
1421  }
1422 
1423  //Skip lower state local quanta
1424  {
1425  Index ell;
1426  extract(ell,line,9);
1427  }
1428  // Accuracy index for intensity
1429  {
1430  Index di0;
1431  // Extract MYTRAN value:
1432  extract(di0,line,1);
1433  //convert to ARTS units (%)
1434  convMytranIER(mdi0,di0);
1435  }
1436 
1437  // Accuracy index for AGAM
1438  {
1439  Index dgam;
1440  // Extract MYTRAN value:
1441  extract(dgam,line,1);
1442  //convert to ARTS units (%)
1443  convMytranIER(mdagam,dgam);
1444  }
1445 
1446  // Accuracy index for NAIR
1447  {
1448  Index dnair;
1449  // Extract MYTRAN value:
1450  extract(dnair,line,1);
1451  //convert to ARTS units (%);
1452  convMytranIER(mdnair,dnair);
1453  }
1454 
1455 
1456  // These were all the parameters that we can extract from
1457  // MYTRAN. However, we still have to set the reference temperatures
1458  // to the appropriate value:
1459 
1460  // Reference temperature for Intensity in K.
1461  // (This is fix for MYTRAN2)
1462  mti0 = 296.0;
1463 
1464  // It is important that you intialize here all the new parameters that
1465  // you added to the line record. (This applies to all the reading
1466  // functions, also for ARTS, JPL, and HITRAN format.) Parameters
1467  // should be either set from the catalogue, or set to -1.)
1468 
1469  // That's it!
1470  return false;
1471 }
1472 
1473 
1475 {
1476  // Global species lookup data:
1478 
1479  // We need a species index sorted by JPL tag. Keep this in a
1480  // static variable, so that we have to do this only once. The ARTS
1481  // species index is JplMap[<JPL tag>]. We need Index in this map,
1482  // because the tag array within the species data is an Index array.
1483  static map<Index, SpecIsoMap> JplMap;
1484 
1485  // Remeber if this stuff has already been initialized:
1486  static bool hinit = false;
1487 
1488  if ( !hinit )
1489  {
1490 
1491  out3 << " JPL index table:\n";
1492 
1493  for ( Index i=0; i<species_data.nelem(); ++i )
1494  {
1495  const SpeciesRecord& sr = species_data[i];
1496 
1497 
1498  for ( Index j=0; j<sr.Isotope().nelem(); ++j)
1499  {
1500 
1501  for (Index k=0; k<sr.Isotope()[j].JplTags().nelem(); ++k)
1502  {
1503 
1504  SpecIsoMap indicies(i,j);
1505 
1506  JplMap[sr.Isotope()[j].JplTags()[k]] = indicies;
1507 
1508  // Print the generated data structures (for debugging):
1509  // The explicit conversion of Name to a c-String is
1510  // necessary, because setw does not work correctly for
1511  // stl Strings.
1512  const Index& i1 = JplMap[sr.Isotope()[j].JplTags()[k]].Speciesindex();
1513  const Index& i2 = JplMap[sr.Isotope()[j].JplTags()[k]].Isotopeindex();
1514 
1515  out3 << " JPL TAG = " << sr.Isotope()[j].JplTags()[k] << " Species = "
1516  << setw(10) << setiosflags(ios::left)
1517  << species_data[i1].Name().c_str()
1518  << "iso = "
1519  << species_data[i1].Isotope()[i2].Name().c_str()
1520  << "\n";
1521  }
1522 
1523  }
1524  }
1525  hinit = true;
1526  }
1527 
1528 
1529  // This contains the rest of the line to parse. At the beginning the
1530  // entire line. Line gets shorter and shorter as we continue to
1531  // extract stuff from the beginning.
1532  String line;
1533 
1534 
1535  // Look for more comments?
1536  bool comment = true;
1537 
1538  while (comment)
1539  {
1540  // Return true if eof is reached:
1541  if (is.eof()) return true;
1542 
1543  // Throw runtime_error if stream is bad:
1544  if (!is) throw runtime_error ("Stream bad.");
1545 
1546  // Read line from file into linebuffer:
1547  getline(is,line);
1548 
1549  // It is possible that we were exactly at the end of the file before
1550  // calling getline. In that case the previous eof() was still false
1551  // because eof() evaluates only to true if one tries to read after the
1552  // end of the file. The following check catches this.
1553  if (line.nelem() == 0 && is.eof()) return true;
1554 
1555  // Because of the fixed FORTRAN format, we need to break up the line
1556  // explicitly in apropriate pieces. Not elegant, but works!
1557 
1558  // Extract center frequency:
1559  // Initialization of v is important, because v stays the same
1560  // if line is empty.
1561  // JPL position in MHz:
1562  Numeric v = 0.0;
1563 
1564  // Extract JPL position:
1565  extract(v,line,13);
1566 
1567  // check for empty line
1568  if (v != 0.0)
1569  {
1570  // ARTS position in Hz:
1571  mf = v * 1E6;
1572 
1573  comment = false;
1574  }
1575  }
1576 
1577  // Accuracy for line position
1578  {
1579  Numeric df;
1580  extract(df,line,8);
1581  //convert to ARTS units (Hz)
1582  mdf = df * 1E6;
1583  }
1584 
1585  // Intensity.
1586  {
1587  // JPL has log (10) of intensity in nm2 MHz at 300 Kelvin.
1588  //
1589  // We need to do the following:
1590  // 1. take 10^intensity
1591  // 2. convert to cm-1/(molecule * cm-2): devide by c * 1e10
1592  // 3. Convert frequency from wavenumber to Hz (factor 1e2 * c)
1593  // 4. Convert [molec * cm-2] to [molec * m-2] (factor 1e-4)
1594 
1595  Numeric s;
1596 
1597  // Extract JPL intensity:
1598  extract(s,line,8);
1599 
1600  // remove log
1601  s = pow(10,s);
1602 
1603  // Convert to ARTS units (Hz / (molec * m-2) ), or shorter: Hz*m^2
1604  mi0 = s / 1E12;
1605  }
1606 
1607  // Degrees of freedom
1608  {
1609  Index dr;
1610 
1611  // Extract degrees of freedom
1612  extract(dr,line,2);
1613  }
1614 
1615  // Lower state energy.
1616  {
1617  // JPL parameter is in wavenumbers (cm^-1).
1618  // We have to convert this to the ARTS unit Joule.
1619 
1620  // Extract from Catalogue line
1621  extract(melow,line,10);
1622 
1623  // Convert to Joule:
1625  }
1626 
1627  // Upper state degeneracy
1628  {
1629  Index gup;
1630 
1631  // Extract upper state degeneracy
1632  extract(gup,line,3);
1633  }
1634 
1635  // Tag number
1636  Index tag;
1637  {
1638  // Extract Tag number
1639  extract(tag,line,7);
1640 
1641  // make sure tag is not negative (damned jpl cat):
1642  tag = tag > 0 ? tag : -tag;
1643  }
1644 
1645  // ok, now for the cool index map:
1646 
1647  // is this tag valid?
1648  const map<Index, SpecIsoMap>::const_iterator i = JplMap.find(tag);
1649  if ( i == JplMap.end() )
1650  {
1651  ostringstream os;
1652  os << "JPL Tag: " << tag << " is unknown.";
1653  throw runtime_error(os.str());
1654  }
1655 
1656  SpecIsoMap id = i->second;
1657 
1658 
1659  // Set mspecies:
1660  mspecies = id.Speciesindex();
1661 
1662  // Set misotope:
1663  misotope = id.Isotopeindex();
1664 
1665  // Air broadening parameters: unknown to jpl, use old iup forward
1666  // model default values, which is mostly set to 0.0025 GHz/hPa, even
1667  // though for some lines the pressure broadening is given explicitly
1668  // in the program code. The explicitly given values are ignored and
1669  // only the default value is set. Self broadening was in general not
1670  // considered in the old forward model.
1671  {
1672  // ARTS parameter in Hz/Pa:
1673  magam = 2.5E4;
1674 
1675  // ARTS parameter in Hz/Pa:
1676  msgam = magam;
1677  }
1678 
1679 
1680  // Temperature coefficient of broadening parameters. Was set to 0.75
1681  // in old forward model, even though for some lines the parameter is
1682  // given explicitly in the program code. The explicitly given values
1683  // are ignored and only the default value is set. Self broadening
1684  // not considered.
1685  {
1686  mnair = 0.75;
1687  mnself = 0.0;
1688  }
1689 
1690 
1691  // Reference temperature for broadening parameter in K, was
1692  // generally set to 300 K in old forward model, with the exceptions
1693  // as already mentioned above:
1694  {
1695  mtgam = 300.0;
1696  }
1697 
1698 
1699  // Pressure shift: not given in JPL, set to 0
1700  {
1701  mpsf = 0.0;
1702  }
1703 
1704 
1705  // These were all the parameters that we can extract from
1706  // JPL. However, we still have to set the reference temperatures
1707  // to the appropriate value:
1708 
1709  // Reference temperature for Intensity in K.
1710  // (This is fix for JPL)
1711  mti0 = 300.0;
1712 
1713  // That's it!
1714  return false;
1715 }
1716 
1717 
1719 {
1720  // Global species lookup data:
1722 
1723  // We need a species index sorted by Arts identifier. Keep this in a
1724  // static variable, so that we have to do this only once. The ARTS
1725  // species index is ArtsMap[<Arts String>].
1726  static map<String, SpecIsoMap> ArtsMap;
1727 
1728  // Remeber if this stuff has already been initialized:
1729  static bool hinit = false;
1730 
1731  if ( !hinit )
1732  {
1733 
1734  out3 << " ARTS index table:\n";
1735 
1736  for ( Index i=0; i<species_data.nelem(); ++i )
1737  {
1738  const SpeciesRecord& sr = species_data[i];
1739 
1740 
1741  for ( Index j=0; j<sr.Isotope().nelem(); ++j)
1742  {
1743 
1744  SpecIsoMap indicies(i,j);
1745  String buf = sr.Name()+"-"+sr.Isotope()[j].Name();
1746 
1747  ArtsMap[buf] = indicies;
1748 
1749  // Print the generated data structures (for debugging):
1750  // The explicit conversion of Name to a c-String is
1751  // necessary, because setw does not work correctly for
1752  // stl Strings.
1753  const Index& i1 = ArtsMap[buf].Speciesindex();
1754  const Index& i2 = ArtsMap[buf].Isotopeindex();
1755 
1756  out3 << " Arts Identifier = " << buf << " Species = "
1757  << setw(10) << setiosflags(ios::left)
1758  << species_data[i1].Name().c_str()
1759  << "iso = "
1760  << species_data[i1].Isotope()[i2].Name().c_str()
1761  << "\n";
1762 
1763  }
1764  }
1765  hinit = true;
1766  }
1767 
1768 
1769  // This always contains the rest of the line to parse. At the
1770  // beginning the entire line. Line gets shorter and shorter as we
1771  // continue to extract stuff from the beginning.
1772  String line;
1773 
1774  // Look for more comments?
1775  bool comment = true;
1776 
1777  while (comment)
1778  {
1779  // Return true if eof is reached:
1780  if (is.eof()) return true;
1781 
1782  // Throw runtime_error if stream is bad:
1783  if (!is) throw runtime_error ("Stream bad.");
1784 
1785  // Read line from file into linebuffer:
1786  getline(is,line);
1787 
1788  // It is possible that we were exactly at the end of the file before
1789  // calling getline. In that case the previous eof() was still false
1790  // because eof() evaluates only to true if one tries to read after the
1791  // end of the file. The following check catches this.
1792  if (line.nelem() == 0 && is.eof()) return true;
1793 
1794  // @ as first character marks catalogue entry
1795  char c;
1796  extract(c,line,1);
1797 
1798  // check for empty line
1799  if (c == '@')
1800  {
1801  comment = false;
1802  }
1803  }
1804 
1805 
1806  // read the arts identifier String
1807  istringstream icecream(line);
1808 
1809  String artsid;
1810  icecream >> artsid;
1811 
1812  if (artsid.length() != 0)
1813  {
1814 
1815  // ok, now for the cool index map:
1816  // is this arts identifier valid?
1817  const map<String, SpecIsoMap>::const_iterator i = ArtsMap.find(artsid);
1818  if ( i == ArtsMap.end() )
1819  {
1820  ostringstream os;
1821  os << "ARTS Tag: " << artsid << " is unknown.";
1822  throw runtime_error(os.str());
1823  }
1824 
1825  SpecIsoMap id = i->second;
1826 
1827 
1828  // Set mspecies:
1829  mspecies = id.Speciesindex();
1830 
1831  // Set misotope:
1832  misotope = id.Isotopeindex();
1833 
1834 
1835  // Extract center frequency:
1836  icecream >> mf;
1837 
1838 
1839  // Extract pressure shift:
1840  icecream >> mpsf;
1841 
1842  // Extract intensity:
1843  icecream >> mi0;
1844 
1845 
1846  // Extract reference temperature for Intensity in K:
1847  icecream >> mti0;
1848 
1849 
1850  // Extract lower state energy:
1851  icecream >> melow;
1852 
1853 
1854  // Extract air broadening parameters:
1855  icecream >> magam;
1856  icecream >> msgam;
1857 
1858  // Extract temperature coefficient of broadening parameters:
1859  icecream >> mnair;
1860  icecream >> mnself;
1861 
1862 
1863  // Extract reference temperature for broadening parameter in K:
1864  icecream >> mtgam;
1865 
1866  // Extract the aux parameters:
1867  Index naux;
1868  icecream >> naux;
1869 
1870  // resize the aux array and read it
1871  maux.resize(naux);
1872 
1873  for (Index i = 0; i<naux; i++)
1874  {
1875  icecream >> maux[i];
1876  //cout << "maux" << i << " = " << maux[i] << "\n";
1877  }
1878 
1879  // Extract accuracies:
1880  try
1881  {
1882  icecream >> mdf;
1883  icecream >> mdi0;
1884  icecream >> mdagam;
1885  icecream >> mdsgam;
1886  icecream >> mdnair;
1887  icecream >> mdnself;
1888  icecream >> mdpsf;
1889  }
1890  catch (runtime_error x)
1891  {
1892  // Nothing to do here, the accuracies are optional, so we
1893  // just set them to -1 and continue reading the next line of
1894  // the catalogue
1895  mdf = -1;
1896  mdi0 = -1;
1897  mdagam = -1;
1898  mdsgam = -1;
1899  mdnair = -1;
1900  mdnself = -1;
1901  mdpsf = -1;
1902  }
1903  }
1904 
1905  // That's it!
1906  return false;
1907 }
1908 
1909 
1910 
1912 {
1913  // Species lookup data:
1914  extern const Array<SpeciesRecord> species_data;
1915  // The species map. This is used to find the species id.
1916  extern std::map<String, Index> SpeciesMap;
1917  // Name of species and isotope (aux variables):
1918  String name, isoname;
1919  // Aux index:
1920  Index n;
1921 
1922 
1923  // Set frequency limits to default values (no limits):
1924  mlf = -1;
1925  muf = -1;
1926 
1927  // We cannot set a default value for the isotope, because the
1928  // default should be `ALL' and the value for `ALL' depends on the
1929  // species.
1930 
1931 
1932  // Extract the species name:
1933  n = (Index)def.find('-'); // find the '-'
1934  if (n != def.npos )
1935  {
1936  name = def.substr(0,n); // Extract before '-'
1937  def.erase(0,n+1); // Remove from def
1938  }
1939  else
1940  {
1941  // n==def.npos means that def does not contain a '-'. In that
1942  // case we assume that it contains just a species name and
1943  // nothing else
1944  name = def;
1945  def = "";
1946  }
1947 
1948 
1949 // cout << "name / def = " << name << " / " << def << endl;
1950 
1951  // Look for species name in species map:
1952  map<String, Index>::const_iterator mi = SpeciesMap.find(name);
1953  if ( mi != SpeciesMap.end() )
1954  {
1955  // Ok, we've found the species. Set mspecies.
1956  mspecies = mi->second;
1957  }
1958  else
1959  {
1960  ostringstream os;
1961  os << "Species " << name << " is not a valid species.";
1962  throw runtime_error(os.str());
1963  }
1964 
1965  // Get a reference to the relevant Species Record:
1966  const SpeciesRecord& spr = species_data[mspecies];
1967 
1968  if ( 0 == def.nelem() )
1969  {
1970  // This means that there is nothing else to parse. Apparently
1971  // the user wants all isotopes and no frequency limits.
1972  // Frequency defaults are already set. Set isotope defaults:
1973  misotope = spr.Isotope().nelem();
1974  // This means all isotopes.
1975 
1976  return;
1977  }
1978 
1979  // Extract the isotope name:
1980  n = (Index)def.find('-'); // find the '-'
1981  if (n != def.npos )
1982  {
1983  isoname = def.substr(0,n); // Extract before '-'
1984  def.erase(0,n+1); // Remove from def
1985  }
1986  else
1987  {
1988  // n==def.npos means that def does not contain a '-'. In that
1989  // case we assume that it contains just the isotope name and
1990  // nothing else.
1991  isoname = def;
1992  def = "";
1993  }
1994 
1995  // Check for joker:
1996  if ( "*" == isoname )
1997  {
1998  // The user wants all isotopes. Set this accordingly:
1999  misotope = spr.Isotope().nelem();
2000  }
2001  else
2002  {
2003  // Make an array containing the isotope names:
2004  ArrayOfString ins;
2005  for ( Index i=0; i<spr.Isotope().nelem(); ++i )
2006  ins.push_back( spr.Isotope()[i].Name() );
2007 
2008  // Use the find algorithm from the STL to find the isotope ID. It
2009  // returns an iterator, so to get the index we take the
2010  // difference to the begin() iterator.
2011  misotope = (Index)(find( ins.begin(),
2012  ins.end(),
2013  isoname ) - ins.begin());
2014 
2015  // Check if we found a matching isotope:
2016  if ( misotope >= ins.nelem() )
2017  {
2018  ostringstream os;
2019  os << "Isotope " << isoname << " is not a valid isotope for "
2020  << "species " << name << ".\n"
2021  << "Valid isotopes are:";
2022  for ( Index i=0; i<ins.nelem(); ++i )
2023  os << " " << ins[i];
2024  throw runtime_error(os.str());
2025  }
2026  }
2027 
2028  if ( 0 == def.nelem() )
2029  {
2030  // This means that there is nothing else to parse. Apparently
2031  // the user wants no frequency limits. Frequency defaults are
2032  // already set, so we can return directly.
2033 
2034  return;
2035  }
2036 
2037 
2038  // Look for the two frequency limits:
2039 
2040  // Extract first frequency
2041  n = (Index)def.find('-'); // find the '-'
2042  if (n != def.npos )
2043  {
2044  // Frequency as a String:
2045  String fname;
2046  fname = def.substr(0,n); // Extract before '-'
2047  def.erase(0,n+1); // Remove from def
2048 
2049  // Check for joker:
2050  if ( "*" == fname )
2051  {
2052  // The default for mlf is already -1, meaning `ALL'.
2053  // So there is nothing to do here.
2054  }
2055  else
2056  {
2057  // Convert to Numeric:
2058  istringstream is(fname);
2059  is >> mlf;
2060  }
2061  }
2062  else
2063  {
2064  // n==def.npos means that def does not contain a '-'. In this
2065  // case that is not allowed!
2066  throw runtime_error("You must either speciefy both frequency limits\n"
2067  "(at least with jokers), or none.");
2068  }
2069 
2070 
2071  // Now there should only be the upper frequency left in def.
2072  // Check for joker:
2073  if ( "*" == def )
2074  {
2075  // The default for muf is already -1, meaning `ALL'.
2076  // So there is nothing to do here.
2077  }
2078  else
2079  {
2080  // Convert to Numeric:
2081  istringstream is(def);
2082  is >> muf;
2083  }
2084 }
2085 
2086 
2088 {
2089  // Species lookup data:
2090  extern const Array<SpeciesRecord> species_data;
2091  // A reference to the relevant record of the species data:
2092  const SpeciesRecord& spr = species_data[mspecies];
2093  // For return value:
2094  ostringstream os;
2095 
2096  // First the species name:
2097  os << spr.Name() << "-";
2098 
2099  // Now the isotope. Can be a single isotope or ALL.
2100  if ( misotope == spr.Isotope().nelem() )
2101  {
2102  // One larger than allowed means all isotopes!
2103  os << "*-";
2104  }
2105  else
2106  {
2107  os << spr.Isotope()[misotope].Name() << "-";
2108  }
2109 
2110  // Now the frequency limits, if there are any. For this we first
2111  // need to determine the floating point precision.
2112 
2113  // Determine the precision, depending on whether Numeric is double
2114  // or float:
2115  Index precision;
2116  switch (sizeof(Numeric)) {
2117  case sizeof(float) : precision = FLT_DIG; break;
2118  case sizeof(double) : precision = DBL_DIG; break;
2119  default: out0 << "Numeric must be double or float\n"; exit(1);
2120  }
2121 
2122  if ( 0 > mlf )
2123  {
2124  // mlf < 0 means no lower limit.
2125  os << "*-";
2126  }
2127  else
2128  {
2129  os << setprecision((int)precision);
2130  os << mlf << "-";
2131  }
2132 
2133  if ( 0 > muf )
2134  {
2135  // muf < 0 means no upper limit.
2136  os << "*";
2137  }
2138  else
2139  {
2140  os << setprecision((int)precision);
2141  os << muf;
2142  }
2143 
2144  return os.str();
2145 }
2146 
2147 ostream& operator << (ostream& os, const OneTag& ot)
2148 {
2149  return os << ot.Name();
2150 }
2151 
2152 
2153 
2177  ArrayOfIndex& tags1_index,
2178  const TagGroups& tags1,
2179  const ArrayOfString& tags2_Strings )
2180 {
2181  const Index n1 = tags1.nelem();
2182  const Index n2 = tags2_Strings.nelem();
2183  TagGroups tags2; // Internal tag names for tag_Strings
2184  Index i1, i2, nj, j, found, ok;
2185 
2186  tags1_index.resize(n2);
2187  // cout << "tags2_Strings: " << tags2_Strings << "\n";
2188  tgsDefine( tags2, tags2_Strings );
2189 
2190  for ( i2=0; i2<n2; i2++ )
2191  {
2192  found = 0;
2193  for ( i1=0; (i1<n1) && !found; i1++ )
2194  {
2195  nj = tags2[i2].nelem();
2196  if ( nj == tags1[i1].nelem() )
2197  {
2198  ok = 1;
2199  for ( j=0; j<nj; j++ )
2200  {
2201  if ( tags2[i2][j].Name() != tags1[i1][j].Name() )
2202  ok = 0;
2203  }
2204  if ( ok )
2205  {
2206  found = 1;
2207  tags1_index[i2] = i1;
2208  }
2209  }
2210  }
2211  if ( !found )
2212  {
2213  ostringstream os;
2214  os << "The tag String \"" << tags2_Strings[i2] <<
2215  "\" does not match any of the given tags.\n";
2216  throw runtime_error(os.str());
2217  }
2218  }
2219 }
2220 
2221 
2234  const TagGroups& tgs1,
2235  const Array<OneTag>& tg2 )
2236 {
2237  bool found = false;
2238 
2239  for ( Index i=0;
2240  i<tgs1.nelem() && !found;
2241  ++i )
2242  {
2243  // Is at least the size correct?
2244  if ( tg2.nelem() == tgs1[i].nelem() )
2245  {
2246  bool ok = true;
2247 
2248  for ( Index j=0; j<tg2.nelem(); ++j )
2249  {
2250  if ( tg2[j].Name() != tgs1[i][j].Name() )
2251  ok = false;
2252  }
2253 
2254  if ( ok )
2255  {
2256  found = true;
2257  tgs1_index = i;
2258  }
2259  }
2260  }
2261 
2262  if ( !found )
2263  {
2264  ostringstream os;
2265  os << "The tag String \"" << tg2 <<
2266  "\" does not match any of the given tags.\n";
2267  throw runtime_error(os.str());
2268  }
2269 }
2270 
2271 
2287 {
2288  String name;
2289  Index i;
2290 
2291  for ( i=0; i<tg.nelem()-1; ++i )
2292  {
2293  name += tg[i].Name() + ", ";
2294  }
2295  name += tg[i].Name();
2296 
2297  return name;
2298 }
2299 
2300 
2309 void write_lines_to_stream(ostream& os,
2310  const ArrayOfLineRecord& lines)
2311 {
2312  // We need this dummy line record, so that we can get the catalogue
2313  // version tag from dummy.Version, even if the line list is empty.
2314  LineRecord dummy;
2315 
2316  // Output version String first, followed by number of lines:
2317  os << dummy.Version() << " " << lines.nelem() << "\n";
2318 
2319  // Now output the line data:
2320  for ( Index i=0; i<lines.nelem(); ++i )
2321  {
2322  // out3 << lines[i] << "\n";
2323  os << lines[i] << "\n";
2324  }
2325 }
2326 
2327 
2329 
2341 {
2342  if( x.nelem() > 1 )
2343  {
2344  for( Index i=1; i<x.nelem(); i++ )
2345  {
2346  if( x[i] < x[i-1] )
2347  return false;
2348  }
2349  }
2350  return true;
2351 }
2352 
2353 
2380  ConstVectorView f_mono,
2381  ConstVectorView p_abs,
2382  ConstVectorView t_abs,
2383  ConstVectorView h2o_abs,
2384  ConstVectorView vmr,
2385  const ArrayOfLineRecord& lines,
2386  const Index ind_ls,
2387  const Index ind_lsn,
2388  const Numeric cutoff)
2389 {
2390  // Make lineshape and species lookup data visible:
2393 
2394  // speed of light constant
2395  extern const Numeric SPEED_OF_LIGHT;
2396 
2397  // Boltzmann constant
2398  extern const Numeric BOLTZMAN_CONST;
2399 
2400  // Avogadros constant
2401  extern const Numeric AVOGADROS_NUMB;
2402 
2403  // Planck constant
2404  extern const Numeric PLANCK_CONST;
2405 
2406  // sqrt(ln(2))
2407  // extern const Numeric SQRT_NAT_LOG_2;
2408 
2409  // Constant within the Doppler Broadening calculation:
2410  static const Numeric doppler_const = sqrt( 2.0 * BOLTZMAN_CONST *
2412 
2413  // dimension of f_mono, lines
2414  Index nf = f_mono.nelem();
2415  Index nl = lines.nelem();
2416 
2417  // Define the vector for the line shape function and the
2418  // normalization factor of the lineshape here, so that we don't need
2419  // so many free store allocations. the last element is used to
2420  // calculate the value at the cutoff frequency
2421  Vector ls(nf+1);
2422  Vector fac(nf+1);
2423 
2424  bool cut = (cutoff != -1) ? true : false;
2425 
2426  // Check that the frequency grid is sorted in the case of lineshape
2427  // with cutoff. Duplicate frequency values are allowed.
2428  if (cut)
2429  {
2430  if ( ! is_sorted( f_mono ) )
2431  {
2432  ostringstream os;
2433  os << "If you use a lineshape function with cutoff, your\n"
2434  << "frequency grid *f_mono* must be sorted.\n"
2435  << "(Duplicate values are allowed.)";
2436  throw runtime_error(os.str());
2437  }
2438  }
2439 
2440  // Check that all temperatures are non-negative
2441  bool negative = false;
2442 
2443  for (Index i = 0; !negative && i < t_abs.nelem (); i++)
2444  {
2445  if (t_abs[i] < 0.)
2446  negative = true;
2447  }
2448 
2449  if (negative)
2450  {
2451  ostringstream os;
2452  os << "t_abs contains at least one negative temperature value.\n"
2453  << "This is not allowed.";
2454  throw runtime_error(os.str());
2455  }
2456 
2457  // We need a local copy of f_mono which is 1 element longer, because
2458  // we append a possible cutoff to it.
2459  // The initialization of this has to be inside the line loop!
2460  Vector f_local( nf + 1 );
2461 
2462  // Voigt generally needs a different frequency grid. If we allocate
2463  // that in the outer loop, instead of in voigt, we don't have the
2464  // free store allocation at each lineshape call. Calculation is
2465  // still done in the voigt routine itself, this is just an auxillary
2466  // parameter, passed to lineshape. For selected lineshapes (e.g.,
2467  // Rosenkranz) it is used additionally to pass parameters needed in
2468  // the lineshape (e.g., overlap, ...). Consequently we have to
2469  // assure that aux has a dimension not less then the number of
2470  // parameters passed.
2471  Index ii = (nf+1 < 10) ? 10 : nf+1;
2472  Vector aux(ii);
2473 
2474  // Check that p_abs, t_abs, and h2o_abs all have the same
2475  // dimension. This could be a user error, so we throw a
2476  // runtime_error. FIXME: why do we do this for each tag? wouldn't a
2477  // check in abscalc be sufficient?
2478 
2479  if ( t_abs.nelem() != p_abs.nelem() )
2480  {
2481  ostringstream os;
2482  os << "Variable t_abs must have the same dimension as p_abs.\n"
2483  << "t_abs.nelem() = " << t_abs.nelem() << '\n'
2484  << "p_abs.nelem() = " << p_abs.nelem();
2485  throw runtime_error(os.str());
2486  }
2487 
2488  if ( vmr.nelem() != p_abs.nelem() )
2489  {
2490  ostringstream os;
2491  os << "Variable vmr must have the same dimension as p_abs.\n"
2492  << "vmr.nelem() = " << vmr.nelem() << '\n'
2493  << "p_abs.nelem() = " << p_abs.nelem();
2494  throw runtime_error(os.str());
2495  }
2496 
2497  if ( h2o_abs.nelem() != p_abs.nelem() )
2498  {
2499  ostringstream os;
2500  os << "Variable h2o_abs must have the same dimension as p_abs.\n"
2501  << "h2o_abs.nelem() = " << h2o_abs.nelem() << '\n'
2502  << "p_abs.nelem() = " << p_abs.nelem();
2503  throw runtime_error(os.str());
2504  }
2505 
2506 
2507  // Check that the dimension of xsec is indeed [f_mono.nelem(),
2508  // p_abs.nelem()]:
2509  if ( xsec.nrows() != nf || xsec.ncols() != p_abs.nelem() )
2510  {
2511  ostringstream os;
2512  os << "Variable xsec must have dimensions [f_mono.nelem(),p_abs.nelem()].\n"
2513  << "[xsec.nrows(),xsec.ncols()] = [" << xsec.nrows()
2514  << ", " << xsec.ncols() << "]\n"
2515  << "f_mono.nelem() = " << nf << '\n'
2516  << "p_abs.nelem() = " << p_abs.nelem();
2517  throw runtime_error(os.str());
2518  }
2519 
2520 
2521  // Loop all pressures:
2522  for ( Index i=0; i<p_abs.nelem(); ++i )
2523  {
2524 
2525  // store variables p_abs[i] and t_abs[i],
2526  // this is slightly faster
2527  Numeric p_i = p_abs[i];
2528  Numeric t_i = t_abs[i];
2529 
2530 
2531  //out3 << " p = " << p_i << " Pa\n";
2532 
2533  // Calculate total number density from pressure and temperature.
2534  // n = n0*T0/p0 * p/T or n = p/kB/t, ideal gas law
2535  Numeric n = p_i / BOLTZMAN_CONST / t_i;
2536 
2537  // For the pressure broadening, we also need the partial pressure:
2538  const Numeric p_partial = p_i * vmr[i];
2539 
2540 
2541  // Loop all lines:
2542  for ( Index l=0; l< nl; ++l )
2543  {
2544  // Copy f_mono to the beginning of f_local. There is one
2545  // element left at the end of f_local.
2546  // THIS HAS TO BE INSIDE THE LINE LOOP, BECAUSE THE CUTOFF
2547  // FREQUENCY IS ALWAYS PUT IN A DIFFERENT PLACE!
2548  f_local[Range(0,nf)] = f_mono;
2549 
2550  // This will hold the actual number of frequencies to add to
2551  // xsec later on:
2552  Index nfl = nf;
2553 
2554  // This will hold the actual number of frequencies for the
2555  // call to the lineshape functions later on:
2556  Index nfls = nf;
2557 
2558  // The baseline to substract for cutoff frequency
2559  Numeric base=0.0;
2560 
2561  // lines[l] is used several times, this construct should be
2562  // faster (Oliver Lemke)
2563  LineRecord l_l = lines[l]; // which line are we dealing with
2564  // Center frequency in vacuum:
2565  Numeric F0 = l_l.F();
2566 
2567  // Intensity is already in the right units (Hz*m^2). It also
2568  // includes already the isotopic ratio. Needs only to be
2569  // coverted to the actual temperature and multiplied by total
2570  // number density and lineshape.
2571  Numeric intensity = l_l.I0();
2572 
2573  // Lower state energy is already in the right unit (Joule).
2574  Numeric e_lower = l_l.Elow();
2575 
2576  // Upper state energy:
2577  Numeric e_upper = e_lower + F0 * PLANCK_CONST;
2578 
2579  // Get the ratio of the partition function.
2580  // This will throw a runtime error if no data exists.
2581  // Important: This function needs both the reference
2582  // temperature and the actual temperature, because the
2583  // reference temperature can be different for each line,
2584  // even of the same species.
2585  Numeric part_fct_ratio =
2587  t_i );
2588 
2589  // Boltzmann factors
2590  Numeric nom = exp(- e_lower / ( BOLTZMAN_CONST * t_i ) ) -
2591  exp(- e_upper / ( BOLTZMAN_CONST * t_i ) );
2592 
2593  Numeric denom = exp(- e_lower / ( BOLTZMAN_CONST * l_l.Ti0() ) ) -
2594  exp(- e_upper / ( BOLTZMAN_CONST * l_l.Ti0() ) );
2595 
2596 
2597  // intensity at temperature
2598  // (calculate the line intensity according to the standard
2599  // expression given in HITRAN)
2600  intensity *= part_fct_ratio * nom / denom;
2601 
2602  if (lineshape_norm_data[ind_lsn].Name() == "quadratic")
2603  {
2604  // in case of the quadratic normalization factor use the
2605  // so called 'microwave approximation' of the line intensity
2606  // given by
2607  // P. W. Rosenkranz, Chapter 2, Eq.(2.16), in M. A. Janssen,
2608  // Atmospheric Remote Sensing by Microwave Radiometry,
2609  // John Wiley & Sons, Inc., 1993
2610  Numeric mafac = (PLANCK_CONST * F0) / (2.000e0 * BOLTZMAN_CONST * t_i);
2611  intensity = intensity * mafac / sinh(mafac);
2612  }
2613 
2614 
2615  // 2. Get pressure broadened line width:
2616  // (Agam is in Hz/Pa, p_abs is in Pa, gamma is in Hz)
2617  const Numeric theta = l_l.Tgam() / t_i;
2618  const Numeric theta_Nair = pow(theta, l_l.Nair());
2619 
2620  Numeric gamma
2621  = l_l.Agam() * theta_Nair * (p_i - p_partial)
2622  + l_l.Sgam() * pow(theta, l_l.Nself()) * p_partial;
2623 
2624  // 3. Doppler broadening without the sqrt(ln(2)) factor, which
2625  // seems to be redundant FIXME: verify .
2626  Numeric sigma = F0 * doppler_const *
2627  sqrt( t_i / l_l.IsotopeData().Mass());
2628 
2629  // 3.a. Put in pressure shift.
2630  // The T dependence is connected to that of agam by:
2631  // n_shift = .25 + 1.5 * n_agam
2632  // Theta has been initialized above.
2633  F0 += l_l.Psf() * p_i *
2634  std::pow( theta , (Numeric).25 + (Numeric)1.5*l_l.Nair() );
2635 
2636  // 4. the rosenkranz lineshape for oxygen calculates the
2637  // pressure broadening, overlap, ... differently. Therefore
2638  // all required parameters are passed in the aux array, the
2639  // given order must agree with how they are accessed in the
2640  // used lineshape (currently only the Rosenkranz routines are
2641  // using this method of passing parameters). Parameters are
2642  // always passed, because the Rosenkranz lineshape function
2643  // can be used without overlap correction, which is then just
2644  // set to 0. I know, this is not very nice, but I suggested to
2645  // put the oxygen absorption into a different workspace
2646  // method, but nobody listened to me, Schnief.
2647  aux[0] = theta;
2648  aux[1] = theta_Nair;
2649  aux[2] = p_i;
2650  aux[3] = vmr[i];
2651  aux[4] = h2o_abs[i];
2652  aux[5] = l_l.Agam();
2653  aux[6] = l_l.Nair();
2654  // is overlap available, otherwise pass zero
2655  if (l_l.Naux() > 1)
2656  {
2657  aux[7] = l_l.Aux()[0];
2658  aux[8] = l_l.Aux()[1];
2659  // cout << "aux7, aux8: " << aux[7] << " " << aux[8] << "\n";
2660  }
2661  else
2662  {
2663  aux[7] = 0.0;
2664  aux[8] = 0.0;
2665  }
2666 
2667 
2668  // Indices pointing at begin/end frequencies of f_mono or at
2669  // the elements that have to be calculated in case of cutoff
2670  Index i_f_min = 0;
2671  Index i_f_max = nf-1;
2672 
2673  // cutoff ?
2674  if ( cut )
2675  {
2676  // Check whether we have elements in ls that can be
2677  // ignored at lower frequencies of f_mono.
2678  //
2679  // Loop through all frequencies, finding min value and
2680  // set all values to zero on that way.
2681  while ( i_f_min < nf && (F0 - cutoff) > f_mono[i_f_min] )
2682  {
2683  // ls[i_f_min] = 0;
2684  ++i_f_min;
2685  }
2686 
2687 
2688  // Check whether we have elements in ls that can be
2689  // ignored at higher frequencies of f_mono.
2690  //
2691  // Loop through all frequencies, finding max value and
2692  // set all values to zero on that way.
2693  while ( i_f_max >= 0 && (F0 + cutoff) < f_mono[i_f_max] )
2694  {
2695  // ls[i_f_max] = 0;
2696  --i_f_max;
2697  }
2698 
2699  // Append the cutoff frequency to f_local:
2700  ++i_f_max;
2701  f_local[i_f_max] = F0 + cutoff;
2702 
2703  // Number of frequencies to calculate:
2704  nfls = i_f_max - i_f_min + 1; // Add one because indices
2705  // are pointing to first and
2706  // last valid element. This
2707  // is for the lineshape
2708  // calls.
2709  nfl = nfls -1; // This is for xsec.
2710  }
2711  else
2712  {
2713  // Nothing to do here. Note that nfl and nfls are both still set to nf.
2714  }
2715 
2716  // cout << "nf, nfl, nfls = " << nf << ", " << nfl << ", " << nfls << ".\n";
2717 
2718  // Maybe there are no frequencies left to compute? Note that
2719  // the number that counts here is nfl, since only these are
2720  // the `real' frequencies, for which xsec is changed. nfls
2721  // will always be at least one, because it contains the cutoff.
2722  if ( nfl > 0 )
2723  {
2724  // Calculate the line shape:
2725  lineshape_data[ind_ls].Function()(ls,
2726  aux,F0,gamma,sigma,
2727  f_local[Range(i_f_min,nfls)],
2728  nfls);
2729 
2730  // Calculate the chosen normalization factor:
2731  lineshape_norm_data[ind_lsn].Function()(fac,F0,
2732  f_local[Range(i_f_min,nfls)],
2733  t_i,nfls);
2734 
2735  // Get a handle on the range of xsec that we want to change.
2736  // We use nfl here, which could be one less than nfls.
2737  VectorView this_xsec = xsec(Range(i_f_min,nfl),i);
2738 
2739  // Get handles on the range of ls and fac that we need.
2740  VectorView this_ls = ls[Range(0,nfl)];
2741  VectorView this_fac = fac[Range(0,nfl)];
2742 
2743  // cutoff ?
2744  if ( cut )
2745  {
2746  // The index nfls-1 should be exactly the index pointing
2747  // to the value at the cutoff frequency.
2748  base = ls[nfls-1];
2749 
2750  // Subtract baseline from xsec.
2751  // this_xsec -= base;
2752  this_ls -= base;
2753  }
2754 
2755  // Add line to xsec.
2756  {
2757  // To make the loop a bit faster, precompute all constant
2758  // factors. These are:
2759  // 1. Total number density of the air.
2760  // 2. Line intensity.
2761  // 3. Isotopic ratio.
2762  //
2763  // The isotopic ratio must be applied here, since we are
2764  // summing up lines belonging to different isotopes.
2765 
2766  const Numeric factors = n * intensity * l_l.IsotopeData().Abundance();
2767 
2768  // We have to do:
2769  // xsec(j,i) += factors * ls[j1] * fac[j1];
2770  //
2771  // We use ls as a dummy to compute the product, then add it
2772  // to this_xsec.
2773 
2774  // FIXME: Maybe try if this is faster with a good
2775  // old-fashioned simple loop.
2776 
2777  this_ls *= this_fac;
2778  this_ls *= factors;
2779  this_xsec += this_ls;
2780  }
2781  }
2782  }
2783  }
2784 }
2785 
2786 
2787 
2788 
2789 //======================================================================
2790 // Functions related to refraction
2791 //======================================================================
2792 
2793 
2795 
2811  Vector& refr_index,
2812  ConstVectorView p_abs,
2813  ConstVectorView t_abs )
2814 {
2815  const Index n = p_abs.nelem();
2816  refr_index.resize( n );
2817 
2818  assert ( n == t_abs.nelem() );
2819 
2820  // N = 77.593e-2 * p / t ppm
2821  for ( Index i=0; i<n; i++ )
2822  refr_index[i] = 1.0 + 77.593e-8 * p_abs[i] / t_abs[i];
2823 }
2824 
2825 
2826 
2828 
2843  Vector& refr_index,
2844  ConstVectorView p_abs,
2845  ConstVectorView t_abs,
2846  ConstVectorView h2o_abs )
2847 {
2848  const Index n = p_abs.nelem();
2849  refr_index.resize( n );
2850 
2851  assert ( n == t_abs.nelem() );
2852  assert ( n == h2o_abs.nelem() );
2853 
2854  Numeric e; // Partial pressure of water in Pa
2855  Numeric p; // Partial pressure of the dry air: p = p_tot - e
2856 
2857  for ( Index i=0; i<n; i++ )
2858  {
2859  e = p_abs[i] * h2o_abs[i];
2860  p = p_abs[i] - e;
2861 
2862  refr_index[i] = 1.0 + 77.593e-8 * p / t_abs[i] +
2863  72e-8 * e / t_abs[i] +
2864  3.754e-3 * e / (t_abs[i]*t_abs[i]);
2865  }
2866 }
2867 
2880 {
2881  // Planck constant [Js]
2882  extern const Numeric PLANCK_CONST;
2883 
2884  // Speed of light [m/s]
2885  extern const Numeric SPEED_OF_LIGHT;
2886 
2887  // Constant to convert lower state energy from cm^-1 to J
2888  static const Numeric lower_energy_const = PLANCK_CONST * SPEED_OF_LIGHT * 1E2;
2889 
2890  return e*lower_energy_const;
2891 }
2892 
2893 //======================================================================
2894 // Functions to convert the accuracy index to ARTS units
2895 //======================================================================
2896 
2897 // ********* for HITRAN database *************
2898 //convert HITRAN index for line position accuracy to ARTS
2899 //units (Hz).
2900 
2902  Numeric& mdf,
2903  const Index& df
2904  )
2905 {
2906  switch ( df )
2907  {
2908  case 0:
2909  {
2910  mdf = -1;
2911  break;
2912  }
2913  case 1:
2914  {
2915  mdf = 30*1E9;
2916  break;
2917  }
2918  case 2:
2919  {
2920  mdf = 3*1E9;
2921  break;
2922  }
2923  case 3:
2924  {
2925  mdf = 300*1E6;
2926  break;
2927  }
2928  case 4:
2929  {
2930  mdf = 30*1E6;
2931  break;
2932  }
2933  case 5:
2934  {
2935  mdf = 3*1E6;
2936  break;
2937  }
2938  case 6:
2939  {
2940  mdf = 0.3*1E6;
2941  break;
2942  }
2943  }
2944 }
2945 
2946 //convert HITRAN index for intensity and halfwidth accuracy to ARTS
2947 //units (relative difference).
2949  Numeric& mdh,
2950  const Index& dh
2951  )
2952 {
2953  switch ( dh )
2954  {
2955  case 0:
2956  {
2957  mdh = -1;
2958  break;
2959  }
2960  case 1:
2961  {
2962  mdh = -1;
2963  break;
2964  }
2965  case 2:
2966  {
2967  mdh = -1;
2968  break;
2969  }
2970  case 3:
2971  {
2972  mdh = 30;
2973  break;
2974  }
2975  case 4:
2976  {
2977  mdh = 20;
2978  break;
2979  }
2980  case 5:
2981  {
2982  mdh = 10;
2983  break;
2984  }
2985  case 6:
2986  {
2987  mdh =5;
2988  break;
2989  }
2990  case 7:
2991  {
2992  mdh =2;
2993  break;
2994  }
2995  case 8:
2996  {
2997  mdh =1;
2998  break;
2999  }
3000  }
3001  mdh=mdh/100;
3002 }
3003 
3004 // ********* for MYTRAN database *************
3005 //convert MYTRAN index for intensity and halfwidth accuracy to ARTS
3006 //units (relative difference).
3008  Numeric& mdh,
3009  const Index & dh
3010  )
3011 {
3012  switch ( dh )
3013  {
3014  case 0:
3015  {
3016  mdh = 200;
3017  break;
3018  }
3019  case 1:
3020  {
3021  mdh = 100;
3022  break;
3023  }
3024  case 2:
3025  {
3026  mdh = 50;
3027  break;
3028  }
3029  case 3:
3030  {
3031  mdh = 30;
3032  break;
3033  }
3034  case 4:
3035  {
3036  mdh = 20;
3037  break;
3038  }
3039  case 5:
3040  {
3041  mdh = 10;
3042  break;
3043  }
3044  case 6:
3045  {
3046  mdh =5;
3047  break;
3048  }
3049  case 7:
3050  {
3051  mdh = 2;
3052  break;
3053  }
3054  case 8:
3055  {
3056  mdh = 1;
3057  break;
3058  }
3059  case 9:
3060  {
3061  mdh = 0.5;
3062  break;
3063  }
3064  }
3065  mdh=mdh/100;
3066 }
3067 
extract
void extract(T &x, String &line, Index n)
Extract something from a catalogue line.
Definition: absorption.cc:127
LineRecord::mdnair
Numeric mdnair
Definition: absorption.h:1001
LineRecord::IsotopeData
const IsotopeRecord & IsotopeData() const
The matching IsotopeRecord from species_data.
Definition: absorption.h:573
LineRecord::Ti0
Numeric Ti0() const
Reference temperature for I0 in K:
Definition: absorption.h:609
IsotopeRecord::CalculatePartitionFctAtTemp
Numeric CalculatePartitionFctAtTemp(Numeric temperature) const
Definition: absorption.cc:48
LineRecord::dNair
Numeric dNair() const
Accuracy for AGAM temperature exponent in relative value :
Definition: absorption.h:665
refr_index_Boudouris
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.
Definition: absorption.cc:2842
MatrixView
The MatrixView class.
Definition: matpackI.h:476
PLANCK_CONST
const Numeric PLANCK_CONST
LineRecord::Sgam
Numeric Sgam() const
Self broadened width in Hz/Pa:
Definition: absorption.h:621
auto_md.h
SpecIsoMap
Definition: absorption.h:1009
AVOGADROS_NUMB
const Numeric AVOGADROS_NUMB
Global constant, the Avogadro's number [molec/kg].
refr_index_BoudourisDryAir
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.
Definition: absorption.cc:2810
LineRecord::melow
Numeric melow
Definition: absorption.h:976
absorption.h
Declarations required for the calculation of absorption coefficients.
OneTag::OneTag
OneTag()
Default constructor.
Definition: absorption.h:1064
LineRecord::Agam
Numeric Agam() const
Air broadened width in Hz/Pa:
Definition: absorption.h:615
LineRecord::mi0
Numeric mi0
Definition: absorption.h:972
LineRecord::Nself
Numeric Nself() const
SGAM temperature exponent (dimensionless):
Definition: absorption.h:633
get_tag_group_name
String get_tag_group_name(const Array< OneTag > &tg)
Print the name of a tag group.
Definition: absorption.cc:2286
tgsDefine
void tgsDefine(TagGroups &, const ArrayOfString &tags)
Definition: m_abs.cc:901
convMytranIER
void convMytranIER(Numeric &mdh, const Index &dh)
Definition: absorption.cc:3007
out0
Out0 out0
Level 0 output stream.
Definition: messages.cc:50
convHitranIERSH
void convHitranIERSH(Numeric &mdh, const Index &dh)
Definition: absorption.cc:2948
operator<<
ostream & operator<<(ostream &os, const LineRecord &lr)
Output operator for LineRecord.
Definition: absorption.cc:76
LineRecord::Version
String Version() const
Return the version String.
Definition: absorption.h:521
Vector::resize
void resize(Index n)
Resize function.
Definition: matpackI.h:1467
SpeciesMap
std::map< String, Index > SpeciesMap
The map associated with species_data.
Definition: absorption.cc:42
LineRecord::mdpsf
Numeric mdpsf
Definition: absorption.h:1005
ConstMatrixView::nrows
Index nrows() const
Returns the number of rows.
Definition: matpackI.h:1491
LineRecord::mpsf
Numeric mpsf
Definition: absorption.h:970
IsotopeRecord::Abundance
const Numeric & Abundance() const
Normal abundance ( = isotopic ratio).
Definition: absorption.h:216
get_tag_group_index_for_tag_group
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.
Definition: absorption.cc:2233
LineRecord::ReadFromHitran2004Stream
bool ReadFromHitran2004Stream(istream &is)
Read one line from a stream associated with a HITRAN 2004 file.
Definition: absorption.cc:609
xsec_species
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.
Definition: absorption.cc:2379
LineRecord::mnself
Numeric mnself
Definition: absorption.h:984
OneTag::Name
String Name() const
Return the full name of this tag according to Bredbeck convention.
Definition: absorption.cc:2087
LineRecord::ReadFromJplStream
bool ReadFromJplStream(istream &is)
Read one line from a stream associated with a JPL file.
Definition: absorption.cc:1474
Array
This can be used to make arrays out of anything.
Definition: array.h:48
LineRecord::Elow
Numeric Elow() const
Lower state energy in cm^-1:
Definition: absorption.h:612
IsotopeRecord::mqcoeff
ArrayOfNumeric mqcoeff
Definition: absorption.h:277
LineRecord::Psf
Numeric Psf() const
The pressure shift parameter in Hz/Pa.
Definition: absorption.h:586
SPEED_OF_LIGHT
const Numeric SPEED_OF_LIGHT
species_data
Array< SpeciesRecord > species_data
Definition: species_data.cc:36
messages.h
Declarations having to do with the four output streams.
LineRecord::ReadFromMytran2Stream
bool ReadFromMytran2Stream(istream &is)
Read one line from a stream associated with a MYTRAN2 file.
Definition: absorption.cc:1075
LineRecord::Name
String Name() const
The full name of the species and isotope.
Definition: absorption.h:540
write_lines_to_stream
void write_lines_to_stream(ostream &os, const ArrayOfLineRecord &lines)
A helper function that writes lines in a line list to a stream.
Definition: absorption.cc:2309
ConstMatrixView::ncols
Index ncols() const
Returns the number of columns.
Definition: matpackI.h:1497
TORR2PA
const Numeric TORR2PA
Global constant, converts torr to Pa.
my_basic_string
The implementation for String, the ARTS string class.
Definition: mystring.h:61
LineRecord::ReadFromHitranStream
bool ReadFromHitranStream(istream &is)
Read one line from a stream associated with a HITRAN 1986-2001 file.
Definition: absorption.cc:159
convHitranIERF
void convHitranIERF(Numeric &mdf, const Index &df)
Definition: absorption.cc:2901
LineRecord::mf
Numeric mf
Definition: absorption.h:968
LineRecord::mti0
Numeric mti0
Definition: absorption.h:974
LineRecord::misotope
Index misotope
Definition: absorption.h:966
LineRecord::mdnself
Numeric mdnself
Definition: absorption.h:1003
LineRecord::dNself
Numeric dNself() const
Accuracy for SGAM temperature exponent in relative value:
Definition: absorption.h:668
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: arts.h:147
ConstVectorView::nelem
Index nelem() const
Returns the number of elements.
Definition: matpackI.h:937
VectorView
The VectorView class.
Definition: matpackI.h:281
is_sorted
bool is_sorted(ConstVectorView x)
Checks if a vector is sorted in ascending order.
Definition: absorption.cc:2340
SpeciesRecord::Isotope
const Array< IsotopeRecord > & Isotope() const
Definition: absorption.h:339
SpecIsoMap::Speciesindex
const Index & Speciesindex() const
Definition: absorption.h:1019
define_species_map
void define_species_map()
Define the species data map.
Definition: absorption.cc:64
math_funcs.h
Contains declerations of basic mathematical and vector/matrix functions.
LineRecord
Spectral line catalog data.
Definition: absorption.h:438
LineRecord::mdf
Numeric mdf
Definition: absorption.h:993
ATM2PA
const Numeric ATM2PA
Global constant, converts atm to Pa.
LineRecord::magam
Numeric magam
Definition: absorption.h:978
LineRecord::maux
ArrayOfNumeric maux
Definition: absorption.h:988
SpeciesRecord
Contains the lookup data for one species.
Definition: absorption.h:284
max
#define max(a, b)
Definition: continua.cc:12949
my_basic_string::nelem
Index nelem() const
Number of elements.
Definition: mystring.h:185
LineRecord::mspecies
Index mspecies
Definition: absorption.h:964
OneTag
A tag group can consist of the sum of several of these.
Definition: absorption.h:1061
LineRecord::dSgam
Numeric dSgam() const
Accuracy for self broadened width in relative value :
Definition: absorption.h:662
LineRecord::dAgam
Numeric dAgam() const
Accuracy for air broadened width in relative value :
Definition: absorption.h:659
LineRecord::dF
Numeric dF() const
Accuracy for line position in Hz :
Definition: absorption.h:653
iso
void iso(Array< IsotopeRecord >::iterator &ii, String name, const ArrayOfNumeric &coeff)
Definition: partition_function_data.cc:706
Range
The range class.
Definition: matpackI.h:139
out3
Out3 out3
Level 3 output stream.
Definition: messages.cc:56
IsotopeRecord::Mass
const Numeric & Mass() const
Mass of the isotope.
Definition: absorption.h:219
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:2879
OneTag::misotope
Index misotope
Definition: absorption.h:1104
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: arts.h:153
LineRecord::mtgam
Numeric mtgam
Definition: absorption.h:986
LineRecord::ReadFromArtsStream
bool ReadFromArtsStream(istream &is)
Read one line from a stream associated with an Arts file.
Definition: absorption.cc:1718
SpeciesRecord::Name
const String & Name() const
Definition: absorption.h:337
LineRecord::Naux
Index Naux() const
Number of auxiliary parameters.
Definition: absorption.h:646
OneTag::mspecies
Index mspecies
Definition: absorption.h:1100
lineshape_data
Array< LineshapeRecord > lineshape_data
Definition: lineshapes.cc:2062
LineRecord::msgam
Numeric msgam
Definition: absorption.h:980
get_tagindex_for_Strings
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.
Definition: absorption.cc:2176
OneTag::mlf
Numeric mlf
Definition: absorption.h:1107
lineshape_norm_data
Array< LineshapeNormRecord > lineshape_norm_data
Definition: lineshapes.cc:2159
LineRecord::mdi0
Numeric mdi0
Definition: absorption.h:995
BOLTZMAN_CONST
const Numeric BOLTZMAN_CONST
LineRecord::mdagam
Numeric mdagam
Definition: absorption.h:997
my_basic_string::npos
static const Index npos
Define npos:
Definition: mystring.h:87
LineRecord::Tgam
Numeric Tgam() const
Reference temperature for AGAM and SGAM in K:
Definition: absorption.h:639
Vector
The Vector class.
Definition: matpackI.h:389
LineRecord::F
Numeric F() const
The line center frequency in Hz.
Definition: absorption.h:580
LineRecord::Aux
const ArrayOfNumeric & Aux() const
Auxiliary parameters.
Definition: absorption.h:649
LineRecord::mdsgam
Numeric mdsgam
Definition: absorption.h:999
ConstVectorView
A constant view of a Vector.
Definition: matpackI.h:232
Array::nelem
Index nelem() const
Number of elements.
Definition: array.h:115
LineRecord::mnair
Numeric mnair
Definition: absorption.h:982
LineRecord::I0
Numeric I0() const
The line intensity in m^2*Hz at the reference temperature Ti0.
Definition: absorption.h:603
LineRecord::dPsf
Numeric dPsf() const
Accuracy for pressure shift in relative value :
Definition: absorption.h:671
LineRecord::dI0
Numeric dI0() const
Accuracy for line intensity in relative value :
Definition: absorption.h:656
arts.h
The global header file for ARTS.
LineRecord::Nair
Numeric Nair() const
AGAM temperature exponent (dimensionless):
Definition: absorption.h:627
OneTag::muf
Numeric muf
Definition: absorption.h:1110
IsotopeRecord::CalculatePartitionFctRatio
Numeric CalculatePartitionFctRatio(Numeric reference_temperature, Numeric actual_temperature) const
Calculate partition function ratio.
Definition: absorption.h:247