-
Notifications
You must be signed in to change notification settings - Fork 0
/
initialization.f90
309 lines (196 loc) · 7.78 KB
/
initialization.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
module array
use precision, only: dp
use input, only: level, partition, number_density, xHI, xHII, temperature, pi
use type, only: pyramid, cartesian, data
use solid_angle, only: ss, cc
implicit none
type(pyramid), dimension(1:level) :: pyramid_grid
type(cartesian), dimension(1:level) :: cartesian_grid
type(data), dimension(1:level) :: pos_x_pos_y_pos_z_dir_z
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_creation()
! allocation arrays of volume, length...
integer :: count
do count = 1,level
allocate(pyramid_grid(count)%volume(1:partition(count),1:partition(count)))
allocate(pyramid_grid(count)%cellsize(1:partition(count),1:partition(count)))
allocate(pyramid_grid(count)%distance(1:partition(count),1:partition(count)))
allocate(pyramid_grid(count)%solid_angle(1:partition(count),1:partition(count)))
allocate(pyramid_grid(count)%case(1:partition(count),1:partition(count)))
allocate(pyramid_grid(count)%transform(1:partition(count),1:partition(count)))
end do
end subroutine pyramid_creation
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine cartesian_creation()
! allocation arrays of volume, length...
integer :: count
do count = 1,level
allocate(cartesian_grid(count)%transform(1:level,1:level))
end do
end subroutine cartesian_creation
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_creation()
! allocation arrays of number_density, xHI...
integer :: count
do count = 1,level
allocate(pos_x_pos_y_pos_z_dir_z(count)%number_density(1:partition(count),1:partition(count)))
allocate(pos_x_pos_y_pos_z_dir_z(count)%xHI(1:partition(count),1:partition(count)))
allocate(pos_x_pos_y_pos_z_dir_z(count)%xHII(1:partition(count),1:partition(count)))
allocate(pos_x_pos_y_pos_z_dir_z(count)%temperature(1:partition(count),1:partition(count)))
end do
end subroutine domain_creation
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_volume()
integer :: k, n
real :: volume
! Volume on the same level are the same
do k = 1,level
n = partition(k)
volume = real(3*k*k-3*k+1)/real(3*n*n)
pyramid_grid(k)%volume = volume
end do
end subroutine pyramid_volume
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_cellsize()
integer :: i,j,k,n
!The radial length of the pyramid_grid
do k = 1,level
do i = 1,partition(k)
do j = 1,partition(k)
n = partition(k)
pyramid_grid(k)%cellsize(i,j) = sqrt(real(i*i+j*j+n*n)+0.5)/real(n)
end do
end do
end do
end subroutine pyramid_cellsize
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_distance()
integer :: i,j,k,n
!Distance of the center of pyramid_grid to the source
do k = 1,level
n = partition(k)
do i = 1,partition(k)
do j = 1,partition(k)
pyramid_grid(k)%distance(i,j) = sqrt(real(i*i+j*j+n*n-i-j+0.5))*(2*k-1)/real(2*n)
end do
end do
end do
end subroutine pyramid_distance
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_solid_angle()
integer :: i,j,k,n
n = partition(level)
! solid angle spanned by the pyramid_grid of the toppest level
do i = 1,partition(level)
do j = 1,partition(level)
if (i.ge.j) then
pyramid_grid(level)%solid_angle(i,j) = ss(i-1,j,i-1,n) - ss(i-1,j-1,i-1,n) - &
ss(i,j,i,n) + ss(i,j-1,i,n) + &
cc(j,j,i,n) - cc(j,j,i-1,n) + &
cc(j-1,j-1,i-1,n) - cc(j-1,j-1,i,n)
elseif (i.lt.j) then
pyramid_grid(level)%solid_angle(i,j) = ss(j-1,i,j-1,n) - ss(j-1,i-1,j-1,n) - &
ss(j,i,j,n) + ss(j,i-1,j,n) + &
cc(i,i,j,n) - cc(i,i,j-1,n) + &
cc(i-1,i-1,j-1,n) - cc(i-1,i-1,j,n)
endif
end do
end do
! solid angle spanned by the pyramid_grid of the other level
do k = level-1,1,-1
if (partition(k).eq.partition(k+1)) then
pyramid_grid(k)%solid_angle = pyramid_grid(k+1)%solid_angle
else
do i = 1,partition(k)
do j = 1,partition(k)
pyramid_grid(k)%solid_angle(i,j) = pyramid_grid(k+1)%solid_angle(2*i-1,2*j-1) + &
pyramid_grid(k+1)%solid_angle(2*i-1,2*j) + &
pyramid_grid(k+1)%solid_angle(2*i,2*j-1) + &
pyramid_grid(k+1)%solid_angle(2*i,2*j)
end do
end do
endif
enddo
end subroutine pyramid_solid_angle
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_number_density()
integer :: i,j,k,n
!Distance of the center of pyramid_grid to the source
do k = 1,level
do i = 1,partition(k)
do j = 1,partition(k)
pos_x_pos_y_pos_z_dir_z(k)%number_density(i,j) = number_density
end do
end do
end do
end subroutine domain_number_density
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_xHI()
integer :: i,j,k,n
! Ionization fraction of the pyramid_grids
do k = 1,level
do i = 1,partition(k)
do j = 1,partition(k)
pos_x_pos_y_pos_z_dir_z(k)%xHI(i,j) = xHI
end do
end do
end do
end subroutine domain_xHI
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_xHII()
integer :: i,j,k,n
! Ionization fraction of the pyramid_grids
do k = 1,level
do i = 1,partition(k)
do j = 1,partition(k)
pos_x_pos_y_pos_z_dir_z(k)%xHII(i,j) = xHII
end do
end do
end do
end subroutine domain_xHII
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_temperature()
integer :: i,j,k,n
! Temperature of the pyramid_grids
do k = 1,level
do i = 1,partition(k)
do j = 1,partition(k)
pos_x_pos_y_pos_z_dir_z(k)%temperature(i,j) = temperature
end do
end do
end do
end subroutine domain_temperature
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine pyramid_destruction()
integer :: k
! To deallocate of the arrays
do k = 1,level
deallocate(pyramid_grid(k)%volume)
deallocate(pyramid_grid(k)%cellsize)
deallocate(pyramid_grid(k)%distance)
deallocate(pyramid_grid(k)%solid_angle)
deallocate(pyramid_grid(k)%case)
deallocate(pyramid_grid(k)%transform)
end do
end subroutine pyramid_destruction
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine cartiesian_destruction()
integer :: k
! To deallocate of the arrays
do k = 1,level
deallocate(cartesian_grid(k)%transform)
end do
end subroutine cartiesian_destruction
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine domain_destruction()
integer :: k
! To deallocate of the arrays
do k = 1,level
deallocate(pos_x_pos_y_pos_z_dir_z(k)%number_density)
deallocate(pos_x_pos_y_pos_z_dir_z(k)%xHI)
deallocate(pos_x_pos_y_pos_z_dir_z(k)%xHII)
deallocate(pos_x_pos_y_pos_z_dir_z(k)%temperature)
end do
end subroutine domain_destruction
end module array