-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmean-2d.f90
301 lines (247 loc) · 9.63 KB
/
mean-2d.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
!.... Use a LNS3D field file for nonparallel stability
!============================================================================!
module meanflow
!============================================================================!
integer :: nxm, nym, nzm, ndofm
!.... flow variables and derivatives
real, allocatable :: vm(:,:,:)
real, allocatable :: g1vm(:,:,:), g2vm(:,:,:)
real, allocatable :: g11vm(:,:,:), g12vm(:,:,:), g22vm(:,:,:)
real, allocatable :: g1vl(:), g2vl(:), g11vl(:), g12vl(:), g22vl(:)
!.... mesh and metrics
real :: dxi, deta
real, allocatable :: xm(:,:), ym(:,:)
real, allocatable :: m1(:,:), m2(:,:), n1(:,:), n2(:,:), &
m11(:,:), m12(:,:), m22(:,:), &
n11(:,:), n12(:,:), n22(:,:)
!.... mean splines
real, allocatable :: vms(:,:), g1vms(:,:), g2vms(:,:), &
g11vms(:,:), g12vms(:,:), g22vms(:,:)
end module meanflow
!============================================================================!
subroutine initmean(imean)
!============================================================================!
use meanflow
implicit none
integer :: imean
!.... derivative parameters
integer :: optx=-1, opty=-1
logical :: xper=.false., yper=.false.
logical :: lsym=.false. , rsym=.false., bsym=.false., tsym=.false.
!.... local
integer :: lstep
real :: time, Ma, Re, Pr, gamma, cv
integer :: i, j, k, idof
real :: tmp
logical :: useIJ=.true.
!============================================================================!
!.... read the grid file
open(imean,file='lstx.dat',form='unformatted',status='old',err=100)
read(imean) nxm, nym, nzm
allocate( xm(nym,nxm), ym(nym,nxm) )
read(imean) (((xm(j,i), i = 1, nxm), j = 1, nym), k = 1, nzm), &
(((ym(j,i), i = 1, nxm), j = 1, nym), k = 1, nzm), &
((( tmp, i = 1, nxm), j = 1, nym), k = 1, nzm)
close(imean)
dxi = 1.0 / float(nxm-1)
deta = 1.0 / float(nym-1)
!.... read the metrics
allocate (m1(nym,nxm), m2(nym,nxm), n1(nym,nxm), n2(nym,nxm), &
m11(nym,nxm), m12(nym,nxm), m22(nym,nxm), &
n11(nym,nxm), n12(nym,nxm), n22(nym,nxm) )
open (imean,file='lstm.dat',form='unformatted', status='old',err=200)
if (useIJ) then
read(imean) ((m1(j,i), i=1,nxm), j=1,nym), &
((m2(j,i), i=1,nxm), j=1,nym), &
((n1(j,i), i=1,nxm), j=1,nym), &
((n2(j,i), i=1,nxm), j=1,nym), &
((m11(j,i), i=1,nxm), j=1,nym),&
((m12(j,i), i=1,nxm), j=1,nym),&
((m22(j,i), i=1,nxm), j=1,nym),&
((n11(j,i), i=1,nxm), j=1,nym),&
((n12(j,i), i=1,nxm), j=1,nym),&
((n22(j,i), i=1,nxm), j=1,nym)
else
read(imean) m1, m2, n1, n2, m11, m12, m22, n11, n12, n22
endif
close(imean)
!.... read the datafile
open(imean, file='lstq.dat', form='unformatted', status='old',err=300)
read(imean) lstep, time, nxm, nym, nzm, ndofm, Re, Ma, Pr, gamma, cv
allocate( vm(nym,nxm,ndofm))
if (useIJ) then
read(imean) (((vm(j,i,k), k=1,ndofm), i=1,nxm), j=1,nym)
else
read(imean) vm
endif
close(imean)
!.... Compute first derivatives of field in the mapped space
allocate( g1vm(nym,nxm,ndofm), g2vm(nym,nxm,ndofm) )
call grad(ndofm, nxm, nym, vm, g1vm, g2vm, dxi, deta, optx, opty, &
xper, yper, lsym, rsym, bsym, tsym)
!.... Compute second derivatives of field
allocate( g11vm(nym,nxm,ndofm), g12vm(nym,nxm,ndofm), &
g22vm(nym,nxm,ndofm) )
call grad2(ndofm, nxm, nym, vm, g1vm, g11vm, g12vm, g22vm, dxi, deta, &
optx, opty, xper, yper, lsym, rsym, bsym, tsym)
!.... transform the gradients to physical space
allocate(g1vl(nym), g2vl(nym), g11vl(nym), g12vl(nym), g22vl(nym))
do idof = 1, ndofm
do i = 1, nxm
g1vl = g1vm(:,i,idof)*m1(:,i) + g2vm(:,i,idof)*n1(:,i)
g2vl = g1vm(:,i,idof)*m2(:,i) + g2vm(:,i,idof)*n2(:,i)
g11vl = g11vm(:,i,idof) * m1(:,i)*m1(:,i) + &
2.0 * g12vm(:,i,idof) * m1(:,i)*n1(:,i) + &
g22vm(:,i,idof) * n1(:,i)*n1(:,i) + &
g1vm(:,i,idof) * m11(:,i) + &
g2vm(:,i,idof) * n11(:,i)
g12vl = g11vm(:,i,idof) * m1(:,i)*m2(:,i) + &
g12vm(:,i,idof) * m1(:,i)*n2(:,i) + &
g12vm(:,i,idof) * m2(:,i)*n1(:,i) + &
g22vm(:,i,idof) * n1(:,i)*n2(:,i) + &
g1vm(:,i,idof) * m12(:,i) + &
g2vm(:,i,idof) * n12(:,i)
g22vl = g11vm(:,i,idof) * m2(:,i)*m2(:,i) + &
2.0 * g12vm(:,i,idof) * m2(:,i)*n2(:,i) + &
g22vm(:,i,idof) * n2(:,i)*n2(:,i) + &
g1vm(:,i,idof) * m22(:,i) + &
g2vm(:,i,idof) * n22(:,i)
g1vm(:,i,idof) = g1vl
g2vm(:,i,idof) = g2vl
g11vm(:,i,idof) = g11vl
g12vm(:,i,idof) = g12vl
g22vm(:,i,idof) = g22vl
end do
end do
deallocate( g1vl, g2vl, g11vl, g12vl, g22vl )
!.... don't need the metrics anymore
deallocate( m1, m2, n1, n2, m11, m12, m22, n11, n12, n22 )
!.... allocate space for the meanflow splines
allocate( vms(nym,ndofm), g1vms(nym,ndofm), g2vms(nym,ndofm), &
g11vms(nym,ndofm), g12vms(nym,ndofm), g22vms(nym,ndofm) )
return
100 write(*,"('ERROR: reading lstx.dat')")
call exit(1)
200 write(*,"('ERROR: reading lstm.dat')")
call exit(1)
300 write(*,"('ERROR: reading lstq.dat')")
call exit(1)
end
!============================================================================!
subroutine mean(imean, i)
!
! Spline the meanflow profile at station i
!
!============================================================================!
use meanflow
implicit none
integer :: imean, i
!.... local variables
integer :: idof
!============================================================================!
do idof = 1, ndofm
call SPLINE(nym, ym(1,i), vm(1,i,idof), vms(1,idof))
call SPLINE(nym, ym(1,i), g1vm(1,i,idof), g1vms(1,idof))
call SPLINE(nym, ym(1,i), g2vm(1,i,idof), g2vms(1,idof))
call SPLINE(nym, ym(1,i), g11vm(1,i,idof), g11vms(1,idof))
call SPLINE(nym, ym(1,i), g12vm(1,i,idof), g12vms(1,idof))
call SPLINE(nym, ym(1,i), g22vm(1,i,idof), g22vms(1,idof))
end do
return
end
!=============================================================================!
subroutine getmean( i, y, v, g1v, g2v, g11v, g12v, g22v )
!=============================================================================!
use meanflow
use global
implicit none
integer :: i
real :: y, v(ndofm), g1v(ndofm), g2v(ndofm)
real :: g11v(ndofm), g12v(ndofm) , g22v(ndofm)
integer :: idof
logical, save :: first_time = .true.
!=============================================================================!
if (useParallel) then
if (first_time) then
first_time = .false.
write(*,"('WARNING: using parallel flow approximation even through')")
write(*,"(' shoot-2d uses a full flow')")
endif
call getmeanp( i, y, v, g2v, g22v )
v(2) = 0.0
g1v = 0.0
g2v(3) = 0.0
g11v = 0.0
g12v = 0.0
g22v(3) = 0.0
else
if ( y .gt. ym(nym,i) ) then
do idof = 1, ndofm
v(idof) = vm(nym,i,idof)
g1v(idof) = g1vm(nym,i,idof)
g2v(idof) = g2vm(nym,i,idof)
g11v(idof) = g11vm(nym,i,idof)
g12v(idof) = g12vm(nym,i,idof)
g22v(idof) = g22vm(nym,i,idof)
end do
return
end if
call MSPEVAL(nym, ym(1,i), ndofm, vm(:,i,:), vms(:,:), y, v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g1vm(:,i,:), g1vms(:,:), y, g1v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g2vm(:,i,:), g2vms(:,:), y, g2v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g11vm(:,i,:), g11vms(:,:), y, g11v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g12vm(:,i,:), g12vms(:,:), y, g12v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g22vm(:,i,:), g22vms(:,:), y, g22v(:))
endif
return
end
!=============================================================================!
subroutine getmeanp( i, y, v, g2v, g22v )
!=============================================================================!
use meanflow
implicit none
integer :: i
real :: y, v(ndofm), g2v(ndofm), g22v(ndofm)
!.... local
integer :: idof
!=============================================================================!
if ( y .gt. ym(nym,i) ) then
do idof = 1, ndofm
v(idof) = vm(nym,i,idof)
g2v(idof) = g2vm(nym,i,idof)
g22v(idof) = g22vm(nym,i,idof)
end do
return
end if
call MSPEVAL(nym, ym(1,i), ndofm, vm(:,i,:), vms(:,:), y, v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g2vm(:,i,:), g2vms(:,:), y, g2v(:))
call MSPEVAL(nym, ym(1,i), ndofm, g22vm(:,i,:), g22vms(:,:), y, g22v(:))
return
end
!=============================================================================!
subroutine mspeval(n, x, ndof, y, ys, xl, yl)
!=============================================================================!
integer :: n, ndof
real :: x(n), y(n,ndof), ys(n,ndof), xl, yl(ndof)
real, parameter :: onesixth = 0.16666666666666666
!=============================================================================!
KLO=1
KHI=N
1 if (KHI-KLO.GT.1) then
K=(KHI+KLO) / 2
if(X(K).GT.XL)then
KHI=K
else
KLO=K
endif
goto 1
endif
DXM = XL - X(KLO)
DXP = X(KHI) - XL
DEL = X(KHI) - X(KLO)
DELINV = 1.0 / DEL
yl(:) = ( YS(KLO,:)*DXP*(DXP*DXP*DELINV - DEL) + &
YS(KHI,:)*DXM*(DXM*DXM*DELINV - DEL) ) * onesixth + &
( Y(KLO,:)*DXP + Y(KHI,:)*DXM ) * DELINV
return
end