forked from PhysicsofFluids/AFiD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
param.F90
177 lines (161 loc) · 7.03 KB
/
param.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
! Declaration of global variables
!***********************************************************
module param
implicit none
!==========================================================
! read from input file bou.in
!==========================================================
integer :: nx, ny, nz
integer :: nsst, nread, ntst, ireset
real :: walltimemax,tout,tmax
real :: alx3,str3
integer :: istr3
real :: ylen,zlen
real :: ray,pra,dt,resid
integer :: inslws,inslwn
integer :: starea,tsta
real :: dtmin,dtmax,limitCFL,limitVel
integer :: nson,idtv
!=================================================
! end of input file
!=================================================
real :: time
!******* Grid parameters**************************
real :: dx,dy,dz,dxq,dyq,dzq
!
real, allocatable, dimension(:) :: xc,xm
real, allocatable, dimension(:) :: yc,ym
real, allocatable, dimension(:) :: zc,zm
real, allocatable, dimension(:) :: g3rc,g3rm
!====================================================
!******* QUANTITIES FOR DERIVATIVES******************
real, allocatable, dimension(:) :: udx3c,udx3m
!==========================================================
!******* Grid indices**************************************
integer, allocatable, dimension(:) :: kmc,kpc,kmv,kpv
!===========================================================
!******* Metric coefficients *******************************
real, allocatable, dimension(:) :: ap3ck,ac3ck,am3ck
real, allocatable, dimension(:) :: ap3sk,ac3sk,am3sk
real, allocatable, dimension(:) :: ap3ssk,ac3ssk,am3ssk
!============================================================
!******* Variables for FFTW and Poisson solver****************
real, allocatable, dimension(:) :: ak2,ap
real, allocatable, dimension(:) :: ak1,ao
real, allocatable, dimension(:) :: amphk,acphk,apphk
!===========================================================
!******* Other variables ***********************************
integer :: nxm, nym, nzm
real :: ren, pec
real :: pi
real :: al,ga,ro
real :: beta
real :: qqmax,qqtot
real :: re
real :: tempmax,tempmin,tempm
integer :: ntime
integer, parameter:: ndv=3
real, dimension(1:ndv) :: vmax
real, dimension(1:3) :: gam,rom,alm
real, allocatable, dimension(:,:) :: tempbp,temptp
logical :: dumpslabs=.false.
logical :: statcal=.false.
logical :: disscal=.false.
logical :: readflow=.false.
logical :: readstats=.false.
logical :: ismaster=.false.
logical :: resetlogstime=.false.
logical :: variabletstep=.true.
integer :: lvlhalo=1
end module param
!************* End of param module******************************
!===============================================================
!******* 2D arrays, dynamically allocated by each process*******
module local_arrays
use param
implicit none
real,allocatable,dimension(:,:,:) :: vx,vy,vz
real,allocatable,dimension(:,:,:) :: pr,temp,rhs
real,allocatable,dimension(:,:,:) :: rux,ruy,ruz,rutemp
real,allocatable,dimension(:,:,:) :: dph,qcap,dq,hro,dphhalo
end module local_arrays
!===============================================================
module stat_arrays
implicit none
real,allocatable, dimension(:) :: vz_me,vz_rms
real,allocatable, dimension(:) :: vy_me,vx_me,vy_rms,vx_rms
real,allocatable, dimension(:) :: temp_me,temp_rms
real, allocatable,dimension(:) :: disste,dissth,tempvx_me
integer :: nstatsamples
end module stat_arrays
!=====================================================
module stat3_param
implicit none
integer :: kslab(1:9)
real :: xslab(1:9)
end module stat3_param
!=====================================================
module mpih
implicit none
include 'mpif.h'
integer :: ierr
integer, parameter :: master=0
integer :: MDP = MPI_DOUBLE_PRECISION
end module mpih
!====================================================
module fftw_params
! use param, only: m2m,m2mh,m1m
use iso_c_binding
type, bind(C) :: fftw_iodim
integer(C_INT) n, is, os
end type fftw_iodim
interface
type(C_PTR) function fftw_plan_guru_dft(rank,dims, &
howmany_rank,howmany_dims,in,out,sign,flags) &
bind(C, name='fftw_plan_guru_dft')
import
integer(C_INT), value :: rank
type(fftw_iodim), dimension(*), intent(in) :: dims
integer(C_INT), value :: howmany_rank
type(fftw_iodim), dimension(*), intent(in) :: howmany_dims
complex(C_DOUBLE_COMPLEX), dimension(*), intent(out) :: in
complex(C_DOUBLE_COMPLEX), dimension(*), intent(out) :: out
integer(C_INT), value :: sign
integer(C_INT), value :: flags
end function fftw_plan_guru_dft
type(C_PTR) function fftw_plan_guru_dft_r2c(rank,dims, &
howmany_rank,howmany_dims,in,out,flags) &
bind(C, name='fftw_plan_guru_dft_r2c')
import
integer(C_INT), value :: rank
type(fftw_iodim), dimension(*), intent(in) :: dims
integer(C_INT), value :: howmany_rank
type(fftw_iodim), dimension(*), intent(in) :: howmany_dims
real(C_DOUBLE), dimension(*), intent(out) :: in
complex(C_DOUBLE_COMPLEX), dimension(*), intent(out) :: out
integer(C_INT), value :: flags
end function fftw_plan_guru_dft_r2c
type(C_PTR) function fftw_plan_guru_dft_c2r(rank,dims, &
howmany_rank,howmany_dims,in,out,flags) &
bind(C, name='fftw_plan_guru_dft_c2r')
import
integer(C_INT), value :: rank
type(fftw_iodim), dimension(*), intent(in) :: dims
integer(C_INT), value :: howmany_rank
type(fftw_iodim), dimension(*), intent(in) :: howmany_dims
complex(C_DOUBLE_COMPLEX), dimension(*), intent(out) :: in
real(C_DOUBLE), dimension(*), intent(out) :: out
integer(C_INT), value :: flags
end function fftw_plan_guru_dft_c2r
end interface
integer FFTW_PATIENT, FFTW_FORWARD, FFTW_BACKWARD,FFTW_ESTIMATE
parameter (FFTW_PATIENT=32)
parameter (FFTW_ESTIMATE=64)
parameter (FFTW_FORWARD=-1)
parameter (FFTW_BACKWARD=1)
type(C_PTR) :: fwd_guruplan_y,bwd_guruplan_y
type(C_PTR) :: fwd_guruplan_z,bwd_guruplan_z
logical :: planned=.false.
real,allocatable,dimension(:,:,:) :: ry1,rz1
complex,allocatable,dimension(:,:,:) :: cy1,cz1,dphc
end module fftw_params