ARTS 2.5.10 (git: 2f1c442c)
m_absorptionlines.cc
Go to the documentation of this file.
1/* Copyright (C) 2019
2 Richard Larsson <ric.larsson@gmail.com>
3
4
5 This program is free software; you can redistribute it and/or modify it
6 under the terms of the GNU General Public License as published by the
7 Free Software Foundation; either version 2, or (at your option) any
8 later version.
9
10
11 This program is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16
17 You should have received a copy of the GNU General Public License
18 along with this program; if not, write to the Free Software
19 Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
20 USA. */
21
30#include "absorptionlines.h"
31#include "array.h"
32#include "arts_omp.h"
33#include "artstime.h"
34#include "auto_md.h"
35#include "debug.h"
36#include "enums.h"
37#include "file.h"
38#include "global_data.h"
39#include "lineshapemodel.h"
40#include "m_xml.h"
41#include "quantum_numbers.h"
42#include <algorithm>
43#include <atomic>
44#include <exception>
45#include <iterator>
46
50
52 const Verbosity&)
53{
54 abs_lines.erase(std::remove_if(abs_lines.begin(), abs_lines.end(),
55 [](auto& x){return x.NumLines() == 0;}),
56 abs_lines.end());
57}
58
60 const Verbosity& verbosity)
61{
62 const Index n = abs_lines.nelem();
63
64 for (Index i=0; i<n; i++) {
65 AbsorptionLines& band = abs_lines[i];
66 if (band.NumLines()) {
67 for (Index j=i+1; j<n; j++) {
68 if (band.Match(abs_lines[j]).first) {
69 for (auto& line: abs_lines[j].lines) {
70 band.AppendSingleLine(line);
71 }
72 abs_lines[j].lines.clear();
73 }
74 }
75 }
76 }
77
78 abs_linesRemoveEmptyBands(abs_lines, verbosity);
79}
80
82 const Verbosity& verbosity)
83{
84 for (auto& abs_lines: abs_lines_per_species) abs_linesFlatten(abs_lines, verbosity);
85}
86
90
91
92
101 const QuantumIdentifier& global_qid) {
102 auto band =
103 std::find_if(abs_lines.begin(), abs_lines.end(), [&](const auto& li) {
104 return li.MatchWithExternal(sline, global_qid);
105 });
106 if (band not_eq abs_lines.end()) {
107 band->AppendSingleLine(sline.line);
108 } else {
109 abs_lines.emplace_back(sline.selfbroadening,
110 sline.bathbroadening,
111 sline.cutoff,
112 sline.mirroring,
113 sline.population,
114 sline.normalization,
115 sline.lineshapetype,
116 sline.T0,
117 sline.cutofffreq,
118 sline.linemixinglimit,
119 global_qid,
120 sline.species,
121 Array<AbsorptionSingleLine>{sline.line});
122 }
123}
124
126
133 const ArrayOfAbsorptionLines& local_lines) {
134 for (auto& band : local_lines) {
135 if (auto ptr = std::find_first_of(
136 abs_lines.begin(),
137 abs_lines.end(),
138 &band,
139 &band,
140 [](const auto& abs_band, const auto& local_band) {
141 return abs_band.Match(local_band).first;
142 });
143 ptr == abs_lines.end())
144 abs_lines.push_back(band);
145 else
146 for (auto& line : band.lines) ptr->lines.push_back(line);
147 }
148}
149
158 for(auto qn: qns) {
159 if (qid.val.has(qn)) {
160 out.val.set(qid.val[qn]);
161 }
162 }
163 return out;
164}
165
172Array<QuantumNumberType> string2vecqn(std::string_view qnstr) {
173 using namespace Quantum::Number;
174
175 if (qnstr == "DEFAULT_GLOBAL") {
176 return Array<Type>(global_types);
177 }
178 if (qnstr == "DEFAULT_LOCAL") {
179 return Array<Type>(local_types);
180 }
181
182 const Index N = count_items(qnstr);
183 Array<Type> nums(N);
184 for (Index i = 0; i < N; i++) {
185 nums[i] = toTypeOrThrow(items(qnstr, i));
186 }
187
188 return nums;
189}
190
191bool check_local(const Array<QuantumNumberType>& local_state) {
192 using namespace Quantum::Number;
193 for (auto qn: local_state) if (common_value_type(common_value_type(qn), ValueType::H) not_eq ValueType::H) return false;
194 return true;
195}
196
197/* Workspace method: Doxygen documentation will be auto-generated */
199 const String& artscat_file,
200 const Numeric& fmin,
201 const Numeric& fmax,
202 const String& globalquantumnumbers,
203 const String& localquantumnumbers,
204 const String& normalization_option,
205 const String& mirroring_option,
206 const String& population_option,
207 const String& lineshapetype_option,
208 const String& cutoff_option,
209 const Numeric& cutoff_value,
210 const Numeric& linemixinglimit_value,
211 const Verbosity& verbosity)
212{
213 abs_lines.resize(0);
214
215 // Global numbers
216 const Array<QuantumNumberType> global_nums = string2vecqn(globalquantumnumbers);
217
218 // Local numbers
219 const Array<QuantumNumberType> local_nums = string2vecqn(localquantumnumbers);
220 ARTS_USER_ERROR_IF(not check_local(local_nums), "Can only have non-string values in the local state")
221
223
224 ArtsXMLTag tag(verbosity);
225 Index nelem;
226
227 // ARTSCAT data
228 shared_ptr<istream> ifs = nullptr;
229 xml_find_and_open_input_file(ifs, artscat_file, verbosity);
230 istream& is_xml = *ifs;
231 auto a = FILE_TYPE_ASCII;
232 auto b = NUMERIC_TYPE_DOUBLE;
233 auto c = ENDIAN_TYPE_LITTLE;
235 verbosity);
236
237 tag.read_from_stream(is_xml);
238 tag.check_name("Array");
239
240 Index num_arrays;
241 tag.get_attribute_value("nelem", num_arrays);
242
243 ArrayOfAbsorptionLines local_bands(0);
244 for (Index i=0; i<num_arrays; i++) {
245 tag.read_from_stream(is_xml);
246 tag.check_name("ArrayOfLineRecord");
247
248 tag.get_attribute_value("nelem", nelem);
249
250 String version;
251 tag.get_attribute_value("version", version);
252
253 Index artscat_version;
254
255 if (version == "3") {
256 artscat_version = 3;
257 } else if (version.substr(0, 8) != "ARTSCAT-") {
259 "The ARTS line file you are trying to read does not contain a valid version tag.\n"
260 "Probably it was created with an older version of ARTS that used different units.")
261 } else {
262 istringstream is(version.substr(8));
263 is >> artscat_version;
264 }
265
266 ARTS_USER_ERROR_IF (artscat_version < 3 or artscat_version > 5,
267 "Unknown ARTS line file version: ", version)
268
269 bool go_on = true;
270 Index n = 0;
271 while (n<nelem) {
272 n++;
273
274 if (go_on) {
276 switch(artscat_version) {
277 case 3:
279 break;
280 case 4:
282 break;
283 case 5:
285 break;
286 default:
287 ARTS_ASSERT (false, "Bad version!");
288 }
289
290 ARTS_USER_ERROR_IF(sline.bad, "Cannot read line ", n)
291 if (sline.line.F0 < fmin) continue;
292 if (sline.line.F0 > fmax) {go_on = false; continue;}
293
295
296 // Get the global quantum number identifier
297 const QuantumIdentifier global_qid = global_quantumidentifier(global_nums, sline.quantumidentity);
298
299 // Get local quantum numbers into the line
300 for(auto qn: local_nums) {
301 if (sline.quantumidentity.val.has(qn)) {
302 sline.line.localquanta.val.set(sline.quantumidentity.val[qn]);
303 }
304 }
305
306 merge_external_line(local_bands, sline, global_qid);
307 if (local_bands.nelem() > merge_local_lines_size) {
308 merge_local_lines(abs_lines, local_bands);
309 local_bands.resize(0);
310 }
311 } else {
312 String line;
313 getline(is_xml, line);
314 }
315 }
316
317 tag.read_from_stream(is_xml);
318 tag.check_name("/ArrayOfLineRecord");
319 }
320
321 merge_local_lines(abs_lines, local_bands);
322
323 abs_linesNormalization(abs_lines, normalization_option, verbosity);
324 abs_linesMirroring(abs_lines, mirroring_option, verbosity);
325 abs_linesPopulation(abs_lines, population_option, verbosity);
326 abs_linesLineShapeType(abs_lines, lineshapetype_option, verbosity);
327 abs_linesCutoff(abs_lines, cutoff_option, cutoff_value, verbosity);
328 abs_linesLinemixingLimit(abs_lines, linemixinglimit_value, verbosity);
329
330 tag.read_from_stream(is_xml);
331 tag.check_name("/Array");
332}
333
334/* Workspace method: Doxygen documentation will be auto-generated */
336 const String& artscat_file,
337 const Numeric& fmin,
338 const Numeric& fmax,
339 const String& globalquantumnumbers,
340 const String& localquantumnumbers,
341 const String& normalization_option,
342 const String& mirroring_option,
343 const String& population_option,
344 const String& lineshapetype_option,
345 const String& cutoff_option,
346 const Numeric& cutoff_value,
347 const Numeric& linemixinglimit_value,
348 const Verbosity& verbosity)
349{
350 abs_lines.resize(0);
351
352 // Global numbers
353 const Array<QuantumNumberType> global_nums = string2vecqn(globalquantumnumbers);
354
355 // Local numbers
356 const Array<QuantumNumberType> local_nums = string2vecqn(localquantumnumbers);
357 ARTS_USER_ERROR_IF(not check_local(local_nums), "Can only have non-string values in the local state")
358
360
361 ArtsXMLTag tag(verbosity);
362 Index nelem;
363
364 // ARTSCAT data
365 shared_ptr<istream> ifs = nullptr;
366 xml_find_and_open_input_file(ifs, artscat_file, verbosity);
367 istream& is_xml = *ifs;
368 auto a = FILE_TYPE_ASCII;
369 auto b = NUMERIC_TYPE_DOUBLE;
370 auto c = ENDIAN_TYPE_LITTLE;
372 verbosity);
373
374 tag.read_from_stream(is_xml);
375 tag.check_name("ArrayOfLineRecord");
376
377 tag.get_attribute_value("nelem", nelem);
378
379 String version;
380 tag.get_attribute_value("version", version);
381
382 Index artscat_version;
383
384 if (version == "3") {
385 artscat_version = 3;
386 } else if (version.substr(0, 8) != "ARTSCAT-") {
388 "The ARTS line file you are trying to read does not contain a valid version tag.\n"
389 "Probably it was created with an older version of ARTS that used different units.")
390 } else {
391 istringstream is(version.substr(8));
392 is >> artscat_version;
393 }
394
395 ARTS_USER_ERROR_IF (artscat_version < 3 or artscat_version > 5,
396 "Unknown ARTS line file version: ", version)
397
398 bool go_on = true;
399 Index n = 0;
400 ArrayOfAbsorptionLines local_bands(0);
401 while (n<nelem) {
402 n++;
403
404 if (go_on) {
406 switch(artscat_version) {
407 case 3:
409 break;
410 case 4:
412 break;
413 case 5:
415 break;
416 default:
417 ARTS_ASSERT (false, "Bad version!");
418 }
419
420 ARTS_USER_ERROR_IF(sline.bad, "Cannot read line ", n)
421 if (sline.line.F0 < fmin) continue;
422 if (sline.line.F0 > fmax) {go_on = false; continue;}
423
425
426 // Get the global quantum number identifier
427 const QuantumIdentifier global_qid = global_quantumidentifier(global_nums, sline.quantumidentity);
428
429 // Get local quantum numbers into the line
430 for(auto qn: local_nums) {
431 if (sline.quantumidentity.val.has(qn)) {
432 sline.line.localquanta.val.set(sline.quantumidentity.val[qn]);
433 }
434 }
435
436 merge_external_line(local_bands, sline, global_qid);
437 if (local_bands.nelem() > merge_local_lines_size) {
438 merge_local_lines(abs_lines, local_bands);
439 local_bands.resize(0);
440 }
441 } else {
442 String line;
443 getline(is_xml, line);
444 }
445 }
446
447 merge_local_lines(abs_lines, local_bands);
448
449 abs_linesNormalization(abs_lines, normalization_option, verbosity);
450 abs_linesMirroring(abs_lines, mirroring_option, verbosity);
451 abs_linesPopulation(abs_lines, population_option, verbosity);
452 abs_linesLineShapeType(abs_lines, lineshapetype_option, verbosity);
453 abs_linesCutoff(abs_lines, cutoff_option, cutoff_value, verbosity);
454 abs_linesLinemixingLimit(abs_lines, linemixinglimit_value, verbosity);
455
456 tag.read_from_stream(is_xml);
457 tag.check_name("/ArrayOfLineRecord");
458}
459
460/* Workspace method: Doxygen documentation will be auto-generated */
462 const ArrayOfArrayOfSpeciesTag& abs_species,
463 const String& basename,
464 const Numeric& fmin,
465 const Numeric& fmax,
466 const String& globalquantumnumbers,
467 const String& localquantumnumbers,
468 const Index& ignore_missing,
469 const String& normalization_option,
470 const String& mirroring_option,
471 const String& population_option,
472 const String& lineshapetype_option,
473 const String& cutoff_option,
474 const Numeric& cutoff_value,
475 const Numeric& linemixinglimit_value,
476 const Verbosity& verbosity)
477{
478 abs_lines.resize(0);
479
480 // Build a set of species indices. Duplicates are ignored.
481 const std::set<Species::Species> unique_species = lbl_species(abs_species);
482
483 String tmpbasename = basename;
484 if (basename.length() && basename[basename.length() - 1] != '/') {
485 tmpbasename += '.';
486 }
487
488 // Read catalogs for each identified species and put them all into
489 // abs_lines.
490 abs_lines.resize(0);
491 for (auto& spec: unique_species) {
492 ArrayOfAbsorptionLines more_abs_lines;
493
494 try {
495 ReadARTSCAT(more_abs_lines,
496 tmpbasename + String(Species::toShortName(spec)) + ".xml",
497 fmin,
498 fmax,
499 globalquantumnumbers,
500 localquantumnumbers,
501 normalization_option,
502 mirroring_option,
503 population_option,
504 lineshapetype_option,
505 cutoff_option,
506 cutoff_value,
507 linemixinglimit_value,
508 verbosity);
509
510 // Either find a line like this in the list of lines or start a new Lines
511 for (auto& newband: more_abs_lines) {
512 bool found = false;
513 for (auto& band: abs_lines) {
514 if (band.Match(newband).first) {
515 for (Index k=0; k<newband.NumLines(); k++) {
516 band.AppendSingleLine(newband.lines[k]);
517 found = true;
518 }
519 }
520 }
521 if (not found) {
522 abs_lines.push_back(newband);
523 }
524 }
525 } catch (const std::exception& e) {
526 ARTS_USER_ERROR_IF (not ignore_missing,
527 "Errors in calls by *propmat_clearskyAddZeeman*:\n",
528 e.what())
529 continue;
530 }
531 }
532
533 for (auto& band: abs_lines)
534 band.sort_by_frequency();
535}
536
537/* Workspace method: Doxygen documentation will be auto-generated */
539 const String& hitran_file,
540 const Numeric& fmin,
541 const Numeric& fmax,
542 const String& globalquantumnumbers,
543 const String& localquantumnumbers,
544 const String& hitran_type,
545 const String& normalization_option,
546 const String& mirroring_option,
547 const String& population_option,
548 const String& lineshapetype_option,
549 const String& cutoff_option,
550 const Numeric& cutoff_value,
551 const Numeric& linemixinglimit_value,
552 const Verbosity& verbosity)
553{
554 abs_lines.resize(0);
555
556 // Global numbers
557 const Array<QuantumNumberType> global_nums = string2vecqn(globalquantumnumbers);
558
559 // Local numbers
560 const Array<QuantumNumberType> local_nums = string2vecqn(localquantumnumbers);
561 ARTS_USER_ERROR_IF(not check_local(local_nums), "Can only have non-string values in the local state")
562
563 // HITRAN type
564 const Options::HitranType hitran_version = Options::toHitranTypeOrThrow(hitran_type);
565
566 // Hitran data
567 ifstream is;
568 open_input_file(is, hitran_file);
569
570 ArrayOfAbsorptionLines local_bands(0);
571 bool go_on = true;
572 while (go_on) {
573
575 switch (hitran_version) {
576 case Options::HitranType::Post2004:
578 break;
579 case Options::HitranType::Pre2004:
581 break;
582 case Options::HitranType::Online:
584 break;
585 default:
586 ARTS_ASSERT(false, "The HitranType enum class has to be fully updated!\n");
587 }
588
589 if (sline.bad) {
590 if (is.eof())
591 break;
592 ARTS_USER_ERROR("Cannot read line ", nelem(abs_lines) + nelem(local_bands) + 1);
593 }
594 if (sline.line.F0 < fmin)
595 continue; // Skip this line
596 if (sline.line.F0 > fmax)
597 break; // We assume sorted so quit here
598
599 // Set Zeeman if implemented
601
602 // Get the global quantum number identifier
603 const QuantumIdentifier global_qid = global_quantumidentifier(global_nums, sline.quantumidentity);
604
605 // Get local quantum numbers into the line
606 for(auto qn: local_nums) {
607 if (sline.quantumidentity.val.has(qn)) {
608 sline.line.localquanta.val.set(sline.quantumidentity.val[qn]);
609 }
610 }
611
612 // Either find a line like this in the list of lines or start a new Lines
613 merge_external_line(local_bands, sline, global_qid);
614 if (local_bands.nelem() > merge_local_lines_size) {
615 merge_local_lines(abs_lines, local_bands);
616 local_bands.resize(0);
617 }
618 }
619
620 merge_local_lines(abs_lines, local_bands);
621
622 abs_linesNormalization(abs_lines, normalization_option, verbosity);
623 abs_linesMirroring(abs_lines, mirroring_option, verbosity);
624 abs_linesPopulation(abs_lines, population_option, verbosity);
625 abs_linesLineShapeType(abs_lines, lineshapetype_option, verbosity);
626 abs_linesCutoff(abs_lines, cutoff_option, cutoff_value, verbosity);
627 abs_linesLinemixingLimit(abs_lines, linemixinglimit_value, verbosity);
628}
629
630/* Workspace method: Doxygen documentation will be auto-generated */
632 const String& lblrtm_file,
633 const Numeric& fmin,
634 const Numeric& fmax,
635 const String& globalquantumnumbers,
636 const String& localquantumnumbers,
637 const String& normalization_option,
638 const String& mirroring_option,
639 const String& population_option,
640 const String& lineshapetype_option,
641 const String& cutoff_option,
642 const Numeric& cutoff_value,
643 const Numeric& linemixinglimit_value,
644 const Verbosity& verbosity)
645{
646 abs_lines.resize(0);
647
648 // Global numbers
649 const Array<QuantumNumberType> global_nums = string2vecqn(globalquantumnumbers);
650
651 // Local numbers
652 const Array<QuantumNumberType> local_nums = string2vecqn(localquantumnumbers);
653 ARTS_USER_ERROR_IF(not check_local(local_nums), "Can only have non-string values in the local state")
654
655 // LBLRTM data
656 ifstream is;
657 open_input_file(is, lblrtm_file);
658
659 std::vector<Absorption::SingleLineExternal> v(0);
660
661 bool go_on = true;
662 while (go_on) {
663 v.push_back(Absorption::ReadFromLBLRTMStream(is));
664
665 if (v.back().bad) {
666 v.pop_back();
667 go_on = false;
668 } else if (v.back().line.F0 < fmin) {
669 v.pop_back();
670 } else if (v.back().line.F0 > fmax) {
671 v.pop_back();
672 go_on = false;
673 }
674 }
675
676 for (auto& x: v)
677 x.line.zeeman = Zeeman::GetAdvancedModel(x.quantumidentity);
678
679 auto x = Absorption::split_list_of_external_lines(v, local_nums, global_nums);
680 abs_lines.resize(0);
681 abs_lines.reserve(x.size());
682 while (x.size()) {
683 abs_lines.push_back(x.back());
684 abs_lines.back().sort_by_frequency();
685 x.pop_back();
686 }
687
688 abs_linesNormalization(abs_lines, normalization_option, verbosity);
689 abs_linesMirroring(abs_lines, mirroring_option, verbosity);
690 abs_linesPopulation(abs_lines, population_option, verbosity);
691 abs_linesLineShapeType(abs_lines, lineshapetype_option, verbosity);
692 abs_linesCutoff(abs_lines, cutoff_option, cutoff_value, verbosity);
693 abs_linesLinemixingLimit(abs_lines, linemixinglimit_value, verbosity);
694}
695
696/* Workspace method: Doxygen documentation will be auto-generated */
698 const String& jpl_file,
699 const Numeric& fmin,
700 const Numeric& fmax,
701 const String& globalquantumnumbers,
702 const String& localquantumnumbers,
703 const String& normalization_option,
704 const String& mirroring_option,
705 const String& population_option,
706 const String& lineshapetype_option,
707 const String& cutoff_option,
708 const Numeric& cutoff_value,
709 const Numeric& linemixinglimit_value,
710 const Verbosity& verbosity)
711{
712 abs_lines.resize(0);
713
714 // Global numbers
715 const Array<QuantumNumberType> global_nums = string2vecqn(globalquantumnumbers);
716
717 // Local numbers
718 const Array<QuantumNumberType> local_nums = string2vecqn(localquantumnumbers);
719 ARTS_USER_ERROR_IF(not check_local(local_nums), "Can only have non-string values in the local state")
720
721 // LBLRTM data
722 ifstream is;
723 open_input_file(is, jpl_file);
724
725 std::vector<Absorption::SingleLineExternal> v(0);
726
727 bool go_on = true;
728 while (go_on) {
729 v.push_back(Absorption::ReadFromJplStream(is));
730
731 if (v.back().bad) {
732 v.pop_back();
733 go_on = false;
734 } else if (v.back().line.F0 < fmin) {
735 v.pop_back();
736 } else if (v.back().line.F0 > fmax) {
737 v.pop_back();
738 go_on = false;
739 }
740 }
741
742 for (auto& x: v)
743 x.line.zeeman = Zeeman::GetAdvancedModel(x.quantumidentity);
744
745 auto x = Absorption::split_list_of_external_lines(v, local_nums, global_nums);
746 abs_lines.resize(0);
747 abs_lines.reserve(x.size());
748 while (x.size()) {
749 abs_lines.push_back(x.back());
750 abs_lines.back().sort_by_frequency();
751 x.pop_back();
752 }
753
754 abs_linesNormalization(abs_lines, normalization_option, verbosity);
755 abs_linesMirroring(abs_lines, mirroring_option, verbosity);
756 abs_linesPopulation(abs_lines, population_option, verbosity);
757 abs_linesLineShapeType(abs_lines, lineshapetype_option, verbosity);
758 abs_linesCutoff(abs_lines, cutoff_option, cutoff_value, verbosity);
759 abs_linesLinemixingLimit(abs_lines, linemixinglimit_value, verbosity);
760}
761
765
766/* Workspace method: Doxygen documentation will be auto-generated */
768 const ArrayOfAbsorptionLines& abs_lines,
769 const String& basename,
770 const Verbosity& verbosity) {
771 // Set the true name of the saving
772 String true_basename = basename;
773 if (not(true_basename.back() == '.' or true_basename.back() == '/'))
774 true_basename += '.';
775
776 // Find all species
777 ArrayOfString specs(0);
778 for (auto& band: abs_lines) {
779 auto specname = band.SpeciesName();
780
781 bool any = false;
782 for (auto& thisname: specs) {
783 if (any) break;
784 if (thisname == specname) any = true;
785 }
786
787 if (not any)
788 specs.push_back(specname);
789 }
790
791 // Make all species into a species tag array
792 Index throwaway;
794 abs_speciesSet(as, throwaway, specs, verbosity);
795
796 // Split lines by species
798 abs_lines_per_speciesCreateFromLines(alps, abs_lines, as, verbosity);
799
800 // Save the arrays
801 for (Index i=0; i<specs.nelem(); i++) {
802 auto& name = specs[i];
803 auto& lines = alps[i];
804
805 WriteXML(output_format, lines,
806 true_basename + name + ".xml",
807 0, "", "", "", verbosity);
808 }
809}
810
811/* Workspace method: Doxygen documentation will be auto-generated */
813 const ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
814 const String& basename,
815 const Verbosity& verbosity) {
816 // Compact to abs_lines
817 ArrayOfAbsorptionLines abs_lines(0);
818 for (auto& lines: abs_lines_per_species) {
819 for (auto& band: lines) {
820 abs_lines.push_back(band);
821 }
822 }
823
824 // Save using the other function
825 abs_linesWriteSpeciesSplitCatalog(output_format, abs_lines, basename, verbosity);
826}
827
828/* Workspace method: Doxygen documentation will be auto-generated */
830 const String& basename,
831 const Index& robust,
832 const Verbosity& verbosity) {
833 abs_lines.resize(0);
834
836 std::size_t bands_found{0};
837
838 String tmpbasename = basename;
839 if (basename.length() && basename[basename.length() - 1] != '/') {
840 tmpbasename += '.';
841 }
842
843 // Read catalogs for each identified species and put them all into
844 // abs_lines
845 for (auto& ir: Species::Isotopologues) {
846 String filename = tmpbasename + ir.FullName() + ".xml";
847 if (find_xml_file_existence(filename)) {
848 ArrayOfAbsorptionLines speclines;
849 xml_read_from_file(filename, speclines, verbosity);
850 for (auto& band: speclines) {
851 abs_lines.push_back(band);
852 bands_found++;
853 }
854 }
855 }
856
857 ARTS_USER_ERROR_IF (not bands_found and not robust,
858 "Cannot find any bands in the directory you are reading");
859 out3 << "Found " << bands_found << " bands\n";
860}
861
862/* Workspace method: Doxygen documentation will be auto-generated */
864 const ArrayOfArrayOfSpeciesTag& abs_species,
865 const String& basename,
866 const Index& robust,
867 const Verbosity& verbosity) {
868 abs_lines_per_species.resize(0);
869
871
872 // Build a set of species indices. Duplicates are ignored.
873 const std::set<Species::Species> unique_species = lbl_species(abs_species);
874
875 String tmpbasename = basename;
876 if (basename.length() && basename[basename.length() - 1] != '/') {
877 tmpbasename += '.';
878 }
879
880 // Read catalogs for each identified species and put them all into
881 // abs_lines in aa parallel manner
883 ArrayOfArrayOfAbsorptionLines par_abs_lines(n);
884 ArrayOfString errors;
885 std::atomic_bool error{false};
886
887 #pragma omp parallel for schedule(dynamic) if (!arts_omp_in_parallel())
888 for (std::size_t ispec=0; ispec<unique_species.size(); ispec++) {
889
890 const auto i = arts_omp_get_thread_num();
891 const auto& spec = *std::next(unique_species.cbegin(), ispec);
892 const auto isots = Species::isotopologues(spec);
893
894 if (not error.load()) {
895 try {
896 for (const auto& isot: isots) {
897 String filename = tmpbasename + isot.FullName() + ".xml";
898 if (find_xml_file_existence(filename)) {
899 ArrayOfAbsorptionLines speclines;
900 xml_read_from_file(filename, speclines, verbosity);
901 for (auto& band: speclines) {
902 par_abs_lines[i].push_back(std::move(band));
903 }
904 }
905 }
906 } catch (std::exception& e) {
907 error.store(true);
908
909 #pragma omp critical
910 errors.push_back(var_string("\n\nError in thread ", i, " for species ", spec, ":\n", e.what(), "\n\n"));
911 }
912 }
913 }
914
915 ARTS_USER_ERROR_IF(error.load(), "Encountered ", errors.nelem(), ':', errors)
916
917 const Index bands_found = nelem(par_abs_lines);
918 ARTS_USER_ERROR_IF (not bands_found and not robust,
919 "Cannot find any bands in the directory you are reading");
920 out3 << "Found " << bands_found << " bands\n";
921
922 ArrayOfAbsorptionLines abs_lines;
923 abs_lines.reserve(bands_found);
924 for (auto& lines: par_abs_lines) {
925 for (auto& band: lines) {
926 abs_lines.push_back(std::move(band));
927 }
928 }
929
930 abs_lines_per_speciesCreateFromLines(abs_lines_per_species, abs_lines, abs_species, verbosity);
931}
932
936
937/* Workspace method: Doxygen documentation will be auto-generated */
939{
940 for (auto& rlines: replacing_lines) {
941 Index number_of_matching_bands = 0;
942 for (auto& tlines: abs_lines) {
943 if (tlines.Match(rlines).first) {
944 number_of_matching_bands++;
945 for (auto& rline: rlines.lines) {
946 Index number_of_matching_single_lines = 0;
947 for (auto& tline: tlines.lines) {
948 if (tline.localquanta.val.check_match(rline.localquanta.val)) {
949 number_of_matching_single_lines++;
950 tline = rline;
951 }
952 }
953
954 ARTS_USER_ERROR_IF (number_of_matching_single_lines not_eq 1,
955 "Error: Did not match to a single single line. "
956 "This means the input data has not been understood. "
957 "This function needs exactly one match.");
958 }
959 tlines.sort_by_frequency();
960 }
961 }
962
963 ARTS_USER_ERROR_IF (number_of_matching_bands not_eq 1,
964 "Error: Did not match to a single set of absorption lines. "
965 "This means the input data has not been understood. "
966 "This function needs exactly one match.");
967 }
968}
969
970/* Workspace method: Doxygen documentation will be auto-generated */
971void abs_linesAppendWithLines(ArrayOfAbsorptionLines& abs_lines, const ArrayOfAbsorptionLines& appending_lines, const Index& safe, const Verbosity&)
972{
973 if (safe) {
974 std::vector<AbsorptionLines> addedlines(0);
975
976 for (auto& alines: appending_lines) {
977 Index number_of_matching_bands = 0;
978 for (auto& tlines: abs_lines) {
979 if (tlines.Match(alines).first) {
980 number_of_matching_bands++;
981 for (auto& aline: alines.lines) {
982 Index number_of_matching_single_lines = 0;
983 for (auto& tline: tlines.lines) {
984 if (tline.localquanta.val.check_match(aline.localquanta.val)) {
985 number_of_matching_single_lines++;
986 }
987 }
988 ARTS_USER_ERROR_IF (number_of_matching_single_lines not_eq 0,
989 "Error: Did match to a single single line. "
990 "This means the input data has not been understood. "
991 "This function needs exactly zero matches.");
992 tlines.AppendSingleLine(aline);
993 }
994 tlines.sort_by_frequency();
995 }
996 }
997
998 if (number_of_matching_bands == 0) addedlines.push_back(alines);
999 ARTS_USER_ERROR_IF (number_of_matching_bands not_eq 1,
1000 "Error: Did not match to a single set of absorption lines. "
1001 "This means the input data has not been understood. "
1002 "This function needs exactly one or zero matches.");
1003 }
1004
1005 for (auto& lines: addedlines) {
1006 abs_lines.push_back(std::move(lines));
1007 }
1008 } else {
1009 for (auto& band: appending_lines)
1010 abs_lines.push_back(band);
1011 }
1012}
1013
1014/* Workspace method: Doxygen documentation will be auto-generated */
1015void abs_linesDeleteBadF0(ArrayOfAbsorptionLines& abs_lines, const Numeric& f0, const Index& lower, const Verbosity&)
1016{
1017 for (auto& lines: abs_lines) {
1018 std::vector<Index> hits;
1019 for (Index i=0; i<lines.NumLines(); i++) {
1020 if (lower and lines.lines[i].F0 < f0)
1021 hits.push_back(i);
1022 else if (not lower and lines.lines[i].F0 > f0)
1023 hits.push_back(i);
1024 }
1025
1026 // Remove the bad values (sort by descending firs)
1027 std::sort(hits.begin(), hits.end());
1028 while(not hits.empty()) {
1029 lines.RemoveLine(hits.back());
1030 hits.pop_back();
1031 }
1032 }
1033}
1034
1035/* Workspace method: Doxygen documentation will be auto-generated */
1037{
1038 for (auto& dlines: deleting_lines) {
1039 for (auto& tlines: abs_lines) {
1040 std::vector<Index> hits(0);
1041
1042 if (tlines.Match(dlines).first) {
1043 for (auto& dline: dlines.lines) {
1044 for (Index i=0; i<tlines.NumLines(); i++) {
1045 if (tlines.lines[i].localquanta.val.check_match(dline.localquanta.val)) {
1046 hits.push_back(i);
1047 }
1048 }
1049 }
1050
1051 // Sort and test the input
1052 std::sort(hits.begin(), hits.end());
1053 auto n = hits.size();
1054 hits.erase(std::unique(hits.begin(), hits.end()), hits.end());
1055 ARTS_USER_ERROR_IF(n not_eq hits.size(),
1056 "Removing the same line more than once is not accepted");
1057
1058 // Remove the bad values
1059 while(not hits.empty()) {
1060 tlines.RemoveLine(hits.back());
1061 hits.pop_back();
1062 }
1063 }
1064 }
1065 }
1066}
1067
1068/* Workspace method: Doxygen documentation will be auto-generated */
1070{
1071 for (auto& band: abs_lines) {
1072 std::array<bool, LineShape::nVars> var_is_empty;
1073
1074 // Species by species can be empty, so loop each species by themselves
1075 for (Index ispec=0; ispec<band.NumBroadeners(); ispec++) {
1076 var_is_empty.fill(true);
1077
1078 // Check if any variable in this band for any line is non-empty
1079 for (Index iline=0; iline<band.NumLines(); iline++) {
1080 for (Index ivar=0; ivar < LineShape::nVars; ivar++) {
1081 if (not LineShape::modelparameterEmpty(band.lines[iline].lineshape.Data()[ispec].Data()[ivar])) {
1082 var_is_empty[ivar] = false;
1083 }
1084 }
1085 }
1086
1087 // Remove empty variables from the writing. This will also speed up some calculations
1088 for (Index iline=0; iline<band.NumLines(); iline++) {
1089 for (Index ivar=0; ivar < LineShape::nVars; ivar++) {
1090 if (var_is_empty[ivar]) {
1091 band.lines[iline].lineshape.Data()[ispec].Data()[ivar].type = LineShape::TemperatureModel::None;
1092 }
1093 }
1094 }
1095 }
1096 }
1097}
1098
1100{
1101 for (auto& band: abs_lines) {
1102 const Quantum::Number::StateMatch lt(qid, band.quantumidentity);
1103 while (lt not_eq Quantum::Number::StateMatchType::Full and band.NumLines()) {
1104 band.RemoveLine(0);
1105 }
1106 }
1107}
1108
1110 const Verbosity&) {
1111 const Index nb = lines.nelem();
1112 for (Index i=0; i<nb; i++) {
1113 for (Index j=i+1; j<nb; j++) {
1115 Quantum::Number::StateMatch(lines[i].quantumidentity, lines[j].quantumidentity) ==
1116 Quantum::Number::StateMatchType::Full,
1117 "Not unique, these bands match:\n", lines[i], "\nand\n", lines[j])
1118 }
1119 }
1120}
1121
1123 const ArrayOfAbsorptionLines& replacing_bands,
1124 const Verbosity&) {
1125 for (auto& replacement: replacing_bands) {
1126 struct {Index band;} pos{-1};
1127
1128 for (Index i=0; i<abs_lines.nelem(); i++) {
1129 auto& band = abs_lines[i];
1130
1131 if (Quantum::Number::StateMatch(band.quantumidentity, replacement.quantumidentity) ==
1132 Quantum::Number::StateMatchType::Full) {
1133 ARTS_USER_ERROR_IF (pos.band not_eq -1, "Duplicate band matches for replacement line:\n",
1134 replacement, "\nThese are for band indexes ", pos.band, " and ", i)
1135
1136 pos.band = i;
1137 }
1138 }
1139
1140 ARTS_USER_ERROR_IF(pos.band == -1,
1141 "There is no match for replacement band:\n", replacement,
1142 "\nYou need to append the entire band")
1143 abs_lines[pos.band] = replacement;
1144 }
1145}
1146
1148 const ArrayOfAbsorptionLines& replacing_lines,
1149 const Verbosity&) {
1150 const Index nb = abs_lines.nelem(); // Evaluate first so new bands are not counted
1151
1152 for (auto& replacement: replacing_lines) {
1153 const Index nl=replacement.NumLines();
1154
1155 struct {Index band; ArrayOfIndex lines;} pos{-1, ArrayOfIndex(nl, -1)};
1156
1157 for (Index i=0; i<nb; i++) {
1158 auto& band = abs_lines[i];
1159
1160 if (Quantum::Number::StateMatch(band.quantumidentity, replacement.quantumidentity) ==
1161 Quantum::Number::StateMatchType::Full) {
1162 ARTS_USER_ERROR_IF (pos.band not_eq -1, "Duplicate band matches for replacement line:\n",
1163 replacement, "\nThese are for band indexes ", pos.band, " and ", i, '\n',
1164 "ID ", i, ": ", band.quantumidentity, '\n',
1165 "ID ", pos.band, ": ", abs_lines[pos.band].quantumidentity, '\n',
1166 "ID target: ", replacement.quantumidentity, '\n')
1167
1168 pos.band = i;
1169
1170 for (Index k=0; k<band.NumLines(); k++) {
1171 for (Index j=0; j<nl; j++) {
1172 if (replacement.lines[j].localquanta.val.check_match(band.lines[k].localquanta.val)) {
1173 // We cannot have multiple entries
1174 ARTS_USER_ERROR_IF(pos.lines[j] not_eq -1 or std::any_of(pos.lines.begin(), pos.lines.end(), [k](auto& a){return a == k;}),
1175 "Found multiple matches of lines in:\n", replacement, "\n\nin mathcing band:\n", band)
1176
1177 pos.lines[j] = k;
1178 }
1179 }
1180 }
1181 }
1182 }
1183
1184 ARTS_USER_ERROR_IF(pos.band == -1 or std::any_of(pos.lines.begin(), pos.lines.end(), [](auto& a){return a == -1;}),
1185 "There is no match for replacement line:\n", replacement,
1186 "\nYou need to append the entire band")
1187
1188 // Add or change the line catalog
1189 auto& band = abs_lines[pos.band];
1190 if (const auto [match, nullable] = band.Match(replacement); nullable) {
1191 for (Index j=0; j<nl; j++) band.lines[pos.lines[j]] = replacement.lines[j];
1192 band.MakeLineShapeModelCommon();
1193 } else {
1194 // Sort to remove from behind
1195 std::sort(pos.lines.begin(), pos.lines.end(), std::greater<>());
1196 for (auto& k: pos.lines) band.RemoveLine(k);
1197 abs_lines.push_back(replacement);
1198 }
1199 }
1200}
1201
1205
1209
1210/* Workspace method: Doxygen documentation will be auto-generated */
1212 const String& type,
1213 const Numeric& x,
1214 const Verbosity&) {
1215 auto t = Absorption::toCutoffTypeOrThrow(type);
1216 for (auto& lines : abs_lines) {
1217 lines.cutoff = t;
1218 lines.cutofffreq = x;
1219 }
1220}
1221
1222/* Workspace method: Doxygen documentation will be auto-generated */
1224 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1225 const String& type,
1226 const Numeric& x,
1227 const Verbosity& v) {
1228 for (auto& abs_lines : abs_lines_per_species)
1229 abs_linesCutoff(abs_lines, type, x, v);
1230}
1231
1232/* Workspace method: Doxygen documentation will be auto-generated */
1234 const String& type,
1235 const Numeric& x,
1236 const QuantumIdentifier& QI,
1237 const Verbosity&) {
1238 auto t = Absorption::toCutoffTypeOrThrow(type);
1239 for (auto& band : abs_lines) {
1240 const Quantum::Number::StateMatch lt(QI, band.quantumidentity);
1241 if (lt == Quantum::Number::StateMatchType::Full) {
1242 band.cutoff = t;
1243 band.cutofffreq = x;
1244 }
1245 }
1246}
1247
1248/* Workspace method: Doxygen documentation will be auto-generated */
1250 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1251 const String& type,
1252 const Numeric& x,
1253 const QuantumIdentifier& QI,
1254 const Verbosity& verbosity) {
1255 for (auto& lines : abs_lines_per_species) {
1256 abs_linesCutoffMatch(lines, type, x, QI, verbosity);
1257 }
1258}
1259
1260/* Workspace method: Doxygen documentation will be auto-generated */
1262 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1263 const ArrayOfArrayOfSpeciesTag& abs_species,
1264 const String& type,
1265 const Numeric& x,
1266 const String& species_tag,
1267 const Verbosity& verbosity) {
1268 Index t1;
1269 ArrayOfArrayOfSpeciesTag target_species;
1270 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1271 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1272 if (std::equal(abs_species[ispec].begin(),
1273 abs_species[ispec].end(),
1274 target_species[0].begin())) {
1275 abs_linesCutoff(abs_lines_per_species[ispec], type, x, verbosity);
1276 }
1277 }
1278}
1279
1283
1284/* Workspace method: Doxygen documentation will be auto-generated */
1286 const String& type,
1287 const Verbosity&) {
1288 auto t = Absorption::toMirroringTypeOrThrow(type);
1289 for (auto& lines : abs_lines) lines.mirroring = t;
1290}
1291
1292/* Workspace method: Doxygen documentation will be auto-generated */
1294 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1295 const String& type,
1296 const Verbosity& v) {
1297 for (auto& abs_lines : abs_lines_per_species)
1298 abs_linesMirroring(abs_lines, type, v);
1299}
1300
1301/* Workspace method: Doxygen documentation will be auto-generated */
1303 const String& type,
1304 const QuantumIdentifier& QI,
1305 const Verbosity&) {
1306 auto t = Absorption::toMirroringTypeOrThrow(type);
1307 for (auto& band : abs_lines) {
1308 const Quantum::Number::StateMatch lt(QI, band.quantumidentity);
1309 if (lt == Quantum::Number::StateMatchType::Full) {
1310 band.mirroring = t;
1311 }
1312 }
1313}
1314
1315/* Workspace method: Doxygen documentation will be auto-generated */
1317 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1318 const String& type,
1319 const QuantumIdentifier& QI,
1320 const Verbosity& v) {
1321 for (auto& abs_lines : abs_lines_per_species)
1322 abs_linesMirroringMatch(abs_lines, type, QI, v);
1323}
1324
1325/* Workspace method: Doxygen documentation will be auto-generated */
1327 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1328 const ArrayOfArrayOfSpeciesTag& abs_species,
1329 const String& type,
1330 const String& species_tag,
1331 const Verbosity& verbosity) {
1332 Index t1;
1333 ArrayOfArrayOfSpeciesTag target_species;
1334 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1335 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1336 if (std::equal(abs_species[ispec].begin(),
1337 abs_species[ispec].end(),
1338 target_species[0].begin())) {
1339 abs_linesMirroring(abs_lines_per_species[ispec], type, verbosity);
1340 }
1341 }
1342}
1343
1344/* Workspace method: Doxygen documentation will be auto-generated */
1346 const Verbosity&) {
1347 const ArrayOfAbsorptionLines abs_lines_copy = abs_lines;
1348 for (AbsorptionLines band : abs_lines_copy) {
1349 band.mirroring = Absorption::MirroringType::Manual;
1350
1353 std::find_if(abs_lines_copy.cbegin(),
1354 abs_lines_copy.cend(),
1355 [&band](const AbsorptionLines& li) {
1356 return band.Match(li).first;
1357 }) not_eq abs_lines_copy.cend(),
1358 "Dual bands with same setup is not allowed for mirroring of band:\n",
1359 band,
1360 '\n')
1361
1362 for (auto& line : band.lines) {
1363 line.F0 *= -1;
1364 }
1365
1366 abs_lines.emplace_back(std::move(band));
1367 }
1368}
1369
1370/* Workspace method: Doxygen documentation will be auto-generated */
1372 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1373 const Verbosity& verbosity) {
1374 for (auto& abs_lines : abs_lines_per_species)
1375 abs_linesManualMirroring(abs_lines, verbosity);
1376}
1377
1378/* Workspace method: Doxygen documentation will be auto-generated */
1380 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1381 const ArrayOfArrayOfSpeciesTag& abs_species,
1382 const ArrayOfSpeciesTag& species,
1383 const Verbosity& verbosity) {
1384 ARTS_USER_ERROR_IF(abs_species.size() not_eq abs_lines_per_species.size(),
1385 "Mismatch abs_species and abs_lines_per_species sizes [",
1386 abs_species.size(),
1387 " vs ",
1388 abs_lines_per_species.size(),
1389 ", respectively]")
1390
1391 if (auto ind = std::distance(
1392 abs_species.cbegin(),
1393 std::find(abs_species.cbegin(), abs_species.cend(), species));
1394 ind not_eq abs_species.nelem()) {
1395 abs_linesManualMirroring(abs_lines_per_species[ind], verbosity);
1396 } else {
1397 ARTS_USER_ERROR("Cannot find species: ",
1398 species,
1399 "\nIn abs_species: [",
1400 abs_species,
1401 ']')
1402 }
1403}
1404
1408
1409/* Workspace method: Doxygen documentation will be auto-generated */
1411 const String& type,
1412 const Verbosity&) {
1413 auto t = Absorption::toPopulationTypeOrThrow(type);
1414 for (auto& lines : abs_lines) lines.population = t;
1415}
1416
1417/* Workspace method: Doxygen documentation will be auto-generated */
1419 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1420 const String& type,
1421 const Verbosity& v) {
1422 for (auto& abs_lines : abs_lines_per_species)
1423 abs_linesPopulation(abs_lines, type, v);
1424}
1425
1426/* Workspace method: Doxygen documentation will be auto-generated */
1428 const String& type,
1429 const QuantumIdentifier& QI,
1430 const Verbosity&) {
1431 auto t = Absorption::toPopulationTypeOrThrow(type);
1432 for (auto& lines : abs_lines) {
1433 const Quantum::Number::StateMatch lt(QI, lines.quantumidentity);
1434 if (lt == Quantum::Number::StateMatchType::Full) {
1435 lines.population = t;
1436 }
1437 }
1438}
1439
1440/* Workspace method: Doxygen documentation will be auto-generated */
1442 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1443 const String& type,
1444 const QuantumIdentifier& QI,
1445 const Verbosity& v) {
1446 for (auto& abs_lines : abs_lines_per_species)
1447 abs_linesPopulationMatch(abs_lines, type, QI, v);
1448}
1449
1450/* Workspace method: Doxygen documentation will be auto-generated */
1452 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1453 const ArrayOfArrayOfSpeciesTag& abs_species,
1454 const String& type,
1455 const String& species_tag,
1456 const Verbosity& verbosity) {
1457 Index t1;
1458 ArrayOfArrayOfSpeciesTag target_species;
1459 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1460 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1461 if (std::equal(abs_species[ispec].begin(),
1462 abs_species[ispec].end(),
1463 target_species[0].begin())) {
1464 abs_linesPopulation(abs_lines_per_species[ispec], type, verbosity);
1465 }
1466 }
1467}
1468
1472
1473/* Workspace method: Doxygen documentation will be auto-generated */
1475 const String& type,
1476 const Verbosity&) {
1477 auto t = Absorption::toNormalizationTypeOrThrow(type);
1478 for (auto& lines : abs_lines) lines.normalization = t;
1479}
1480
1481/* Workspace method: Doxygen documentation will be auto-generated */
1483 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1484 const String& type,
1485 const Verbosity& v) {
1486 for (auto& abs_lines : abs_lines_per_species)
1487 abs_linesNormalization(abs_lines, type, v);
1488}
1489
1490/* Workspace method: Doxygen documentation will be auto-generated */
1492 const String& type,
1493 const QuantumIdentifier& QI,
1494 const Verbosity&) {
1495 auto t = Absorption::toNormalizationTypeOrThrow(type);
1496 for (auto& lines : abs_lines) {
1497 const Quantum::Number::StateMatch lt(QI, lines.quantumidentity);
1498 if (lt == Quantum::Number::StateMatchType::Full) {
1499 lines.normalization = t;
1500 }
1501 }
1502}
1503
1504/* Workspace method: Doxygen documentation will be auto-generated */
1506 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1507 const String& type,
1508 const QuantumIdentifier& QI,
1509 const Verbosity& v) {
1510 for (auto& abs_lines : abs_lines_per_species)
1511 abs_linesNormalizationMatch(abs_lines, type, QI, v);
1512}
1513
1514/* Workspace method: Doxygen documentation will be auto-generated */
1516 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1517 const ArrayOfArrayOfSpeciesTag& abs_species,
1518 const String& type,
1519 const String& species_tag,
1520 const Verbosity& verbosity) {
1521 Index t1;
1522 ArrayOfArrayOfSpeciesTag target_species;
1523 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1524 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1525 if (std::equal(abs_species[ispec].begin(),
1526 abs_species[ispec].end(),
1527 target_species[0].begin())) {
1528 abs_linesNormalization(abs_lines_per_species[ispec], type, verbosity);
1529 }
1530 }
1531}
1532
1536
1537/* Workspace method: Doxygen documentation will be auto-generated */
1539 const String& type,
1540 const Verbosity&) {
1541 auto t = LineShape::toType(type);
1542 check_enum_error(t, "Cannot understand type: ", type);
1543 for (auto& lines : abs_lines) lines.lineshapetype = t;
1544}
1545
1546/* Workspace method: Doxygen documentation will be auto-generated */
1548 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1549 const String& type,
1550 const Verbosity& v) {
1551 for (auto& abs_lines : abs_lines_per_species)
1552 abs_linesLineShapeType(abs_lines, type, v);
1553}
1554
1555/* Workspace method: Doxygen documentation will be auto-generated */
1557 const String& type,
1558 const QuantumIdentifier& QI,
1559 const Verbosity&) {
1560 auto t = LineShape::toTypeOrThrow(type);
1561 for (auto& lines : abs_lines) {
1562 const Quantum::Number::StateMatch lt(QI, lines.quantumidentity);
1563 if (lt == Quantum::Number::StateMatchType::Full) {
1564 lines.lineshapetype = t;
1565 }
1566 }
1567}
1568
1569/* Workspace method: Doxygen documentation will be auto-generated */
1571 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1572 const String& type,
1573 const QuantumIdentifier& QI,
1574 const Verbosity& v) {
1575 for (auto& abs_lines : abs_lines_per_species)
1576 abs_linesLineShapeTypeMatch(abs_lines, type, QI, v);
1577}
1578
1579/* Workspace method: Doxygen documentation will be auto-generated */
1581 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1582 const ArrayOfArrayOfSpeciesTag& abs_species,
1583 const String& type,
1584 const String& species_tag,
1585 const Verbosity& verbosity) {
1586 Index t1;
1587 ArrayOfArrayOfSpeciesTag target_species;
1588 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1589 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1590 if (std::equal(abs_species[ispec].begin(),
1591 abs_species[ispec].end(),
1592 target_species[0].begin())) {
1593 abs_linesLineShapeType(abs_lines_per_species[ispec], type, verbosity);
1594 }
1595 }
1596}
1597
1601
1602/* Workspace method: Doxygen documentation will be auto-generated */
1604 const Numeric& x,
1605 const Verbosity&) {
1606 for (auto& lines : abs_lines) lines.linemixinglimit = x;
1607}
1608
1609/* Workspace method: Doxygen documentation will be auto-generated */
1611 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1612 const Numeric& x,
1613 const Verbosity& v) {
1614 for (auto& abs_lines : abs_lines_per_species)
1615 abs_linesLinemixingLimit(abs_lines, x, v);
1616}
1617
1618/* Workspace method: Doxygen documentation will be auto-generated */
1620 const Numeric& x,
1621 const QuantumIdentifier& QI,
1622 const Verbosity&) {
1623 for (auto& lines : abs_lines) {
1624 const Quantum::Number::StateMatch lt(QI, lines.quantumidentity);
1625 if (lt == Quantum::Number::StateMatchType::Full) {
1626 lines.linemixinglimit = x;
1627 }
1628 }
1629}
1630
1631/* Workspace method: Doxygen documentation will be auto-generated */
1633 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1634 const Numeric& x,
1635 const QuantumIdentifier& QI,
1636 const Verbosity& v) {
1637 for (auto& abs_lines : abs_lines_per_species)
1638 abs_linesLinemixingLimitMatch(abs_lines, x, QI, v);
1639}
1640
1641/* Workspace method: Doxygen documentation will be auto-generated */
1643 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1644 const ArrayOfArrayOfSpeciesTag& abs_species,
1645 const Numeric& x,
1646 const String& species_tag,
1647 const Verbosity& verbosity) {
1648 Index t1;
1649 ArrayOfArrayOfSpeciesTag target_species;
1650 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1651 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1652 if (std::equal(abs_species[ispec].begin(),
1653 abs_species[ispec].end(),
1654 target_species[0].begin())) {
1655 abs_linesLinemixingLimit(abs_lines_per_species[ispec], x, verbosity);
1656 }
1657 }
1658}
1659
1663
1664/* Workspace method: Doxygen documentation will be auto-generated */
1666 const Numeric& x,
1667 const Verbosity&) {
1668 for (auto& lines : abs_lines) lines.T0 = x;
1669}
1670
1671/* Workspace method: Doxygen documentation will be auto-generated */
1673 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1674 const Numeric& x,
1675 const Verbosity& v) {
1676 for (auto& abs_lines : abs_lines_per_species) abs_linesT0(abs_lines, x, v);
1677}
1678
1679/* Workspace method: Doxygen documentation will be auto-generated */
1681 const Numeric& x,
1682 const QuantumIdentifier& QI,
1683 const Verbosity&) {
1684 for (auto& lines : abs_lines) {
1685 const Quantum::Number::StateMatch lt(QI, lines.quantumidentity);
1686 if (lt == Quantum::Number::StateMatchType::Full) {
1687 lines.T0 = x;
1688 }
1689 }
1690}
1691
1692/* Workspace method: Doxygen documentation will be auto-generated */
1694 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1695 const Numeric& x,
1696 const QuantumIdentifier& QI,
1697 const Verbosity& v) {
1698 for (auto& abs_lines : abs_lines_per_species)
1699 abs_linesT0Match(abs_lines, x, QI, v);
1700}
1701
1702/* Workspace method: Doxygen documentation will be auto-generated */
1704 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1705 const ArrayOfArrayOfSpeciesTag& abs_species,
1706 const Numeric& x,
1707 const String& species_tag,
1708 const Verbosity& verbosity) {
1709 Index t1;
1710 ArrayOfArrayOfSpeciesTag target_species;
1711 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1712 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1713 if (std::equal(abs_species[ispec].begin(),
1714 abs_species[ispec].end(),
1715 target_species[0].begin())) {
1716 abs_linesT0(abs_lines_per_species[ispec], x, verbosity);
1717 }
1718 }
1719}
1720
1724
1725/* Workspace method: Doxygen documentation will be auto-generated */
1727 const QuantumIdentifier& QI,
1728 const String& parameter_name,
1729 const Numeric& change,
1730 const Index& relative,
1731 const Verbosity&)
1732{
1733 Index parameter_switch = -1;
1734
1735 ARTS_USER_ERROR_IF (parameter_name.nelem() == 0,
1736 "parameter_name is empty.\n");
1737
1738 if (parameter_name == "Central Frequency" or
1739 parameter_name == "Line Center")
1740 parameter_switch = 0;
1741 else if (parameter_name == "Line Strength")
1742 parameter_switch = 1;
1743 else if (parameter_name == "Lower State Energy")
1744 parameter_switch = 4;
1745 else if (parameter_name == "Einstein Coefficient")
1746 parameter_switch = 5;
1747 else if (parameter_name == "Lower Statistical Weight")
1748 parameter_switch = 6;
1749 else if (parameter_name == "Upper Statistical Weight")
1750 parameter_switch = 7;
1751 else if (parameter_name == "Lower Zeeman Coefficient")
1752 parameter_switch = 8;
1753 else if (parameter_name == "Upper Zeeman Coefficient")
1754 parameter_switch = 9;
1755
1756 for (auto& band: abs_lines) {
1757 for (Index k=0; k<band.NumLines(); k++) {
1758 const Quantum::Number::StateMatch lt(QI, band.lines[k].localquanta, band.quantumidentity);
1759 if (lt == Quantum::Number::StateMatchType::Full) {
1760 switch (parameter_switch) {
1761 case 0: // "Central Frequency":
1762 if (relative == 0)
1763 band.lines[k].F0 += change;
1764 else
1765 band.lines[k].F0 *= 1.0e0 + change;
1766 break;
1767 case 1: // "Line Strength":
1768 if (relative == 0)
1769 band.lines[k].I0 += change;
1770 else
1771 band.lines[k].I0 *= 1.0e0 + change;
1772 break;
1773 case 4: // "Lower State Energy":
1774 if (relative == 0)
1775 band.lines[k].E0 += change;
1776 else
1777 band.lines[k].E0 *= 1.0e0 + change;
1778 break;
1779 case 5: // "Einstein":
1780 if (relative == 0)
1781 band.lines[k].A += change;
1782 else
1783 band.lines[k].A *= 1.0e0 + change;
1784 break;
1785 case 6: // "Lower Statistical Weight":
1786 if (relative == 0)
1787 band.lines[k].glow += change;
1788 else
1789 band.lines[k].glow *= 1.0e0 + change;
1790 break;
1791 case 7: // "Upper Statistical Weight":
1792 if (relative == 0)
1793 band.lines[k].gupp += change;
1794 else
1795 band.lines[k].gupp *= 1.0e0 + change;
1796 break;
1797 case 8: // "Lower Zeeman Coefficient":
1798 if (relative == 0)
1799 band.lines[k].zeeman.gl() += change;
1800 else
1801 band.lines[k].zeeman.gl() *= 1.0e0 + change;
1802 break;
1803 case 9: // "Upper Zeeman Coefficient":
1804 if (relative == 0)
1805 band.lines[k].zeeman.gu() += change;
1806 else
1807 band.lines[k].zeeman.gu() *= 1.0e0 + change;
1808 break;
1809 default: {
1811 "Usupported paramter_name\n", parameter_name,
1812 "\nSee method description for supported parameter names.\n")
1813 }
1814 }
1815 }
1816 }
1817 }
1818}
1819
1820/* Workspace method: Doxygen documentation will be auto-generated */
1822 const QuantumIdentifier& QI,
1823 const String& parameter_name,
1824 const Numeric& change,
1825 const Index& relative,
1826 const Verbosity& verbosity)
1827{
1828 for (auto& lines: abs_lines_per_species)
1829 abs_linesChangeBaseParameterForMatchingLines(lines, QI, parameter_name, change, relative, verbosity);
1830}
1831
1832/* Workspace method: Doxygen documentation will be auto-generated */
1834 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1835 const ArrayOfArrayOfSpeciesTag& abs_species,
1836 const QuantumIdentifier& QI,
1837 const String& parameter_name,
1838 const Numeric& change,
1839 const Index& relative,
1840 const String& species_tag,
1841 const Verbosity& verbosity)
1842{
1843 Index t1;
1844 ArrayOfArrayOfSpeciesTag target_species;
1845 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1846 for (Index ispec=0; ispec<abs_species.nelem(); ispec++) {
1847 if (std::equal(abs_species[ispec].begin(), abs_species[ispec].end(), target_species[0].begin())) {
1848 abs_linesChangeBaseParameterForMatchingLines(abs_lines_per_species[ispec], QI, parameter_name, change, relative, verbosity);
1849 }
1850 }
1851}
1852
1853/* Workspace method: Doxygen documentation will be auto-generated */
1855 const QuantumIdentifier& QI,
1856 const String& parameter_name,
1857 const Numeric& x,
1858 const Verbosity&) {
1859 Index parameter_switch = -1;
1860
1861 ARTS_USER_ERROR_IF(parameter_name.nelem() == 0, "parameter_name is empty.\n");
1862
1863 if (parameter_name == "Central Frequency" or parameter_name == "Line Center")
1864 parameter_switch = 0;
1865 else if (parameter_name == "Line Strength")
1866 parameter_switch = 1;
1867 else if (parameter_name == "Lower State Energy")
1868 parameter_switch = 4;
1869 else if (parameter_name == "Einstein Coefficient")
1870 parameter_switch = 5;
1871 else if (parameter_name == "Lower Statistical Weight")
1872 parameter_switch = 6;
1873 else if (parameter_name == "Upper Statistical Weight")
1874 parameter_switch = 7;
1875 else if (parameter_name == "Lower Zeeman Coefficient")
1876 parameter_switch = 8;
1877 else if (parameter_name == "Upper Zeeman Coefficient")
1878 parameter_switch = 9;
1879
1880 for (auto& band : abs_lines) {
1881 for (Index k = 0; k < band.NumLines(); k++) {
1883 QI, band.lines[k].localquanta, band.quantumidentity);
1884 if (lt == Quantum::Number::StateMatchType::Full) {
1885 switch (parameter_switch) {
1886 case 0: // "Central Frequency":
1887 band.lines[k].F0 = x;
1888 break;
1889 case 1: // "Line Strength":
1890 band.lines[k].I0 = x;
1891 break;
1892 case 4: // "Lower State Energy":
1893 band.lines[k].E0 = x;
1894 break;
1895 case 5: // "Einstein":
1896 band.lines[k].A = x;
1897 break;
1898 case 6: // "Lower Statistical Weight":
1899 band.lines[k].glow = x;
1900 break;
1901 case 7: // "Upper Statistical Weight":
1902 band.lines[k].gupp = x;
1903 break;
1904 case 8:
1905 band.lines[k].zeeman.gl() = x;
1906 break;
1907 case 9:
1908 band.lines[k].zeeman.gu() = x;
1909 break;
1910 default: {
1912 "Usupported paramter_name\n",
1913 parameter_name,
1914 "\nSee method description for supported parameter names.\n")
1915 }
1916 }
1917 }
1918 }
1919 }
1920}
1921
1922/* Workspace method: Doxygen documentation will be auto-generated */
1924 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1925 const QuantumIdentifier& QI,
1926 const String& parameter_name,
1927 const Numeric& change,
1928 const Verbosity& verbosity) {
1929 for (auto& lines : abs_lines_per_species)
1931 lines, QI, parameter_name, change, verbosity);
1932}
1933
1934/* Workspace method: Doxygen documentation will be auto-generated */
1936 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
1937 const ArrayOfArrayOfSpeciesTag& abs_species,
1938 const QuantumIdentifier& QI,
1939 const String& parameter_name,
1940 const Numeric& change,
1941 const String& species_tag,
1942 const Verbosity& verbosity) {
1943 Index t1;
1944 ArrayOfArrayOfSpeciesTag target_species;
1945 abs_speciesSet(target_species, t1, {species_tag}, verbosity);
1946 for (Index ispec = 0; ispec < abs_species.nelem(); ispec++) {
1947 if (std::equal(abs_species[ispec].begin(),
1948 abs_species[ispec].end(),
1949 target_species[0].begin())) {
1951 abs_lines_per_species[ispec], QI, parameter_name, change, verbosity);
1952 }
1953 }
1954}
1955
1956/* Workspace method: Doxygen documentation will be auto-generated */
1958 ArrayOfAbsorptionLines& abs_lines,
1959 const QuantumIdentifier& QI,
1960 const String& parameter,
1961 const String& species,
1962 const String& temperaturemodel,
1963 const Vector& new_values,
1964 const Verbosity& verbosity) {
1966
1967 const bool do_self = species == LineShape::self_broadening;
1968 const bool do_bath = species == LineShape::bath_broadening;
1969
1970 // Set the spec index if possible
1971 const Species::Species spec = do_self ? Species::Species::FINAL
1972 : do_bath ? Species::Species::Bath
1973 : SpeciesTag(species).Spec();
1974
1975 const LineShape::Variable var = LineShape::toVariableOrThrow(parameter);
1976
1978 "Mismatch between input and expected number of variables\n"
1979 "\tInput is: ",
1980 new_values.nelem(),
1981 " long but expects: ",
1983 " values\n")
1984
1986 newdata.type = LineShape::toTemperatureModelOrThrow(temperaturemodel);
1987 newdata.X0 = new_values[0];
1988 newdata.X1 = new_values[1];
1989 newdata.X2 = new_values[2];
1990 newdata.X3 = new_values[3];
1991
1992 for (auto& band : abs_lines) {
1993 for (Index k = 0; k < band.NumLines(); k++) {
1995 QI, band.lines[k].localquanta, band.quantumidentity);
1996 if (lt == Quantum::Number::StateMatchType::Full) {
1997 if (do_self and band.selfbroadening) {
1998 out3 << "Changing self\n";
1999 band.lines[k].lineshape.Data().front().Data()[Index(var)] = newdata;
2000 } else if (do_bath and band.bathbroadening) {
2001 out3 << "Changing bath\n";
2002 band.lines[k].lineshape.Data().back().Data()[Index(var)] = newdata;
2003 } else {
2004 for (Index i = band.selfbroadening;
2005 i < band.broadeningspecies.nelem() - band.bathbroadening;
2006 i++) {
2007 if (spec == band.broadeningspecies[i]) {
2008 out3 << "Changing species: " << Species::toShortName(spec)
2009 << '\n';
2010 band.lines[k].lineshape.Data()[i].Data()[Index(var)] = newdata;
2011 }
2012 }
2013 }
2014 }
2015 }
2016 }
2017}
2018
2019/* Workspace method: Doxygen documentation will be auto-generated */
2021 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2022 const QuantumIdentifier& QI,
2023 const String& parameter,
2024 const String& species,
2025 const String& temperaturemodel,
2026 const Vector& new_values,
2027 const Verbosity& verbosity) {
2028 for (auto& lines : abs_lines_per_species)
2030 lines, QI, parameter, species, temperaturemodel, new_values, verbosity);
2031}
2032
2036
2037/* Workspace method: Doxygen documentation will be auto-generated */
2039 const QuantumIdentifier& QI,
2040 const String& parameter_name,
2041 const Numeric& change,
2042 const Index& relative,
2043 const Verbosity&)
2044{
2045 Index parameter_switch = -1;
2046
2047 ARTS_USER_ERROR_IF (parameter_name.nelem() == 0,
2048 "parameter_name is empty.\n");
2049 if (parameter_name == "Statistical Weight")
2050 parameter_switch = 1;
2051 else if (parameter_name == "Zeeman Coefficient")
2052 parameter_switch = 2;
2053
2054 for (auto& band: abs_lines) {
2055 for (Index k=0; k<band.NumLines(); k++) {
2056 const Quantum::Number::StateMatch lt(QI, band.lines[k].localquanta, band.quantumidentity);
2057 if (lt == Quantum::Number::StateMatchType::Level and lt.low) {
2058 switch (parameter_switch) {
2059 case 1: // "Statistical Weight":
2060 if (relative == 0)
2061 band.lines[k].glow += change;
2062 else
2063 band.lines[k].glow *= 1.0e0 + change;
2064 break;
2065 case 2: // "Zeeman Coefficient":
2066 if (relative == 0)
2067 band.lines[k].zeeman.gl() += change;
2068 else
2069 band.lines[k].zeeman.gl() *= 1.0e0 + change;
2070 break;
2071 default: {
2073 "Usupported paramter_name\n", parameter_name,
2074 "\nSee method description for supported parameter names.\n")
2075 }
2076 }
2077 }
2078
2079 if (lt == Quantum::Number::StateMatchType::Level and lt.upp) {
2080 switch (parameter_switch) {
2081 case 1: // "Statistical Weight":
2082 if (relative == 0)
2083 band.lines[k].gupp += change;
2084 else
2085 band.lines[k].gupp *= 1.0e0 + change;
2086 break;
2087 case 2: // "Zeeman Coefficient":
2088 if (relative == 0)
2089 band.lines[k].zeeman.gu() += change;
2090 else
2091 band.lines[k].zeeman.gu() *= 1.0e0 + change;
2092 break;
2093 default: {
2095 "Usupported paramter_name\n", parameter_name,
2096 "\nSee method description for supported parameter names.\n")
2097 }
2098 }
2099 }
2100 }
2101 }
2102}
2103
2104/* Workspace method: Doxygen documentation will be auto-generated */
2106 const QuantumIdentifier& QI,
2107 const String& parameter_name,
2108 const Numeric& change,
2109 const Index& relative,
2110 const Verbosity& verbosity)
2111{
2112 for (auto& lines: abs_lines_per_species)
2113 abs_linesChangeBaseParameterForMatchingLevel(lines, QI, parameter_name, change, relative, verbosity);
2114}
2115
2116/* Workspace method: Doxygen documentation will be auto-generated */
2118 const ArrayOfQuantumIdentifier& QID,
2119 const String& parameter_name,
2120 const Vector& change,
2121 const Index& relative,
2122 const Verbosity& verbosity)
2123{
2124 ARTS_USER_ERROR_IF (QID.nelem() not_eq change.nelem(),
2125 "Mismatch between QID and change input lengths not allowed");
2126
2127 for (Index iq=0; iq<QID.nelem(); iq++)
2128 abs_linesChangeBaseParameterForMatchingLevel(abs_lines, QID[iq], parameter_name, change[iq], relative, verbosity);
2129}
2130
2131/* Workspace method: Doxygen documentation will be auto-generated */
2133 const ArrayOfQuantumIdentifier& QID,
2134 const String& parameter_name,
2135 const Vector& change,
2136 const Index& relative,
2137 const Verbosity& verbosity)
2138{
2139 ARTS_USER_ERROR_IF (QID.nelem() not_eq change.nelem(),
2140 "Mismatch between QID and change input lengths not allowed");
2141
2142 for (Index iq=0; iq<QID.nelem(); iq++)
2143 for (auto& lines: abs_lines_per_species)
2144 abs_linesChangeBaseParameterForMatchingLevel(lines, QID[iq], parameter_name, change[iq], relative, verbosity);
2145}
2146
2147/* Workspace method: Doxygen documentation will be auto-generated */
2149 const QuantumIdentifier& QI,
2150 const String& parameter_name,
2151 const Numeric& x,
2152 const Verbosity&) {
2153 Index parameter_switch = -1;
2154
2155 ARTS_USER_ERROR_IF(parameter_name.nelem() == 0, "parameter_name is empty.\n");
2156 if (parameter_name == "Statistical Weight")
2157 parameter_switch = 1;
2158 else if (parameter_name == "Zeeman Coefficient")
2159 parameter_switch = 2;
2160
2161 for (auto& band : abs_lines) {
2162 for (Index k = 0; k < band.NumLines(); k++) {
2164 QI, band.lines[k].localquanta, band.quantumidentity);
2165 if (lt == Quantum::Number::StateMatchType::Level and lt.low) {
2166 switch (parameter_switch) {
2167 case 1: // "Statistical Weight":
2168 band.lines[k].glow = x;
2169 break;
2170 case 2: // "Zeeman Coefficient":
2171 band.lines[k].zeeman.gl() = x;
2172 break;
2173 default: {
2175 "Usupported paramter_name\n",
2176 parameter_name,
2177 "\nSee method description for supported parameter names.\n")
2178 }
2179 }
2180 }
2181
2182 if (lt == Quantum::Number::StateMatchType::Level and lt.upp) {
2183 switch (parameter_switch) {
2184 case 1: // "Statistical Weight":
2185 band.lines[k].gupp = x;
2186 break;
2187 case 2: // "Zeeman Coefficient":
2188 band.lines[k].zeeman.gu() = x;
2189 break;
2190 default: {
2192 "Usupported paramter_name\n",
2193 parameter_name,
2194 "\nSee method description for supported parameter names.\n")
2195 }
2196 }
2197 }
2198 }
2199 }
2200}
2201
2202/* Workspace method: Doxygen documentation will be auto-generated */
2204 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2205 const QuantumIdentifier& QI,
2206 const String& parameter_name,
2207 const Numeric& change,
2208 const Verbosity& verbosity) {
2209 for (auto& lines : abs_lines_per_species)
2211 lines, QI, parameter_name, change, verbosity);
2212}
2213
2214/* Workspace method: Doxygen documentation will be auto-generated */
2216 const ArrayOfQuantumIdentifier& QID,
2217 const String& parameter_name,
2218 const Vector& change,
2219 const Verbosity& verbosity) {
2221 QID.nelem() not_eq change.nelem(),
2222 "Mismatch between QID and change input lengths not allowed");
2223
2224 for (Index iq = 0; iq < QID.nelem(); iq++)
2226 abs_lines, QID[iq], parameter_name, change[iq], verbosity);
2227}
2228
2229/* Workspace method: Doxygen documentation will be auto-generated */
2231 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2232 const ArrayOfQuantumIdentifier& QID,
2233 const String& parameter_name,
2234 const Vector& change,
2235 const Verbosity& verbosity) {
2237 QID.nelem() not_eq change.nelem(),
2238 "Mismatch between QID and change input lengths not allowed");
2239
2240 for (Index iq = 0; iq < QID.nelem(); iq++)
2241 for (auto& lines : abs_lines_per_species)
2243 lines, QID[iq], parameter_name, change[iq], verbosity);
2244}
2245
2249
2250/* Workspace method: Doxygen documentation will be auto-generated */
2252 Index& nlte_do,
2253 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2254 const EnergyLevelMap& nlte_field,
2255 const Verbosity&) {
2256 nlte_field.ThrowIfNotOK();
2257
2258 if (nlte_field.value.empty()) {
2259 nlte_do = 0;
2260 return;
2261 }
2262 nlte_do = 1;
2263
2264 const Absorption::PopulationType poptyp =
2265 nlte_field.vib_energy.empty() ? Absorption::PopulationType::NLTE
2266 : Absorption::PopulationType::VibTemps;
2267
2268 for (auto& spec_lines : abs_lines_per_species) {
2269 for (auto& band : spec_lines) {
2270 Index low = 0, upp = 0;
2271 for (auto& id : nlte_field.levels) {
2272 for (auto& line : band.lines) {
2273 const auto lt =
2274 poptyp == Absorption::PopulationType::NLTE
2276 id, line.localquanta, band.quantumidentity)
2277 : Quantum::Number::StateMatch(id, band.quantumidentity);
2278 low += lt.low;
2279 upp += lt.upp;
2280 }
2281 }
2282
2284 not(low == 0 or low == band.NumLines()) or
2285 not(upp == 0 or upp == band.NumLines()),
2286 "The band with metadata:\n",
2287 band.MetaData(),
2288 "\n\nwith lines:\n",
2289 band.lines,
2290 "\n\nThe number of levels don't add upp correctly.\n There are ",
2291 band.NumLines(),
2292 " lines but there are ",
2293 low,
2294 " lower level matches and ",
2295 upp,
2296 " upper level matches.")
2297
2298 if (upp or low) band.population = poptyp;
2299 }
2300 }
2301}
2302
2303
2307
2308/* Workspace method: Doxygen documentation will be auto-generated */
2310 ArrayOfArrayOfAbsorptionLines & abs_lines_per_species,
2311 const ArrayOfArrayOfSpeciesTag& abs_species,
2312 const Verbosity&) {
2313 abs_lines_per_species = ArrayOfArrayOfAbsorptionLines(
2314 abs_species.nelem(), ArrayOfAbsorptionLines(0));
2315}
2316
2317/* Workspace method: Doxygen documentation will be auto-generated */
2319 const Vector& f_grid,
2320 const Verbosity&) {
2321 const Numeric fmax = max(f_grid);
2322 const Numeric fmin = min(f_grid);
2323
2324 for (auto& band : abs_lines) {
2325 for (Index k = band.NumLines() - 1; k >= 0; k--) {
2326 const Numeric fcut_upp = band.CutoffFreq(k);
2327 const Numeric fcut_low = band.CutoffFreqMinus(k);
2328
2329 if (fmax < fcut_low or fmin > fcut_upp) {
2330 band.RemoveLine(k);
2331 }
2332 }
2333 }
2334}
2335
2336/* Workspace method: Doxygen documentation will be auto-generated */
2338 ArrayOfArrayOfAbsorptionLines & abs_lines_per_species,
2339 const Vector& f_grid,
2340 const Verbosity& verbosity) {
2341 for (auto& lines : abs_lines_per_species) {
2342 abs_linesCompact(lines, f_grid, verbosity);
2343 }
2344}
2345
2346/* Workspace method: Doxygen documentation will be auto-generated */
2348 const QuantumIdentifier& qid,
2349 const Verbosity&) {
2350 for (Index i = 0; i < abs_lines.nelem(); i++) {
2351 const Quantum::Number::StateMatch lt(qid, abs_lines[i].quantumidentity);
2352 if (lt == Quantum::Number::StateMatchType::Full) {
2353 abs_lines.erase(abs_lines.begin() + i);
2354 break;
2355 }
2356 }
2357}
2358
2362
2363template <class T>
2364std::vector<T> linspace(
2365 T s, T e, typename std::vector<T>::size_type count) noexcept {
2366 std::vector<T> ls(count);
2367
2368 if (count == 0) return ls;
2369 if (count == 1) {
2370 ls.front() = (e + s) / 2;
2371 return ls;
2372 }
2373 const T step = (e - s) / T(count - 1);
2374 ls.front() = s;
2375 ls.back() = e;
2376 for (typename std::vector<T>::size_type i = 1; i < count - 1; ++i)
2377 ls[i] = s + step * T(i);
2378 return ls;
2379}
2380
2381/* Workspace method: Doxygen documentation will be auto-generated */
2383 Vector & f_grid,
2384 const ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2385 const Numeric& delta_f_low,
2386 const Numeric& delta_f_upp,
2387 const Index& num_freqs,
2388 const Verbosity&) {
2389 const Index n = nelem(abs_lines_per_species);
2390
2391 ARTS_USER_ERROR_IF(delta_f_low >= delta_f_upp,
2392 "The lower frequency delta has to be smaller "
2393 "than the upper frequency delta");
2394 ARTS_USER_ERROR_IF(num_freqs == 0, "Need more than zero frequency points");
2395 ARTS_USER_ERROR_IF(n < 1,
2396 "No lines found. Error? Use *VectorSet* "
2397 "to resize *f_grid*");
2398
2399 std::vector<Numeric> fout(0);
2400 for (auto& lines : abs_lines_per_species) {
2401 for (auto& band : lines) {
2402 for (Index k = 0; k < band.NumLines(); k++) {
2403 if (num_freqs > 1) {
2404 auto ftmp =
2405 linspace<Numeric>(band.lines[k].F0 + delta_f_low,
2406 band.lines[k].F0 + delta_f_upp,
2407 std::vector<Numeric>::size_type(num_freqs));
2408 for (auto& f : ftmp) {
2409 if (f > 0) fout.push_back(f);
2410 }
2411 } else {
2412 fout.push_back(band.lines[k].F0);
2413 }
2414 }
2415 }
2416 }
2417
2418 std::sort(fout.begin(), fout.end());
2419 fout.erase(std::unique(fout.begin(), fout.end()), fout.end());
2420 f_grid.resize(fout.size());
2421 for (Index i = 0; i < f_grid.nelem(); i++) f_grid[i] = fout[i];
2422}
2423
2427
2429 const ArrayOfSpeciesTag& species,
2430 const Numeric lower_frequency,
2431 const Numeric upper_frequency,
2432 const Numeric lower_intensity,
2433 const Index safe,
2434 const Index flip_flims,
2435 const Verbosity& verbosity) {
2436 const bool care_about_species = species.nelem();
2437
2438 for (auto& band : abs_lines) {
2439 if (care_about_species and
2440 species[0].Isotopologue() not_eq band.Isotopologue())
2441 continue;
2442
2443 auto& lines = band.lines;
2444
2445 if (not safe) {
2446 std::vector<std::size_t> rem;
2447 if (flip_flims) {
2448 for (std::size_t i=lines.size()-1; i<lines.size(); i--) {
2449 auto& line = lines[i];
2450 if ((line.F0 >= lower_frequency and
2451 line.F0 <= upper_frequency) or
2452 line.I0 < lower_intensity)
2453 rem.push_back(i);
2454 }
2455 } else {
2456 for (std::size_t i = lines.size() - 1; i < lines.size(); i--) {
2457 auto& line = lines[i];
2458 if (line.F0 < lower_frequency or line.F0 > upper_frequency or
2459 line.I0 < lower_intensity)
2460 rem.push_back(i);
2461 }
2462 }
2463 for (auto i : rem) band.RemoveLine(i);
2464 } else {
2465 ARTS_USER_ERROR_IF(flip_flims,
2466 "Not allowed to combine GINs flip_flims and safe.")
2467 const bool all_low = std::all_of(
2468 lines.begin(), lines.end(), [lower_frequency](auto& line) {
2469 return line.F0 < lower_frequency;
2470 });
2471 const bool all_upp = std::all_of(
2472 lines.begin(), lines.end(), [upper_frequency](auto& line) {
2473 return line.F0 > upper_frequency;
2474 });
2475 const bool low_int = std::all_of(
2476 lines.begin(), lines.end(), [lower_intensity](auto& line) {
2477 return line.I0 < lower_intensity;
2478 });
2479 if (all_low or all_upp or low_int) lines.resize(0);
2480 }
2481 }
2482
2483 // Removes empty bands
2484 abs_linesRemoveEmptyBands(abs_lines, verbosity);
2485}
2486
2488 const Numeric& lower_frequency,
2489 const Numeric& upper_frequency,
2490 const Numeric& lower_intensity,
2491 const Index& safe,
2492 const Index& flip_flims,
2493 const Verbosity& verbosity) {
2494 remove_impl(abs_lines,
2495 {},
2496 lower_frequency,
2497 upper_frequency,
2498 lower_intensity,
2499 safe,
2500 flip_flims,
2501 verbosity);
2502}
2503
2505 ArrayOfArrayOfAbsorptionLines & abs_lines_per_species,
2506 const Numeric& lower_frequency,
2507 const Numeric& upper_frequency,
2508 const Numeric& lower_intensity,
2509 const Index& safe,
2510 const Index& flip_flims,
2511 const Verbosity& verbosity) {
2512 for (auto& abs_lines : abs_lines_per_species)
2513 abs_linesRemoveLines(abs_lines,
2514 lower_frequency,
2515 upper_frequency,
2516 lower_intensity,
2517 safe,
2518 flip_flims,
2519 verbosity);
2520}
2521
2523 const ArrayOfSpeciesTag& species,
2524 const Numeric& lower_frequency,
2525 const Numeric& upper_frequency,
2526 const Numeric& lower_intensity,
2527 const Index& safe,
2528 const Index& flip_flims,
2529 const Verbosity& verbosity) {
2531 species.nelem() not_eq 1, "Must have a single species, got: ", species)
2532 ARTS_USER_ERROR_IF(species[0].Isotopologue().joker(),
2533 "Cannot give joker species, got: ",
2534 species)
2535
2536 remove_impl(abs_lines,
2537 species,
2538 lower_frequency,
2539 upper_frequency,
2540 lower_intensity,
2541 safe,
2542 flip_flims,
2543 verbosity);
2544}
2545
2547 ArrayOfArrayOfAbsorptionLines & abs_lines_per_species,
2548 const ArrayOfSpeciesTag& species,
2549 const Numeric& lower_frequency,
2550 const Numeric& upper_frequency,
2551 const Numeric& lower_intensity,
2552 const Index& safe,
2553 const Index& flip_flims,
2554 const Verbosity& verbosity) {
2555 for (auto& abs_lines : abs_lines_per_species)
2557 species,
2558 lower_frequency,
2559 upper_frequency,
2560 lower_intensity,
2561 safe,
2562 flip_flims,
2563 verbosity);
2564}
2565
2567 const String& option,
2568 const Verbosity& verbosity) {
2569 const auto opt = Options::toSortingOptionOrThrow(option);
2570
2571 abs_linesRemoveEmptyBands(abs_lines, verbosity);
2572
2573 switch (opt) {
2574 case Options::SortingOption::ByFrequency:
2575 for (auto& band : abs_lines) band.sort_by_frequency();
2576 std::sort(abs_lines.begin(), abs_lines.end(), [](auto& a, auto& b) {
2577 return a.lines[0].F0 <= b.lines[0].F0;
2578 });
2579 break;
2580 case Options::SortingOption::ByEinstein:
2581 for (auto& band : abs_lines) band.sort_by_einstein();
2582 std::sort(abs_lines.begin(), abs_lines.end(), [](auto& a, auto& b) {
2583 return a.lines[0].A <= b.lines[0].A;
2584 });
2585 break;
2586 case Options::SortingOption::FINAL: { /*leave last*/
2587 }
2588 }
2589}
2590
2592 const Verbosity&) {
2593 for (auto& band : abs_lines) {
2594 for (auto& line : band.lines) {
2595 for (auto& slsm : line.lineshape) {
2596 // Set the line mixing stuff to empty
2597 slsm.Y() = LineShape::ModelParameters();
2598 slsm.G() = LineShape::ModelParameters();
2599 slsm.DV() = LineShape::ModelParameters();
2600 }
2601 }
2602 }
2603}
2604
2606 ArrayOfArrayOfAbsorptionLines& abs_lines_per_species,
2607 const Verbosity& verbosity) {
2608 for (auto& abs_lines : abs_lines_per_species) {
2609 abs_linesTurnOffLineMixing(abs_lines, verbosity);
2610 }
2611}
Contains the absorption namespace.
Array< ArrayOfAbsorptionLines > ArrayOfArrayOfAbsorptionLines
This file contains the definition of Array.
Array< Index > ArrayOfIndex
An array of Index.
Definition: array.h:136
base max(const Array< base > &x)
Max function.
Definition: array.h:145
base min(const Array< base > &x)
Min function.
Definition: array.h:161
int arts_omp_get_max_threads()
Wrapper for omp_get_max_threads.
Definition: arts_omp.cc:46
int arts_omp_get_thread_num()
Wrapper for omp_get_thread_num.
Definition: arts_omp.cc:76
Header file for helper functions for OpenMP.
Stuff related to time in ARTS.
This can be used to make arrays out of anything.
Definition: array.h:48
Index nelem() const ARTS_NOEXCEPT
Definition: array.h:92
The ARTS XML tag class.
Definition: xml_io.h:45
void get_attribute_value(const String &aname, SpeciesTag &value)
Returns value of attribute as type SpeciesTag.
Definition: xml_io.cc:70
bool empty() const noexcept
Definition: matpackIV.h:152
Index nelem() const noexcept
Returns the number of elements.
Definition: matpackI.h:547
bool empty() const noexcept
Returns true if variable size is zero.
Definition: matpackI.h:536
void set(Value v)
Sets the value if it exists or adds it otherwise.
bool has(Types... ts) const ARTS_NOEXCEPT
Returns whether all the Types are part of the list, the types must be sorted.
The Vector class.
Definition: matpackI.h:910
void resize(Index n)
Resize function.
Definition: matpackI.cc:390
void check_name(const String &expected_name)
Check tag name.
Definition: xml_io_base.cc:54
void read_from_stream(istream &is)
Reads next XML tag.
Definition: xml_io_base.cc:201
Index nelem() const
Definition: mystring.h:189
Helper macros for debugging.
#define ARTS_ASSERT(condition,...)
Definition: debug.h:102
#define ARTS_USER_ERROR(...)
Definition: debug.h:169
std::string var_string(Args &&... args)
Definition: debug.h:36
#define ARTS_USER_ERROR_IF(condition,...)
Definition: debug.h:153
constexpr void check_enum_error(EnumType type, Messages... args)
Checks if the enum class type is good and otherwise throws an error message composed by variadic inpu...
Definition: enums.h:90
bool find_xml_file_existence(String &filename)
As find_xml_file but does not throw in the main body.
Definition: file.cc:394
void open_input_file(ifstream &file, const std::string_view name)
Open a file for reading.
Definition: file.cc:128
This file contains basic functions to handle ASCII files.
Contains the line shape namespace.
void abs_lines_per_speciesCreateFromLines(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfAbsorptionLines &abs_lines, const ArrayOfArrayOfSpeciesTag &abs_species, const Verbosity &)
WORKSPACE METHOD: abs_lines_per_speciesCreateFromLines.
Definition: m_abs.cc:106
void abs_speciesSet(ArrayOfArrayOfSpeciesTag &abs_species, Index &propmat_clearsky_agenda_checked, const ArrayOfString &names, const Verbosity &verbosity)
WORKSPACE METHOD: abs_speciesSet.
QuantumIdentifier global_quantumidentifier(const Array< QuantumNumberType > &qns, const QuantumIdentifier &qid)
Selects the global quantum numbers.
void abs_linesCutoff(ArrayOfAbsorptionLines &abs_lines, const String &type, const Numeric &x, const Verbosity &)
WORKSPACE METHOD: abs_linesCutoff.
void abs_lines_per_speciesManualMirroringSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const ArrayOfSpeciesTag &species, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesManualMirroringSpecies.
void abs_lines_per_speciesBaseParameterMatchingLines(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Verbosity &verbosity)
void abs_lines_per_speciesChangeBaseParameterForMatchingLevel(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Index &relative, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesChangeBaseParameterForMatchingLevel.
void abs_lines_per_speciesPopulation(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesPopulation.
void abs_linesPopulationMatch(ArrayOfAbsorptionLines &abs_lines, const String &type, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesPopulationMatch.
void abs_linesDeleteBadF0(ArrayOfAbsorptionLines &abs_lines, const Numeric &f0, const Index &lower, const Verbosity &)
WORKSPACE METHOD: abs_linesDeleteBadF0.
void abs_linesNormalization(ArrayOfAbsorptionLines &abs_lines, const String &type, const Verbosity &)
WORKSPACE METHOD: abs_linesNormalization.
void abs_linesT0(ArrayOfAbsorptionLines &abs_lines, const Numeric &x, const Verbosity &)
WORKSPACE METHOD: abs_linesT0.
void abs_linesFlatten(ArrayOfAbsorptionLines &abs_lines, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesFlatten.
void abs_linesNormalizationMatch(ArrayOfAbsorptionLines &abs_lines, const String &type, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesNormalizationMatch.
void abs_lines_per_speciesLinemixingLimitSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const Numeric &x, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesLinemixingLimitSpecies.
bool check_local(const Array< QuantumNumberType > &local_state)
void abs_lines_per_speciesChangeBaseParameterForMatchingLines(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Index &relative, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesChangeBaseParameterForMatchingLines.
void abs_lines_per_speciesTurnOffLineMixing(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesTurnOffLineMixing.
void abs_linesKeepBand(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &qid, const Verbosity &)
WORKSPACE METHOD: abs_linesKeepBand.
void abs_lines_per_speciesLineShapeModelParametersMatchingLines(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const QuantumIdentifier &QI, const String &parameter, const String &species, const String &temperaturemodel, const Vector &new_values, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesLineShapeModelParametersMatchingLines.
void abs_lines_per_speciesLineShapeTypeMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesLineShapeTypeMatch.
void abs_linesMirroringMatch(ArrayOfAbsorptionLines &abs_lines, const String &type, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesMirroringMatch.
void ReadARTSCAT(ArrayOfAbsorptionLines &abs_lines, const String &artscat_file, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadARTSCAT.
void abs_lines_per_speciesChangeBaseParameterForSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Index &relative, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesChangeBaseParameterForSpecies.
void abs_lines_per_speciesPopulationNlteField(Index &nlte_do, ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const EnergyLevelMap &nlte_field, const Verbosity &)
WORKSPACE METHOD: abs_lines_per_speciesPopulationNlteField.
void abs_linesTurnOffLineMixing(ArrayOfAbsorptionLines &abs_lines, const Verbosity &)
WORKSPACE METHOD: abs_linesTurnOffLineMixing.
void abs_lines_per_speciesCutoff(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Numeric &x, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesCutoff.
void abs_lines_per_speciesLinemixingLimit(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &x, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesLinemixingLimit.
void abs_lines_per_speciesLinemixingLimitMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesLinemixingLimitMatch.
void abs_linesChangeBaseParameterForMatchingLines(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Index &relative, const Verbosity &)
WORKSPACE METHOD: abs_linesChangeBaseParameterForMatchingLines.
void abs_linesEmptyBroadeningParameters(ArrayOfAbsorptionLines &abs_lines, const Verbosity &)
WORKSPACE METHOD: abs_linesEmptyBroadeningParameters.
void abs_lines_per_speciesT0(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &x, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesT0.
void abs_lines_per_speciesMirroringMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesMirroringMatch.
void abs_lines_per_speciesReadSpeciesSplitCatalog(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &basename, const Index &robust, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesReadSpeciesSplitCatalog.
void ReadJPL(ArrayOfAbsorptionLines &abs_lines, const String &jpl_file, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadJPL.
void abs_lines_per_speciesManualMirroring(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesManualMirroring.
void abs_linesReplaceLines(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &replacing_lines, const Verbosity &)
WORKSPACE METHOD: abs_linesReplaceLines.
void abs_linesLineShapeModelParametersMatchingLines(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &QI, const String &parameter, const String &species, const String &temperaturemodel, const Vector &new_values, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesLineShapeModelParametersMatchingLines.
void abs_linesCompact(ArrayOfAbsorptionLines &abs_lines, const Vector &f_grid, const Verbosity &)
WORKSPACE METHOD: abs_linesCompact.
void remove_impl(ArrayOfAbsorptionLines &abs_lines, const ArrayOfSpeciesTag &species, const Numeric lower_frequency, const Numeric upper_frequency, const Numeric lower_intensity, const Index safe, const Index flip_flims, const Verbosity &verbosity)
void abs_linesRemoveBand(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &qid, const Verbosity &)
WORKSPACE METHOD: abs_linesRemoveBand.
constexpr Index merge_local_lines_size
void abs_linesCutoffMatch(ArrayOfAbsorptionLines &abs_lines, const String &type, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesCutoffMatch.
void abs_linesChangeBaseParameterForMatchingLevels(ArrayOfAbsorptionLines &abs_lines, const ArrayOfQuantumIdentifier &QID, const String &parameter_name, const Vector &change, const Index &relative, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesChangeBaseParameterForMatchingLevels.
void abs_linesBaseParameterMatchingLevels(ArrayOfAbsorptionLines &abs_lines, const ArrayOfQuantumIdentifier &QID, const String &parameter_name, const Vector &change, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesBaseParameterMatchingLevels.
void abs_lines_per_speciesRemoveLinesFromSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfSpeciesTag &species, const Numeric &lower_frequency, const Numeric &upper_frequency, const Numeric &lower_intensity, const Index &safe, const Index &flip_flims, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesRemoveLinesFromSpecies.
void abs_linesRemoveEmptyBands(ArrayOfAbsorptionLines &abs_lines, const Verbosity &)
WORKSPACE METHOD: abs_linesRemoveEmptyBands.
void abs_lines_per_speciesNormalizationSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &type, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesNormalizationSpecies.
void ReadArrayOfARTSCAT(ArrayOfAbsorptionLines &abs_lines, const String &artscat_file, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadArrayOfARTSCAT.
void abs_linesDeleteWithLines(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &deleting_lines, const Verbosity &)
void ReadSplitARTSCAT(ArrayOfAbsorptionLines &abs_lines, const ArrayOfArrayOfSpeciesTag &abs_species, const String &basename, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const Index &ignore_missing, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadSplitARTSCAT.
void abs_linesPopulation(ArrayOfAbsorptionLines &abs_lines, const String &type, const Verbosity &)
WORKSPACE METHOD: abs_linesPopulation.
void abs_lines_per_speciesT0Species(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const Numeric &x, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesT0Species.
void abs_lines_per_speciesNormalization(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesNormalization.
void abs_lines_per_speciesRemoveLines(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &lower_frequency, const Numeric &upper_frequency, const Numeric &lower_intensity, const Index &safe, const Index &flip_flims, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesRemoveLines.
void abs_linesAppendWithLines(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &appending_lines, const Index &safe, const Verbosity &)
void abs_lines_per_speciesMirroring(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesMirroring.
void abs_lines_per_speciesCutoffSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &type, const Numeric &x, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesCutoffSpecies.
void abs_linesLineShapeTypeMatch(ArrayOfAbsorptionLines &abs_lines, const String &type, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesLineShapeTypeMatch.
Array< QuantumNumberType > string2vecqn(std::string_view qnstr)
Get a list of quantum numbers from a string.
void abs_linesChangeBaseParameterForMatchingLevel(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Index &relative, const Verbosity &)
WORKSPACE METHOD: abs_linesChangeBaseParameterForMatchingLevel.
void abs_lines_per_speciesNormalizationMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesNormalizationMatch.
void abs_lines_per_speciesPopulationMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesPopulationMatch.
void ReadLBLRTM(ArrayOfAbsorptionLines &abs_lines, const String &lblrtm_file, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadLBLRTM.
void abs_linesBaseParameterMatchingLevel(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &x, const Verbosity &)
WORKSPACE METHOD: abs_linesBaseParameterMatchingLevel.
void abs_linesT0Match(ArrayOfAbsorptionLines &abs_lines, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesT0Match.
void abs_linesRemoveLines(ArrayOfAbsorptionLines &abs_lines, const Numeric &lower_frequency, const Numeric &upper_frequency, const Numeric &lower_intensity, const Index &safe, const Index &flip_flims, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesRemoveLines.
void abs_lines_per_speciesPopulationSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &type, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesPopulationSpecies.
void abs_linesManualMirroring(ArrayOfAbsorptionLines &abs_lines, const Verbosity &)
WORKSPACE METHOD: abs_linesManualMirroring.
void ReadHITRAN(ArrayOfAbsorptionLines &abs_lines, const String &hitran_file, const Numeric &fmin, const Numeric &fmax, const String &globalquantumnumbers, const String &localquantumnumbers, const String &hitran_type, const String &normalization_option, const String &mirroring_option, const String &population_option, const String &lineshapetype_option, const String &cutoff_option, const Numeric &cutoff_value, const Numeric &linemixinglimit_value, const Verbosity &verbosity)
WORKSPACE METHOD: ReadHITRAN.
void abs_lines_per_speciesT0Match(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesT0Match.
void abs_lines_per_speciesSetEmpty(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const Verbosity &)
WORKSPACE METHOD: abs_lines_per_speciesSetEmpty.
void abs_lines_per_speciesBaseParameterMatchingLevel(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &change, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesBaseParameterMatchingLevel.
void merge_external_line(ArrayOfAbsorptionLines &abs_lines, const Absorption::SingleLineExternal &sline, const QuantumIdentifier &global_qid)
Merge an external line to abs_lines.
void abs_lines_per_speciesBaseParameterMatchingLevels(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfQuantumIdentifier &QID, const String &parameter_name, const Vector &change, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesBaseParameterMatchingLevels.
void abs_linesReplaceWithLines(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &replacing_lines, const Verbosity &)
void abs_lines_per_speciesMirroringSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &type, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesMirroringSpecies.
void abs_linesLinemixingLimit(ArrayOfAbsorptionLines &abs_lines, const Numeric &x, const Verbosity &)
WORKSPACE METHOD: abs_linesLinemixingLimit.
void abs_lines_per_speciesCutoffMatch(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesCutoffMatch.
void abs_linesRemoveLinesFromSpecies(ArrayOfAbsorptionLines &abs_lines, const ArrayOfSpeciesTag &species, const Numeric &lower_frequency, const Numeric &upper_frequency, const Numeric &lower_intensity, const Index &safe, const Index &flip_flims, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesRemoveLinesFromSpecies.
void abs_lines_per_speciesFlatten(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesFlatten.
void abs_linesBaseParameterMatchingLines(ArrayOfAbsorptionLines &abs_lines, const QuantumIdentifier &QI, const String &parameter_name, const Numeric &x, const Verbosity &)
WORKSPACE METHOD: abs_linesBaseParameterMatchingLines.
void abs_lines_per_speciesWriteSpeciesSplitCatalog(const String &output_format, const ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &basename, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesWriteSpeciesSplitCatalog.
void abs_linesLineShapeType(ArrayOfAbsorptionLines &abs_lines, const String &type, const Verbosity &)
WORKSPACE METHOD: abs_linesLineShapeType.
void abs_lines_per_speciesChangeBaseParameterForMatchingLevels(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfQuantumIdentifier &QID, const String &parameter_name, const Vector &change, const Index &relative, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesChangeBaseParameterForMatchingLevels.
void abs_linesReplaceBands(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &replacing_bands, const Verbosity &)
WORKSPACE METHOD: abs_linesReplaceBands.
void abs_lines_per_speciesLineShapeType(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const String &type, const Verbosity &v)
WORKSPACE METHOD: abs_lines_per_speciesLineShapeType.
void abs_linesMirroring(ArrayOfAbsorptionLines &abs_lines, const String &type, const Verbosity &)
WORKSPACE METHOD: abs_linesMirroring.
void f_gridFromAbsorptionLines(Vector &f_grid, const ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Numeric &delta_f_low, const Numeric &delta_f_upp, const Index &num_freqs, const Verbosity &)
WORKSPACE METHOD: f_gridFromAbsorptionLines.
void abs_linesLinemixingLimitMatch(ArrayOfAbsorptionLines &abs_lines, const Numeric &x, const QuantumIdentifier &QI, const Verbosity &)
WORKSPACE METHOD: abs_linesLinemixingLimitMatch.
void abs_linesSort(ArrayOfAbsorptionLines &abs_lines, const String &option, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesSort.
void merge_local_lines(ArrayOfAbsorptionLines &abs_lines, const ArrayOfAbsorptionLines &local_lines)
Merge lines to abs_lines.
void abs_linesWriteSpeciesSplitCatalog(const String &output_format, const ArrayOfAbsorptionLines &abs_lines, const String &basename, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesWriteSpeciesSplitCatalog.
void abs_lines_per_speciesLineShapeTypeSpecies(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const ArrayOfArrayOfSpeciesTag &abs_species, const String &type, const String &species_tag, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesLineShapeTypeSpecies.
std::vector< T > linspace(T s, T e, typename std::vector< T >::size_type count) noexcept
void CheckUnique(const ArrayOfAbsorptionLines &lines, const Verbosity &)
WORKSPACE METHOD: CheckUnique.
void abs_linesReadSpeciesSplitCatalog(ArrayOfAbsorptionLines &abs_lines, const String &basename, const Index &robust, const Verbosity &verbosity)
WORKSPACE METHOD: abs_linesReadSpeciesSplitCatalog.
void abs_lines_per_speciesCompact(ArrayOfArrayOfAbsorptionLines &abs_lines_per_species, const Vector &f_grid, const Verbosity &verbosity)
WORKSPACE METHOD: abs_lines_per_speciesCompact.
Workspace methods and template functions for supergeneric XML IO.
void WriteXML(const String &file_format, const T &v, const String &f, const Index &no_clobber, const String &v_name, const String &f_name, const String &no_clobber_name, const Verbosity &verbosity)
WORKSPACE METHOD: WriteXML.
Definition: m_xml.h:118
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:33
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:39
const Joker joker
#define CREATE_OUT3
Definition: messages.h:206
#define CREATE_OUT2
Definition: messages.h:205
my_basic_string< char > String
The String type for ARTS.
Definition: mystring.h:216
SingleLineExternal ReadFromHitranOnlineStream(istream &is)
Read from HITRAN online.
SingleLineExternal ReadFromJplStream(istream &is)
Read from JPL.
SingleLineExternal ReadFromArtscat3Stream(istream &is)
Read from ARTSCAT-3.
SingleLineExternal ReadFromLBLRTMStream(istream &is)
Read from LBLRTM.
SingleLineExternal ReadFromArtscat5Stream(istream &is)
Read from ARTSCAT-5.
SingleLineExternal ReadFromHitran2004Stream(istream &is)
Read from newer HITRAN.
SingleLineExternal ReadFromHitran2001Stream(istream &is)
Read from HITRAN before 2004.
SingleLineExternal ReadFromArtscat4Stream(istream &is)
Read from ARTSCAT-4.
std::vector< Lines > split_list_of_external_lines(std::vector< SingleLineExternal > &external_lines, const std::vector< QuantumNumberType > &localquantas={}, const std::vector< QuantumNumberType > &globalquantas={})
Splits a list of lines into proper Lines.
constexpr std::string_view bath_broadening
Name for bath broadening in printing and reading user input.
constexpr bool modelparameterEmpty(const ModelParameters mp) noexcept
constexpr std::string_view self_broadening
Name for self broadening in printing and reading user input.
constexpr Index nVars
Current max number of line shape variables.
constexpr std::array Isotopologues
A list of all ARTS isotopologues, note how the species enum class input HAS to be sorted.
Definition: isotopologues.h:57
consteval std::array< IsotopeRecord, count_isotopologues< spec >()> isotopologues() noexcept
Model GetAdvancedModel(const QuantumIdentifier &qid) ARTS_NOEXCEPT
Returns an advanced Zeeman model.
Definition: zeemandata.cc:131
#define N
Definition: rng.cc:164
std::set< Species::Species > lbl_species(const ArrayOfArrayOfSpeciesTag &abs_species) noexcept
Species::Tag SpeciesTag
Definition: species_tags.h:100
Index NumLines() const noexcept
Number of lines.
std::pair< bool, bool > Match(const Lines &l) const noexcept
Checks if another line list matches this structure.
void AppendSingleLine(SingleLine &&sl)
Appends a single line to the absorption lines.
Single line reading output.
Numeric F0
Central frequency.
Quantum::Number::LocalState localquanta
Local quantum numbers.
Zeeman::Model zeeman
Zeeman model.
ArrayOfQuantumIdentifier levels
void ThrowIfNotOK() const ARTS_NOEXCEPT
Coefficients and temperature model for SingleSpeciesModel.
static constexpr Index N
A logical struct for global quantum numbers with species identifiers.
Species::IsotopeRecord Isotopologue() const noexcept
StateMatchType operates so that a check less than a level should be 'better', bar None.
#define v
#define a
#define c
#define b
void xml_find_and_open_input_file(std::shared_ptr< istream > &ifs, const String &filename, const Verbosity &verbosity)
Open plain or zipped xml file.
Definition: xml_io.cc:136
void xml_read_from_file(const String &filename, T &type, const Verbosity &verbosity)
Reads data from XML file.
Definition: xml_io.h:151
void xml_read_header_from_stream(istream &is, FileType &ftype, NumericType &ntype, EndianType &etype, const Verbosity &verbosity)
Reads XML header and root tag.
Definition: xml_io_base.cc:574
@ FILE_TYPE_ASCII
Definition: xml_io_base.h:43
@ NUMERIC_TYPE_DOUBLE
Definition: xml_io_base.h:48
@ ENDIAN_TYPE_LITTLE
Definition: xml_io_base.h:49