forked from cpraveen/hit3d
-
Notifications
You must be signed in to change notification settings - Fork 0
/
m_stats.f90
494 lines (341 loc) · 13.3 KB
/
m_stats.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
module m_stats
use m_parameters
real*8, allocatable :: e_spec(:), e_spec1(:), moments(:,:), sc_diss(:), sc_min(:), sc_max(:)
integer*8, allocatable :: hits(:), hits1(:)
real*8 :: energy, eps_v, eta, etakmax, enstrophy, re_lambda, uvar, x_length
real*8 :: lambda, re_lambda1, tau_e
real*8 :: sctmp
contains
!================================================================================
! This is a small module so the arrays get allocated in all parts of the code
!================================================================================
subroutine m_stats_init
implicit none
allocate(e_spec(kmax), e_spec1(kmax), moments(3+n_scalars,4), hits(kmax), hits1(kmax),&
sc_diss(n_scalars), sc_min(n_scalars), sc_max(n_scalars),stat=ierr)
write(out,*) "Stat allocated.", ierr
e_spec = zip
e_spec1 = zip
hits = 0
hits1 = 0
sc_diss = zip
sc_min = zip
sc_max = zip
return
end subroutine m_stats_init
!================================================================================
!================================================================================
subroutine stat_main
use m_parameters
use m_fields
implicit none
integer :: n
logical :: there2
call stat_velocity
if (int_scalars) call stat_scalars
if (int_scalars .and. myid.eq.0) then
! now outputting the scalar statistics
do n = 1, n_scalars
write(fname,"('sc',i2.2,'.gp')") n
inquire(file=fname, exist=there, opened=there2)
if (.not.there) then
open(100+n,file=fname,form='formatted')
write(100+n,'(A)') '# 1.itime 2.time 3.sc.diss 4. mean 5.variance 6.min 7.max'
end if
if(there.and..not.there2) then
open(100+n,file=fname,position='append')
end if
write(100+n,"(i7,10e15.6)") itime, time, sc_diss(n), moments(3+n,1:2), sc_min(n), sc_max(n)
call flush(100+n)
end do
end if
if (task_split) then
write(out,9000) ITIME, TIME
call flush(out)
end if
return
9000 format('ITIME=',i7,' TIME=',e15.7,' Stat files are written.')
end subroutine stat_main
!================================================================================
subroutine stat_velocity
implicit none
logical :: there2
integer :: k, n
! getting the enstrophy
call get_gradient_statistics
! getting the energy spectrum e_spec to the main process
call get_e_spec
! outputting the statistics into files
if (myid.eq.0) then
! getting the total energy
energy = sum(e_spec(1:kmax)) + half * e_spec(1)
! finding dissipation spectrum and total dissipation
do k = 1,kmax
e_spec1(k) = e_spec(k) * real(k**2,8) * two * nu
end do
eps_v = sum(e_spec1(1:kmax)) + half*e_spec1(1)
! finding Kolmogorov scale
eta = (nu**3/eps_v)**0.25
etakmax = eta * real(kmax,8)
! variance
uvar = two/three*energy
! integral length scale
sctmp = zip
do k = 1, kmax
sctmp = sctmp + e_spec(k) / real(k,8)
end do
sctmp = sctmp + half*e_spec(1)
x_length = PI / two * sctmp / uvar
! Taylor microscale
lambda = sqrt(15.d0 * uvar * nu / eps_v)
! Taylor-Reynolds number
re_lambda = uvar*sqrt(15.d0/eps_v*RE)
re_lambda1 = sqrt(uvar)*lambda / nu
! Eddy turnover time
tau_e = x_length / sqrt(uvar)
! outputting all this in the stat1 file
inquire(file='stat1.gp', exist=there, opened=there2)
if (.not.there) then
open(69,file='stat1.gp',form='formatted')
write(69,*) '# 1.itime 2.time 3.energy 4.diss 5.eta 6.enstrophy 7.R_lambda'
end if
if(there.and..not.there2) then
open(69,file='stat1.gp',position='append')
end if
write(69,"(i8,20e15.6)") itime, time, energy, eps_v, eta, enstrophy, re_lambda
call flush(69)
! outputting all this in the stat2 file
inquire(file='stat2.gp', exist=there, opened=there2)
if (.not.there) then
open(70,file='stat2.gp',form='formatted')
write(70,'(A)') '# 1.itime 2.time 3.int LS 4. lambda 5.R_lambda1 6.tau_e 7.etakmax'
end if
if(there.and..not.there2) then
open(70,file='stat2.gp',position='append')
end if
write(70,"(i8,20e15.6)") itime, time, x_length, lambda, re_lambda1, tau_e, etakmax
call flush(70)
! outputting the energy spectrum
open(900,file='es.gp',position='append')
write(900,"()")
write(900,"()")
write(900,"('# ITIME=',i7,' TIME=',e17.8)") ITIME, TIME
do k = 1,kmax !min(kmax,nx/3)
write(900,"(i4,2e15.6)") k,e_spec(k), e_spec1(k)!, hits(k)
end do
close(900)
end if
return
end subroutine stat_velocity
!================================================================================
!================================================================================
!================================================================================
subroutine get_e_spec
use m_io
use m_fields
use x_fftw
implicit none
real*8 :: sc_rad1, sc_rad2, fac
integer :: i, j, k, n_shell
! need this normalization factor because the FFT is unnormalized
fac = one / real(nx*ny*nz_all)**2
e_spec1 = zip
e_spec = zip
hits = 0
hits1 = 0
! assembling the total energy in each shell and number of hits in each shell
do k = 1,nz
do j = 1,ny
do i = 1,nx
n_shell = nint(sqrt(real(akx(i)**2 + aky(k)**2 + akz(j)**2, 4)))
if (n_shell .gt. 0 .and. n_shell .le. kmax) then
if (akx(i).le.nx/3 .and. aky(k).le.nx/3 .and. akz(j).le.nx/3) then
hits1(n_shell) = hits1(n_shell) + 1
e_spec1(n_shell) = e_spec1(n_shell) + &
fac * (fields(i,j,k,1)**2 + fields(i,j,k,2)**2 + fields(i,j,k,3)**2)
end if
end if
end do
end do
end do
! reducing the number of hits and energy to two arrays on master node
count = kmax
call MPI_REDUCE(hits1,hits,count,MPI_INTEGER8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
call MPI_REDUCE(e_spec1,e_spec,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
! now the master node counts the energy density in each shell
if (myid.eq.0) then
! 4/3 PI is prefactor for the shell volume, and 1/2 is the
! factor from the energy (u^2+v^2+w^2)/2, omitted in the above cycle
fac = four/three * PI / two
do k = 1,kmax
sc_rad1 = real(k,8) + half
sc_rad2 = real(k,8) - half
if (k.eq.1) sc_rad2 = 0.d0
if (hits(k).gt.0) then
e_spec(k) = e_spec(k) / hits(k) * fac * (sc_rad1**3 - sc_rad2**3)
else
e_spec(k) = zip
end if
end do
end if
return
end subroutine get_e_spec
!================================================================================-
!================================================================================-
subroutine get_gradient_statistics
use m_fields
use m_work
use x_fftw
implicit none
integer :: i
real*8 :: fac
! normalization factor
fac = one / real(nx*ny*nz_all,8)
do i = 1,3
wrk(:,:,:,i) = fields(:,:,:,i)
end do
! Taking derivatives
call x_derivative(3,'y',6)
call x_derivative(3,'x',5)
call x_derivative(2,'z',4)
call x_derivative(2,'x',3)
call x_derivative(1,'z',2)
call x_derivative(1,'y',1)
!------------------------------------------------------------
! getting vorticity and enstrophy
!------------------------------------------------------------
wrk(:,:,:,3) = wrk(:,:,:,3) - wrk(:,:,:,1) ! omega_3 = v_x - u_y
wrk(:,:,:,2) = wrk(:,:,:,2) - wrk(:,:,:,5) ! omega_2 = u_z - w_x
wrk(:,:,:,1) = wrk(:,:,:,6) - wrk(:,:,:,4) ! omega_1 = w_y - v_z
call xFFT3d(-1,1)
call xFFT3d(-1,2)
call xFFT3d(-1,3)
! getting mean enstrophy
wrk(:,:,:,0) = wrk(:,:,:,1)**2 + wrk(:,:,:,2)**2 + wrk(:,:,:,3)**2
sctmp = sum(wrk(1:nx,:,:,0)) * fac
count = 1
call MPI_REDUCE(sctmp,enstrophy,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
!------------------------------------------------------------
return
end subroutine get_gradient_statistics
!================================================================================
!================================================================================
!================================================================================
subroutine stat_scalars
use m_openmpi
use m_fields
use m_work
use x_fftw
implicit none
integer :: i, j, k, n
real*8 :: q1, q2, fac
if (.not. int_scalars) return
! getting the spectra of the scalar variances
call get_scalar_spectra
! scaling factor
fac = one / real(nx*ny*nz_all)
! --- Calculating moments of scalars
do n = 1, n_scalars
! putting the scalar in wrk0
wrk(:,:,:,0) = fields(:,:,:,3+n)
! taking derivatives
call x_derivative(0,'x',1)
call x_derivative(0,'y',2)
call x_derivative(0,'z',3)
! converting the derivatives to X-space
call xFFT3d(-1,1)
call xFFT3d(-1,2)
call xFFT3d(-1,3)
! getting the dissipation rate of the variance
wrk(:,:,:,4) = wrk(:,:,:,1)**2 + wrk(:,:,:,2)**2 + wrk(:,:,:,3)**2
q1 = two * pe(n) * sum(wrk(1:nx,:,:,4)) * fac
count = 1
call MPI_REDUCE(q1,q2,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
if (myid.eq.0) sc_diss(n) = q2
! converting the scalar itself to X-space
call xFFT3d(-1,0)
! First moment - mean
q1 = sum(wrk(1:nx,:,:,0)) * fac
count = 1
call MPI_REDUCE(q1,q2,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
if (myid.eq.0) moments(3+n,1) = q2
! Second moment - variance
q1 = sum(wrk(1:nx,:,:,0)**2) * fac
count = 1
call MPI_REDUCE(q1,q2,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
if (myid.eq.0) moments(3+n,2) = q2 - moments(3+n,1)**2
! Min and max of the scalar
q1 = minval(wrk(1:nx,:,:,0))
q2 = maxval(wrk(1:nx,:,:,0))
count = 1
call MPI_REDUCE(q1,sc_min(n),count,MPI_REAL8,MPI_MIN,0,MPI_COMM_TASK,mpi_err)
call MPI_REDUCE(q2,sc_max(n),count,MPI_REAL8,MPI_MAX,0,MPI_COMM_TASK,mpi_err)
end do
return
end subroutine stat_scalars
!================================================================================
!================================================================================
subroutine get_scalar_spectra
use m_io
use m_fields
use x_fftw
implicit none
real*8 :: sc_rad1, sc_rad2, fac
integer :: i, j, k, n, n_shell
! cycle over the scalars
do n = 1,n_scalars
! need this normalization factor because the FFT is unnormalized
fac = one / real(nx*ny*nz_all)**2
! using the neergy spectra arrays to keep the scalar spectra
e_spec1 = zip
e_spec = zip
hits = 0
hits1 = 0
! assembling the total scalar energy in each shell and number of hits in each shell
do k = 1,nz
do j = 1,ny
do i = 1,nx
n_shell = nint(sqrt(real(akx(i)**2 + aky(k)**2 + akz(j)**2, 4)))
if (n_shell .gt. 0 .and. n_shell .le. kmax) then
hits1(n_shell) = hits1(n_shell) + 1
e_spec1(n_shell) = e_spec1(n_shell) + fac * fields(i,j,k,3+n)**2
end if
end do
end do
end do
! reducing the number of hits and energy to two arrays on master node
count = kmax
call MPI_REDUCE(hits1,hits,count,MPI_INTEGER8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
call MPI_REDUCE(e_spec1,e_spec,count,MPI_REAL8,MPI_SUM,0,MPI_COMM_TASK,mpi_err)
! now the master node counts the energy density in each shell
master_node: if (myid.eq.0) then
! 4/3 PI is prefactor for the shell volume, and 1/2 is the
! factor from the energy (u^2+v^2+w^2)/2, omitted in the above cycle
fac = four/three * PI / two
do k = 1,kmax
sc_rad1 = real(k,8) + half
sc_rad2 = real(k,8) - half
if (k.eq.1) sc_rad2 = 0.d0
if (hits(k).gt.0) then
e_spec(k) = e_spec(k) / hits(k) * fac * (sc_rad1**3 - sc_rad2**3)
else
e_spec(k) = zip
end if
end do
! now the master node puts the scalar energy in the file es_sc##.gp
write(fname,"('es_',i2.2,'.gp')") n
open(900,file=fname,position='append')
write(900,"()")
write(900,"()")
write(900,"('# ITIME=',i7,' TIME=',e17.8)") ITIME, TIME
do k = 1,kmax
write(900,"(i4,4e15.6)") k,e_spec(k)
end do
close(900)
end if master_node
end do
return
end subroutine get_scalar_spectra
!================================================================================
!================================================================================
end module m_stats