-
Notifications
You must be signed in to change notification settings - Fork 8
/
ratecoeff.cc
1415 lines (1193 loc) · 60.7 KB
/
ratecoeff.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include "ratecoeff.h"
#include <mpi.h>
#if !USE_SIMPSON_INTEGRATOR
#include <gsl/gsl_integration.h>
#endif
#include <array>
#include <cmath>
#include <cstring>
#include <tuple>
// #define D_POSIX_SOURCE
#include <gsl/gsl_errno.h>
#include <cstdio>
#include <cstdlib>
#include <ctime>
#include "artisoptions.h"
#include "atomic.h"
#include "constants.h"
#include "globals.h"
#include "grid.h"
#include "ltepop.h"
#include "macroatom.h"
#include "md5.h"
#include "radfield.h"
#include "rpkt.h"
#include "sn3d.h"
namespace {
double *spontrecombcoeffs{};
// for USE_LUT_PHOTOION = true
double *corrphotoioncoeffs{};
double *bfcooling_coeffs{};
struct GSLIntegralParasGammaCorr {
double nu_edge;
double departure_ratio;
const float *photoion_xs;
float T_e;
int nonemptymgi;
};
char adatafile_hash[33];
char compositionfile_hash[33];
std::array<char[33], 3> phixsfile_hash;
auto read_ratecoeff_dat(FILE *ratecoeff_file) -> bool
// Try to read in the precalculated rate coefficients from file
// return true if successful or false otherwise
{
// Check whether current atomic data and temperature range match
// the precalculated rate coefficients
char adatafile_hash_in[33] = "UNKNOWN";
if (fscanf(ratecoeff_file, "%32s\n", adatafile_hash_in) != 1) {
return false;
}
printout("ratecoeff.dat: MD5 adata.txt = %s ", adatafile_hash_in);
if (strcmp(adatafile_hash, adatafile_hash_in) == 0) {
printout("(pass)\n");
} else {
printout("MISMATCH: MD5 adata.txt = %s\n", adatafile_hash);
return false;
}
char compositionfile_hash_in[33] = "UNKNOWN";
if (fscanf(ratecoeff_file, "%32s\n", compositionfile_hash_in) != 1) {
return false;
}
printout("ratecoeff.dat: MD5 compositiondata.txt %s ", compositionfile_hash_in);
if (strcmp(compositionfile_hash, compositionfile_hash_in) == 0) {
printout("(pass)\n");
} else {
printout("\nMISMATCH: MD5 compositiondata.txt = %s\n", compositionfile_hash);
return false;
}
for (int phixsver = 1; phixsver <= 2; phixsver++) {
if (phixs_file_version_exists[phixsver]) {
char phixsfile_hash_in[33] = "UNKNOWN";
if (fscanf(ratecoeff_file, "%32s\n", phixsfile_hash_in) != 1) {
return false;
}
printout("ratecoeff.dat: MD5 %s = %s ", phixsdata_filenames[phixsver], phixsfile_hash_in);
if (strcmp(phixsfile_hash[phixsver], phixsfile_hash_in) == 0) {
printout("(pass)\n");
} else {
printout("\nMISMATCH: MD5 %s = %s\n", phixsdata_filenames[phixsver], phixsfile_hash[phixsver]);
return false;
}
}
}
double in_T_min = -1.;
double in_T_max = -1.;
int in_tablesize = -1;
int in_nlines = -1;
int in_nbfcontinua = -1;
double in_ratecoeff_integral_accuracy = -1.;
const int items_read = fscanf(ratecoeff_file, "%la %la %d %d %d %la\n", &in_T_min, &in_T_max, &in_tablesize,
&in_nlines, &in_nbfcontinua, &in_ratecoeff_integral_accuracy);
if (items_read != 6) {
printout("\nMISMATCH: error reading header line\n");
return false;
}
printout("ratecoeff.dat: Tmin %g Tmax %g TABLESIZE %d nlines %d nbfcontinua %d in_ratecoeff_integral_accuracy %g ",
in_T_min, in_T_max, in_tablesize, in_nlines, in_nbfcontinua, in_ratecoeff_integral_accuracy);
if (in_T_min != MINTEMP) {
printout("\nMISMATCH: this simulation has MINTEMP %g\n", MINTEMP);
return false;
}
if (in_T_max != MAXTEMP) {
printout("\nMISMATCH: this simulation has MAXTEMP %g\n", MAXTEMP);
return false;
}
if (in_tablesize != TABLESIZE) {
printout("\nMISMATCH: this simulation has TABLESIZE %d\n", TABLESIZE);
return false;
}
if (in_nlines != globals::nlines) {
printout("\nMISMATCH: this simulation has nlines %d\n", globals::nlines);
return false;
}
if (in_nbfcontinua != globals::nbfcontinua) {
printout("\nMISMATCH: this simulation has nbfcontinua %d\n", globals::nbfcontinua);
return false;
}
if (in_ratecoeff_integral_accuracy != RATECOEFF_INTEGRAL_ACCURACY) {
printout("\nMISMATCH: this simulation has RATECOEFF_INTEGRAL_ACCURACY %g\n", RATECOEFF_INTEGRAL_ACCURACY);
return false;
}
printout("(pass)\n");
// this is redundant if the adata and composition data matches, consider removing
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element);
for (int ion = 0; ion < nions; ion++) {
int in_element = 0;
int in_ionstage = 0;
int in_levels = 0;
int in_ionisinglevels = 0;
assert_always(
fscanf(ratecoeff_file, "%d %d %d %d\n", &in_element, &in_ionstage, &in_levels, &in_ionisinglevels) == 4);
const int nlevels = get_nlevels(element, ion);
const int ionisinglevels = get_nlevels_ionising(element, ion);
if (get_atomicnumber(element) != in_element || get_ionstage(element, ion) != in_ionstage ||
nlevels != in_levels || ionisinglevels != in_ionisinglevels) {
printout(
"Levels or ionising levels count mismatch! element %d %d ionstage %d %d nlevels %d %d ionisinglevels "
"%d %d\n",
get_atomicnumber(element), in_element, get_ionstage(element, ion), in_ionstage, nlevels, in_levels,
ionisinglevels, in_ionisinglevels);
return false;
}
}
}
printout("Existing ratecoeff.dat is valid. Reading this file...\n");
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element) - 1;
for (int ion = 0; ion < nions; ion++) {
// nlevels = get_nlevels(element,ion);
const int nlevels = get_nlevels_ionising(element, ion); // number of ionising levels associated with current ion
for (int level = 0; level < nlevels; level++) {
// Loop over the phixs target states
const int nphixstargets = get_nphixstargets(element, ion, level);
for (int phixstargetindex = 0; phixstargetindex < nphixstargets; phixstargetindex++) {
// Loop over the temperature grid
for (int iter = 0; iter < TABLESIZE; iter++) {
double in_alpha_sp{NAN};
double in_bfcooling_coeff{NAN};
double in_corrphotoioncoeff{NAN};
double in_bfheating_coeff{NAN};
assert_always(fscanf(ratecoeff_file, "%la %la %la %la\n", &in_alpha_sp, &in_bfcooling_coeff,
&in_corrphotoioncoeff, &in_bfheating_coeff) == 4);
// assert_always(std::isfinite(alpha_sp) && alpha_sp >= 0);
spontrecombcoeffs[get_bflutindex(iter, element, ion, level, phixstargetindex)] = in_alpha_sp;
// assert_always(std::isfinite(bfcooling_coeff) && bfcooling_coeff >= 0);
bfcooling_coeffs[get_bflutindex(iter, element, ion, level, phixstargetindex)] = in_bfcooling_coeff;
if constexpr (USE_LUT_PHOTOION) {
if (in_corrphotoioncoeff >= 0) {
corrphotoioncoeffs[get_bflutindex(iter, element, ion, level, phixstargetindex)] = in_corrphotoioncoeff;
} else {
printout(
"ERROR: USE_LUT_PHOTOION is on, but there are no corrphotoioncoeff values in ratecoeff file\n");
std::abort();
}
}
if constexpr (USE_LUT_BFHEATING) {
if (in_bfheating_coeff >= 0) {
globals::bfheating_coeff[get_bflutindex(iter, element, ion, level, phixstargetindex)] =
in_bfheating_coeff;
} else {
printout(
"ERROR: USE_LUT_BFHEATING is on, but there are no bfheating_coeff values in the ratecoeff "
"file\n");
std::abort();
}
}
}
}
}
}
}
return true;
}
void write_ratecoeff_dat() {
FILE *ratecoeff_file = fopen_required("ratecoeff.dat", "w");
fprintf(ratecoeff_file, "%32s\n", adatafile_hash);
fprintf(ratecoeff_file, "%32s\n", compositionfile_hash);
for (int phixsver = 1; phixsver <= 2; phixsver++) {
if (phixs_file_version_exists[phixsver]) {
fprintf(ratecoeff_file, "%32s\n", phixsfile_hash[phixsver]);
}
}
fprintf(ratecoeff_file, "%la %la %d %d %d %la\n", MINTEMP, MAXTEMP, TABLESIZE, globals::nlines, globals::nbfcontinua,
RATECOEFF_INTEGRAL_ACCURACY);
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element);
for (int ion = 0; ion < nions; ion++) {
fprintf(ratecoeff_file, "%d %d %d %d\n", get_atomicnumber(element), get_ionstage(element, ion),
get_nlevels(element, ion), get_nlevels_ionising(element, ion));
}
}
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element) - 1;
for (int ion = 0; ion < nions; ion++) {
// nlevels = get_nlevels(element,ion);
const int nlevels = get_nlevels_ionising(element, ion);
for (int level = 0; level < nlevels; level++) {
// Loop over the phixs targets
const auto nphixstargets = get_nphixstargets(element, ion, level);
for (int phixstargetindex = 0; phixstargetindex < nphixstargets; phixstargetindex++) {
// Loop over the temperature grid
for (int iter = 0; iter < TABLESIZE; iter++) {
const int bflutindex = get_bflutindex(iter, element, ion, level, phixstargetindex);
fprintf(ratecoeff_file, "%la %la %la %la\n", spontrecombcoeffs[bflutindex], bfcooling_coeffs[bflutindex],
USE_LUT_PHOTOION ? corrphotoioncoeffs[bflutindex] : -1,
USE_LUT_BFHEATING ? globals::bfheating_coeff[bflutindex] : -1);
}
}
}
}
}
fclose(ratecoeff_file);
}
// Integrand to calculate the rate coefficient for spontaneous recombination
auto alpha_sp_integrand_gsl(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegrationParas *>(voidparas);
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, params->nu_edge, nu);
const double x = TWOOVERCLIGHTSQUARED * sigma_bf * pow(nu, 2) * exp(-HOVERKB * nu / params->T);
// in formula this looks like
// x = sigma_bf/H/nu * 2*H*pow(nu,3)/pow(CLIGHT,2) * exp(-H*nu/KB/T);
// set contributions from Lyman continuum artificially to zero to overcome it's large opacity
return x;
}
// Integrand to calculate the rate coefficient for spontaneous recombination
auto alpha_sp_E_integrand_gsl(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegrationParas *>(voidparas);
const float T = params->T;
const double nu_edge = params->nu_edge;
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, nu_edge, nu);
const double x = TWOOVERCLIGHTSQUARED * sigma_bf * pow(nu, 3) / nu_edge * exp(-HOVERKB * nu / T);
// in formula this looks like
// x = sigma_bf/H/nu * 2*H*pow(nu,3)/pow(CLIGHT,2) * exp(-H*nu/KB/T);
// set contributions from Lyman continuum artificially to zero to overcome it's large opacity
return x;
}
// Integrand to calculate the rate coefficient for photoionization corrected for stimulated recombination.
auto gammacorr_integrand_gsl(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegrationParas *>(voidparas);
const float T = params->T;
const double nu_edge = params->nu_edge;
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, nu_edge, nu);
// Dependence on dilution factor W is linear. This allows to set it here to
// 1. and scale to its actual value later on.
// Assumption T_e = T_R makes n_kappa/n_i * (n_i/n_kappa)* = 1
return sigma_bf * ONEOVERH / nu * radfield::dbb(nu, T, 1) * (1 - exp(-HOVERKB * nu / T));
}
// Integrand to precalculate the bound-free heating ratecoefficient in an approximative way
// on a temperature grid using the assumption that T_e=T_R and W=1 in the ionisation
// formula. The radiation fields dependence on W is taken into account by multiplying
// the resulting expression with the correct W later on.
auto approx_bfheating_integrand_gsl(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegrationParas *>(voidparas);
const float T = params->T;
const double nu_edge = params->nu_edge;
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, nu_edge, nu);
// Precalculation for T_e=T_R and W=1
const double x = sigma_bf * (1 - nu_edge / nu) * radfield::dbb(nu, T, 1) * (1 - exp(-HOVERKB * nu / T));
return x;
}
// Integrand to precalculate the bound-free heating ratecoefficient in an approximative way
// on a temperature grid using the assumption that T_e=T_R and W=1 in the ionisation
// formula. The radiation fields dependence on W is taken into account by multiplying
// the resulting expression with the correct W later on.
auto bfcooling_integrand_gsl(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegrationParas *>(voidparas);
const float T = params->T;
const double nu_edge = params->nu_edge;
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, nu_edge, nu);
// return sigma_bf * (1-nu_edge/nu) * TWOHOVERCLIGHTSQUARED * pow(nu,3) * exp(-HOVERKB*nu/T);
return sigma_bf * (nu - nu_edge) * TWOHOVERCLIGHTSQUARED * nu * nu * exp(-HOVERKB * nu / T);
}
void precalculate_rate_coefficient_integrals() {
// target fractional accuracy of the integrator //=1e-5 took 8 hours with Fe I to V!
const double epsrelwarning = 1e-2; // fractional error to emit a warning
// Calculate the rate coefficients for each level of each ion of each element
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element) - 1;
#ifdef _OPENMP
#pragma omp parallel for
#endif
for (int ion = 0; ion < nions; ion++) {
const int atomic_number = get_atomicnumber(element);
const int ionstage = get_ionstage(element, ion);
const int nlevels = get_nlevels_ionising(element, ion);
printout("Performing rate integrals for Z = %d, ionstage %d...\n", atomic_number, ionstage);
gsl_error_handler_t *previous_handler = gsl_set_error_handler(gsl_error_handler_printout);
for (int level = 0; level < nlevels; level++) {
if ((level > 0) && (level % 50 == 0)) {
printout(" completed up to level %d of %d\n", level, nlevels);
}
// coefficients are stored in node shared memory, so divide up the work on the node
if ((level % globals::node_nprocs) != globals::rank_in_node) {
continue;
}
const int nphixstargets = get_nphixstargets(element, ion, level);
for (int phixstargetindex = 0; phixstargetindex < nphixstargets; phixstargetindex++) {
const int upperlevel = get_phixsupperlevel(element, ion, level, phixstargetindex);
const double phixstargetprobability = get_phixsprobability(element, ion, level, phixstargetindex);
// const double E_threshold = epsilon(element,ion+1,upperlevel) - epsilon(element,ion,level);
const double E_threshold = get_phixs_threshold(element, ion, level, phixstargetindex);
const double nu_threshold = E_threshold / H;
const double nu_max_phixs =
nu_threshold * last_phixs_nuovernuedge; // nu of the uppermost point in the phixs table
// Loop over the temperature grid
for (int iter = 0; iter < TABLESIZE; iter++) {
const int bflutindex = get_bflutindex(iter, element, ion, level, phixstargetindex);
double error{NAN};
int status = 0;
const float T_e = MINTEMP * exp(iter * T_step_log);
const double sfac = calculate_sahafact(element, ion, level, upperlevel, T_e, E_threshold);
assert_always(get_phixs_table(element, ion, level) != nullptr);
// the threshold of the first target gives nu of the first phixstable point
const GSLIntegrationParas intparas = {
.nu_edge = nu_threshold, .T = T_e, .photoion_xs = get_phixs_table(element, ion, level)};
// Spontaneous recombination and bf-cooling coefficient don't depend on the radiation field
double alpha_sp = 0.;
status =
integrator<alpha_sp_integrand_gsl>(intparas, nu_threshold, nu_max_phixs, 0, RATECOEFF_INTEGRAL_ACCURACY,
GSL_INTEG_GAUSS61, &alpha_sp, &error);
if (status != 0 && (status != 18 || (error / alpha_sp) > epsrelwarning)) {
printout("alpha_sp integrator status %d. Integral value %9.3e +/- %9.3e\n", status, alpha_sp, error);
}
alpha_sp *= FOURPI * sfac * phixstargetprobability;
if (!std::isfinite(alpha_sp) || alpha_sp < 0) {
printout(
"WARNING: alpha_sp was negative or non-finite for level %d Te %g. alpha_sp %g sfac %g "
"phixstargetindex %d phixstargetprobability %g\n",
level, T_e, alpha_sp, sfac, phixstargetindex, phixstargetprobability);
alpha_sp = 0;
}
spontrecombcoeffs[bflutindex] = alpha_sp;
if constexpr (USE_LUT_PHOTOION) {
double gammacorr = 0.;
status = integrator<gammacorr_integrand_gsl>(intparas, nu_threshold, nu_max_phixs, 0,
RATECOEFF_INTEGRAL_ACCURACY, GSL_INTEG_GAUSS61, &gammacorr,
&error);
if (status != 0 && (status != 18 || (error / gammacorr) > epsrelwarning)) {
printout("gammacorr integrator status %d. Integral value %9.3e +/- %9.3e\n", status, gammacorr, error);
}
gammacorr *= FOURPI * phixstargetprobability;
assert_always(gammacorr >= 0);
if (gammacorr < 0) {
printout("WARNING: gammacorr was negative for level %d\n", level);
gammacorr = 0;
}
corrphotoioncoeffs[bflutindex] = gammacorr;
}
if constexpr (USE_LUT_BFHEATING) {
double this_bfheating_coeff = 0.;
status = integrator<approx_bfheating_integrand_gsl>(intparas, nu_threshold, nu_max_phixs, 0,
RATECOEFF_INTEGRAL_ACCURACY, GSL_INTEG_GAUSS61,
&this_bfheating_coeff, &error);
if (status != 0 && (status != 18 || (error / this_bfheating_coeff) > epsrelwarning)) {
printout("bfheating_coeff integrator status %d. Integral value %9.3e +/- %9.3e\n", status,
this_bfheating_coeff, error);
}
this_bfheating_coeff *= FOURPI * phixstargetprobability;
if (this_bfheating_coeff < 0) {
printout("WARNING: bfheating_coeff was negative for level %d\n", level);
this_bfheating_coeff = 0;
}
globals::bfheating_coeff[bflutindex] = this_bfheating_coeff;
}
double this_bfcooling_coeff = 0.;
status = integrator<bfcooling_integrand_gsl>(intparas, nu_threshold, nu_max_phixs, 0,
RATECOEFF_INTEGRAL_ACCURACY, GSL_INTEG_GAUSS61,
&this_bfcooling_coeff, &error);
if (status != 0 && (status != 18 || (error / this_bfcooling_coeff) > epsrelwarning)) {
printout("bfcooling_coeff integrator status %d. Integral value %9.3e +/- %9.3e\n", status,
this_bfcooling_coeff, error);
}
this_bfcooling_coeff *= FOURPI * sfac * phixstargetprobability;
if (!std::isfinite(this_bfcooling_coeff) || this_bfcooling_coeff < 0) {
printout(
"WARNING: bfcooling_coeff was negative or non-finite for level %d Te %g. bfcooling_coeff %g sfac %g "
"phixstargetindex %d phixstargetprobability %g\n",
level, T_e, this_bfcooling_coeff, sfac, phixstargetindex, phixstargetprobability);
this_bfcooling_coeff = 0;
}
bfcooling_coeffs[bflutindex] = this_bfcooling_coeff;
}
}
}
gsl_set_error_handler(previous_handler);
}
}
}
// multiply the cross sections associated with a level by some factor and
// also update the quantities integrated from (and proportional to) the cross sections
void scale_level_phixs(const int element, const int ion, const int level, const double factor) {
// if we store the data in node shared memory, then only one rank should update it
if (globals::rank_in_node == 0) {
const int nphixstargets = get_nphixstargets(element, ion, level);
if (nphixstargets == 0) {
return;
}
auto *phixstable = get_phixs_table(element, ion, level);
for (int n = 0; n < globals::NPHIXSPOINTS; n++) {
phixstable[n] *= factor;
}
for (int phixstargetindex = 0; phixstargetindex < nphixstargets; phixstargetindex++) {
for (int iter = 0; iter < TABLESIZE; iter++) {
const auto bflutindex = get_bflutindex(iter, element, ion, level, phixstargetindex);
spontrecombcoeffs[bflutindex] *= factor;
if constexpr (USE_LUT_PHOTOION) {
corrphotoioncoeffs[bflutindex] *= factor;
}
if constexpr (USE_LUT_BFHEATING) {
globals::bfheating_coeff[bflutindex] *= factor;
}
bfcooling_coeffs[bflutindex] *= factor;
}
}
}
}
// calibrate the recombination rates to tabulated values by scaling the photoionisation cross sections
void read_recombrate_file() {
use_cellcache = false;
FILE *recombrate_file = fopen("recombrates.txt", "r");
if (recombrate_file == nullptr) {
printout("No recombrates.txt file found. Skipping recombination rate scaling...\n");
return;
}
printout("Reading recombination rate file (recombrates.txt)...\n");
const double Te_estimate = RECOMBCALIBRATION_T_ELEC;
const double log_Te_estimate = log10(Te_estimate);
printout("Calibrating recombination rates for a temperature of %.1f K\n", Te_estimate);
struct RRCRow {
double log_Te;
double rrc_low_n;
double rrc_total;
};
int atomicnumber = 0;
int upperionstage = 0;
int tablerows = 0;
while (fscanf(recombrate_file, "%d %d %d\n", &atomicnumber, &upperionstage, &tablerows) > 0) {
RRCRow T_highestbelow = {.log_Te = 0, .rrc_low_n = 0, .rrc_total = 0};
RRCRow T_lowestabove = {.log_Te = 0, .rrc_low_n = 0, .rrc_total = 0};
T_highestbelow.log_Te = -1;
T_lowestabove.log_Te = -1;
for (int i = 0; i < tablerows; i++) {
RRCRow row{};
assert_always(fscanf(recombrate_file, "%lg %lg %lg\n", &row.log_Te, &row.rrc_low_n, &row.rrc_total) == 3);
if (row.log_Te < log_Te_estimate && row.log_Te > T_highestbelow.log_Te) {
T_highestbelow.log_Te = row.log_Te;
T_highestbelow.rrc_low_n = row.rrc_low_n;
T_highestbelow.rrc_total = row.rrc_total;
}
if (row.log_Te > log_Te_estimate && (row.log_Te < T_lowestabove.log_Te || T_lowestabove.log_Te < 0)) {
T_lowestabove.log_Te = row.log_Te;
T_lowestabove.rrc_low_n = row.rrc_low_n;
T_lowestabove.rrc_total = row.rrc_total;
}
}
const int element = get_elementindex(atomicnumber);
if (element >= 0) {
const int ion = upperionstage - get_ionstage(element, 0); // the index of the upper ion
if (ion > 0 && ion < get_nions(element)) {
printout("Z=%d ionstage %d->%d\n", atomicnumber, upperionstage, upperionstage - 1);
assert_always(T_highestbelow.log_Te > 0);
assert_always(T_lowestabove.log_Te > 0);
const int nlevels = get_nlevels_ionising(element, ion - 1);
const double x = (log_Te_estimate - T_highestbelow.log_Te) / (T_lowestabove.log_Te - T_highestbelow.log_Te);
const double input_rrc_low_n = (x * T_highestbelow.rrc_low_n) + ((1 - x) * T_lowestabove.rrc_low_n);
const double input_rrc_total = (x * T_highestbelow.rrc_total) + ((1 - x) * T_lowestabove.rrc_total);
const bool assume_lte = true;
const bool printdebug = false;
const bool per_groundmultipletpop = true;
double rrc = calculate_ionrecombcoeff(-1, Te_estimate, element, ion, assume_lte, false, printdebug, false,
per_groundmultipletpop, false);
printout(" rrc: %10.3e\n", rrc);
if (input_rrc_low_n >= 0) // if it's < 0, ignore it
{
printout(" input_rrc_low_n: %10.3e\n", input_rrc_low_n);
const double phixs_multiplier = input_rrc_low_n / rrc;
if (phixs_multiplier < 0.05 || phixs_multiplier >= 2.0) {
printout(" Not scaling phixs of all levels by %.3f (because < 0.05 or >= 2.0)\n", phixs_multiplier);
} else {
printout(" scaling phixs of all levels by %.3f\n", phixs_multiplier);
for (int level = 0; level < nlevels; level++) {
scale_level_phixs(element, ion - 1, level, phixs_multiplier);
}
rrc = calculate_ionrecombcoeff(-1, Te_estimate, element, ion, assume_lte, false, printdebug, false,
per_groundmultipletpop, false);
printout(" rrc: %10.3e\n", rrc);
}
}
// hopefully the RRC now matches the low_n value well, if it was defined
// Next, use the superlevel recombination rates to make up the excess needed to reach the total RRC
printout(" input_rrc_total: %10.3e\n", input_rrc_total);
if (rrc < input_rrc_total) {
const double rrc_superlevel = calculate_ionrecombcoeff(-1, Te_estimate, element, ion, assume_lte, false,
printdebug, true, per_groundmultipletpop, false);
printout(" rrc(superlevel): %10.3e\n", rrc_superlevel);
if (rrc_superlevel > 0) {
const double phixs_multiplier_superlevel = 1.0 + ((input_rrc_total - rrc) / rrc_superlevel);
printout(" scaling phixs of levels in the superlevel by %.3f\n", phixs_multiplier_superlevel);
assert_always(phixs_multiplier_superlevel >= 0);
const int first_superlevel_level = get_nlevels_nlte(element, ion - 1) + 1;
for (int level = first_superlevel_level; level < nlevels; level++) {
scale_level_phixs(element, ion - 1, level, phixs_multiplier_superlevel);
}
} else {
printout("There is no superlevel recombination, so multiplying all levels instead\n");
const double phixs_multiplier = input_rrc_total / rrc;
printout(" scaling phixs of all levels by %.3f\n", phixs_multiplier);
assert_always(phixs_multiplier >= 0);
for (int level = 0; level < nlevels; level++) {
scale_level_phixs(element, ion - 1, level, phixs_multiplier);
}
}
} else {
printout("rrc >= input_rrc_total!\n");
const double phixs_multiplier = input_rrc_total / rrc;
printout(" scaling phixs of all levels by %.3f\n", phixs_multiplier);
assert_always(phixs_multiplier >= 0);
for (int level = 0; level < nlevels; level++) {
scale_level_phixs(element, ion - 1, level, phixs_multiplier);
}
}
rrc = calculate_ionrecombcoeff(-1, Te_estimate, element, ion, assume_lte, false, printdebug, false,
per_groundmultipletpop, false);
printout(" rrc: %10.3e\n", rrc);
}
}
}
fclose(recombrate_file);
}
void precalculate_ion_alpha_sp() {
globals::ion_alpha_sp.resize(get_includedions() * TABLESIZE);
for (int iter = 0; iter < TABLESIZE; iter++) {
const float T_e = MINTEMP * exp(iter * T_step_log);
for (int element = 0; element < get_nelements(); element++) {
const int nions = get_nions(element) - 1;
for (int ion = 0; ion < nions; ion++) {
const auto uniqueionindex = get_uniqueionindex(element, ion);
const int nionisinglevels = get_nlevels_ionising(element, ion);
double zeta = 0.;
for (int level = 0; level < nionisinglevels; level++) {
const auto nphixstargets = get_nphixstargets(element, ion, level);
for (int phixstargetindex = 0; phixstargetindex < nphixstargets; phixstargetindex++) {
const double zeta_level = get_spontrecombcoeff(element, ion, level, phixstargetindex, T_e);
zeta += zeta_level;
}
}
globals::ion_alpha_sp[(uniqueionindex * TABLESIZE) + iter] = zeta;
}
}
}
}
auto integrand_stimrecombination_custom_radfield(const double nu, void *const voidparas) -> double {
const auto *const params = static_cast<const GSLIntegralParasGammaCorr *>(voidparas);
const int nonemptymgi = params->nonemptymgi;
const float T_e = params->T_e;
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, params->nu_edge, nu);
const double Jnu = radfield::radfield(nu, nonemptymgi);
// TODO: MK thesis page 41, use population ratios and Te?
return ONEOVERH * sigma_bf / nu * Jnu * exp(-HOVERKB * nu / T_e);
}
auto calculate_stimrecombcoeff_integral(const int element, const int lowerion, const int level,
const int phixstargetindex, const int nonemptymgi) -> double {
const double epsrel = 1e-3;
const double epsabs = 0.;
const double E_threshold = get_phixs_threshold(element, lowerion, level, phixstargetindex);
const double nu_threshold = ONEOVERH * E_threshold;
const double nu_max_phixs = nu_threshold * last_phixs_nuovernuedge; // nu of the uppermost point in the phixs table
const auto T_e = grid::get_Te(nonemptymgi);
const auto intparas = GSLIntegralParasGammaCorr{
.nu_edge = nu_threshold,
.departure_ratio = 1., // not used, but must be set to something
.photoion_xs = get_phixs_table(element, lowerion, level),
.T_e = T_e,
.nonemptymgi = nonemptymgi,
};
const int upperionlevel = get_phixsupperlevel(element, lowerion, level, phixstargetindex);
const double sf = calculate_sahafact(element, lowerion, level, upperionlevel, T_e, H * nu_threshold);
double error = 0.;
gsl_error_handler_t *previous_handler = gsl_set_error_handler(gsl_error_handler_printout);
double stimrecombcoeff = 0.;
// const int status =
integrator<integrand_stimrecombination_custom_radfield>(intparas, nu_threshold, nu_max_phixs, epsabs, epsrel,
GSL_INTEG_GAUSS61, &stimrecombcoeff, &error);
gsl_set_error_handler(previous_handler);
stimrecombcoeff *= FOURPI * sf * get_phixsprobability(element, lowerion, level, phixstargetindex);
// if (status != 0)
// {
// error *= FOURPI * get_phixsprobability(element, ion, level, phixstargetindex);
// printout("stimrecombcoeff gsl integrator warning %d. modelgridindex %d Z=%d ionstage %d lower %d
// phixstargetindex %d gamma %g error %g\n",
// status, modelgridindex, get_atomicnumber(element), get_ionstage(element, ion), level,
// phixstargetindex, gammacorr, error);
// }
return stimrecombcoeff;
}
auto integrand_corrphotoioncoeff_custom_radfield(const double nu, void *const voidparas) -> double
// Integrand to calculate the rate coefficient for photoionization
// using gsl integrators. Corrected for stimulated recombination.
{
const GSLIntegralParasGammaCorr *const params = static_cast<GSLIntegralParasGammaCorr *>(voidparas);
const int nonemptymgi = params->nonemptymgi;
#if (SEPARATE_STIMRECOMB)
const double corrfactor = 1.;
#else
const float T_e = params->T_e;
double corrfactor = 1. - (params->departure_ratio * exp(-HOVERKB * nu / T_e));
if (corrfactor < 0) {
corrfactor = 0.;
}
#endif
const float sigma_bf = photoionization_crosssection_fromtable(params->photoion_xs, params->nu_edge, nu);
const double Jnu = radfield::radfield(nu, nonemptymgi);
// TODO: MK thesis page 41, use population ratios and Te?
return ONEOVERH * sigma_bf / nu * Jnu * corrfactor;
}
auto calculate_corrphotoioncoeff_integral(int element, const int ion, const int level, const int phixstargetindex,
int nonemptymgi) -> double {
constexpr double epsrel = 1e-3;
constexpr double epsrelwarning = 1e-1;
constexpr double epsabs = 0.;
const double E_threshold = get_phixs_threshold(element, ion, level, phixstargetindex);
const double nu_threshold = ONEOVERH * E_threshold;
const double nu_max_phixs = nu_threshold * last_phixs_nuovernuedge; // nu of the uppermost point in the phixs table
const auto T_e = grid::get_Te(nonemptymgi);
#if SEPARATE_STIMRECOMB
const double departure_ratio = 0.; // zero the stimulated recomb contribution
#else
// stimulated recombination is negative photoionisation
const double nnlevel = get_levelpop(nonemptymgi, element, ion, level);
const double nne = grid::get_nne(nonemptymgi);
const int upperionlevel = get_phixsupperlevel(element, ion, level, phixstargetindex);
const double sf = calculate_sahafact(element, ion, level, upperionlevel, T_e, H * nu_threshold);
const double nnupperionlevel = get_levelpop(nonemptymgi, element, ion + 1, upperionlevel);
double departure_ratio = nnlevel > 0. ? nnupperionlevel / nnlevel * nne * sf : 1.; // put that to phixslist
if (!std::isfinite(departure_ratio)) {
departure_ratio = 0.;
}
#endif
const auto intparas = GSLIntegralParasGammaCorr{
.nu_edge = nu_threshold,
.departure_ratio = departure_ratio,
.photoion_xs = get_phixs_table(element, ion, level),
.T_e = T_e,
.nonemptymgi = nonemptymgi,
};
double error = 0.;
gsl_error_handler_t *previous_handler = gsl_set_error_handler(gsl_error_handler_printout);
double gammacorr = 0.;
const int status = integrator<integrand_corrphotoioncoeff_custom_radfield>(
intparas, nu_threshold, nu_max_phixs, epsabs, epsrel, GSL_INTEG_GAUSS61, &gammacorr, &error);
gsl_set_error_handler(previous_handler);
if (status != 0 && (status != 18 || (error / gammacorr) > epsrelwarning)) {
printout(
"corrphotoioncoeff gsl integrator warning %d. modelgridindex %d Z=%d ionstage %d lower %d phixstargetindex %d "
"integral %g error %g\n",
status, grid::get_mgi_of_nonemptymgi(nonemptymgi), get_atomicnumber(element), get_ionstage(element, ion), level,
phixstargetindex, gammacorr, error);
if (!std::isfinite(gammacorr)) {
gammacorr = 0.;
}
}
gammacorr *= FOURPI * get_phixsprobability(element, ion, level, phixstargetindex);
return gammacorr;
}
// get the number of levels that make up a fraction of the ion population
// of at least IONGAMMA_POPFRAC_LEVELS_INCLUDED
auto get_nlevels_important(const int nonemptymgi, const int element, const int ion, const bool assume_lte,
const float T_e) -> std::tuple<int, double> {
// get the stored ion population for comparison with the cumulative sum of level pops
const double nnion_real = get_nnion(nonemptymgi, element, ion);
if (IONGAMMA_POPFRAC_LEVELS_INCLUDED >= 1.) {
return {get_nlevels(element, ion), nnion_real};
}
double nnlevelsum = 0.;
int nlevels_important = get_nlevels_ionising(element, ion); // levels needed to get majority of ion pop
// debug: treat all ionising levels as important
// *nnlevelsum_out = nnion_real;
// return nlevels_important;
for (int lower = 0; lower < get_nlevels_ionising(element, ion); lower++) {
double nnlowerlevel{NAN};
if (assume_lte) {
const double T_exc = T_e; // remember, other parts of the code in LTE mode use TJ, not T_e
const double E_level = epsilon(element, ion, lower);
const double E_ground = epsilon(element, ion, 0);
const double nnground = get_groundlevelpop(nonemptymgi, element, ion);
nnlowerlevel = (nnground * stat_weight(element, ion, lower) / stat_weight(element, ion, 0) *
exp(-(E_level - E_ground) / KB / T_exc));
} else {
nnlowerlevel = get_levelpop(nonemptymgi, element, ion, lower);
}
nnlevelsum += nnlowerlevel;
nlevels_important = lower + 1;
if ((nnlevelsum / nnion_real) >= IONGAMMA_POPFRAC_LEVELS_INCLUDED) {
break;
}
}
assert_always(nlevels_important <= get_nlevels(element, ion));
return {nlevels_important, nnlevelsum};
}
} // anonymous namespace
void setup_photoion_luts() {
assert_always(globals::nbfcontinua > 0);
size_t mem_usage_photoionluts = 2 * TABLESIZE * globals::nbfcontinua * sizeof(double);
spontrecombcoeffs = MPI_shared_malloc<double>(TABLESIZE * globals::nbfcontinua);
assert_always(spontrecombcoeffs != nullptr);
if constexpr (USE_LUT_PHOTOION) {
corrphotoioncoeffs = MPI_shared_malloc<double>(TABLESIZE * globals::nbfcontinua);
assert_always(corrphotoioncoeffs != nullptr);
mem_usage_photoionluts += TABLESIZE * globals::nbfcontinua * sizeof(double);
}
if constexpr (USE_LUT_BFHEATING) {
globals::bfheating_coeff = MPI_shared_malloc<double>(TABLESIZE * globals::nbfcontinua);
assert_always(globals::bfheating_coeff != nullptr);
mem_usage_photoionluts += TABLESIZE * globals::nbfcontinua * sizeof(double);
}
bfcooling_coeffs = MPI_shared_malloc<double>(TABLESIZE * globals::nbfcontinua);
printout(
"[info] mem_usage: lookup tables derived from photoionisation (spontrecombcoeff, bfcooling and "
"corrphotoioncoeff/bfheating if enabled) occupy %.3f MB\n",
mem_usage_photoionluts / 1024. / 1024.);
}
__host__ __device__ auto select_continuum_nu(int element, const int lowerion, const int lower, const int upperionlevel,
float T_e) -> double {
const int phixstargetindex = get_phixtargetindex(element, lowerion, lower, upperionlevel);
const double E_threshold = get_phixs_threshold(element, lowerion, lower, phixstargetindex);
const double nu_threshold = ONEOVERH * E_threshold;
const double nu_max_phixs = nu_threshold * last_phixs_nuovernuedge; // nu of the uppermost point in the phixs table
const int npieces = globals::NPHIXSPOINTS;
const GSLIntegrationParas intparas = {
.nu_edge = nu_threshold, .T = T_e, .photoion_xs = get_phixs_table(element, lowerion, lower)};
const double zrand = 1. - rng_uniform(); // Make sure that 0 < zrand <= 1
const double deltanu = (nu_max_phixs - nu_threshold) / npieces;
double error{NAN};
#if !USE_SIMPSON_INTEGRATOR
gsl_error_handler_t *previous_handler = gsl_set_error_handler(gsl_error_handler_printout);
#endif
double total_alpha_sp = 0.;
integrator<alpha_sp_E_integrand_gsl>(intparas, nu_threshold, nu_max_phixs, 0, CONTINUUM_NU_INTEGRAL_ACCURACY,
GSL_INTEG_GAUSS31, &total_alpha_sp, &error);
double alpha_sp_old = total_alpha_sp;
double alpha_sp = total_alpha_sp;
int i = 1;
for (i = 1; i < npieces; i++) {
alpha_sp_old = alpha_sp;
const double xlow = nu_threshold + (i * deltanu);
// Spontaneous recombination and bf-cooling coefficient don't depend on the radiation field
integrator<alpha_sp_E_integrand_gsl>(intparas, xlow, nu_max_phixs, 0, CONTINUUM_NU_INTEGRAL_ACCURACY,
GSL_INTEG_GAUSS31, &alpha_sp, &error);
if (zrand >= alpha_sp / total_alpha_sp) {
break;
}
}
#if !USE_SIMPSON_INTEGRATOR
gsl_set_error_handler(previous_handler);
#endif
const double nuoffset =
(alpha_sp != alpha_sp_old) ? (total_alpha_sp * zrand - alpha_sp_old) / (alpha_sp - alpha_sp_old) * deltanu : 0.;
const double nu_lower = nu_threshold + ((i - 1) * deltanu) + nuoffset;
assert_testmodeonly(std::isfinite(nu_lower));
return nu_lower;
}
// Return the rate coefficient for spontaneous recombination.
__host__ __device__ auto get_spontrecombcoeff(int element, const int ion, const int level, const int phixstargetindex,
float T_e) -> double {
double Alpha_sp{NAN};
const int lowerindex = floor(log(T_e / MINTEMP) / T_step_log);
assert_always(lowerindex >= 0);
if (lowerindex < TABLESIZE - 1) {
const int upperindex = lowerindex + 1;
const double T_lower = MINTEMP * exp(lowerindex * T_step_log);
const double T_upper = MINTEMP * exp(upperindex * T_step_log);
const double f_upper = spontrecombcoeffs[get_bflutindex(upperindex, element, ion, level, phixstargetindex)];
const double f_lower = spontrecombcoeffs[get_bflutindex(lowerindex, element, ion, level, phixstargetindex)];
Alpha_sp = (f_lower + (f_upper - f_lower) / (T_upper - T_lower) * (T_e - T_lower));
} else {
Alpha_sp = spontrecombcoeffs[get_bflutindex(TABLESIZE - 1, element, ion, level, phixstargetindex)];
}
return Alpha_sp;
}
// multiply by upper ion population (or ground population if per_groundmultipletpop is true) and nne to get a rate
auto calculate_ionrecombcoeff(const int modelgridindex, const float T_e, const int element, const int upperion,
const bool assume_lte, const bool collisional_not_radiative, const bool printdebug,
const bool lower_superlevel_only, const bool per_groundmultipletpop, const bool stimonly)
-> double {
const int lowerion = upperion - 1;
if (lowerion < 0) {
return 0.;
}
double alpha = 0.;
if (lowerion < get_nions(element) - 1) {
const auto nonemptymgi = (modelgridindex >= 0) ? grid::get_nonemptymgi_of_mgi(modelgridindex) : -1;
// this gets divided and cancelled out in the radiative case anyway
const double nne = (modelgridindex >= 0) ? grid::get_nne(nonemptymgi) : 1.;
double nnupperion = 0;
// nnupperion = get_groundmultiplet_pop(modelgridindex, T_e, element, upperion, assume_lte);
int upper_nlevels = 0;
if (per_groundmultipletpop) {
// assume that photoionisation of the ion below is only to the ground multiplet levels of the current ion
// const int nphixstargets = get_nphixstargets(element, lowerion, 0);
// upper_nlevels = get_phixsupperlevel(element, lowerion, 0, nphixstargets - 1) + 1;
upper_nlevels = get_nlevels_groundterm(element, lowerion + 1);
} else {
upper_nlevels = get_nlevels(element, lowerion + 1);
}
for (int upper = 0; upper < upper_nlevels; upper++) {
double nnupperlevel{NAN};
if (assume_lte) {
const double T_exc = T_e;
const double E_level = epsilon(element, lowerion + 1, upper);
const double E_ground = epsilon(element, lowerion + 1, 0);
const double nnground = (modelgridindex >= 0) ? get_groundlevelpop(nonemptymgi, element, lowerion + 1) : 1.;
nnupperlevel = (nnground * stat_weight(element, lowerion + 1, upper) / stat_weight(element, lowerion + 1, 0) *
exp(-(E_level - E_ground) / KB / T_exc));
} else {
nnupperlevel = get_levelpop(nonemptymgi, element, lowerion + 1, upper);
}