-
Notifications
You must be signed in to change notification settings - Fork 1
/
equiSources.f90
5060 lines (4259 loc) · 175 KB
/
equiSources.f90
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
!--- Alexei Razoumov, 2005-
!--- Main driver for the Fully Threaded Transport Engine: diffuse and point-source solvers.
!--- Diffuse solver: please cite Razoumov A. and Cardall C.Y., 2005, MNRAS, 362, 1413.
!--- Point-source solver: please cite Razoumov A. and Sommer-Larsen J., 2006, ApJ, 651, L89.
module localDefinitions
use definitions
integer, parameter :: nradius = 7, maxPixelLevel = 6
real(kind=RealKind), dimension(nradius) :: outputRadius = (/ 0.1, 0.3, 1., 3., 10., 30., 100. /) ! [kpc]
real(kind=RealKind), dimension(nradius) :: ndotRemaining, ndotBoundary, fraction
real(kind=RealKind) :: ndotDust, maxDepth
real(kind=RealKind), dimension(:), pointer :: lymanFraction
end module localDefinitions
program pointTransfer
use definitions
use localDefinitions
use stellarPopulationModule, only: stellarPopulation
use utilities, only: sortStellarParticles
use dust
use rotateIndicesModule
use transportRoutinesModule
implicit none
integer :: nlevels, icell, level, i, j, k, nx, ny, nz, itmp, i0, j0, k0, &
nside, iAngularLevel, imax, ir, itime, iStar, jStar, iSpectrum, iWavelength, &
iMetal, restart, &
nStars, nStarsSpecificAge, iradius, ienergy, ilambda, iBadSource, &
iRedshift, iDensity, iSpecificEnergy, iref, mode, nvariables, counter, &
reionizationModel, nxtransfer, nytransfer, nztransfer, imean
integer*8 :: longIntegerArgument, iray
real(kind=RealKind) :: tmp, tmp1, tmp2, tmp3, x0, y0, z0, xnew, ynew, znew, &
xa, xb, ya, yb, za, zb, phi, theta, weight, totalIntegral, &
ndot1, xneu, nh, nhe, coefSpectrum, coefMetal, &
depth1, depth2, depth3, depthDust, timeReadTable, &
upperAgeLimit, phiLarge, thetaLarge, phi1, theta1, &
xmean, ymean, zmean, localTemperature, maxOxygenAbundance, &
totalOxygen, rho, uvbCoefficient, oldNeutralFraction, actualRate, rateCoef
real(kind=RealKind) :: Jmean1, Jmean2, Jmean3, cellSizeAbsoluteUnits, Iin1, Iin2, Iin3, &
dpath, tau1, tau2, tau3, tmpemi1, tmpemi2, tmpemi3, nemi1, nemi2, nemi3, &
tmpabs1, tmpabs2, tmpabs3
type(readLevelType), dimension(:), pointer :: readLevel
real*4, dimension(maxNumberReadVariables) :: readArray
type(zoneType), pointer :: currentCell
type(nodeType), pointer :: currentNode, newNode
type(pixelType), target :: sphere
type(pixelType), pointer :: currentPixel, tmpPixel, leafPixel
type(starType), dimension(:), pointer :: star
type(starType), pointer :: currentStar, nextStar
type(pointType) :: startingPoint
type(segmentType), pointer :: newSegment
type(nodeType), target :: box
integer :: sfstart, sffinfo, sfselect, sfginfo, sfrdata, sfendacc, sfend, sfcreate, sfwdata
integer, dimension(3) :: baseGridSize
logical :: readNewSpectrum
character(60) :: synthesisFile, synthfile(5), sphDir, synthesisDir, grid, sources
character(90) :: scratchString, restartCellArrayName, filename
real(kind=RealKind), parameter :: alphaQuasar = 1.8
real(kind=RealKind), parameter :: alphaStellar = 5.
real(kind=RealKind), parameter :: contributionQuasar = 1.
real(kind=RealKind), parameter :: contributionStellar = 1.
integer, parameter :: stellarTransferThinUVB = 1, plotPDFs = 2, initialConfiguration = 3, &
printNumberOfCells = 4, noStarsThinUVB = 6, clumpingFactor = 7, &
bothStellarUVBTransfer = 8, UVBTransferOnly = 9
logical :: runStellarTransfer, runUVBTransfer
integer :: nxmap, nymap, kstart, kend, jcell, kcell
integer(1) :: izone
integer, dimension(2,2,2) :: is, js, ks
real*4, dimension(:,:), pointer :: mapArray
real(kind=RealKind), dimension(2) :: startPoint, endPoint
real(kind=RealKind), dimension(3) :: centerPoint
real(kind=RealKind) :: zoomFactor, zstart, zend
character(60) :: dirname, mapname
real(kind=RealKind), dimension(:,:,:), pointer :: unigrid, tmpgrid
real(kind=RealKind), dimension(17) :: z, rate
type(angleType) :: angles
real(kind=RealKind), dimension(3) :: alpha
type(patternType), dimension(:), pointer :: pattern
synthfile(1) = 'model41-salpeter-burst34/spectrum.out'
synthfile(2) = 'model42-salpeter-burst34/spectrum.out'
synthfile(3) = 'model43-salpeter-burst34/spectrum.out'
synthfile(4) = 'model44-salpeter-burst34/spectrum.out'
synthfile(5) = 'model45-salpeter-burst34/spectrum.out'
dustApproximation = 0
selfShieldingThreshold = 1.*kpc
currentRedshift = 3.
massStellarParticle = 1 ! normal(8x) = 1, hiRes(64x) = 2, superHiRes(512x) = 3, massive = 10, hiResHeavy(64x_heavy) = 4, crazyHiRes(512x8x) = 5
grid = ''
sources = ''
restartCellArrayName = ''
mode = 1
upperAgeLimit = 10.*Myr ! 34.*Myr 10.*Myr 3.4*Myr
sphDir = ''
synthesisDir = ''
open(unit=10, file='inputParameters', action="read", status="old")
status = 0
restart = 0
uvbCoefficient = 1.
reionizationModel = 0
do while (status.eq.0)
read(unit=10, fmt="(a90)", IOSTAT=status) scratchString
if (scratchString(1:20).eq.'dustApproximation = ') read (scratchString(21:90),*) dustApproximation
if (scratchString(1:25).eq.'selfShieldingThreshold = ') then
read (scratchString(26:90),*) selfShieldingThreshold
selfShieldingThreshold = selfShieldingThreshold*kpc
endif
if (scratchString(1:18).eq.'currentRedshift = ') read (scratchString(19:90),*) currentRedshift
if (scratchString(1:17).eq.'uvbCoefficient = ') read (scratchString(18:90),*) uvbCoefficient
if (scratchString(1:22).eq.'massStellarParticle = ') read (scratchString(23:90),*) massStellarParticle
if (scratchString(1:7).eq.'grid = ') read (scratchString(8:90),*) grid
if (scratchString(1:10).eq.'sources = ') read (scratchString(11:90),*) sources
if (scratchString(1:7).eq.'mode = ') read (scratchString(8:90),*) mode
if (scratchString(1:15).eq.'upperAgeLimit = ') then
read (scratchString(16:90),*) upperAgeLimit
upperAgeLimit = upperAgeLimit * Myr
endif
if (scratchString(1:9).eq.'sphDir = ') read (scratchString(10:90),*) sphDir
if (scratchString(1:15).eq.'synthesisDir = ') read (scratchString(16:90),*) synthesisDir
if (scratchString(1:10).eq.'restart = ') read (scratchString(11:90),*) restart
if (scratchString(1:23).eq.'restartCellArrayName = ') read (scratchString(24:90),*) restartCellArrayName
if (scratchString(1:20).eq.'reionizationModel = ') read (scratchString(21:90),*) reionizationModel
enddo
close(10)
print*, 'dustApproximation =', dustApproximation
print*, 'selfShieldingThreshold =', selfShieldingThreshold/kpc
print*, 'currentRedshift =', currentRedshift
print*, 'massStellarParticle =', massStellarParticle
print*, 'grid = ', trim(grid)
print*, 'sources = ', trim(sources)
print*, 'mode = ', mode
print*, 'upperAgeLimit = ', upperAgeLimit/Myr
print*, 'sphDir = ', trim(sphDir)
print*, 'synthesisDir = ', trim(synthesisDir)
print*, 'restart = ', restart
print*, 'restartCellArrayName = ', trim(restartCellArrayName)
print*, 'uvbCoefficient = ', uvbCoefficient
print*, 'reionizationModel =', reionizationModel
readKinematics = .false.
itmp = len(trim(grid))
do i = 1, itmp
if (grid(i:i).eq.'v'.and.i.lt.itmp-1) then
if (grid(i:i+2).eq.'vel') readKinematics = .true.
endif
enddo
readMetals = .false.
itmp = len(trim(grid))
do i = 1, itmp
if (grid(i:i).eq.'m'.and.i.lt.itmp-1) then
if (grid(i:i+2).eq.'met') readMetals = .true.
endif
enddo
if (mode.eq.stellarTransferThinUVB .or. mode.eq.bothStellarUVBTransfer) then
runStellarTransfer = .true.
else
runStellarTransfer = .false.
endif
if (mode.eq.UVBTransferOnly .or. mode.eq.bothStellarUVBTransfer) then
runUVBTransfer = .true.
else
runUVBTransfer = .false.
endif
! initialize collisional/recombination rate data
logtem0 = log(temstart)
logtem9 = log(temend)
dlogtem = (log(temend) - log(temstart))/real(nratec-1)
call calc_rates(nratec, temstart, temend, &
ceHIa, ceHeIa, ceHeIIa, ciHIa, ciHeIa, &
ciHeISa, ciHeIIa, reHIIa, reHeII1a, &
reHeII2a, reHeIIIa, brema, lineHIa, compa, &
hyd01ka, h2k01a, vibha, rotha, rotla, &
gpldl, gphdl, hdltea, hdlowa, &
k1a, k2a, k3a, k4a, k5a, k6a, k7a, k8a, k9a, k10a, &
k11a, k12a, k13a, k13dda, k14a, k15a, k16a, k17a, &
k18a, k19a, k20a, k21a, k22a, &
k50a, k51a, k52a, k53a, k54a, k55a, k56a, &
recombinationType)
! compute radiative rates for the uniform fiducial background
call uniformTable(nfbins,frequencyBinWidth,alphaQuasar,alphaStellar)
! set background radiation
! Abel & Haehnelt 99 stellar component
stellar99 = 1./(1.+(7./(1.+currentRedshift))**4)*exp(-(currentRedshift/4.)**3)
! Abel & Haehnelt 99 quasar component
quasar99 = 10./(1.+(7./(1.+currentRedshift))**4)*exp(-(currentRedshift/2.5)**3)
! Paschos 02 total
pascal02 = 0.0188*exp(-(currentRedshift-0.5)**2/(1.+0.0625*(currentRedshift+2.09)**2.075))*(1.+currentRedshift)**3.35
! Razoumov 02 stellar component
component1 = stellar99
component2 = pascal02
step = 0.5*(tanh((currentRedshift-4.2)*1.5)+1.)
stellar02 = (1.-step)*component1 + step*component2
! Razoumov 02 quasar component
quasar02 = 10./(1.+(7./(1.+currentRedshift))**4)*exp(-(currentRedshift/2.5)**3)
gaussian = exp(-((currentRedshift-4.5)/2.)**2)*0.3
newQuasar02 = gaussian*stellar02 + (1.-gaussian)*quasar02
newStellar02 = (1.-gaussian)*stellar02 + gaussian*quasar02
step = 0.5*(tanh((currentRedshift-14.)*0.5)+1.)
newStellar02 = step*0. + (1.-step)*newStellar02
! Razoumov 02 total
total02 = newStellar02 + newQuasar02
uniformQuasar = newQuasar02 * 1.d-21 * contributionQuasar * uvbCoefficient
uniformStellar = newStellar02 * 1.d-21 * contributionStellar * uvbCoefficient
if (runUVBTransfer) then
uvbStellar1 = newStellar02 * 1.d-21 * contributionStellar * uvbCoefficient
uvbStellar2 = uvbStellar1 * (nu2/nu1)**(-alphaStellar)
uvbStellar3 = uvbStellar2 * (nu3/nu2)**(-alphaStellar)
uvbQuasar1 = newQuasar02 * 1.d-21 * contributionQuasar * uvbCoefficient
uvbQuasar2 = uvbQuasar1 * (nu2/nu1)**(-alphaQuasar)
uvbQuasar3 = uvbQuasar2 * (nu3/nu2)**(-alphaQuasar)
call powerSpectrumIndex(uvbStellar1,alphaStellar,uvbQuasar1,alphaQuasar,uvb1,alpha(1),nu1,nu2,.true.)
call powerSpectrumIndex(uvbStellar2,alphaStellar,uvbQuasar2,alphaQuasar,uvb2,alpha(2),nu2,nu3,.true.)
call powerSpectrumIndex(uvbStellar3,alphaStellar,uvbQuasar3,alphaQuasar,uvb3,alpha(3),nu3,nu3,.false.)
write(*,'("band1: ",3f15.10)') uvbStellar1/1.d-21, uvbQuasar1/1.d-21, uvb1/1.d-21
write(*,'("band2: ",3f15.10)') uvbStellar2/1.d-21, uvbQuasar2/1.d-21, uvb2/1.d-21
write(*,'("band3: ",3f15.10)') uvbStellar3/1.d-21, uvbQuasar3/1.d-21, uvb3/1.d-21
print*, 'alpha =', alpha
call uvbBetaTable(nfbins,frequencyBinWidth,alpha)
endif
! hydrogen ionization rates for models z_re=6 and z_re=10
if (reionizationModel.ne.0) then
print*, 'assuming reionization at z=', reionizationModel
select case (reionizationModel)
case (6)
z = (/ 0., 0.316, 0.697, 1.187, 1.513, 2.343, 2.547, 2.765, &
3.024, 3.296, 3.772, 4.316, 4.657, 4.997, 5.302, 5.609, 100. /)
rate = (/ 0.0045, 0.0100, 0.0248, 0.0585, 0.0968, 0.1594, 0.1621, 0.1564, &
0.1403, 0.1159, 0.0683, 0.0248, 0.0112, 0.0058, 0.0017, 0.0004, 0.0000 /) * 1.d-11
case (10)
z = (/ 0., 0.316, 0.697, 1.187, 1.513, 2.343, 2.547, 2.972, &
3.432, 3.976, 5.065, 6.221, 6.902, 7.650, 8.331, 9.419, 100. /)
rate = (/ 0.0045, 0.0100, 0.0248, 0.0585, 0.0968, 0.1594, 0.1621, 0.1570, &
0.1444, 0.1240, 0.0710, 0.0262, 0.0128, 0.0058, 0.0014, 0.0003, 0.0000 /) * 1.d-11
case default
print*, 'select proper reionization model'
stop
endselect
i = 1
do while (currentRedshift.gt.z(i))
i = i + 1
enddo
actualRate = (currentRedshift-z(i-1)) / (z(i)-z(i-1)) * (rate(i)-rate(i-1)) + rate(i-1)
rateCoef = actualRate / (4.*pi * (uniformQuasar*quasar%ksi24 + uniformStellar*stellar%ksi24))
uniformQuasar = uniformQuasar * rateCoef
uniformStellar = uniformStellar * rateCoef
if (runUVBTransfer) then
uvb1 = uvb1 * rateCoef
uvb2 = uvb2 * rateCoef
uvb3 = uvb3 * rateCoef
endif
endif
call dustInitialize ! initialize dust data
rmax(1) = 1.984
rmax(2) = 3.264
rmax(3) = 5.739
rmax(4) = 10.65
rmax(5) = 20.46
rmax(6) = 40.05
rmax(7) = 79.25
rmax(8) = 157.6
rmax(9) = 314.4
rmax(10) = 627.9
do ir = 1, nrmax
rmax(ir) = sqrt(3.)*(sqrt(0.5*4.**(ir-1)-1./12.)+0.5)
if (mode.ne.printNumberOfCells) print*, ir, rmax(ir)
enddo
rmax = rmax/2.
print*, 'maxPixelLevel =', maxPixelLevel
sphere%level = 0
sphere%refined = .false.
sd_id = sfstart(trim(sphDir)//trim(grid)//'.h4', dfacc_read)
status = sffinfo(sd_id, n_datasets, n_file_attrs)
start = 0
edges(1) = 1
stride = 1
sds_id = sfselect(sd_id, 0)
status = sfrdata(sds_id, start, stride, edges, nlevels)
status = sfendacc(sds_id)
nvariables = 4
if (readKinematics) nvariables = nvariables + 1
if (readMetals) nvariables = nvariables + 1
if (nlevels.ne.(n_datasets-1)/nvariables) then
print*, trim(sphDir)//trim(grid)//'.h4'
print*, 'error in nlevels', nlevels, n_datasets
stop
endif
allocate(readLevel(nlevels))
if (readMetals) then
maxOxygenAbundance = 0.
totalOxygen = 0.
endif
do level = 1, nlevels
sds_id = sfselect(sd_id, nvariables*(level-1)+1)
status = sfginfo(sds_id,sds_name,rank,dim_sizes,data_type,n_attrs)
! if (mode.ne.printNumberOfCells) print*, trim(sds_name), rank, dim_sizes(1), dim_sizes(2)
if (status.ne.0 .or. rank.ne.2) then
write(*,*) 'error reading ...', level, status, rank, nlevels
stop
endif
readLevel(level)%ncell = dim_sizes(1)
allocate(readLevel(level)%pos(readLevel(level)%ncell,3))
allocate(readLevel(level)%lT(readLevel(level)%ncell))
allocate(readLevel(level)%lnH(readLevel(level)%ncell))
allocate(readLevel(level)%lx(readLevel(level)%ncell))
if (readKinematics) allocate(readLevel(level)%vel(readLevel(level)%ncell,3))
if (readMetals) allocate(readLevel(level)%abun(readLevel(level)%ncell,4))
! print*, 'level =', level, ' ncells =', readLevel(level)%ncell
edges(1:2) = (/ readLevel(level)%ncell, 3 /)
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%pos)
status = sfendacc(sds_id)
sds_id = sfselect(sd_id, nvariables*(level-1)+2)
edges(1) = readLevel(level)%ncell
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%lT)
status = sfendacc(sds_id)
sds_id = sfselect(sd_id, nvariables*(level-1)+3)
edges(1) = readLevel(level)%ncell
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%lnH)
status = sfendacc(sds_id)
if (mode.eq.printNumberOfCells) then
tmp = 0.
do icell = 1, readLevel(level)%ncell
tmp = max(tmp,readLevel(level)%lnH(icell))
enddo
write(*,*) 'level =', level, readLevel(level)%ncell, tmp
endif
sds_id = sfselect(sd_id, nvariables*(level-1)+4)
edges(1) = readLevel(level)%ncell
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%lx)
status = sfendacc(sds_id)
counter = 5
if (readMetals) then
sds_id = sfselect(sd_id, nvariables*(level-1)+counter)
edges(1:2) = (/ readLevel(level)%ncell, 4 /)
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%abun)
status = sfendacc(sds_id)
counter = counter + 1
endif
if (readKinematics) then
sds_id = sfselect(sd_id, nvariables*(level-1)+counter)
edges(1:2) = (/ readLevel(level)%ncell, 3 /)
status = sfrdata(sds_id, start, stride, edges, readLevel(level)%vel)
status = sfendacc(sds_id)
counter = counter + 1
endif
if (readMetals) then
do icell = 1, readLevel(level)%ncell
rho = 10.**readLevel(level)%lnH(icell)
if (readLevel(level)%abun(icell,2).gt.maxOxygenAbundance) &
maxOxygenAbundance = readLevel(level)%abun(icell,2)
totalOxygen = totalOxygen + readLevel(level)%abun(icell,2)*rho/(2.**(3*(level-1)))
enddo
endif
enddo
if (readMetals) print*, 'maximum oxygen abundance =', maxOxygenAbundance, totalOxygen
status = sfend(sd_id)
if (mode.eq.printNumberOfCells) stop
i = 0
itmp = 0
do while (itmp.lt.readLevel(1)%ncell)
i = i + 1
itmp = i**3
enddo
if (itmp.ne.readLevel(1)%ncell) then
write(*,*) 'base grid needs to be of size n^3', itmp, readLevel(1)%ncell
stop
endif
nx = i
ny = i
nz = i
write(*,*) 'grid:', nx,'^3'
print*, 'nlevels =', nlevels
do i = 1, nlevels
print*, 'level ', i, ' has', readLevel(i)%ncell, ' cells'
enddo
level = 1
xa = 1.d10
xb = - 1.d10
ya = 1.d10
yb = - 1.d10
za = 1.d10
zb = - 1.d10
do icell = 1, readLevel(level)%ncell
xa = min(xa,dble(readLevel(level)%pos(icell,1)))
xb = max(xb,dble(readLevel(level)%pos(icell,1)))
ya = min(ya,dble(readLevel(level)%pos(icell,2)))
yb = max(yb,dble(readLevel(level)%pos(icell,2)))
za = min(za,dble(readLevel(level)%pos(icell,3)))
zb = max(zb,dble(readLevel(level)%pos(icell,3)))
enddo
tmp1 = 0.5*(xa+xb)
tmp2 = 0.5*(xb-xa)*float(nx)/float(nx-1)
xa = tmp1 - tmp2
xb = tmp1 + tmp2
tmp1 = 0.5*(ya+yb)
tmp2 = 0.5*(yb-ya)*float(ny)/float(ny-1)
ya = tmp1 - tmp2
yb = tmp1 + tmp2
tmp1 = 0.5*(za+zb)
tmp2 = 0.5*(zb-za)*float(nz)/float(nz-1)
za = tmp1 - tmp2
zb = tmp1 + tmp2
write(*,*) 'edges of computational grid in kpc:'
write(*,*) xa, xb
write(*,*) ya, yb
write(*,*) za, zb
physicalBoxSize = abs(xa-xb)*kpc
do level = 1, nlevels
do icell = 1, readLevel(level)%ncell
readLevel(level)%pos(icell,1) = (readLevel(level)%pos(icell,1)-xa)/(xb-xa)
readLevel(level)%pos(icell,2) = (readLevel(level)%pos(icell,2)-ya)/(yb-ya)
readLevel(level)%pos(icell,3) = (readLevel(level)%pos(icell,3)-za)/(zb-za)
enddo
enddo
print*, 'setting up a fully nested 3D grid'
allocate(baseGrid%cell(nx,ny,nz))
do i = 1, nx
do j = 1, ny
do k = 1, nz
baseGrid%cell(i,j,k)%refined = .false.
baseGrid%cell(i,j,k)%tgas = 0.d0
baseGrid%cell(i,j,k)%rho = 0.d0
baseGrid%cell(i,j,k)%HI = 0.d0
baseGrid%cell(i,j,k)%HeI = 0.d0
baseGrid%cell(i,j,k)%HeII = 0.d0
baseGrid%cell(i,j,k)%rhoCoef = 1.
baseGrid%cell(i,j,k)%velx = 0.d0
baseGrid%cell(i,j,k)%vely = 0.d0
baseGrid%cell(i,j,k)%velz = 0.d0
baseGrid%cell(i,j,k)%krate24 = 0.d0
baseGrid%cell(i,j,k)%krate25 = 0.d0
baseGrid%cell(i,j,k)%krate26 = 0.d0
baseGrid%cell(i,j,k)%crate24 = 0.d0
baseGrid%cell(i,j,k)%crate25 = 0.d0
baseGrid%cell(i,j,k)%crate26 = 0.d0
baseGrid%cell(i,j,k)%level = 0
baseGrid%cell(i,j,k)%parent => baseGrid
nullify(baseGrid%cell(i,j,k)%cell)
nullify(baseGrid%cell(i,j,k)%segment)
nullify(baseGrid%cell(i,j,k)%lastSegment)
enddo
enddo
enddo
! smoothing level 1 abundancies
level = 1
allocate(unigrid(nx,ny,nz))
allocate(tmpgrid(nx,ny,nz))
do icell = 1, readLevel(level)%ncell
i0 = int(readLevel(level)%pos(icell,1)*nx) + 1
j0 = int(readLevel(level)%pos(icell,2)*ny) + 1
k0 = int(readLevel(level)%pos(icell,3)*nz) + 1
unigrid(i0,j0,k0) = readLevel(level)%abun(icell,2)
enddo
do itmp = 1, 2
tmpgrid = 0.
do i = 1, nx
do j = 1, ny
do k = 1, nz
tmpgrid(i,j,k) = tmpgrid(i,j,k) + 0.5*unigrid(i,j,k)
if (i.gt.1) tmpgrid(i-1,j,k) = tmpgrid(i-1,j,k) + 0.25*unigrid(i,j,k)
if (i.lt.nx) tmpgrid(i+1,j,k) = tmpgrid(i+1,j,k) + 0.25*unigrid(i,j,k)
enddo
enddo
enddo
unigrid = tmpgrid
tmpgrid = 0.
do i = 1, nx
do j = 1, ny
do k = 1, nz
tmpgrid(i,j,k) = tmpgrid(i,j,k) + 0.5*unigrid(i,j,k)
if (j.gt.1) tmpgrid(i,j-1,k) = tmpgrid(i,j-1,k) + 0.25*unigrid(i,j,k)
if (j.lt.ny) tmpgrid(i,j+1,k) = tmpgrid(i,j+1,k) + 0.25*unigrid(i,j,k)
enddo
enddo
enddo
unigrid = tmpgrid
tmpgrid = 0.
do i = 1, nx
do j = 1, ny
do k = 1, nz
tmpgrid(i,j,k) = tmpgrid(i,j,k) + 0.5*unigrid(i,j,k)
if (k.gt.1) tmpgrid(i,j,k-1) = tmpgrid(i,j,k-1) + 0.25*unigrid(i,j,k)
if (k.lt.nz) tmpgrid(i,j,k+1) = tmpgrid(i,j,k+1) + 0.25*unigrid(i,j,k)
enddo
enddo
enddo
unigrid = tmpgrid
enddo
do icell = 1, readLevel(level)%ncell
i0 = int(readLevel(level)%pos(icell,1)*nx) + 1
j0 = int(readLevel(level)%pos(icell,2)*ny) + 1
k0 = int(readLevel(level)%pos(icell,3)*nz) + 1
readLevel(level)%abun(icell,2) = unigrid(i0,j0,k0)
enddo
deallocate(unigrid,tmpgrid)
do level = 1, nlevels
do icell = 1, readLevel(level)%ncell
x0 = readLevel(level)%pos(icell,1)
y0 = readLevel(level)%pos(icell,2)
z0 = readLevel(level)%pos(icell,3)
i0 = int(x0*nx) + 1
j0 = int(y0*ny) + 1
k0 = int(z0*nz) + 1
xnew = x0*float(nx) - float(i0-1)
ynew = y0*float(ny) - float(j0-1)
znew = z0*float(nz) - float(k0-1)
readArray(1) = readLevel(level)%lT(icell)
readArray(2) = readLevel(level)%lnH(icell)
readArray(3) = readLevel(level)%lx(icell)
nvariables = 3
if (readKinematics) then
readArray(nvariables+1) = readLevel(level)%vel(icell,1)
readArray(nvariables+2) = readLevel(level)%vel(icell,2)
readArray(nvariables+3) = readLevel(level)%vel(icell,3)
nvariables = nvariables + 3
endif
if (readMetals) then
readArray(nvariables+1) = readLevel(level)%abun(icell,1)
readArray(nvariables+2) = readLevel(level)%abun(icell,2)
readArray(nvariables+3) = readLevel(level)%abun(icell,3)
readArray(nvariables+4) = readLevel(level)%abun(icell,4)
nvariables = nvariables + 4
endif
call placeCellProjectWithVelocity(baseGrid%cell(i0,j0,k0),level,xnew,ynew,znew,readArray)
enddo
enddo
do level = 1, nlevels
deallocate(readLevel(level)%pos)
deallocate(readLevel(level)%lT)
deallocate(readLevel(level)%lnH)
deallocate(readLevel(level)%lx)
if (readKinematics) deallocate(readLevel(level)%vel)
if (readMetals) deallocate(readLevel(level)%abun)
enddo
deallocate(readLevel)
! densestCell => baseGrid%cell(1,1,1)
! do i = 1, nx
! do j = 1, ny
! do k = 1, nz
! call printCell(baseGrid%cell(i,j,k),0,(/i,j,k/))
! enddo
! enddo
! enddo
! densestCell => baseGrid%cell(1,1,1)
! write(*,*) densestCell%nh, densestCell%level
! do i = 1, nx
! do j = 1, ny
! do k = 1, nz
! call countCells(baseGrid%cell(i,j,k))
! enddo
! enddo
! enddo
! write(*,*) 'found', icosmic, ' cells'
! write(*,*) densestCell%nh, densestCell%level
! currentCell => baseGrid%cell(1,1,1)
! if (associated(currentCell%lastSegment)) then
! allocate(currentCell%lastSegment%segment)
! currentCell%lastSegment => currentCell%lastSegment%segment
! else
! allocate(currentCell%segment)
! currentCell%lastSegment => currentCell%segment
! endif
! newSegment => currentCell%lastSegment
if (mode.eq.clumpingFactor) then
volSum = 0.
volSumDensity = 0.
volSumDensity2 = 0.
do i = 1, nx
do j = 1, ny
do k = 1, nz
call computeClumping(baseGrid%cell(i,j,k))
enddo
enddo
enddo
volSumDensity = volSumDensity / volSum
volSumDensity2 = volSumDensity2 / volSum
print*, 'clumping =', volSumDensity2 / volSumDensity**2, volSum
stop
endif
if (mode.eq.initialConfiguration) then
izone = 3
nxmap = 1024
nymap = nxmap
centerPoint = (/ 0.5, 0.5, 0.5 /)
zoomFactor = 2.
mapname = 'abun-S33-06-64-10.4.h4'
startPoint = (/ max(centerPoint(1) - 0.5/zoomFactor,0.), &
max(centerPoint(2) - 0.5/zoomFactor,0.) /)
zstart = max(centerPoint(3) - 0.5/zoomFactor,0.)
zend = min(centerPoint(3) + 0.5/zoomFactor,1.)
endPoint(1:2) = min(startPoint(1:2)+1./zoomFactor,1.)
do i = 1, 2
do j = 1, 2
do k = 1, 2
call rotateIndices(i,j,k,2,2,2,izone,is(i,j,k),js(i,j,k),ks(i,j,k))
enddo
enddo
enddo
allocate(mapArray(nxmap,nymap))
endPoint(1:2) = min(startPoint(1:2)+1./zoomFactor,1.)
kstart = max(int(zstart*nz)+1,1)
kend = min(int(zend*nz)+1,nz)
do i = 1, nxmap
print*, 'doing', i, ' out of', nxmap
x0 = (endPoint(1)-startPoint(1))*(float(i)-0.5)/float(nxmap) + startPoint(1)
i0 = int(x0*nx) + 1
do j = 1, nymap
y0 = (endPoint(2)-startPoint(2))*(float(j)-0.5)/float(nymap) + startPoint(2)
j0 = int(y0*ny) + 1
xnew = x0*float(nx) - float(i0-1)
ynew = y0*float(ny) - float(j0-1)
mapArray(i,j) = 0.
totalMass = 0.
do k0 = kstart, kend
call rotateIndices(i0,j0,k0,nx,ny,nz,izone,icell,jcell,kcell)
call projectVariableToMap(baseGrid%cell(icell,jcell,kcell),mapArray(i,j), &
xnew,ynew,physicalBoxSize/dfloat(nx),is,js,ks)
enddo
mapArray(i,j) = mapArray(i,j) / totalMass
enddo
enddo
sd_id = sfstart(trim(mapname),dfacc_create)
edges(1) = nxmap
edges(2) = nymap
sds_id = sfcreate(sd_id,'map',dfnt_float32,2,edges)
start = 0
stride = 1
status = sfwdata(sds_id,start,stride,edges,mapArray)
status = sfendacc(sds_id)
status = sfend(sd_id)
deallocate(mapArray)
stop
endif
print*, 'read positions and ages of all stars'
open(unit=15,file=trim(sphDir)//trim(sources),action='read',status='old')
status = 0
iStar = 0
do while (status.eq.0)
read(15,*,iostat=status) j
iStar = iStar + 1
enddo
nStars = iStar - 1
close(15)
allocate(star(nStars))
open(unit=15,file=trim(sphDir)//trim(sources),action='read',status='old')
maxOxygenAbundance = 0.
do iStar = 1, nStars
currentStar => star(iStar)
read(15,*,iostat=status) currentStar%level, tmp1, tmp2, tmp3, currentStar%age
currentStar%age = currentStar%age * Myr
xbase = (tmp1-xa)/(xb-xa)
ybase = (tmp2-ya)/(yb-ya)
zbase = (tmp3-za)/(zb-za)
call localizeSplitContinuationCell(xbase,ybase,zbase,0,nx,ny,nz)
! print*, iStar, currentStar%level-1, splitContinuationCell%level
currentStar%level = splitContinuationCell%level ! overwrite the levels
allocate(currentStar%position(3*currentStar%level+3))
currentStar%position = splitContinuationCallSequence(1:3*currentStar%level+3)
startingPoint = splitContinuation
! print*, '--------'
! print*, 'source #', iStar
! print*, 'cube =', xbase, ybase, zbase
! print*, 'level =', splitContinuationCell%level, ' rho =', splitContinuationCell%rho
! print*, 'stellar metallicity =', splitContinuationCell%abun2
! print*, iStar, ' level =', currentStar%level
! print*, 'position =', currentStar%position
! print*, 'point =', startingPoint
if (splitContinuationCell%abun2.gt.maxOxygenAbundance) &
maxOxygenAbundance = splitContinuationCell%abun2
enddo
print*, 'maximum stellar oxygen abundance =', maxOxygenAbundance
close(15)
nStarsSpecificAge = 0
do iStar = 1, nStars
currentStar => star(iStar)
if (currentStar%age.le.upperAgeLimit) then
nStarsSpecificAge = nStarsSpecificAge + 1
currentStar%weight = 1
! print*, iStar, currentStar%level, currentStar%age/Myr
else
currentStar%weight = 0
endif
enddo
if (mode.eq.plotPDFs) then
print*, 'computing stellar and gas PDF'
meanHostDensity = 0.
pdfStar = 0
pdfStarOutside = 0
do iStar = 1, nStars
currentStar => star(iStar)
if (currentStar%weight.gt.0) then
startingPoint%x = 0.5
startingPoint%y = 0.5
startingPoint%z = 0.5
i = currentStar%position(1)
j = currentStar%position(2)
k = currentStar%position(3)
call localizeCellFromStar(baseGrid%cell(i,j,k),currentStar)
currentCell => currentStar%hostCell
write(*,2015) iStar, currentCell%level, psi*currentCell%rho/mh
meanHostDensity = meanHostDensity + currentCell%HI
tmp = dlog10(currentCell%rho/msun*pc**3)
if (tmp.gt.apdf .and. tmp.lt.bpdf) then
itmp = int((tmp-apdf)/(bpdf-apdf)*float(npdf)) + 1
pdfStar(itmp) = pdfStar(itmp) + 1
else
pdfStarOutside = pdfStarOutside + 1
endif
2015 format(2i5,es18.8)
endif
enddo
print*, 'meanHostDensity =', meanHostDensity/float(nStarsSpecificAge), nStarsSpecificAge
pdfGas = 0.
pdfGasOutside = 0.
do i = 1, nx
do j = 1, ny
do k = 1, nz
call computeGasPDF(baseGrid%cell(i,j,k))
enddo
enddo
enddo
itmp = pdfStarOutside
tmp = pdfGasOutside
do i = 1, npdf
itmp = itmp + pdfStar(i)
tmp = tmp + pdfGas(i)
enddo
do i = 1, npdf
tmp1 = max(1.d-38,pdfGas(i)/tmp)
tmp2 = max(1.d-38,float(pdfStar(i))/float(itmp))
write(*,2025) (float(i)-0.5)/float(npdf)*(bpdf-apdf)+apdf, dlog10(tmp1), dlog10(tmp2)
enddo
2025 format(3f15.6)
stop
endif
! initialize stellar population synthesis data
if (nMetallicity.ne.5) then
print*, 'please update the list of metallicities below'
stop
endif
metallicity(1:5) = (/ 0.0004, 0.004, 0.008, 0.020, 0.050 /)
metallicity = dlog10(metallicity)
do iMetal = 1, nMetallicity
filename = synthfile(iMetal)
open (unit=11,file=trim(synthesisDir)//trim(filename),action='read',status='old')
status = 0
i = 0
readNewSpectrum = .false.
iSpectrum = 0
do while (status.eq.0)
i = i + 1
read(unit=11, fmt='(a65)', IOSTAT=status) scratchString
if (status.eq.0) then
! print*, i, scratchString(1:65)
if (scratchString(2:10).eq.'TIME [YR]') then
read(unit=11, fmt='(a65)', IOSTAT=status) scratchString
read(unit=11, fmt='(a65)', IOSTAT=status) scratchString
readNewSpectrum = .true.
iWavelength = 0
else
if (readNewSpectrum .and. scratchString(2:6).ne.'MODEL') then
if (iWavelength.eq.0) iSpectrum = iSpectrum + 1
iWavelength = iWavelength + 1
read (scratchString(2:13),*) spectrumTime(iSpectrum)
spectrumTime(iSpectrum) = spectrumTime(iSpectrum) * yr
read (scratchString(14:28),*) wavelength(iWavelength)
wavelength(iWavelength) = wavelength(iWavelength) * angstrom
read (scratchString(29:41),*) specificLuminosity(iMetal,iSpectrum,iWavelength)
else
readNewSpectrum = .false.
endif
endif
endif
enddo
close(11)
enddo
! specificLuminosity was computed for the SF rate of 11.6 msun/yr distributed among 34 star particles
! specificLuminosity * 34 -- luminosity of 11.6 msun/yr constant SF rate
! nStars / 347 * 11.6 Msun/yr -- SF rate in the volume
! nStars / 347 * 10.**specificLuminosity * 34 -- luminosity of the volume
! nStars / 347 * 10.**specificLuminosity * 34 / nStarsSpecificAge -- luminosity per star particle
specificLuminosity = specificLuminosity + log10(float(nStars) / 347. * 34. / float(nStarsSpecificAge)) ! per star particle
select case (massStellarParticle)
case (hiRes)
print*, 'high-res model: stellar sources are 1/8th the mass of normal model'
specificLuminosity = specificLuminosity - log10(8.)
case (superHiRes)
print*, 'high-res model: stellar sources are 1/64th the mass of normal model'
specificLuminosity = specificLuminosity - log10(64.)
case (crazyHiRes)
print*, 'high-res model: stellar sources are 1/512th the mass of normal model'
specificLuminosity = specificLuminosity - log10(512.)
case (massive)
print*, 'massive star particles: stellar sources are 2.7818X the mass of normal model'
specificLuminosity = specificLuminosity + log10(2.7818)
case (hiResHeavy)
print*, 'high-res model: stellar sources are (5.832/8)X the mass of normal model'
specificLuminosity = specificLuminosity + log10(5.832/8.)
case (light)
print*, 'high-res model: stellar sources are 0.6^3/512th the mass of normal model'
specificLuminosity = specificLuminosity + 3.*log10(0.6) - log10(512.)
case (lyAlpha)
print*, 'lyAlpha model: stellar sources are 0.9286/8th the mass of normal model'
specificLuminosity = specificLuminosity + log10(65./(70.*8.))
endselect
! ! specify sources by physical coordinates
! nStars = 4
! allocate(star(nStars))
! star(1)%data = (/ 0.9765625, 0.0000000, 0.0000000, 151. /)
! star(2)%data = (/ -4.882812, 22.46094, 30.27344, 21. /)
! star(3)%data = (/ -2.929688, -0.9765625, -2.929688, 35. /)
! star(4)%data = (/ -0.9765625, 0.0000000, -0.9765625, 31. /)
! do iStar = 1, nStars
! currentStar => star(iStar)
! xbase = (currentStar%data(1)-xa)/(xb-xa)
! ybase = (currentStar%data(2)-ya)/(yb-ya)
! zbase = (currentStar%data(3)-za)/(zb-za)
! call localizeSplitContinuationCell(xbase,ybase,zbase,0,nx,ny,nz)
! currentStar%level = splitContinuationCell%level
! allocate(currentStar%position(3*currentStar%level+3))
! currentStar%position = splitContinuationCallSequence(1:3*currentStar%level+3)
! startingPoint = splitContinuation
! currentStar%luminosity = currentStar%data(4) * 1.d54
! print*, '--------'
! print*, 'source #', iStar
! print*, 'cube =', xbase, ybase, zbase
! print*, 'level =', currentStar%level
! print*, 'position =', currentStar%position
! print*, 'point =', startingPoint
! print*, 'luminosity =', currentStar%luminosity
! enddo
! currentStar%position = (/65,64,65,1,2,1,1,2,1,2,2,2,1,2,2,1,2,1/)
! ! specify sources by host cells
! allocate(star(nStars))
! do iStar = 1, nStars
! currentStar => star(iStar)
! select case (iStar)
! case(1)
! currentStar%level = 5
! allocate(currentStar%position(3*currentStar%level+3))
! currentStar%position = (/65,64,64, 1,2,2, 1,2,2, 1,1,2, 2,1,2, 2,2,2/) ! rho = 3.95296746662E-22
! currentStar%luminosity = 5.d55 ! total photon # luminosity in all bands [1/s]
! currentStar%mass = 151. * starParticleMass
! case(2)
! currentStar%level = 5
! allocate(currentStar%position(3*currentStar%level+3))
! currentStar%position = (/63,64,63, 2,1,2, 2,2,2, 1,2,2, 1,1,1, 1,1,1/) ! rho = 8.34276624451E-22
! currentStar%luminosity = 1.d30 ! total photon # luminosity in all bands [1/s]
! currentStar%mass = 35. * starParticleMass
! case(3)
! currentStar%level = 4
! allocate(currentStar%position(3*currentStar%level+3))
! currentStar%position = (/64,64,64, 2,2,2, 2,2,1, 2,2,2, 2,1,2/) ! rho = 1.39107558854E-22
! currentStar%luminosity = 1.d30 ! total photon # luminosity in all bands [1/s]
! currentStar%mass = 31. * starParticleMass
! case(4)
! currentStar%level = 4
! allocate(currentStar%position(3*currentStar%level+3))
! currentStar%position = (/62,76,80, 1,2,1, 1,1,2, 2,1,1, 1,1,1/) ! rho = 7.53792642535E-23
! currentStar%luminosity = 1.d30 ! total photon # luminosity in all bands [1/s]
! currentStar%mass = 21. * starParticleMass
! end select
! startingPoint%x = 0.5
! startingPoint%y = 0.5
! startingPoint%z = 0.5
! call absoluteCoordinates(currentStar%level,currentStar%position,startingPoint,nx,ny,nz)
! print*, '--------'
! print*, 'source #', iStar
! print*, 'cube =', xbase, ybase, zbase
! print*, 'level =', currentStar%level
! print*, 'position =', currentStar%position
! print*, 'point =', startingPoint
! ! print*, 'luminosity =', currentStar%luminosity
! call localizeSplitContinuationCell(xbase,ybase,zbase,0,nx,ny,nz)
! print*, 'level,rho =', splitContinuationCell%level, splitContinuationCell%rho
! enddo
print*, 'counting cells'
icosmic = 0
do i = 1, nx
do j = 1, ny