forked from PhysicsofFluids/AFiD
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathCorrectVelocity.F90
49 lines (44 loc) · 1.53 KB
/
CorrectVelocity.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: CorrectVelocity.F90 !
! CONTAINS: subroutine CorrectVelocity !
! !
! PURPOSE: Update velocities with the pressure !
! correction to enforce incompresibility !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine CorrectVelocity
use param
use local_arrays, only: vy,vx,dphhalo,vz
use decomp_2d, only: xstart,xend
implicit none
integer :: jc,jm,kc,km,ic,im
real :: usukm,udy,udz,locdph
udy = al*dt*dy
udz = al*dt*dz
!$OMP PARALLEL DO &
!$OMP DEFAULT(none) &
!$OMP SHARED(vz,vy,vx,dphhalo,udz,udy,udx3c) &
!$OMP SHARED(xstart,xend,nxm,kmv,dt,al) &
!$OMP PRIVATE(ic,jc,kc) &
!$OMP PRIVATE(im,jm,km,usukm,locdph)
do ic=xstart(3),xend(3)
im=ic-1
do jc=xstart(2),xend(2)
jm=jc-1
do kc=1,nxm
km=kmv(kc)
usukm = al*dt*udx3c(kc)
locdph=dphhalo(kc,jc,ic)
vx(kc,jc,ic)=vx(kc,jc,ic)- &
(locdph-dphhalo(km,jc,ic))*usukm
vy(kc,jc,ic)=vy(kc,jc,ic)- &
(locdph-dphhalo(kc,jm,ic))*udy
vz(kc,jc,ic)=vz(kc,jc,ic)- &
(locdph-dphhalo(kc,jc,im))*udz
enddo
enddo
enddo
!$OMP END PARALLEL DO
return
end