-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathbou.f90
executable file
·460 lines (381 loc) · 10 KB
/
bou.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
module mod_bou
use mod_global
implicit none
private
public bounduvw, boundp, boundc,boundT,updthalos, updthalosBig,boundPFM,boundChem,boundloadd,&
boundBigQ,boundphib
contains
!
subroutine bounduvw(u,v,w,PFM_phiB,rholB,vislB)
use mod_cons
implicit none
integer :: i,j,k
real, dimension(0:i1,0:j1,0:k1), intent(inout) :: u,v,w
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(in):: PFM_phiB
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(in):: rholB,vislB
real:: DphiDx,DphiDy,DphiDz,Phi_grad_abs
real :: first_term ,second_term_a,second_term_b,second_term
real :: mu_bound,dvdz_wall
real::g_prime_phi
real :: phi_wall
real:: phi_jp, phi_jm,phi_kp,phi_km
select case(top)
case ('wall')
!this is the wall boundary
do j=0,j1
do i=0,i1
u(i,j,0) = -u(i,j,1)
w(i,j,0) = 0.0
u(i,j,k1) = -u(i,j,kmax)
w(i,j,kmax) = 0.0
w(i,j,k1) = w(i,j,kmax-1)
select case (PhaseField)
case ('Droplet')
v(i,j,k1) = -v(i,j,kmax)
end select
enddo
enddo
case ('outlet')
!here I use outlet boundary condition
do j=0,j1
do i=0,i1
u(i,j,0) = -u(i,j,1)
w(i,j,0) = 0.0
u(i,j,k1) = u(i,j,kmax)
w(i,j,k1) = w(i,j,kmax)
v(i,j,k1) = v(i,j,kmax)
enddo
enddo
end select
call updthalos(u,1)
call updthalos(v,1)
call updthalos(w,1)
call updthalos(u,2)
call updthalos(v,2)
call updthalos(w,2)
!bottom wall
do j=0,j1
do i=0,i1
v(i,j,0) = 0
enddo
enddo
! communicate data in x direction (periodic b.c.'s incorporated)
call updthalos(u,1)
call updthalos(v,1)
call updthalos(w,1)
! communicate data in y direction (periodic b.c.'s incorporated)
call updthalos(u,2)
call updthalos(v,2)
call updthalos(w,2)
return
end subroutine bounduvw
!
subroutine boundp(p)
use mod_cons
implicit none
integer :: i,j
real, dimension(0:,0:,0:) :: p
!
select case(top)
case ('wall')
do j=0,j1
do i=0,i1
p(i,j,0) = p(i,j,1) ! Newmann (consistent with no/free-slip)
p(i,j,k1) = p(i,j,kmax) ! Newmann (consistent with no/free-slip)
enddo
enddo
case ('outlet')
do j=0,j1
do i=0,i1
p(i,j,0) = p(i,j,1) ! Newmann (consistent with no/free-slip)
p(i,j,k1) = 0.0 ! top wall
enddo
enddo
end select
call updthalos(p,1)
!
call updthalos(p,2)
!
return
end subroutine boundp
subroutine boundBigQ(BigQ)
use mod_cons
implicit none
integer :: i,j
real, dimension(0:,0:,0:) :: BigQ
!
do j=0,j1
do i=0,i1
BigQ(i,j,0) = BigQ(i,j,1) ! Newmann (consistent with no/free-slip)
BigQ(i,j,k1) = BigQ(i,j,kmax) ! Newmann (consistent with no/free-slip)
enddo
enddo
call updthalos(BigQ,1)
!
call updthalos(BigQ,2)
!
return
end subroutine boundBigQ
subroutine boundc(c)
use mod_cons
implicit none
integer :: i,j
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2),intent(inout) :: c
!
do j=-2,j1+2
do i=-2,i1+2
! I am not sure about the wall condition
c(i,j,0) = c(i,j,1) ! Newmann (consistent with no/free-slip)
c(i,j,k1) = c(i,j,kmax) ! Newmann (consistent with no/free-slip)
enddo
enddo
do j=0,j1
do i=0,i1
c(i,j,-1) = c(i,j,0)
c(i,j,-2) = c(i,j,0)
c(i,j,k1+1) = c(i,j,k1)
c(i,j,k1+2) = c(i,j,k1)
enddo
enddo
!
call updthalosBig(c,1)
!
call updthalosBig(c,2)
!
return
end subroutine boundc
subroutine boundT(T)
use mod_cons
implicit none
integer :: i,j
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout) :: T
!
do j=-2,j1+2
do i=-2,i1+2
!T(i,j,0) = 270 ! fix value
T(i,j,0) = T(i,j,1) ! fix value
T(i,j,k1) = T(i,j,kmax) ! Newmann, adiabetic
enddo
enddo
!do j=0,j1
! do i=0,i1
do j=-2,j1+2
do i=-2,i1+2
T(i,j,-1) = T(i,j,0)
T(i,j,-2) = T(i,j,0)
T(i,j,k1+1) = T(i,j,k1)
T(i,j,k1+2) = T(i,j,k1)
enddo
enddo
call updthalosBig(T,1)
call updthalosBig(T,2)
!
return
!return
end subroutine boundT
subroutine boundPFM(PFM_phiB,dPFM_boundd,dPFM_boundd_top,u,v,w)
use mod_cons
use mod_chem
implicit none
integer :: i,j,k
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout) :: PFM_phiB
real, dimension(0:i1,0:j1,0:k1), intent(in) :: u,v,w
real, dimension(0:i1,0:j1), intent(out) ::dPFM_boundd,dPFM_boundd_top
real:: first_term,second_term,third_term,DphiDx,DphiDy,DphiDz,Phi_grad_abs
real:: g_prime_phi,v_wall
real::phi_jp,phi_jm,phi_wall
real:: DphiDzii
real::PFM_mu_f_1
select case (PhaseField)
case ('Droplet')
do j=-2,j1+2
do i=-2,i1+2
PFM_phiB(i,j,k1:k1+2) = PFM_phiB(i,j,kmax)
enddo
enddo
dPFM_boundd_top(:,:) = 0.
end select
!Bottom wall
do j=1,jmax
do i=1,imax
PFM_phiB(i,j,0)=PFM_phiB(i,j,1)
enddo
enddo
do j=0,j1
do i=0,i1
PFM_phiB(i,j,-1) = PFM_phiB(i,j,0)
PFM_phiB(i,j,-2) = PFM_phiB(i,j,0)
PFM_phiB(i,j,k1+1) = PFM_phiB(i,j,k1)
PFM_phiB(i,j,k1+2) = PFM_phiB(i,j,k1)
enddo
enddo
call updthalosBig(PFM_phiB,1)
call updthalosBig(PFM_phiB,2)
!
return
return
end subroutine boundPFM
subroutine boundphib(phi_b)
use mod_cons
implicit none
integer :: i,j,k
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout) :: phi_b
real:: g_prime_phi
real::phi_wall
!top surface
do j=-2,j1+2
do i=-2,i1+2
phi_b(i,j,k1:k1+2) = phi_b(i,j,kmax)
enddo
enddo
!Bottom wall
do j=1,jmax
do i=1,imax
phi_wall = 0.5*(phi_b(i,j,0)+phi_b(i,j,1))
g_prime_phi = 1.5*phi_wall*(1-phi_wall)
phi_b(i,j,0) = phi_b(i,j,1)+ (dz/PFM_l)*g_prime_phi*cos(PFM_thetta)*(2.*sqrt(2.)/3)
enddo
enddo
do j=0,j1
do i=0,i1
phi_b(i,j,-1) = phi_b(i,j,0)
phi_b(i,j,-2) = phi_b(i,j,0)
phi_b(i,j,k1+1) = phi_b(i,j,k1)
phi_b(i,j,k1+2) = phi_b(i,j,k1)
enddo
enddo
call updthalosBig(phi_b,1)
call updthalosBig(phi_b,2)
!
return
return
end subroutine boundphib
subroutine boundChem(chem_potB)
use mod_cons
implicit none
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout):: chem_potB
integer:: i,j,k
do j=-2,j1+2
do i=-2,i1+2
chem_PotB(i,j,0) = Chem_potB(i,j,1)
chem_PotB(i,j,-1) = Chem_potB(i,j,1)
chem_PotB(i,j,-2) = Chem_potB(i,j,1)
chem_PotB(i,j,k1+2) = Chem_potB(i,j,kmax)
chem_PotB(i,j,k1+1) = Chem_potB(i,j,kmax)
chem_PotB(i,j,k1) = Chem_potB(i,j,kmax)
enddo
enddo
!
call updthalosBig(chem_PotB,1)
call updthalosBig(chem_PotB,2)
end subroutine boundChem
subroutine boundloadd(u,v,w,p,PFM_phiB)
use mod_cons
implicit none
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout) :: PFM_phiB
real, dimension(0:i1,0:j1,0:k1), intent(inout) :: u,v,w,p
integer::i,j
do j=0,j1
do i=0,i1 !I changed here!
PFM_phiB(i,j,-1) = PFM_phiB(i,j,0)
PFM_phiB(i,j,-2) = PFM_phiB(i,j,0)
PFM_phiB(i,j,k1+1) = PFM_phiB(i,j,k1)
PFM_phiB(i,j,k1+2) = PFM_phiB(i,j,k1)
enddo
enddo
call updthalos(u,1)
call updthalos(v,1)
call updthalos(w,1)
call updthalos(u,2)
call updthalos(v,2)
call updthalos(w,2)
call updthalos(p,1)
call updthalos(p,2)
call updthalosBig(PFM_phiB,1)
call updthalosBig(PFM_phiB,2)
end subroutine boundloadd
subroutine updthalos(var,dir)
use mpi
use mod_cons
use mod_global
implicit none
real, dimension(0:,0:,0:), intent(inout) :: var
integer, intent(in) :: dir
integer :: requests(4), statuses(MPI_STATUS_SIZE,4)
!
!
select case(dir)
case(1) ! x direction
call MPI_SENDRECV(var(1,0,0),1,xhalo,left,0, &
var(i1,0,0),1,xhalo,right,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(imax,0,0),1,xhalo,right,0, &
var(0,0,0),1,xhalo,left,0, &
comm_cart,status,error)
case(2) ! y direction
call MPI_SENDRECV(var(0,1,0),1,yhalo,front,0, & !jmax+1=1
var(0,j1,0),1,yhalo,back,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(0,jmax,0),1,yhalo,back,0, & !jmax=0
var(0,0,0),1,yhalo,front,0, &
comm_cart,status,error)
end select
!
return
end subroutine updthalos
!
!
subroutine updthalosBig(var,dir)
use mpi
use mod_cons
use mod_global
implicit none
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2), intent(inout) :: var
integer, intent(in) :: dir
!
!
select case(dir)
case(1) ! x direction
call MPI_SENDRECV(var(1,-2,-2),1,xhalo3,left,0, &
var(i1,-2,-2),1,xhalo3,right,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(imax,-2,-2),1,xhalo3,right,0, &
var(0,-2,-2),1,xhalo3,left,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(1+1,-2,-2),1,xhalo3,left,0, &
var(i1+1,-2,-2),1,xhalo3,right,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(imax-1,-2,-2),1,xhalo3,right,0, &
var(0-1,-2,-2),1,xhalo3,left,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(1+2,-2,-2),1,xhalo3,left,0, &
var(i1+2,-2,-2),1,xhalo3,right,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(imax-2,-2,-2),1,xhalo3,right,0, &
var(0-2,-2,-2),1,xhalo3,left,0, &
comm_cart,status,error)
case(2) ! y direction
call MPI_SENDRECV(var(-2,1,-2),1,yhalo3,front,0, &
var(-2,j1,-2),1,yhalo3,back,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(-2,jmax,-2),1,yhalo3,back,0, &
var(-2,0,-2),1,yhalo3,front,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(-2,1+1,-2),1,yhalo3,front,0, &
var(-2,j1+1,-2),1,yhalo3,back,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(-2,jmax-1,-2),1,yhalo3,back,0, &
var(-2,0-1,-2),1,yhalo3,front,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(-2,1+2,-2),1,yhalo3,front,0, &
var(-2,j1+2,-2),1,yhalo3,back,0, &
comm_cart,status,error)
call MPI_SENDRECV(var(-2,jmax-2,-2),1,yhalo3,back,0, &
var(-2,0-2,-2),1,yhalo3,front,0, &
comm_cart,status,error)
end select
!
return
end subroutine updthalosBig
!
end module mod_bou