forked from PhysicsofFluids/AFiD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSlabDumpRoutines.F90
133 lines (108 loc) · 3.62 KB
/
SlabDumpRoutines.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: SlabDumper.F90 !
! CONTAINS: subroutines SlabDumper, InitializeSlabDump !
! DumpSingleSlab !
! !
! PURPOSE: Auxiliary routines used for memory allocs !
! and memory freeing !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine SlabDumper
use param
use local_arrays, only: temp,vx,vy,vz
use stat3_param
use decomp_2d, only: xstart,xend
implicit none
integer :: i,j,m
real,dimension(xstart(2):xend(2),xstart(3):xend(3)) :: &
& vxcc,vycc,vzcc,tempcc
character*70 :: filnam
character*1 :: charm
!EP Slabs
!EP cell center only vx
do m=1,9
!$OMP PARALLEL DO DEFAULT(SHARED) &
!$OMP PRIVATE(i,j)
do i=xstart(3),xend(3)
do j=xstart(2),xend(2)
vxcc(j,i) = (vx(kslab(m),j,i)+vx(kslab(m)+1,j,i))*0.5
vycc(j,i) = vy(kslab(m),j,i)
vzcc(j,i) = vz(kslab(m),j,i)
tempcc(j,i) = temp(kslab(m),j,i)
enddo
enddo
!$OMP END PARALLEL DO
write(charm,28) m
28 format(i1.1)
filnam='slab'//charm//'vx_'
call DumpSingleSlab(vxcc,filnam)
filnam='slab'//charm//'vy_'
call DumpSingleSlab(vycc,filnam)
filnam='slab'//charm//'vz_'
call DumpSingleSlab(vzcc,filnam)
filnam='slab'//charm//'temp_'
call DumpSingleSlab(tempcc,filnam)
enddo
return
end subroutine SlabDumper
!===========================================================================
subroutine InitializeSlabDump
use param
use stat3_param
implicit none
integer :: i,k,j
real :: xmloc
character(len=4) :: dummy
!EP Read from stst3.in
open(unit=19,file='stst3.in',status='old')
read(19,301) dummy
read(19,*) (xslab(i),i=2,9)
301 format(a4)
close(19)
!EP Compute which kslab corresponds to which xslab
kslab = 2
do k=2,nxm
xmloc=xm(k)
do j=2,9
if(xm(k).gt.xslab(j).and.xm(k-1).lt.xslab(j)) then
kslab(j) = k
endif
enddo
enddo
!EP Write probe and slab locations
if (ismaster) then
open(unit=23,file='stst3locs.out',status='unknown')
rewind(23)
write(23,*) (kslab(i),i=1,9)
close(23)
endif
return
end subroutine InitializeSlabDump
!==================================================================
subroutine DumpSingleSlab(var,filnam)
USE param
use mpih
USE hdf5
use decomp_2d, only: xstart,xend
IMPLICIT none
real, intent(in) :: var(xstart(2):xend(2) &
& ,xstart(3):xend(3))
real :: tprfi
integer :: itime
character*70,intent(in) :: filnam
character*70 :: namfile,dsetname
character*8 :: ipfi
tprfi = 1/tout
itime=nint(time*tprfi)
write(ipfi,82)itime
82 format(i8.8)
namfile=trim('./stst3/'//trim(filnam)//trim(ipfi)//'.h5')
dsetname = trim('var')
call HdfWriteReal2D(dsetname,namfile,var)
if(ismaster) then
dsetname = trim('time')
call HdfSerialWriteRealScalar(dsetname,namfile,time)
endif
return
end subroutine DumpSingleSlab