-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathprojection_method_csayres.f90
541 lines (464 loc) · 18.8 KB
/
projection_method_csayres.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
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
! ==========================================
! Program: /Users/mandli/Documents/School/grad/ns_course
! File: projection_method
! Created: 2010-05-02
! Author: Kyle Mandli
! =========================================
! Copyright (C) 2010-05-02 Kyle Mandli <[email protected]>
!
! Distributed under the terms of the Berkeley Software Distribution (BSD)
! license
! http://www.opensource.org/licenses/
! ========================================
module grid_module
implicit none
! Grid parameters
integer :: N_x,N_y, N_leadIn
integer, parameter :: mbc = 1
double precision, parameter :: gamma = 0.66d0
double precision :: L_x,L_y,dx,dzeta
! Grid arrays
! i-1,j | i,j | i+1,j u <-- x_edge, matches p indices
! ---v---|---v---|---v---
! | |
! i-1,j | i,j | i+1,j
! p u p u p <--- x_center, y_center, zeta_center, F_center
! | |
! i-1,j-1| i,j-1 | i+1,j-1
! ---v---|---v---|---v--- <--- y_edge, zeta_edge, F_edge
! | |
double precision, allocatable :: x_edge(:),x_center(:)
double precision, allocatable :: zeta_edge(:),zeta_center(:)
double precision, allocatable :: y_edge(:),y_center(:)
double precision, allocatable :: F_edge(:),F_center(:)
contains
! Setup and allocate all grid quantities
subroutine setup_grid()
implicit none
! Local
integer :: i
! Grid spacing
dx = (L_x - 0.d0) / (N_x)
dzeta = (L_y - 0.d0) / (N_y)
! Grid array allocation
allocate(x_edge(0:N_x+1),x_center(0:N_x+1))
allocate(zeta_edge(0:N_y+1),zeta_center(0:N_y+1))
allocate(y_edge(0:N_y+1),y_center(0:N_y+1))
allocate(F_edge(0:N_y+1),F_center(0:N_y+1))
! Calculate grid arrays
forall (i=1-mbc:N_x+mbc)
x_edge(i) = 0.d0 + i * dx
x_center(i) = 0.d0 + (i-0.5d0) * dx
end forall
forall (i=1-mbc:N_y+mbc)
zeta_edge(i) = 0.d0 + i * dzeta
zeta_center(i) = 0.d0 + (i-0.5d0) * dzeta
F_edge(i) = F(zeta_edge(i))
F_center(i) = F(zeta_center(i))
end forall
! Really don't need these in this implementation
y_edge = L_y * (1.d0 - tanh(gamma * (L_y - zeta_edge)) / tanh(gamma*L_y))
y_center = L_y * (1.d0 - tanh(gamma * (L_y - zeta_center)) / tanh(gamma*L_y))
end subroutine setup_grid
! Grid stretching function
double precision pure function F(zeta)
implicit none
double precision, intent(in) :: zeta
F = tanh(gamma*L_y) / (gamma*L_y) * cosh(gamma*(L_y - zeta))**2
!F = L_y
end function F
! Output grid
subroutine output_grid(frame,t,u,v,p)
implicit none
! Arguments
integer, intent(in) :: frame
double precision, intent(in) :: t
double precision, dimension(0:N_x+1,0:N_y+1) :: u,v,p
! Locals
integer :: i,j
! Open output file and write file header:
if (t==0) then
open(unit=70,file='_output/UVP.dat',access='sequential',status='unknown')
write(70,*) ' VARIABLES= "x", "y", "u", "v", "p"'
write(70,100) t,N_X,N_Y
endif
100 FORMAT('ZONE T="t = ',e26.16,'"',' F=POINT, I=',I5,' J=', I5)
if (t>0) then
write(70,100) t,N_X,N_Y
endif
! Write out data
do j=1,N_y
do i=1,N_x
! Reset to zero if exponent is too large
if (abs(u(i,j)) < 1d-99) u(i,j) = 0.d0
if (abs(v(i,j)) < 1d-99) v(i,j) = 0.d0
if (abs(p(i,j)) < 1d-99) p(i,j) = 0.d0
write(70,"(5e26.16)") x_center(i),y_center(j),u(i,j),v(i,j),p(i,j)
enddo
enddo
end subroutine output_grid
end module grid_module
! =========================================
! Program: /Users/mandli/Documents/School/grad/ns_course
! File: main
! Created: 2010-05-02
! Author: Kyle Mandli
! ========================================
! Copyright (C) 2010-05-02 Kyle Mandli <[email protected]>
!
! Distributed under the terms of the Berkeley Software Distribution (BSD)
! license
! http://www.opensource.org/licenses/
! =========================================
program main
use grid_module
implicit none
! ====================================
! Solver parameters
integer, parameter :: MAX_ITERATIONS = 4000
double precision, parameter :: TOLERANCE = 1.d-4, CFL = 0.8d0
logical, parameter :: write_star = .false.
integer :: n_steps
! ===================================
! Physics
double precision :: U_inf = 1.d0
double precision :: rho,nu,Re !rho = 1.d0, nu=1.d-3
! ===================================
! Velocity and pressures
double precision, allocatable :: u(:,:),v(:,:),p(:,:),u_star(:,:),v_star(:,:)
double precision, allocatable :: u_old(:,:), v_old(:,:)
! ===================================
! Locals
character*20 :: arg
integer :: i,j,n,m,frame,i_R,j_R
double precision :: R,t,dt,a
double precision, allocatable :: Flux_ux(:,:),Flux_uy(:,:),Flux_vy(:,:)
double precision :: uu_x,uv_y,uv_x,vv_y,u_xx,u_yy,v_xx,v_yy
double precision, allocatable :: Q(:,:),b(:),cp(:),cm(:)
! ===================================
! Get command line arguments
! if (iargc() /= 6) then
! print *,"Wrong number of command line arguments, expected 6."
! print *," blasius Nx Ny Lx Ly Nsteps Re"
! print *," Nx - Number of grid points in x"
! print *," Ny - Number of grid points in y"
! print *," Lx - Length of domain in x"
! print *," Ly - Length of domain in y"
! print *," Nsteps - Number of steps in between output"
! print *," Re - Reynolds number of flow"
! stop
! else
! call getarg(1,arg)
! read (arg,'(I10)') N_x
! call getarg(2,arg)
! read (arg,'(I10)') N_y
! call getarg(3,arg)
! read (arg,'(e16.8)') L_x
! call getarg(4,arg)
! read (arg,'(e16.8)') L_y
! call getarg(5,arg)
! read (arg,'(I10)') n_steps
! call getarg(6,arg)
! read (arg,'(e16.8)') Re
! endif
!!!!!!!!!!!!!!!!!!
!! Parameters: !!
!!!!!!!!!!!!!!!!!!
N_x=20 !Number of grid points in x-direction
N_y=20 !Number of grid points in y-direction
N_leadIn = 5 !Number of grid points before wall boundary
L_x=10.0 !Length of box in x-direction
L_y=4.0 !Length of box in y-direction
n_steps=500 !Interval that u,v and p are printed to UVP.dat
Re=100.0 !Reynolds number
print *,"Running blasius with following parameters: "
print "(a,i3,a,i3,a)"," (N_x,N_y) = (",N_x,",",N_y,")"
print "(a,e16.8,a,e16.8,a)"," (L_x,L_y) = (",L_x,",",L_y,")"
print "(a,i4,a)"," Output every ",n_steps," steps"
print "(a,e16.8)"," Reynold's number = ",Re
! ===================================
! Setup grid and arrays
call setup_grid()
allocate(Flux_ux(1:N_x+1,0:N_y+1))
allocate(Flux_uy(0:N_x,0:N_y))
allocate(Flux_vy(0:N_x+1,0:N_y))
allocate(Q(1:N_x,1:N_y))
allocate(b(1:N_y),cp(1:N_y),cm(1:N_y))
! Calculate matrix coefficients for Poisson solve
a = 1.d0 / dx**2
forall (j=1:N_y)
b(j) = - (2.d0 / dx**2 + F_center(j) / dzeta**2 * (F_edge(j) + F_edge(j-1)))
cp(j) = (F_center(j) * F_edge(j)) / dzeta**2
cm(j) = (F_center(j) * F_edge(j-1)) / dzeta**2
end forall
! Velocity and pressure arrays
allocate(u(1-mbc:N_x+mbc,1-mbc:N_y+mbc),u_star(1-mbc:N_x+mbc,1-mbc:N_y+mbc))
allocate(v(1-mbc:N_x+mbc,1-mbc:N_y+mbc),v_star(1-mbc:N_x+mbc,1-mbc:N_y+mbc))
allocate(p(1-mbc:N_x+mbc,1-mbc:N_y+mbc))
allocate(u_old(1:N_x,1:N_y),v_old(1:N_x,1:N_y))
! Inital conditions
u = U_inf
v = 0.d0
p = 0.d0
dt = CFL * dx / (Re * U_inf)
t = 0.d0
frame = 0
nu = 1.d0/Re !1.d-3 !1.d0/Re
rho = 1.d0
! Output inital condition
call output_grid(frame,t,u,v,p)
print "(a,i3,a,i4,a,e16.8)","Writing frame ",frame," during step n=",0," t=",t
! Open up file to store residual information in
open(unit=13, file='residual.dat', status="unknown", action="write")
! ===================================
! Main algorithm loop
do n=1,MAX_ITERATIONS
! Store old step for convergence test
u_old = u(1:N_x,1:N_y)
v_old = v(1:N_x,1:N_y)
! ===================================
! Apply BC V = x, U = o, P = +
call bc(u,v,U_inf)
! ===================================
! Step 1: Update velocity to intermediate step
! Calculate fluxes at each boundary
!
! U-Fluxes
! | u |
! FUY_i,j i,j| i,j+1 |i+1,j
! -------|---o---|------- -------v-------v-------
! | | i-1,j | i,j | i+1,j
! FUX_i,j FUX_i+1,j u | u | u
! | | | |
! -------|---o---|------- -------v-------v-------
! FUY_i,j-1 i,j-1| i,j-1 | i+1,j-1
! | u |
!
! V-Fluxes
! | v |
! FVY_i,j i-1,j+1| i,j+1 |i,j+1
! -------|---o---|------- -------u---o---u-------
! | | i-1,j | i,j | i+1,j
! FUY_i-1,j FUY_i,j v o v o v
! | | | |
! -------|---o---|------- -------u---o---u-------
! FVY_i,j-1 i-1,j| i,j-1 |i,j
! | v |
forall (i=1:N_x+1,j=0:N_y+1) Flux_ux(i,j) = (u(i-1,j)+u(i,j))**2 / 4.d0
forall (i=0:N_x,j=0:N_y) Flux_uy(i,j) = (u(i,j)+u(i,j+1)) * (v(i+1,j)+v(i,j)) / 4.d0
forall (i=0:N_x+1,j=0:N_y) Flux_vy(i,j) = (v(i,j+1) + v(i,j))**2 / 4.d0
do j=1,N_y
do i=1,N_x
!********************************************
!*** ADD RHS TO THE FOLLOWING FORMULAS: ***
!********************************************
! Advective terms, see Lecture 4 notes, page 12
! Grid stretching function defined above F(zeta), F_center
uu_x = (Flux_ux(i+1,j) - Flux_ux(i,j)) / dx
uv_y = F_center(j) * (Flux_uy(i,j) - Flux_uy(i,j-1)) / dzeta !
uv_x = (Flux_uy(i,j) - Flux_uy(i-1,j)) / dx
vv_y = F_edge(j) * (Flux_vy(i,j) - Flux_vy(i,j-1)) / dzeta
! Diffusive terms
u_xx = (u(i+1,j) - 2*u(i,j) + u(i-1,j)) / (dx**2)
u_yy = (F_center(j) / dzeta**2) * (F_edge(j) * (u(i,j+1) - u(i,j)) - F_edge(j-1) * (u(i,j) - u(i,j-1))) !edge leads center by 0.5, verify by looking at poisson solver
v_xx = (v(i+1,j) - 2.d0*v(i,j) + v(i-1,j)) / (dx**2)
v_yy = (F_edge(j) / dzeta**2) * (F_center(j+1) * (v(i,j+1) - v(i,j)) - F_center(j) * (v(i,j) - v(i,j-1)))
! Update to u* and v* values
u_star(i,j) = u(i,j) + dt * (-1*(uu_x + uv_y) + nu*(u_xx + u_yy))
v_star(i,j) = v(i,j) + dt * (-1*(uv_x+vv_y) + nu*(v_xx + v_yy))
end do
end do
! Debug, save out u_star,v_star,p
if (write_star) then
frame = frame + 1
print "(a,i3,a,i4,a,e16.8)","Writing frame ",frame," during step n=",n," t=",t
call output_grid(frame,t,u_star,v_star,p)
endif
! ===================================
! Step 2: Solve projection poisson problem
call bc(u_star,v_star,U_inf)
forall(i=1:N_x,j=1:N_y)
! page 16 lecture 4 notes
Q(i,j) = 1.d0/dt * ((u_star(i,j)-u_star(i-1,j)) / dx + (v_star(i,j)-v_star(i,j-1)) / dzeta * F_center(j))
end forall
! Solve poisson problem
call solve_poisson(p,Q,a,b,cm,cp)
! ===================================
! Step 3: Update velocities to end time
!**********************************************
!*** ADD RHS TO VELOCITY UPDATE FORLUMAS: ***
!**********************************************
! see work on page 18 lecture 4
forall (i=1:N_x,j=1:N_y)
! appears rho is 1.d0 from Q matrix
u(i,j) = u_star(i,j) - dt * (p(i+1,j) - p(i,j)) / dx
v(i,j) = v_star(i,j) - dt * (p(i,j+1) - p(i,j)) * F_edge(j) / dzeta
end forall
! ===================================
! Step 4: Check convergence
R = 0.d0
do j=1,N_y
do i=1,N_x
if (R < abs(u(i,j)-u_old(i,j)) .or. R < abs(v(i,j)-v_old(i,j))) then
R = max(R,abs(u(i,j)-u_old(i,j)),abs(v(i,j)-v_old(i,j)))
i_R = i
j_R = j
endif
enddo
enddo
! Finish up loop
!print "(a,i4,a,i3,a,i3,a,e16.8)","Loop ",n,": (",i_R,",",j_R,") - R = ",R
write (13,"(i4,i4,i4,e16.8)") n,i_R,j_R,R
! Write out u,v,p every n_steps
if (mod(n,n_steps) == 0) then
frame = frame + 1
call output_grid(frame,t,u,v,p)
print "(a,i3,a,i4,a,e16.8)","Writing frame ",frame," during step n=",n," t=",t
endif
! Check tolerance
if (R < TOLERANCE) then
print *, "Convergence reached, R = ",R
call output_grid(frame,t,u,v,p)
print "(a,i3,a,i4,a,e16.8)","Writing frame ",frame," during step n=",n," t=",t
exit
endif
! We did not reach our tolerance, iterate again
t = t + dt
enddo
if (R > TOLERANCE) then
print "(a,e16.8)","Convergence was never reached, R = ", R
endif
call output_grid(frame,t,u,v,p) ! ouput last grid?
! print "(a,i3,a,i4,a,e16.8)","Tolerance reaced!, Writing frame ",frame," during step n=",n," t=",t
close(13)
end program main
! ===================================
! Boundary Conditions
! ===================================
subroutine bc(u,v,U_inf)
use grid_module
implicit none
! Input/Output arguments
double precision, intent(in) :: U_inf
double precision, intent(inout) :: u(1-mbc:N_x+mbc,1-mbc:N_y+mbc)
double precision, intent(inout) :: v(1-mbc:N_x+mbc,1-mbc:N_y+mbc)
! Locals
integer :: i,j
! Lower Boundary Upper Boundary
! i,1 i,1 |
! x + x i,1 | o |
! | | |i,N_y+1|
! -------o-------- x + x
! | i,0 | | |
! x + x -|---o---|-
! | i,0 i,0 | i,N_y |
! | | x + x
! i,N_y
! Wall Free shear
! u(i,0) = -u(i,1) u(i,N_y+1) = u(i,N_y)
! v(i,0) = 0.d0 v(i,N_y+1) = v(i,N_y)
! p(i,0) = p(i,1) p(i,N_y+1) = 0.d0
forall(i=0:N_x+1)
u(i,0) = -u(i,1)
v(i,0) = 0.d0
u(i,N_y+1) = u(i,N_y)
v(i,N_y+1) = v(i,N_y)
end forall
! overwrite the lead-in
! boundary to have a slip condition
forall(i=0:N_leadIn)
u(i,0) = u(i,1)
v(i,0) = v(i,1)
end forall
! Left Boundaries
! x 0,j | x 1,j
! 0,j | u(0,j) = U_inf
! + o 0,j + 1,j o 1,j v(0,j) = 0.d0
! | p(0,j) = p(1,j) (P_x = 0)
! Right Boundaries (outflow)
! x N_x,j | x N_x+1,j u(N_x+1,j) = 2 u(N_x,j) - u(N_x-1,j)
! | v(N_x+1,j) = 2 v(N_x,j) - v(N_x-1,j)
! + N_x,j o N_x,j + N_x+1,j o N_x+1,j p(N_x+1,j) = p(N_x,j)
! |
forall(j=1:N_y+1)
u(0,j) = U_inf
v(0,j) = 0.d0
! u(N_x+1,j) = 2.d0*u(N_x,j) - u(N_x-1,j)
u(N_x+1,j) = u(N_x,j)
v(N_x+1,j) = v(N_x,j)
! v(N_x+1,j) = 2.d0*v(N_x,j) - v(N_x-1,j)
end forall
end subroutine bc
! ===================================
! Solve the poisson problem
! laplace(P)_ij = Q_ij
! ===================================
subroutine solve_poisson(P,Q,a,b,cm,cp)
use grid_module
implicit none
! Input
double precision, intent(in) :: Q(1:N_x,1:N_y),a
double precision, intent(in) :: b(1:N_y),cm(1:N_y),cp(1:N_y)
double precision, intent(inout) :: P(1-mbc:N_x+mbc,1-mbc:N_y+mbc)
! Solver parameters
logical, parameter :: verbose = .false.
integer, parameter :: MAX_ITERATIONS = 10000
double precision, parameter :: TOLERANCE = 10.d-4
double precision, parameter :: w = 1.6d0 ! 1 (GS) < w < 2
! Local variables
integer :: i,j,n
double precision :: R,P_old
do n=1,MAX_ITERATIONS
R = 0.d0
! Boundary conditions
forall (j=0:N_y+1)
P(0,j) = P(1,j) ! Left
P(N_x+1,j) = P(N_x,j) ! Right
end forall
forall (i=0:N_x+1)
P(i,0) = P(i,1) ! Bottom wall
P(i,N_y+1) = 0.d0 ! Free stream
end forall
! add in the lead-in bit
forall (i=0:N_leadIn)
P(i,0) = 0.d0 !Bottom wall acts like free stream
end forall
do j=1,N_y
do i=1,N_x
! Update p
P_old = P(i,j)
P(i,j) = (Q(i,j) - (a * (P(i+1,j) + P(i-1,j)) + cp(j)*P(i,j+1) + cm(j)*P(i,j-1))) / b(j)
! relaxation
P(i,j) = w * P(i,j) + (1-w) * P_old
! Calculate Residual
R = max(R,abs(P(i,j)-P_old))
enddo
enddo
!print *, "R: ", R
! Print out convergence and exit
if (verbose) then
print "(a,i5,a,e16.8,a)", "(",n,",",R,")"
endif
if (R < TOLERANCE) then
exit
endif
enddo
! Check to see if we converged and quit if we did not
if (n > MAX_ITERATIONS) then
print *,"Poisson solver did not converge, exiting!"
print "(a,e16.8)","R = ",R
print *,"iteration: ", n
stop
else
print "(a,i6,a,e16.8)", "Solver converged in ",n," steps: R = ",R
! Boundary conditions
forall (j=0:N_y+1)
P(0,j) = P(1,j) ! Left
P(N_x+1,j) = P(N_x,j) ! Right
end forall
forall (i=0:N_x+1)
P(i,0) = P(i,1) ! Bottom wall
P(i,N_y+1) = 0.d0 ! Free stream
end forall
endif
end subroutine solve_poisson