forked from garrelt/C2-Ray3Dm
-
Notifications
You must be signed in to change notification settings - Fork 0
/
sourceprops_cubep3m_compr.F90
259 lines (224 loc) · 9.28 KB
/
sourceprops_cubep3m_compr.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
!>
!! \brief This module contains data and routines for handling the sources
!! of reionization.
!!
!!
!! \b Author: Garrelt Mellema, Ilian Iliev
!!
!! \b Date: 30-Jan-2008
!!
!! \b Version: cubep3m Nbody simulation, with source suppression and compressed ionization fractions
module sourceprops
use precision, only: dp
use my_mpi
use file_admin, only: logf
use cgsconstants, only: m_p
use astroconstants, only: M_SOLAR
use cosmology_parameters, only: Omega_B, Omega0
use nbody, only: id_str, M_grid, dir_src
use material, only: xh_compr,ionized_from_compr
use grid, only: x,y,z
use c2ray_parameters, only: phot_per_atom1, phot_per_atom2, lifetime, &
S_star_nominal, StillNeutral
integer :: NumSrc
integer,dimension(:,:),allocatable :: srcpos
real(kind=dp),dimension(:,:),allocatable :: rsrcpos
real(kind=dp),dimension(:),allocatable :: srcMass
real(kind=dp),dimension(:),allocatable :: NormFlux
integer,dimension(:),allocatable :: srcSeries
integer,private :: NumSrc0=0
integer,dimension(3),private :: srcpos0
real(kind=dp),private :: srcMass00,srcMass01
character(len=6) :: z_str
integer,private :: NumMassiveSrc,NumSupprbleSrc,NumSupprsdSrc
contains
! =======================================================================
subroutine source_properties(zred_now,nz,lifetime2,restart)
! Input routine: establish the source properties
! Authors: Garrelt Mellema, Ilian Iliev
! Update: 30-Jan-2008 (20-Sep-2006 (3-jan-2005, 15-Apr-2004))
! For random permutation of sources
use m_ctrper
real(kind=dp),intent(in) :: zred_now ! current redshift
real(kind=dp),intent(in) :: lifetime2 ! time step
integer,intent(in) :: nz
integer,intent(in) :: restart
character(len=512) :: sourcelistfile,sourcelistfilesuppress
integer :: ns,ns0
#ifdef MPI
integer :: mympierror
#endif
#ifdef MPILOG
write(logf,*) "Check sourceprops: ",zred_now,nz,lifetime2,restart
#endif
! Deallocate source arrays
if (allocated(srcpos)) deallocate(srcpos)
if (allocated(rsrcpos)) deallocate(rsrcpos)
if (allocated(srcMass)) deallocate(srcMass)
if (allocated(NormFlux)) deallocate(NormFlux)
if (allocated(srcSeries)) deallocate(srcSeries)
! Ask for input
if (rank == 0) then
! Sources are read from file
write(z_str,'(f6.3)') zred_now
! Construct the file name
sourcelistfile=trim(adjustl(dir_src))//&
trim(adjustl(z_str))//"-"//trim(adjustl(id_str))//"_sources.dat"
sourcelistfilesuppress=trim(adjustl(dir_src))//&
trim(adjustl(z_str))//"-"//trim(adjustl(id_str))// &
"_sources_used_wfgamma.dat"
if (restart == 0 .or. restart == 1) then
open(unit=50,file=sourcelistfile,status='old')
! Number of sources
read(50,*) NumSrc0
! Report
write(logf,*) "Total number of source locations, no suppression: ", &
NumSrc0
! Read in source positions and mass to establish number
! of non-suppressed sources
NumSrc = 0
NumMassiveSrc = 0
NumSupprbleSrc = 0
NumSupprsdSrc = 0
do ns0=1,NumSrc0
read(50,*) srcpos0(1),srcpos0(2),srcpos0(3),SrcMass00,SrcMass01
! the cell is still neutral, no suppression
if (SrcMass00 /= 0.0 .or. &
ionized_from_compr(xh_compr(srcpos0(1),srcpos0(2),srcpos0(3))) < StillNeutral) &
NumSrc=NumSrc+1
! Count different types of sources
if (SrcMass00 /= 0.0) NumMassiveSrc=NumMassiveSrc+1
if (SrcMass01 /= 0.0) NumSupprbleSrc=NumSupprbleSrc+1
! How many suppressed?
if (SrcMass01 /= 0.0 .and. &
ionized_from_compr(xh_compr(srcpos0(1),srcpos0(2),srcpos0(3))) >= StillNeutral) &
NumSupprsdSrc=NumSupprsdSrc+1
enddo
close(50)
write(logf,*) "Number of suppressable sources: ",NumSupprbleSrc
write(logf,*) "Number of suppressed sources: ",NumSupprsdSrc
write(logf,*) "Number of massive sources: ",NumMassiveSrc
if (NumSupprbleSrc > 0) write(logf,*) "Suppressed fraction: ", &
real(NumSupprsdSrc)/real(NumSupprbleSrc)
else
! Upon restart use the previously calculated suppressed source
! list
open(unit=49,file=sourcelistfilesuppress,status='unknown')
! Number of sources
read(49,*) NumSrc
close(49)
endif
write(logf,*) "Number of sources, with suppression: ",NumSrc
endif ! end of rank 0 test
#ifdef MPI
! Distribute source number to all other nodes
call MPI_BCAST(NumSrc,1,MPI_INTEGER,0,MPI_COMM_NEW,mympierror)
#endif
#ifdef MPILOG
if (rank /=0) write(logf,*) "Number of sources, with suppression: ",NumSrc
#endif
! Allocate arrays for this NumSrc
allocate(srcpos(3,NumSrc))
allocate(rsrcpos(3,NumSrc))
allocate(SrcMass(NumSrc))
allocate(NormFlux(NumSrc))
allocate(SrcSeries(NumSrc))
if (rank == 0) then
if (restart == 0 .or. restart == 1) then
open(unit=50,file=sourcelistfile,status='old')
! Number of sources
read(50,*) NumSrc0
! Read in source positions and mass
ns=0
do ns0=1,NumSrc0
read(50,*) srcpos0(1),srcpos0(2),srcpos0(3), &
SrcMass00,SrcMass01
if (ionized_from_compr(xh_compr(srcpos0(1),srcpos0(2),srcpos0(3))) < StillNeutral) then
! the cell is still neutral, no suppression
ns=ns+1
! Source positions in file start at 1!
srcpos(1,ns)=srcpos0(1)
srcpos(2,ns)=srcpos0(2)
srcpos(3,ns)=srcpos0(3)
! Source is always at cell centre!!
rsrcpos(1,ns)=x(srcpos(1,ns))
rsrcpos(2,ns)=y(srcpos(2,ns))
rsrcpos(3,ns)=z(srcpos(3,ns))
SrcMass(ns)=SrcMass00*phot_per_atom1 & !massive sources
+SrcMass01*phot_per_atom2 !small sources
elseif (SrcMass00 > 0.0d0) then
!the cell is ionized but source is massive enough to survive
!and is assumed Pop. II
ns=ns+1
! Source positions in file start at 1!
srcpos(1,ns)=srcpos0(1)
srcpos(2,ns)=srcpos0(2)
srcpos(3,ns)=srcpos0(3)
! Source is always at cell centre!!
rsrcpos(1,ns)=x(srcpos(1,ns))
rsrcpos(2,ns)=y(srcpos(2,ns))
rsrcpos(3,ns)=z(srcpos(3,ns))
SrcMass(ns)=SrcMass00*phot_per_atom1
!else
! ! Report
! write(logf,*) 'Source dropped: ', &
! srcpos0(1),srcpos0(2),srcpos0(3)
! write(logf,*) 'mass= ',SrcMass01*(M_grid/M_SOLAR), &
! xh(srcpos0(1),srcpos0(2),srcpos0(3),0)
endif
enddo
! Save new source list, without the suppressed ones
open(unit=49,file=sourcelistfilesuppress,status='unknown')
write(49,*) NumSrc
do ns0=1,NumSrc
write(49,*) srcpos(1,ns0),srcpos(2,ns0),srcpos(3,ns0), &
SrcMass(ns0)
enddo
close(49)
else ! of restart test
! Read source list from file saved previously
open(unit=49,file=sourcelistfilesuppress,status="old")
write(logf,*) "Reading ",NumSrc," sources from ", &
trim(adjustl(sourcelistfilesuppress))
read(49,*) NumSrc
do ns0=1,NumSrc
read(49,*) srcpos(1,ns0),srcpos(2,ns0),srcpos(3,ns0), &
SrcMass(ns0)
! Source is always at cell centre!!
rsrcpos(1,ns0)=x(srcpos(1,ns0))
rsrcpos(2,ns0)=y(srcpos(2,ns0))
rsrcpos(3,ns0)=z(srcpos(3,ns0))
enddo
close(49)
endif ! of restart test
endif ! of rank 0 test
#ifdef MPI
! Distribute the source parameters to the other nodes
call MPI_BCAST(srcpos,3*NumSrc,MPI_INTEGER,0,MPI_COMM_NEW,mympierror)
call MPI_BCAST(rsrcpos,3*NumSrc,MPI_DOUBLE_PRECISION,0,MPI_COMM_NEW,mympierror)
call MPI_BCAST(SrcMass,NumSrc,MPI_DOUBLE_PRECISION,0,MPI_COMM_NEW,mympierror)
#endif
! Turn masses into luminosities
do ns=1,NumSrc
NormFlux(ns)=SrcMass(ns)*M_grid* &!note that now photons/atom are included in SrcMass
Omega_B/(Omega0*m_p)/S_star_nominal
!NormFlux(ns)=NormFlux(ns)/lifetime
NormFlux(ns)=NormFlux(ns)/lifetime2
enddo
if (rank == 0) then
write(logf,*) 'Source lifetime=', lifetime2/3.1536e13
!write(logf,*) 'Source lifetime=', lifetime/3.1536e13
write(logf,*) 'Total flux= ',sum(NormFlux)
! Create array of source numbers for generating random order
do ns=1,NumSrc
SrcSeries(ns)=ns
enddo
! Make a random order
call ctrper(SrcSeries(1:NumSrc),1.0)
endif
#ifdef MPI
! Distribute the source series to the other nodes
call MPI_BCAST(SrcSeries,NumSrc,MPI_INTEGER,0,MPI_COMM_NEW,mympierror)
#endif
end subroutine source_properties
end module sourceprops