-
Notifications
You must be signed in to change notification settings - Fork 7
/
cooling_metal_H2.c
2911 lines (2503 loc) · 112 KB
/
cooling_metal_H2.c
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 "mode.h"
#ifdef GASOLINE
#ifndef NOCOOLING
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <math.h>
#include <assert.h>
/* #include <rpc/xdr.h> */
/* Usage/Interfaces:
Functions starting with
Cool are public: intended to be used by outside callers
cl are private: only intended for using by cooling code
*/
/* General accuracy target */
#define EPS 1e-5
#define MAXABUNDITERATIONS 20
/* Accuracy for Temperature estimation for E,rho assuming eqm abundances */
#define EPSTEMP 1e-5
#define ETAMINTIMESTEP 1e-4
/* This min is to prevent Y_e = 0 which shuts off all collisional processes (cooling, recomb or ionization)
and prevents unphysical negative Y_e (resulting from roundoff)
This is important for cases with no cosmic rays or UV
*/
#ifndef Y_EMIN
#define Y_EMIN 1e-7
#endif
#define TESTRATE 1e-3
/* #define CUBICTABLEINTERP */
/* #define USETABLE */
#define CLRATES( _cl, _Rate, _T, _rho, _ZMetal, _correL, _Rate_Phot_H2_stellar) clRates( _cl, _Rate, _T, _rho, _ZMetal, _correL, _Rate_Phot_H2_stellar)
#define CLEDOTINSTANT( _cl, _Y, _Rate, _rho, _ZMetal, _Heat, _Cool ) clEdotInstant( _cl, _Y, _Rate, _rho, _ZMetal, _Heat, _Cool)
#define TABLEFACTOR 1
#define TABLEINTERP( _rname ) (wTln0*RT0->_rname+wTln1*RT1->_rname)
#include "stiff.h"
#if defined(COOLDEBUG) || defined(STARFORM) || defined(SIMPLESF)
#include "pkd.h"
#else
#include "cooling.h"
#endif
#include "outtype.h"
/* Integrator constants */
/* Accuracy target for integrators */
#define EPSINTEG 1e-3
#define WARNINTEGITS 100
#define MAXINTEGITS 10000
/* Debugging Information */
#define PARTICLEIORD 100
#define CL_eH2 (4.476*CL_eV_erg) /*Energy lost during H2 dissociation, Shapira & Kang (1987) through Abel 1997, page 194 */
#define CL_eHI (13.60*CL_eV_erg)
#define CL_eHeI (24.59*CL_eV_erg)
#define CL_eHeII (54.42*CL_eV_erg)
#define CL_E2HeII (3.0*13.6*CL_eV_erg)
#define EMUL (1.000001)
#define M_H 1.672e-24
#define SUM_Metal 0.0184707
COOL *CoolInit( )
{
COOL *cl;
cl = (COOL *) malloc(sizeof(COOL));
assert(cl!=NULL);
cl->nUV = 0;
cl->UV = NULL;
cl->nTable = 0;
cl->RT = NULL;
cl->MetalCoolln = NULL;
cl->MetalHeatln = NULL;
cl->nTableRead = 0; /* Internal Tables read from Files */
cl->DerivsData = malloc(sizeof(clDerivsData));
assert(cl->DerivsData != NULL);
((clDerivsData *) (cl->DerivsData))->IntegratorContext =
StiffInit( EPSINTEG, 5, cl->DerivsData, clDerivs); /*Change array length to 5 to include H2*/
return cl;
}
void CoolFinalize(COOL *cl )
{
StiffFinalize( ((clDerivsData *) (cl->DerivsData))->IntegratorContext );
free(cl->DerivsData);
if (cl->UV != NULL) free(cl->UV);
if (cl->RT != NULL) free(cl->RT);
if (cl->MetalCoolln != NULL) free(cl->MetalCoolln);
if (cl->MetalHeatln != NULL) free (cl->MetalHeatln);
free(cl);
}
void clInitConstants( COOL *cl, double dGmPerCcUnit, double dComovingGmPerCcUnit,
double dErgPerGmUnit, double dSecUnit, double dKpcUnit, double dMsolUnit, double dInitStarMass, COOLPARAM CoolParam)
{
assert(cl!=NULL);
cl->dGmPerCcUnit = dGmPerCcUnit;
cl->dComovingGmPerCcUnit = dComovingGmPerCcUnit;
cl->dErgPerGmUnit = dErgPerGmUnit;
cl->dSecUnit = dSecUnit;
cl->dErgPerGmPerSecUnit = cl->dErgPerGmUnit / cl->dSecUnit;
cl->diErgPerGmUnit = 1./dErgPerGmUnit;
cl->dKpcUnit = dKpcUnit;
cl->dMsolUnit = dMsolUnit;
cl->dMassFracHelium = CoolParam.dMassFracHelium;
cl->dInitStarMass = dInitStarMass;
cl->bUV = CoolParam.bUV;
cl->bUVTableUsesTime = CoolParam.bUVTableUsesTime;
cl->bUVTableLinear = CoolParam.bUVTableUsesTime; /* Linear if using time */
cl->bLowTCool = CoolParam.bLowTCool;
cl->bSelfShield = CoolParam.bSelfShield;
cl->bMetal = CoolParam.bMetal;
cl->bShieldHI = CoolParam.bShieldHI;
cl->dClump = CoolParam.dClump;
{
clDerivsData *Data = cl->DerivsData;
assert(Data != NULL);
Data->cl = cl;
Data->dlnE = (log(EMUL)-log(1/EMUL));
}
}
/* caculate the number fraction YH, YHe as a function of metalicity. Cosmic
production rate of helium relative to metals (in mass)
delta Y/delta Z = 2.1 and primordial He Yp = 0.236 (Jimenez et al. 2003,Science
299, 5612.
piecewise linear
Y = Yp + dY/dZ*ZMetal up to ZMetal = 0.1, then linear decrease to 0 at Z=1)
SUM_Metal = sum(ni/nH *mi),it is a fixed number for cloudy abundance.
Massfraction fmetal = Z*SUM_metal/(1 + 4*nHe/nH + Z*SUM_metal) (1)
4*nHe/nH = mHe/mH = fHe/fH
also fH + fHe + fMetal = 1 (2)
if fHe specified, combining the 2 eq above will solve for
fH and fMetal
a = CoolParam.dMassFracHelium;
c = SUM_Metal*ZMetal;
fH = (-(2.0*a + a*c -1) + sqrt(pow((2.*a+a*c-1),2.)-4.0*a*(1+c)*(a-1)))/(2.*(1+c));
fMetal = 1.0 - fH - CoolParam.dMassFracHelium;
*/
void clSetAbundanceTotals(COOL *cl, double ZMetal, double *pY_H, double *pY_He, double *pY_eMax) {
double Y_H, Y_He;
if (ZMetal <= 0.1) {
Y_He = (0.236 + 2.1*ZMetal)/4.0;
}
else {
Y_He = (-0.446*(ZMetal - 0.1)/0.9 + 0.446)/4.0;
}
Y_H = 1.0 - Y_He*4.0 - ZMetal;
*pY_H = Y_H;
*pY_He = Y_He;
*pY_eMax = Y_H+ Y_He*2; /* Ignoring any electrons from metals */
}
void CoolPARTICLEtoPERBARYON(COOL *cl, PERBARYON *Y, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
Y->HI = cp->f_HI*Y_H;
Y->H2 = cp->f_H2/2.0*Y_H; /*Number of H2 molecules total from fraction of H atoms in H2*/
Y->HII = Y_H - Y->HI - 2.0*Y->H2;
Y->HeI = cp->f_HeI*Y_He;
Y->HeII = cp->f_HeII*Y_He;
Y->HeIII = Y_He - Y->HeI - Y->HeII;
Y->e = Y->HII + Y->HeII + 2*Y->HeIII;
Y->Total = Y->e + Y_H + Y_He + ZMetal/MU_METAL;
}
void CoolPERBARYONtoPARTICLE(COOL *cl, PERBARYON *Y, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
cp->f_HI = Y->HI/Y_H;
cp->f_H2 = 2.0*(Y->H2)/Y_H; /*Fraction of H atoms that are in H2 from number of H2 molecules CC*/
cp->f_HeI = Y->HeI/Y_He;
cp->f_HeII = Y->HeII/Y_He;
}
/*Returns baryonic fraction for a given species*/
FLOAT COOL_ARRAY0(COOL *cl, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
return (cp->f_HI*Y_H);
}
FLOAT COOL_ARRAY1(COOL *cl, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
return (cp->f_HeI*Y_He);
}
FLOAT COOL_ARRAY2(COOL *cl, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
return (cp->f_HeII*Y_He);
}
FLOAT COOL_ARRAY3(COOL *cl, COOLPARTICLE *cp, double ZMetal) {
double Y_H, Y_He, Y_eMax;
clSetAbundanceTotals(cl,ZMetal,&Y_H,&Y_He,&Y_eMax);
return (cp->f_H2/2.0*Y_H);
}
#define COOL_ARRAY4_EXT "dummy"
#define COOL_ARRAY4(x,y,z) 0
#define COOL_IN_ARRAY4(w,x,y,z)
#define COOL_ARRAY5_EXT "dummy"
#define COOL_ARRAY5(x,y,z) 0
#define COOL_IN_ARRAY5(w,x,y,z)
#define COOL_ARRAY6_EXT "dummy"
#define COOL_ARRAY6(x,y,z) 0
#define COOL_IN_ARRAY6(w,x,y,z)
#define COOL_ARRAY7_EXT "dummy"
#define COOL_ARRAY7(x,y,z) 0
#define COOL_IN_ARRAY7(w,x,y,z)
#define COOL_ARRAY8_EXT "dummy"
#define COOL_ARRAY8(x,y,z) 0
#define COOL_IN_ARRAY8(w,x,y,z)
#define COOL_ARRAY9_EXT "dummy"
#define COOL_ARRAY9(x,y,z) 0
#define COOL_IN_ARRAY9(w,x,y,z)
#define COOL_ARRAY10_EXT "dummy"
#define COOL_ARRAY10(x,y,z) 0
#define COOL_IN_ARRAY10(w,x,y,z)
#define COOL_ARRAY11_EXT "dummy"
#define COOL_ARRAY11(x,y,z) 0
#define COOL_IN_ARRAY11(w,x,y,z)
#define COOL_ARRAY12_EXT "dummy"
#define COOL_ARRAY12(x,y,z) 0
#define COOL_IN_ARRAY12(w,x,y,z)
#define COOL_ARRAY13_EXT "dummy"
#define COOL_ARRAY13(x,y,z) 0
#define COOL_IN_ARRAY13(w,x,y,z)
#define COOL_ARRAY14_EXT "dummy"
#define COOL_ARRAY14(x,y,z) 0
#define COOL_IN_ARRAY14(w,x,y,z)
#define COOL_ARRAY15_EXT "dummy"
#define COOL_ARRAY15(x,y,z) 0
#define COOL_IN_ARRAY15(w,x,y,z)
void clReadMetalTable(COOL *cl, COOLPARAM clParam)
{
int i,j,k,l, nt, nnH, nz;
double dt, dnH, dz, tminlog, tmaxlog, nHminlog, nHmaxlog,zmin,zmax;
FILE *fp;
XDR xdrs;
fp = fopen(clParam.CoolInFile, "r");
// fp = fopen("cooltable_xdr", "r");
assert(fp != NULL);
fscanf(fp, "%d %lf %lf %lf \n",&nz, &zmin, &zmax, &dz);
fscanf(fp, "%d %lf %lf %lf \n",&nnH, &nHminlog ,&nHmaxlog,&dnH);
fscanf(fp, "%d %lf %lf %lf \n",&nt, &tminlog, &tmaxlog, &dt);
cl->bMetal = clParam.bMetal;
cl->nzMetalTable = nz;
cl->nnHMetalTable = nnH;
cl->nTMetalTable = nt;
cl->MetalTMin = pow(10.0, tminlog);
cl->MetalTMax = pow(10.0, tmaxlog);
cl->MetalnHMin = pow(10.0, nHminlog); /* Note: here nH means mass density */
cl->MetalnHMax = pow(10.0, nHmaxlog);
cl->MetalzMin = zmin;
cl->MetalzMax = zmax;
cl->MetalTlogMin = tminlog;
cl->MetalTlogMax = tmaxlog;
cl->MetalnHlogMin = nHminlog;
cl->MetalnHlogMax = nHmaxlog;
cl->rDeltaTlog = 1./dt;
cl->rDeltanHlog = 1./dnH;
cl->rDeltaz = 1./dz;
cl->nUV = nz -1;
cl->UV = (UVSPECTRUM *)malloc(sizeof(UVSPECTRUM)*(cl->nUV));
/* 3D table, f = f(T, rho, z), in redshift the order is from high z-->low z,
density and temperature are in the order low --> high */
cl->MetalCoolln = (float ***)malloc(nz*sizeof(float ***));
for (i=0; i<=nz-1;i++){
cl->MetalCoolln[i] = (float **) malloc(nnH*sizeof(float **));
}
for (i=0; i<=nz-1;i++){
for(j=0; j<=nnH-1; j++){
cl->MetalCoolln[i][j] = (float *) malloc(nt*sizeof(float *));
}
}
cl->MetalHeatln = (float ***)malloc(nz*sizeof(float ***));
for (i=0; i<=nz-1;i++){
cl->MetalHeatln[i] = (float **) malloc(nnH*sizeof(float **));
}
for (i=0; i<=nz-1;i++){
for(j=0; j<=nnH-1; j++){
cl->MetalHeatln[i][j] = (float *) malloc(nt*sizeof(float *));
}
}
/* for(;;){
fscanf(fp, "%s", comments);
if (comments[0] eq "$") break;
} */
xdrstdio_create(&xdrs, fp, XDR_DECODE);
for(j=0; j<=nnH-1; j++){
for(k=0; k<=nt-1; k++){
xdr_float(&xdrs, &cl->MetalCoolln[0][j][k]); /*redshift: i, density: j, temperature: k */
xdr_float(&xdrs, &cl->MetalHeatln[0][j][k]);
}
}
for(i=nz-2, l=0; i>=0; i--,l++){
xdr_double(&xdrs, &cl->UV[l].zTime);
xdr_double(&xdrs, &cl->UV[l].Rate_Phot_HI);
xdr_double(&xdrs, &cl->UV[l].Rate_Phot_HeI);
xdr_double(&xdrs, &cl->UV[l].Rate_Phot_HeII);
xdr_double(&xdrs, &cl->UV[l].Heat_Phot_HI);
xdr_double(&xdrs, &cl->UV[l].Heat_Phot_HeI);
xdr_double(&xdrs, &cl->UV[l].Heat_Phot_HeII);
for(j=0; j<=nnH-1; j++){
for(k=0; k<=nt-1; k++){
xdr_float(&xdrs, &cl->MetalCoolln[l+1][j][k]);
xdr_float(&xdrs, &cl->MetalHeatln[l+1][j][k]);
}
}
}
fclose(fp);
return;
}
void clInitUV(COOL *cl, int nTableColumns, int nTableRows, double *dTableData )
{
int i;
printf("call clInit UV...\n");
assert(cl!=NULL);
assert(cl->UV == NULL);
assert(nTableColumns == 7);
cl->nUV = nTableRows;
cl->UV = (UVSPECTRUM *) malloc(nTableRows*sizeof(UVSPECTRUM));
assert(cl->UV!=NULL);
for (i=0;i<nTableRows;i++) {
(cl->UV)[i].zTime = dTableData[i*nTableColumns];
(cl->UV)[i].Rate_Phot_HI = dTableData[i*nTableColumns+1];
(cl->UV)[i].Rate_Phot_HeI = dTableData[i*nTableColumns+2];
(cl->UV)[i].Rate_Phot_HeII = dTableData[i*nTableColumns+3];
(cl->UV)[i].Rate_Phot_H2_cosmo = dTableData[i*nTableColumns+4];
(cl->UV)[i].Heat_Phot_HI = dTableData[i*nTableColumns+5];
(cl->UV)[i].Heat_Phot_HeI = dTableData[i*nTableColumns+6];
(cl->UV)[i].Heat_Phot_HeII = dTableData[i*nTableColumns+7];
(cl->UV)[i].Heat_Phot_H2 = dTableData[i*nTableColumns+8];
/* Make sure the heating is in units of ergs per ionization */
assert( (cl->UV)[i].Heat_Phot_HI>1e-15 && (cl->UV)[i].Heat_Phot_HI<1e-10);
if (i) assert (((cl->UV)[i-1].zTime > (cl->UV)[i].zTime
&& !cl->bUVTableUsesTime) ||
((cl->UV)[i-1].zTime < (cl->UV)[i].zTime
&& cl->bUVTableUsesTime));
}
}
void clInitRatesTable( COOL *cl, double TMin, double TMax, int nTable ) {
/*-----------------------------------------------------------------
* Function of ln(T) tables:
* cl assumed to be predefined (allocated)
* A table spacing of 1.23e-3 in log(T) eg. nTable=15001 10->1e9 K
* Maximum 1% errors except at discontinuities in the functions.
* Storage for such a table is (15001*15*4|8b) => 0.9|1.7 Mb
*-----------------------------------------------------------------*/
int i,j;
double DeltaTln, Tln, T;
assert(cl!=NULL);
assert(cl->RT == NULL);
cl->R.Cool_Coll_HI = CL_eHI*CL_B_gm;
cl->R.Cool_Coll_HeI = CL_eHeI*CL_B_gm;
cl->R.Cool_Coll_HeII = CL_eHeII*CL_B_gm;
cl->R.Cool_Diel_HeII = (CL_E2HeII+CL_eHeI)*CL_B_gm;
cl->R.Cool_Coll_H2 = CL_eH2*CL_B_gm;
cl->nTable = nTable;
cl->TMin = TMin;
cl->TMax = TMax;
cl->TlnMin = log( TMin );
cl->TlnMax = log( TMax );
DeltaTln = ( cl->TlnMax-cl->TlnMin )/( nTable - 1 );
cl->rDeltaTln = 1./DeltaTln;
cl->RT = (RATES_T *) malloc( nTable * sizeof(RATES_T) * TABLEFACTOR );
assert(cl->RT != NULL);
for ( j=0; j<nTable; j++ ) {
Tln = cl->TlnMin + DeltaTln*j;
T=exp(Tln);
i = j*TABLEFACTOR;
(cl->RT+i)->Rate_Coll_HI = clRateCollHI( T );
if ( (cl->RT+i)->Rate_Coll_HI < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_HI = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_HeI = clRateCollHeI( T );
if ( (cl->RT+i)->Rate_Coll_HeI < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_HeI = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_HeII = clRateCollHeII( T );
if ( (cl->RT+i)->Rate_Coll_HeII < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_HeII = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_e_H2 = clRateColl_e_H2( T );
if ( (cl->RT+i)->Rate_Coll_e_H2 < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_e_H2 = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_HI_H2 = clRateColl_HI_H2( T );
if ( (cl->RT+i)->Rate_Coll_HI_H2 < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_HI_H2 = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_H2_H2 = clRateColl_H2_H2( T );
if ( (cl->RT+i)->Rate_Coll_H2_H2 < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_H2_H2 = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_Hm_e = clRateColl_Hm_e( T );/*gas phase form of H2 */
if ( (cl->RT+i)->Rate_Coll_Hm_e < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_Hm_e = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_Hm_HII = clRateColl_Hm_HII( T );/*gas phase form of H2 */
if ( (cl->RT+i)->Rate_Coll_Hm_HII < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_Hm_HII = CL_RT_MIN;
(cl->RT+i)->Rate_Coll_HII_H2 = clRateColl_HII_H2( T );/*gas phase form of H2 */
if ( (cl->RT+i)->Rate_Coll_HII_H2 < CL_RT_MIN ) (cl->RT+i)->Rate_Coll_HII_H2 = CL_RT_MIN;
(cl->RT+i)->Rate_HI_e = clRateHI_e( T );/*gas phase form of H2 */
if ( (cl->RT+i)->Rate_HI_e < CL_RT_MIN ) (cl->RT+i)->Rate_HI_e = CL_RT_MIN;
(cl->RT+i)->Rate_HI_Hm = clRateHI_Hm( T );/*gas phase form of H2 */
if ( (cl->RT+i)->Rate_HI_Hm < CL_RT_MIN ) (cl->RT+i)->Rate_HI_Hm = CL_RT_MIN;
(cl->RT+i)->Rate_Radr_HII = clRateRadrHII( T );
(cl->RT+i)->Rate_Radr_HeII = clRateRadrHeII( T );
(cl->RT+i)->Rate_Radr_HeIII = clRateRadrHeIII( T );
(cl->RT+i)->Rate_Diel_HeII = clRateDielHeII( T );
(cl->RT+i)->Rate_Chtr_HeII = clRateChtrHeII( T );
(cl->RT+i)->Cool_Brem_1 = clCoolBrem1( T );
(cl->RT+i)->Cool_Brem_2 = clCoolBrem2( T );
(cl->RT+i)->Cool_Radr_HII = clCoolRadrHII( T );
(cl->RT+i)->Cool_Radr_HeII = clCoolRadrHeII( T );
(cl->RT+i)->Cool_Radr_HeIII = clCoolRadrHeIII( T );
(cl->RT+i)->Cool_Line_H2_H = clCoolLineH2_HI( T ); /*H2 Rot-Vib transitions out of collisionally-induced excited states. H2-H collisions*/
(cl->RT+i)->Cool_Line_H2_H2 = clCoolLineH2_H2( T ); /*H2 Rot-Vib transitions out of collisionally-induced excited states. H2-H2 collisions*/
(cl->RT+i)->Cool_Line_H2_He = clCoolLineH2_He( T ); /*H2 Rot-Vib transitions out of collisionally-induced excited states. H2-He collisions*/
(cl->RT+i)->Cool_Line_H2_e = clCoolLineH2_e( T ); /*H2 Rot-Vib transitions out of collisionally-induced excited states. H2-e collisions*/
(cl->RT+i)->Cool_Line_H2_HII = clCoolLineH2_HII( T ); /*H2 Rot-Vib transitions out of collisionally-induced excited states. H2-p collisions*/
(cl->RT+i)->Cool_Line_HI = clCoolLineHI( T );
(cl->RT+i)->Cool_Line_HeI = clCoolLineHeI( T );
(cl->RT+i)->Cool_Line_HeII = clCoolLineHeII( T );
(cl->RT+i)->Cool_LowT = clCoolLowT( T );
}
}
void clRatesTableError( COOL *cl ) {
/* This estimates the error for a table of half the size of
* the current one.
* The collisional, dielectric and line cooling rates can all go to
* zero (minimum) and will even for double precision (exp(-1e5/T))
* the critical values that must never get down to zero
* are the radiative recombination rates: Check this before using
* CL_RT_FLOAT float for the table.
*/
int i,j,ierr[15];
double min[15],max[15];
double err,maxerr[15];
CL_RT_FLOAT *p,*p1,*p2;
for (j=0;j<15;j++) {
maxerr[j]=0.0;
min[j]=1e300;
max[j]=-1e300;
}
for ( i=1; i<cl->nTable-1; i+=2 ) {
p = (CL_RT_FLOAT *) &((cl->RT+i)->Rate_Coll_HI);
p1 = (CL_RT_FLOAT *) &((cl->RT+i-1)->Rate_Coll_HI);
p2 = (CL_RT_FLOAT *) &((cl->RT+i+1)->Rate_Coll_HI);
if (i==10001) {
printf(" Comp %i %i %e %e\n",i,(int) sizeof(CL_RT_FLOAT),p[1],(cl->RT+i)->Rate_Coll_HeI);
}
for (j=0;j<15;j++) {
if (p[j] < min[j]) min[j]=p[j];
if (p[j] > max[j]) max[j]=p[j];
err= fabs((0.5*(p1[j]+p2[j])-p[j])/(p[j]+1e-30));
if (err > maxerr[j] ) {
maxerr[j]=err;
ierr[j]=i;
}
}
}
for (j=0;j<15;j++) {
printf("Col %i Max Error: %e at T=%e dlnT=%e (min %e max %e)\n",j,
maxerr[j],exp(cl->TlnMin+ierr[j]/cl->rDeltaTln),2/cl->rDeltaTln,min[j],max[j]);
}
}
#define CL_Ccomp0 0.565e-9
#define CL_Tcmb0 2.735
#define CL_Ccomp (CL_Ccomp0*CL_Tcmb0)
void clRatesRedshift( COOL *cl, double zIn, double dTimeIn ) {
int i;
double xx;
double zTime;
UVSPECTRUM *UV,*UV0;
double Y_H, Y_He, Y_eMax;
double ten = 10.0,output, expon;
/* printf("Redshift: %f \n", zIn); */
/* cl->z = 0.0; */
cl->z = zIn;
cl->dTime = dTimeIn;
cl->dComovingGmPerCcUnit = cl->dGmPerCcUnit*pow(1.+zIn,3.);
cl->dExpand = 1.0/(1.0+zIn);
cl->R.Cool_Comp = pow((1+zIn)*CL_Ccomp,4.0)*CL_B_gm;
cl->R.Tcmb = CL_Tcmb0*(1+zIn);
clSetAbundanceTotals(cl,0.0,&Y_H,&Y_He,&Y_eMax); /* Hack to estimate Y_H */
cl->R.Cool_LowTFactor = (cl->bLowTCool ? CL_B_gm*Y_H*Y_H/0.001 : 0 );
/* Photo-Ionization rates */
UV = cl->UV;
if (cl->bUV) {
assert( UV != NULL );
if (cl->bUVTableUsesTime) {
/*
** Table in order of increasing time
*/
zTime = dTimeIn;
for ( i=0; i < cl->nUV && zTime >= UV->zTime ; i++,UV++ );
}
else {
/*
** Table in order of high to low redshift
*/
zTime = zIn;
for ( i=0; i < cl->nUV && zTime <= UV->zTime ; i++,UV++ );
}
}
if (!cl->bUV || i==0) {
cl->R.Rate_Phot_HI = CL_RT_MIN;
cl->R.Rate_Phot_HeI = CL_RT_MIN;
cl->R.Rate_Phot_HeII = CL_RT_MIN;
cl->R.Rate_Phot_H2_cosmo = CL_RT_MIN;
cl->R.Heat_Phot_HI = 0.0;
cl->R.Heat_Phot_HeI = 0.0;
cl->R.Heat_Phot_HeII = 0.0;
cl->R.Heat_Phot_H2 = 0.0;
return;
}
UV0=UV-1;
if (i == cl->nUV ) {
cl->R.Rate_Phot_HI = UV0->Rate_Phot_HI;
cl->R.Rate_Phot_HeI = UV0->Rate_Phot_HeI;
cl->R.Rate_Phot_HeII = UV0->Rate_Phot_HeII;
expon = 0.90632725*zTime - 0.16790918*zTime*zTime + 0.010241484*zTime*zTime*zTime - 12.518825;
output = pow(ten,expon); /*haardt_madau gal+quasar, z = 0*/
cl->R.Rate_Phot_H2_cosmo = output;
cl->R.Heat_Phot_HI = UV0->Heat_Phot_HI*CL_B_gm;
cl->R.Heat_Phot_HeI = UV0->Heat_Phot_HeI*CL_B_gm;
cl->R.Heat_Phot_HeII = UV0->Heat_Phot_HeII*CL_B_gm;
cl->R.Heat_Phot_H2 = 6.4e-13*CL_B_gm;
}
else {
if (cl->bUVTableLinear) { /* use Linear interpolation */
xx = (zTime - UV0->zTime)/(UV->zTime - UV0->zTime);
cl->R.Rate_Phot_HI = UV0->Rate_Phot_HI*(1-xx)+UV->Rate_Phot_HI*xx;
cl->R.Rate_Phot_HeI = UV0->Rate_Phot_HeI*(1-xx)+UV->Rate_Phot_HeI*xx;
cl->R.Rate_Phot_HeII = UV0->Rate_Phot_HeII*(1-xx)+UV->Rate_Phot_HeII*xx;
expon = 0.90632725*zTime - 0.16790918*zTime*zTime + 0.010241484*zTime*zTime*zTime - 12.518825;
cl->R.Rate_Phot_H2_cosmo = pow(ten,expon); /*haardt_madau gal+quasar, z = 0*/
cl->R.Heat_Phot_HI = (UV0->Heat_Phot_HI*(1-xx)+UV->Heat_Phot_HI*xx)*CL_B_gm;
cl->R.Heat_Phot_HeI = (UV0->Heat_Phot_HeI*(1-xx)+UV->Heat_Phot_HeI*xx)*CL_B_gm;
cl->R.Heat_Phot_HeII = (UV0->Heat_Phot_HeII*(1-xx)+UV->Heat_Phot_HeII*xx)*CL_B_gm;
cl->R.Heat_Phot_H2 = 6.4e-13*CL_B_gm;
}
else { /* use Log interpolation with 1+zTime */
xx = log((1+zTime)/(1+UV0->zTime))/log((1+UV->zTime)/(1+UV0->zTime));
cl->R.Rate_Phot_HI = pow(UV0->Rate_Phot_HI,1-xx)*pow(UV->Rate_Phot_HI,xx);
cl->R.Rate_Phot_HeI = pow(UV0->Rate_Phot_HeI,1-xx)*pow(UV->Rate_Phot_HeI,xx);
cl->R.Rate_Phot_HeII = pow(UV0->Rate_Phot_HeII,1-xx)*pow(UV->Rate_Phot_HeII,xx);
expon = 0.90632725*zTime - 0.16790918*zTime*zTime + 0.010241484*zTime*zTime*zTime - 12.518825;
cl->R.Rate_Phot_H2_cosmo = pow(ten,expon); /*haardt_madau gal+quasar, z = 0*/
cl->R.Heat_Phot_HI = pow(UV0->Heat_Phot_HI,1-xx)*pow(UV->Heat_Phot_HI,xx)*CL_B_gm;
cl->R.Heat_Phot_HeI = pow(UV0->Heat_Phot_HeI,1-xx)*pow(UV->Heat_Phot_HeI,xx)*CL_B_gm;
cl->R.Heat_Phot_HeII = pow(UV0->Heat_Phot_HeII,1-xx)*pow(UV->Heat_Phot_HeII,xx)*CL_B_gm;
cl->R.Heat_Phot_H2 = 6.4e-13*CL_B_gm;
}
}
if (cl->R.Rate_Phot_HI < CL_RT_MIN) cl->R.Rate_Phot_HI = CL_RT_MIN;
if (cl->R.Rate_Phot_HeI < CL_RT_MIN) cl->R.Rate_Phot_HeI = CL_RT_MIN;
if (cl->R.Rate_Phot_HeII < CL_RT_MIN) cl->R.Rate_Phot_HeII = CL_RT_MIN;
if (cl->R.Rate_Phot_H2_cosmo < CL_RT_MIN) cl->R.Rate_Phot_H2_cosmo = CL_RT_MIN;
return;
}
double AP_log_den_mp_percm3[] = { -10.25, -9.75, -9.25, -8.75, -8.25, -7.75, -7.25,
-6.75, -6.25, -5.75, -5.25, -4.75, -4.25, -3.75,
-3.25, -2.75, -2.25,
-1.75, -1.25, -0.75, -0.25, 0.25, 0.75, 1.25, 1.75, 2.25 };
double AP_Gamma_HI_factor[] = { 0.99805271764596307, 0.99877911567687988, 0.99589340865612034,
0.99562060764857702, 0.99165170359332663, 0.9900889877822455,
0.98483276828954668, 0.97387675312245325, 0.97885673164000397,
0.98356305803821331, 0.96655786672182487, 0.9634906824933207,
0.95031917373653985, 0.87967606627349137, 0.79917533618355074,
0.61276011113763151, 0.16315185162187529, 0.02493663181368239,
0.0044013580765645335, 0.00024172553511936628, 1.9576102058649783e-10,
0.0, 0.0, 0.0, 0.0, 0.0 };
double AP_Gamma_HeI_factor[] = { 0.99284882980782224, 0.9946618686265097, 0.98641914356740497,
0.98867015777574848, 0.96519214493135597, 0.97188336387980656,
0.97529866247535113 , 0.97412477991428936, 0.97904139838765991,
0.98368372768570034, 0.96677432842215549, 0.96392622083382651,
0.95145730833093178, 0.88213871255482879, 0.80512823597731886,
0.62474472739578646, 0.17222786134467002, 0.025861959933038869,
0.0045265030237581529, 0.00024724339438128221, 1.3040144591221284e-08,
0.0, 0.0, 0.0, 0.0, 0.0};
double AP_Gamma_HeII_factor[] = { 0.97990208047216765, 0.98606251822654412,
0.97657215632444849, 0.97274858503068629, 0.97416108746560681,
0.97716929017896703, 0.97743607605974214, 0.97555305775319012,
0.97874250764784809, 0.97849791914637996, 0.95135572977973504,
0.92948461312852582, 0.89242272355549912, 0.79325512242742746 ,
0.6683745597121028, 0.51605924897038324, 0.1840253816147828,
0.035905775349044489, 0.0045537756654992923, 0.00035933897136804514,
1.2294426136470751e-6, 0.0, 0.0, 0.0, 0.0, 0.0 };
double AP_Gamma_H2_factor[] = {1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0,
1.0, 1.0, 1.0, 1.0, 1.0};
/*This relates the density of the surroundings to the likelihood of photodissociation
we have another way to do shielding for H2 so this is filled in only to prevent a crash
It's also untested with H2
*/
void clRates( COOL *cl, RATE *Rate, double T, double rho, double ZMetal, double correL, double Rate_Phot_H2_stellar) {
double Tln;
if (T >= cl->TMax) T=cl->TMax*(1.0 - EPS);
if (T < cl->TMin) T=cl->TMin;
Tln = log(T); /* Deprecated but log's, sqrt's etc... used in raw rate functions */
Rate->T = T;
Rate->Tln = Tln;
Rate->Coll_HI = clRateCollHI( T );
Rate->Coll_HeI = clRateCollHeI( T );
Rate->Coll_HeII = clRateCollHeII( T );
Rate->Coll_e_H2 = clRateColl_e_H2( T );
Rate->Coll_HI_H2 = clRateColl_HI_H2( T );
Rate->Coll_H2_H2 = clRateColl_H2_H2( T );
Rate->Coll_Hm_e = clRateColl_Hm_e(T); /*gas phase formation of H2 */
Rate->Coll_Hm_HII = clRateColl_Hm_HII(T); /*gas phase formation of H2 */
Rate->Coll_HII_H2 = clRateColl_HII_H2(T); /*gas phase formation of H2 */
Rate->HI_e = clRateHI_e(T); /*gas phase formation of H2 */
Rate->HI_Hm = clRateHI_Hm(T); /*gas phase formation of H2 */
Rate->Radr_HII = clRateRadrHII( T );
Rate->Radr_HeII = clRateRadrHeII( T );
Rate->Diel_HeII = clRateDielHeII( T );
Rate->Chtr_HeII = clRateChtrHeII( T );
Rate->Totr_HeII = Rate->Radr_HeII + Rate->Diel_HeII + Rate->Chtr_HeII;
Rate->Radr_HeIII = clRateRadrHeIII( T );
Rate->DustForm_H2 = clRateDustFormH2( ZMetal, cl->dClump ); /*H2 Formation on dust*/
Rate->Phot_HI = cl->R.Rate_Phot_HI;
Rate->Phot_HeI = cl->R.Rate_Phot_HeI;
Rate->Phot_HeII = cl->R.Rate_Phot_HeII;
Rate->Phot_H2 = cl->R.Rate_Phot_H2_cosmo
;
Rate->CorreLength = correL * cl->dKpcUnit * 3.08568025e21 * cl->dExpand; /* From system units to cm*/
if (cl->bSelfShield) {
double logen_B;
logen_B = log10(rho*CL_B_gm);
if (logen_B > 2.2499) {
Rate->Phot_HI = 0;
Rate->Phot_HeI = 0;
Rate->Phot_HeII = 0;
}
else if (logen_B > -10.25) {
double x = (logen_B+10.25)*2.0;
int ix;
ix = floor(x);
x -= ix;
Rate->Phot_HI *= (AP_Gamma_HI_factor[ix]*(1-x)+AP_Gamma_HI_factor[ix+1]*x);
Rate->Phot_HeI *= (AP_Gamma_HeI_factor[ix]*(1-x)+AP_Gamma_HeI_factor[ix+1]*x);
Rate->Phot_HeII *= (AP_Gamma_HeII_factor[ix]*(1-x)+AP_Gamma_HeII_factor[ix+1]*x);
}
}
}
#define TABLEINTERPLIN( _rname ) (wTln0*RT0->_rname+wTln1*RT1->_rname)
void clRates_Table_Lin( COOL *cl, RATE *Rate, double T, double rho, double ZMetal, double correL, double Rate_Phot_H2_stellar) {
double Tln;
double xTln,wTln0,wTln1;/*,wTln0d,wTln1d;*/
RATES_T *RT0,*RT1; /**RT0d,*RT1d;*/
int iTln;
if (T >= cl->TMax) T=cl->TMax*(1.0 - EPS);
if (T < cl->TMin) T=cl->TMin;
Tln = log(T);
Rate->T = T;
Rate->Tln = Tln;
xTln = (Tln-cl->TlnMin)*cl->rDeltaTln;
iTln = xTln;
RT0 = (cl->RT+iTln*TABLEFACTOR);
RT1 = RT0+TABLEFACTOR;
xTln = xTln-iTln;
wTln1 = xTln;
wTln0 = 1-xTln;
Rate->Coll_HI = TABLEINTERPLIN( Rate_Coll_HI );
Rate->Coll_e_H2 = TABLEINTERPLIN( Rate_Coll_e_H2 );
Rate->Coll_HI_H2 = TABLEINTERPLIN( Rate_Coll_HI_H2 );
Rate->Coll_H2_H2 = TABLEINTERPLIN( Rate_Coll_H2_H2 );
Rate->Coll_HeI = TABLEINTERPLIN( Rate_Coll_HeI );
Rate->Coll_HeII = TABLEINTERPLIN( Rate_Coll_HeII );
Rate->Coll_Hm_e = TABLEINTERPLIN(Rate_Coll_Hm_e); /*gas phase form of H2 */
Rate->Coll_Hm_HII = TABLEINTERPLIN(Rate_Coll_Hm_HII); /*gas phase form of H2 */
Rate->Coll_HII_H2 = TABLEINTERPLIN(Rate_Coll_HII_H2);
Rate->HI_e = TABLEINTERPLIN(Rate_HI_e); /*gas phase form of H2 */
Rate->HI_Hm = TABLEINTERPLIN(Rate_HI_Hm); /*gas phase form of H2 */
Rate->Radr_HII = TABLEINTERPLIN( Rate_Radr_HII );
Rate->Radr_HeII = TABLEINTERPLIN( Rate_Radr_HeII );
Rate->Diel_HeII = TABLEINTERPLIN( Rate_Diel_HeII );
Rate->Chtr_HeII = TABLEINTERPLIN( Rate_Chtr_HeII );
Rate->Totr_HeII = Rate->Radr_HeII + Rate->Diel_HeII + Rate->Chtr_HeII;
Rate->Radr_HeIII = TABLEINTERPLIN( Rate_Radr_HeIII );
Rate->DustForm_H2 = clRateDustFormH2( ZMetal, cl->dClump );
Rate->Phot_HI = cl->R.Rate_Phot_HI;
Rate->Phot_HeI = cl->R.Rate_Phot_HeI;
Rate->Phot_HeII = cl->R.Rate_Phot_HeII;
Rate->Phot_H2 = cl->R.Rate_Phot_H2_cosmo
;
Rate->CorreLength = correL * cl->dKpcUnit * 3.08568025e21 * cl->dExpand; /* From system units to cm*/
if (cl->bSelfShield) {
double logen_B;
logen_B = log10(rho*CL_B_gm);
if (logen_B > 2.2499) {
Rate->Phot_HI = 0;
Rate->Phot_HeI = 0;
Rate->Phot_HeII = 0;
Rate->Phot_H2 = 0;
}
else if (logen_B > -10.25) {
double x = (logen_B+10.25)*2.0;
int ix;
ix = floor(x);
x -= ix;
Rate->Phot_HI *= (AP_Gamma_HI_factor[ix]*(1-x)+AP_Gamma_HI_factor[ix+1]*x);
Rate->Phot_HeI *= (AP_Gamma_HeI_factor[ix]*(1-x)+AP_Gamma_HeI_factor[ix+1]*x);
Rate->Phot_HeII *= (AP_Gamma_HeII_factor[ix]*(1-x)+AP_Gamma_HeII_factor[ix+1]*x);
Rate->Phot_H2 *= (AP_Gamma_H2_factor[ix]*(1-x)+AP_Gamma_H2_factor[ix+1]*x);
}
}
}
void clRates_Table( COOL *cl, RATE *Rate, double T, double rho, double ZMetal, double correL, double Rate_Phot_H2_stellar) {
double Tln;
double xTln,wTln0,wTln1;/*,wTln0d,wTln1d;*/
RATES_T *RT0,*RT1;/*,*RT0d,*RT1d;*/
int iTln;
#ifdef TESTRATE
RATE test;
clRates( cl, &test, T, rho, ZMetal, correL, Rate_Phot_H2_stellar);
#endif
if (T >= cl->TMax) T=cl->TMax*(1.0 - EPS);
if (T < cl->TMin) T=cl->TMin;
Tln = log(T);
Rate->T = T;
Rate->Tln = Tln;
xTln = (Tln-cl->TlnMin)*cl->rDeltaTln;
iTln = xTln;
RT0 = (cl->RT+iTln*TABLEFACTOR);
RT1 = RT0+TABLEFACTOR;
xTln = xTln-iTln;
wTln1 = xTln;
wTln0 = 1-xTln;
Rate->Coll_HI = TABLEINTERP( Rate_Coll_HI );
Rate->Coll_HeI = TABLEINTERP( Rate_Coll_HeI );
Rate->Coll_HeII = TABLEINTERP( Rate_Coll_HeII );
Rate->Coll_e_H2 = TABLEINTERP( Rate_Coll_e_H2 );
Rate->Coll_HI_H2 = TABLEINTERP( Rate_Coll_HI_H2 );
Rate->Coll_H2_H2 = TABLEINTERP( Rate_Coll_H2_H2 );
Rate->Coll_HII_H2 = TABLEINTERPLIN(Rate_Coll_HII_H2);
Rate->Coll_Hm_e = TABLEINTERPLIN(Rate_Coll_Hm_e); /*gas phase form of H2 */
Rate->Coll_Hm_HII = TABLEINTERPLIN(Rate_Coll_Hm_HII); /*gas phase form of H2 */
Rate->HI_e = TABLEINTERPLIN(Rate_HI_e); /*gas phase form of H2 */
Rate->HI_Hm = TABLEINTERPLIN(Rate_HI_Hm); /*gas phase form of H2 */
Rate->Radr_HII = TABLEINTERP( Rate_Radr_HII );
Rate->Radr_HeII = TABLEINTERP( Rate_Radr_HeII );
Rate->Diel_HeII = TABLEINTERP( Rate_Diel_HeII );
Rate->Chtr_HeII = TABLEINTERP( Rate_Chtr_HeII );
Rate->Totr_HeII = Rate->Radr_HeII + Rate->Diel_HeII + Rate->Chtr_HeII;
Rate->Radr_HeIII = TABLEINTERP( Rate_Radr_HeIII );
Rate->DustForm_H2 = clRateDustFormH2( ZMetal, cl->dClump );
Rate->Phot_HI = cl->R.Rate_Phot_HI;
Rate->Phot_HeI = cl->R.Rate_Phot_HeI;
Rate->Phot_HeII = cl->R.Rate_Phot_HeII;
Rate->Phot_H2 = cl->R.Rate_Phot_H2_cosmo
;
Rate->CorreLength = correL * cl->dKpcUnit*3.08568025e21 * cl->dExpand; /* From system units to cm*/
if (cl->bSelfShield) {
double logen_B;
logen_B = log10(rho*CL_B_gm);
if (logen_B > 2.2499) {
Rate->Phot_HI = 0;
Rate->Phot_HeI = 0;
Rate->Phot_HeII = 0;
Rate->Phot_H2 = 0;
}
else if (logen_B > -10.25) {
double x = (logen_B+10.25)*2.0;
int ix;
ix = floor(x);
x -= ix;
Rate->Phot_HI *= (AP_Gamma_HI_factor[ix]*(1-x)+AP_Gamma_HI_factor[ix+1]*x);
Rate->Phot_HeI *= (AP_Gamma_HeI_factor[ix]*(1-x)+AP_Gamma_HeI_factor[ix+1]*x);
Rate->Phot_HeII *= (AP_Gamma_HeII_factor[ix]*(1-x)+AP_Gamma_HeII_factor[ix+1]*x);
Rate->Phot_H2 *= (AP_Gamma_H2_factor[ix]*(1-x)+AP_Gamma_H2_factor[ix+1]*x); /* Untested */
}
}
#ifdef TESTRATE
#define RATEVAR( _name ) (fabs((test._name - Rate->_name)/(test._name)) > TESTRATE)
if ( (T <.5e6 && rho > 1e3*1e-31) && ((RATEVAR(Coll_HI)) || RATEVAR(Coll_HeI) || RATEVAR(Coll_HeII) ||
RATEVAR(Radr_HII) || RATEVAR(Radr_HeII) || RATEVAR(Diel_HeII) || RATEVAR(Totr_HeII) || RATEVAR(Radr_HeIII) ||
RATEVAR(Phot_HI) || RATEVAR(Phot_HeI) || RATEVAR(Phot_HeII) || RATEVAR(Phot_H2) )) {
printf("Bad interpolated rates at T=%e rho=%e\n",T,rho);
printf("%12f %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %e12 %e12 %e12\n",
T, test.Coll_HI, test.Coll_e_H2, test.Coll_HI_H2, test.Coll_H2_H2, test.Coll_HeI, test.Coll_HeII, test.Radr_HII, test.Radr_HeII, test.Diel_HeII, test.Totr_HeII, test.Radr_HeIII, test.Phot_HI, test.Phot_H2, test.Phot_HeI, test.Phot_HeII);
printf("%12f %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e \n",
T, Rate->Coll_HI, Rate->Coll_e_H2, Rate->Coll_HI_H2, Rate->Coll_H2_H2, Rate->Coll_HeI, Rate->Coll_HeII, Rate->Radr_HII, Rate->Radr_HeII, Rate->Diel_HeII, Rate->Chtr_HeII, Rate->Totr_HeII, Rate->Radr_HeIII, Rate->Phot_HI, Rate->Phot_HI, Rate->Phot_HeI, Rate->Phot_HeII );
clRates_Table_Lin( cl, &test, T, rho, ZMetal, correL, Rate_Phot_H2_stellar);
printf("%12f %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %e12 %12e %e12\n", T, test.Coll_HI, test.Coll_e_H2, test.Coll_HI_H2, test.Coll_H2_H2, test.Coll_HeI, test.Coll_HeII, test.Radr_HII, test.Radr_HeII, test.Diel_HeII, test.Totr_HeII, test.Radr_HeIII, test.Phot_HI, test.Phot_H2, test.Phot_HeI, test.Phot_HeII );
clRates( cl, &test, T*1.00123, rho, ZMetal, correL, Rate_Phot_H2_stellar);
printf("%12f %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %e12 %12e %e12\n", T*1.00123, test.Coll_HI, test.Coll_e_H2, test.Coll_HI_H2, test.Coll_H2_H2, test.Coll_HeI, test.Coll_HeII, test.Radr_HII, test.Radr_HeII, test.Diel_HeII, test.Totr_HeII, test.Radr_HeIII, test.Phot_HI, test.Phot_H2, test.Phot_HeI, test.Phot_HeII );
clRates( cl, &test, T/1.00123, rho, ZMetal, correL, Rate_Phot_H2_stellar);
printf("%12f %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %12e %e12 %12e %12e\n", T/1.00123, test.Coll_HI, test.Coll_e_H2, test.Coll_HI_H2, test.Coll_H2_H2, test.Coll_HeI, test.Coll_HeII, test.Radr_HII, test.Radr_HeII, test.Diel_HeII, test.Totr_HeII, test.Radr_HeIII, test.Phot_HI, test.Phot_H2, test.Phot_HeI, test.Phot_HeII );
assert(0);
}
#endif
}
void clRateMetalTable(COOL *cl, RATE *Rate, double T, double rho, double Y_H, double ZMetal)
{
double tempT, tempnH,tempz, nH;
double Tlog, nHlog;
double xTlog, wTlog0, wTlog1, xz, wz0, wz1, xnHlog, wnHlog0, wnHlog1;
int iTlog, iz, inHlog;
double Cool000, Cool010, Cool100, Cool110, Cool001, Cool011, Cool101, Cool111;
double Cool00, Cool01, Cool10, Cool11, Cool0, Cool1, Cool;
double Heat000, Heat010, Heat100, Heat110, Heat001, Heat011, Heat101, Heat111;
double Heat00, Heat01, Heat10, Heat11, Heat0, Heat1, Heat;
if(!cl->bMetal) {
Rate->Cool_Metal = 0.0;
Rate->Heat_Metal = 0.0;
return;
}
nH = rho*Y_H/M_H;
tempT = T;
tempnH = nH;
tempz = cl->z;
if (T >= cl->MetalTMax) tempT = cl->MetalTMax*(1.0-EPS);
if (T < cl->MetalTMin) tempT=cl->MetalTMin;
if (nH >= cl->MetalnHMax) tempnH = cl->MetalnHMax*(1.0-EPS);
if (nH < cl->MetalnHMin) tempnH = cl->MetalnHMin;
if (cl->z <= cl->MetalzMin) tempz = cl->MetalzMin+EPS;
/* if redshift is too high or no UV, use the no UV metal cooling table*/
if (cl->z > cl->MetalzMax || !cl->bUV) tempz = cl->MetalzMax;
Tlog = log10(tempT);
nHlog = log10(tempnH);
xz = cl->nzMetalTable -1 - (tempz - cl->MetalzMin)*cl->rDeltaz;
iz = xz;
xTlog = (Tlog - cl->MetalTlogMin)*cl->rDeltaTlog;
assert(xTlog >= 0.0);
iTlog = xTlog;
xnHlog = (nHlog - cl->MetalnHlogMin)*cl->rDeltanHlog;
inHlog = xnHlog;
if (inHlog == cl->nnHMetalTable - 1) inHlog == cl->nnHMetalTable - 2; /*CC; To prevent running over the table. Should not be used*/
Cool000 = cl->MetalCoolln[iz][inHlog][iTlog];
Cool001 = cl->MetalCoolln[iz][inHlog][iTlog+1];
Cool010 = cl->MetalCoolln[iz][inHlog+1][iTlog];
Cool011 = cl->MetalCoolln[iz][inHlog+1][iTlog+1];
Cool100 = cl->MetalCoolln[iz+1][inHlog][iTlog];
Cool101 = cl->MetalCoolln[iz+1][inHlog][iTlog+1];
Cool110 = cl->MetalCoolln[iz+1][inHlog+1][iTlog];
Cool111 = cl->MetalCoolln[iz+1][inHlog+1][iTlog+1];
Heat000 = cl->MetalHeatln[iz][inHlog][iTlog];
Heat001 = cl->MetalHeatln[iz][inHlog][iTlog+1];
Heat010 = cl->MetalHeatln[iz][inHlog+1][iTlog];
Heat011 = cl->MetalHeatln[iz][inHlog+1][iTlog+1];
Heat100 = cl->MetalHeatln[iz+1][inHlog][iTlog];
Heat101 = cl->MetalHeatln[iz+1][inHlog][iTlog+1];
Heat110 = cl->MetalHeatln[iz+1][inHlog+1][iTlog];
Heat111 = cl->MetalHeatln[iz+1][inHlog+1][iTlog+1];
xz = xz - iz;
wz1 = xz;
wz0 = 1-xz;
Cool00 = wz0*Cool000 + wz1*Cool100;
Cool01 = wz0*Cool001 + wz1*Cool101;
Cool10 = wz0*Cool010 + wz1*Cool110;
Cool11 = wz0*Cool011 + wz1*Cool111;
Heat00 = wz0*Heat000 + wz1*Heat100;
Heat01 = wz0*Heat001 + wz1*Heat101;
Heat10 = wz0*Heat010 + wz1*Heat110;
Heat11 = wz0*Heat011 + wz1*Heat111;
xnHlog = xnHlog - inHlog;
wnHlog1 = xnHlog;
wnHlog0 = 1-xnHlog;
Cool0 = wnHlog0*Cool00 + wnHlog1*Cool10;