-
Notifications
You must be signed in to change notification settings - Fork 1
Expand file tree
/
Copy pathvariables.f90
More file actions
160 lines (139 loc) · 6.16 KB
/
variables.f90
File metadata and controls
160 lines (139 loc) · 6.16 KB
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
!################################################################################
!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 variables
USE param
! define all major arrays here
real(mytype), save, allocatable, dimension(:,:,:) :: ux1, ux2, ux3,po3,dv3,pp3
real(mytype), save, allocatable, dimension(:,:,:) :: uy1, uy2, uy3
real(mytype), save, allocatable, dimension(:,:,:) :: uz1, uz2, uz3
real(mytype), save, allocatable, dimension(:,:,:) :: phi1, phi2, phi3
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
!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,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.)
#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(di1);call alloc_x(ep1)
allocate(sx(xsize(2),xsize(3)),vx(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)))
!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)
allocate(sy(ysize(1),ysize(3)),vy(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)
allocate(sz(zsize(1),zsize(2)),vz(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.)
return
end subroutine init_variables
end module var