forked from ImperialCollegeLondon/WInc3D
-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathvariables.f90
277 lines (247 loc) · 12.5 KB
/
variables.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
!################################################################################
!This file is part of Incompact3d.
!
!Incompact3d
!Copyright (c) 2012 Eric Lamballais and Sylvain Laizet
!eric.lamballais@univ-poitiers.fr / sylvain.laizet@gmail.com
!
! Incompact3d is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation.
!
! Incompact3d is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with the code. If not, see <http://www.gnu.org/licenses/>.
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------
! We kindly request that you cite Incompact3d in your publications and
! presentations. The following citations are suggested:
!
! 1-Laizet S. & Lamballais E., 2009, High-order compact schemes for
! incompressible flows: a simple and efficient method with the quasi-spectral
! accuracy, J. Comp. Phys., vol 228 (15), pp 5989-6015
!
! 2-Laizet S. & Li N., 2011, Incompact3d: a powerful tool to tackle turbulence
! problems with up to 0(10^5) computational cores, Int. J. of Numerical
! Methods in Fluids, vol 67 (11), pp 1735-1757
!################################################################################
module var
use decomp_2d
USE param
! define all major arrays here
real(mytype), save, allocatable, dimension(:,:,:) :: ux1, ux2, ux3,po3,dv3,pp3,ucx1,ucy1,ucz1,uxf1,uxf2,uxf3
real(mytype), save, allocatable, dimension(:,:,:) :: uy1, uy2, uy3,uyf1,uyf2,uyf3
real(mytype), save, allocatable, dimension(:,:,:) :: uz1, uz2, uz3,uzf1,uzf2,uzf3
real(mytype), save, allocatable, dimension(:,:,:) :: phi1, phi2, phi3, phif1,phif2,phif3
real(mytype), save, allocatable, dimension(:,:,:) :: gx1, gy1, gz1, hx1, hy1, hz1, phis1,phiss1
real(mytype), save, allocatable, dimension(:,:,:) :: px1, py1, pz1
real(mytype), save, allocatable, dimension(:,:,:) :: ep1
real(mytype), save, allocatable, dimension(:,:,:) :: nut1, shrt_coeff
real(mytype), save, allocatable, dimension(:,:,:,:) :: uxt,uyt,uzt
! Define inflow/outflow file
real(mytype), save, allocatable, dimension(:,:,:) :: ux_recOutflow, uy_recOutflow, uz_recOutflow
real(mytype), save, allocatable, dimension(:,:,:) :: ux_inflow, uy_inflow, uz_inflow, phi_inflow
! define the Momentum Source arrays for the three directions
real(mytype), save, allocatable, dimension(:,:,:) :: FTx, FTy, FTz
real(mytype), save, allocatable, dimension(:,:,:) :: Fdiscx, Fdiscy, Fdiscz, GammaDisc
real(mytype), save, allocatable, dimension(:,:,:) :: Ftripx, randomtrip
! define Probe location and velocity/scalar quantities
real(mytype), save, allocatable, dimension(:) :: xprobe,yprobe,zprobe,uprobe,vprobe,wprobe
real(mytype), save, allocatable, dimension(:) :: uprobe_part,vprobe_part,wprobe_part
!arrays for statistic collection
real(mytype), save, allocatable, dimension(:,:,:) :: umean,vmean,wmean,uumean,vvmean,wwmean,uvmean,uwmean,vwmean,tmean
real(mytype), save, allocatable, dimension(:,:,:) :: phimean, phiphimean
real(mytype), save, allocatable, dimension(:,:,:) :: vortxmean,vortymean,vortzmean,vortxxmean,vortyymean, vortzzmean
real(mytype), save, allocatable, dimension(:,:,:) :: vortmean,qcritmean
real(mytype), save, allocatable, dimension(:,:,:) :: pmean, upmean,vpmean, wpmean, ppmean, uvisu1
real(mytype), save, allocatable, dimension(:,:,:) :: sgszmean,sgsxmean,sgsymean,eadvxmean,eadvymean,eadvzmean
real(mytype), save, allocatable, dimension(:,:,:) :: divdiva,curldiva,tauxymean
!arrays for visualization
real(mytype), save, allocatable, dimension(:,:,:) :: uvisu
! define all work arrays here
real(mytype), save, allocatable, dimension(:,:,:) :: ta1,tb1,tc1,td1,&
te1,tf1,tg1,th1,ti1,tj1,di1
real(mytype), save, allocatable, dimension(:,:,:) :: ta2,tb2,tc2,td2,&
te2,tf2,tg2,th2,ti2,tj2,di2
real(mytype), save, allocatable, dimension(:,:,:) :: ta3,tb3,tc3,td3,&
te3,tf3,tg3,th3,ti3,di3
integer, save :: nxmsize, nymsize, nzmsize
contains
subroutine init_variables
TYPE(DECOMP_INFO), save :: ph ! decomposition object
if (nclx==0) then
nxmsize = xsize(1)
else
nxmsize = xsize(1) -1
endif
if (ncly==0) then
nymsize = ysize(2)
else
nymsize = ysize(2) -1
endif
if (nclz==0) then
nzmsize = zsize(3)
else
nzmsize = zsize(3) -1
endif
call decomp_info_init(nxmsize, nymsize, nzmsize, ph)
!X PENCILS
call alloc_x(ux1, opt_global=.true.)
call alloc_x(uy1, opt_global=.true.)
#ifndef TWOD
call alloc_x(uz1, opt_global=.true.)
call alloc_x(pz1, opt_global=.true.)
! Allocate the Momentum Source term in the x direction
#else
allocate (uz1(1,1,1))
allocate (pz1(1,1,1))
#endif
call alloc_x(px1, opt_global=.true.)
call alloc_x(py1, opt_global=.true.)
call alloc_x(phi1, opt_global=.true.)
call alloc_x(gx1);call alloc_x(gy1);call alloc_x(gz1);call alloc_x(phis1)
call alloc_x(hx1);call alloc_x(hy1);call alloc_x(hz1);call alloc_x(phiss1)
call alloc_x(ta1);call alloc_x(tb1);call alloc_x(tc1)
call alloc_x(td1);call alloc_x(te1);call alloc_x(tf1)
call alloc_x(tg1);call alloc_x(th1);call alloc_x(ti1)
call alloc_x(tj1);call alloc_x(di1);call alloc_x(ep1)
call alloc_x(nut1);call alloc_x(ucx1);call alloc_x(ucy1);call alloc_x(ucz1);
call alloc_x(uxf1);call alloc_x(uyf1);call alloc_x(uzf1);call alloc_x(phif1);
call alloc_x(shrt_coeff);
allocate(sx(xsize(2),xsize(3)),vx(xsize(2),xsize(3)))
allocate(fisx(xsize(2),xsize(3)),fivx(xsize(2),xsize(3)))
!inflow/ouflow 2d arrays
allocate(bxx1(xsize(2),xsize(3)),bxy1(xsize(2),xsize(3)))
allocate(bxz1(xsize(2),xsize(3)),bxxn(xsize(2),xsize(3)))
allocate(bxyn(xsize(2),xsize(3)),bxzn(xsize(2),xsize(3)))
allocate(bxo(xsize(2),xsize(3)),byo(xsize(2),xsize(3)))
allocate(bzo(xsize(2),xsize(3)))
allocate(byx1(xsize(1),xsize(3)),byy1(xsize(1),xsize(3)))
allocate(byz1(xsize(1),xsize(3)),byxn(xsize(1),xsize(3)))
allocate(byyn(xsize(1),xsize(3)),byzn(xsize(1),xsize(3)))
allocate(bzx1(xsize(1),xsize(2)),bzy1(xsize(1),xsize(2)))
allocate(bzz1(xsize(1),xsize(2)),bzxn(xsize(1),xsize(2)))
allocate(bzyn(xsize(1),xsize(2)),bzzn(xsize(1),xsize(2)))
!pre_correc 2d array
allocate(dpdyx1(xsize(2),xsize(3)),dpdyxn(xsize(2),xsize(3)))
allocate(dpdzx1(xsize(2),xsize(3)),dpdzxn(xsize(2),xsize(3)))
allocate(dpdxy1(xsize(1),xsize(3)),dpdxyn(xsize(1),xsize(3)))
allocate(dpdzy1(xsize(1),xsize(3)),dpdzyn(xsize(1),xsize(3)))
allocate(dpdxz1(xsize(1),xsize(2)),dpdxzn(xsize(1),xsize(2)))
allocate(dpdyz1(xsize(1),xsize(2)),dpdyzn(xsize(1),xsize(2)))
if (iin==3) then
allocate(ux_inflow(NTimeSteps,xsize(2),xsize(3)))
allocate(uy_inflow(NTimeSteps,xsize(2),xsize(3)))
allocate(uz_inflow(NTimeSteps,xsize(2),xsize(3)))
endif
if (ioutflow==1) then
allocate(ux_recOutflow(NTimeSteps,xsize(2),xsize(3)))
allocate(uy_recOutflow(NTimeSteps,xsize(2),xsize(3)))
allocate(uz_recOutflow(NTimeSteps,xsize(2),xsize(3)))
endif
! Allocate Momentum Source Terms
if (ialm==1) then
allocate(FTx(xsize(1),xsize(2),xsize(3)))
allocate(FTy(xsize(1),xsize(2),xsize(3)))
allocate(FTz(xsize(1),xsize(2),xsize(3)))
endif
if (iadm==1) then
allocate(Fdiscx(xsize(1),xsize(2),xsize(3)))
allocate(Fdiscy(xsize(1),xsize(2),xsize(3)))
allocate(Fdiscz(xsize(1),xsize(2),xsize(3)))
allocate(Gammadisc(xsize(1),xsize(2),xsize(3)))
endif
if(itripping==1.or.itripping==2) then
allocate(randomtrip(xsize(1),xsize(2),xsize(3)))
allocate(Ftripx(xsize(1),xsize(2),xsize(3)))
endif
!arrays for statistic collection!pay attention to the size!
allocate (umean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (wmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (uumean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vvmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (wwmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (uvmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (uwmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vwmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (tmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (tauxymean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
if (iscalar==1) then
allocate (phimean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (phiphimean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
else
allocate (phimean(1,1,1))
allocate (phiphimean(1,1,1))
endif
allocate (vortxmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortymean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortzmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortxxmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortyymean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortzzmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vortmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (qcritmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (pmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (upmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (vpmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (wpmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (ppmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (uvisu1(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
!arrays for LES Models
if(jLES==7 .OR. jLES==8) then
allocate (uxt(xstart(1):xend(1),xstart(2):xend(2),xstart(3):xend(3),25))
allocate (uyt(xstart(1):xend(1),xstart(2):xend(2),xstart(3):xend(3),25))
allocate (uzt(xstart(1):xend(1),xstart(2):xend(2),xstart(3):xend(3),25))
endif
if(jLES .ge. 2) then
allocate (sgsxmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (sgsymean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (sgszmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
if(jLES .gt. 3) then
allocate (divdiva(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (curldiva(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
endif
if(jADV==1) then
allocate (eadvxmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (eadvymean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
allocate (eadvzmean(xstS(1):xenS(1),xstS(2):xenS(2),xstS(3):xenS(3)))
endif
endif
!arrays for visualization!pay attention to the size!
allocate (uvisu(xstV(1):xenV(1),xstV(2):xenV(2),xstV(3):xenV(3)))
!Y PENCILS
call alloc_y(ux2);call alloc_y(uy2);call alloc_y(uz2)
call alloc_y(ta2);call alloc_y(tb2);call alloc_y(tc2)
call alloc_y(td2);call alloc_y(te2);call alloc_y(tf2)
call alloc_y(tg2);call alloc_y(th2);call alloc_y(ti2)
call alloc_y(tj2);call alloc_y(di2);call alloc_y(phi2)
call alloc_y(uxf2);call alloc_y(uyf2);call alloc_y(uzf2); call alloc_y(phif2)
allocate(sy(ysize(1),ysize(3)),vy(ysize(1),ysize(3)))
allocate(fisy(ysize(1),ysize(3)),fivy(ysize(1),ysize(3)))
!Z PENCILS
call alloc_z(ux3);call alloc_z(uy3);call alloc_z(uz3)
call alloc_z(ta3);call alloc_z(tb3);call alloc_z(tc3)
call alloc_z(td3);call alloc_z(te3);call alloc_z(tf3)
call alloc_z(tg3);call alloc_z(th3);call alloc_z(ti3)
call alloc_z(di3);call alloc_z(phi3);
call alloc_z(uxf3);call alloc_z(uyf3);call alloc_z(uzf3); call alloc_z(phif3)
allocate(sz(zsize(1),zsize(2)),vz(zsize(1),zsize(2)))
allocate(fisz(zsize(1),zsize(2)),fivz(zsize(1),zsize(2)))
! if all periodic
! allocate (pp3(ph%zst(1):ph%zen(1),ph%zst(2):ph%zen(2),ph%zst(3):ph%zen(3)))
! allocate (dv3(ph%zst(1):ph%zen(1),ph%zst(2):ph%zen(2),ph%zst(3):ph%zen(3)))
! allocate (po3(ph%zst(1):ph%zen(1),ph%zst(2):ph%zen(2),ph%zst(3):ph%zen(3)))
call alloc_z(pp3,ph,.true.)
call alloc_z(dv3,ph,.true.)
call alloc_z(po3,ph,.true.)
!z_modes=int(zlz/zs_tr)
!allocate(h_coeff(z_modes))
!allocate(h_nxt(xsize(3)),h_i(xsize(3)))
return
end subroutine init_variables
end module var