-
Notifications
You must be signed in to change notification settings - Fork 39
/
CorrectPressure.F90
55 lines (52 loc) · 1.75 KB
/
CorrectPressure.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: CorrectPressure.F90 !
! CONTAINS: subroutine CorrectPressure !
! !
! PURPOSE: Apply the pressure correction to the !
! pressure !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine CorrectPressure
use param
use local_arrays, only: pr,dphhalo
use decomp_2d, only: xstart,xend
implicit none
integer :: kp,km,jm,jp,jc,kc,ic,ip,im
real :: be,amm,acc,app
be=al*beta
!$OMP PARALLEL DO &
!$OMP DEFAULT(none) &
!$OMP SHARED(pr,dphhalo,be,amphk,acphk,apphk) &
!$OMP SHARED(xstart,xend,nxm,kmv,kpv,dzq,dyq) &
!$OMP PRIVATE(ic,jc,kc) &
!$OMP PRIVATE(im,jm,km,ip,jp,kp) &
!$OMP PRIVATE(amm,acc,app)
do ic=xstart(3),xend(3)
im=ic-1
ip=ic+1
do jc=xstart(2),xend(2)
jm=jc-1
jp=jc+1
do kc=1,nxm
kp=kpv(kc)
km=kmv(kc)
amm=amphk(kc)
acc=acphk(kc)
app=apphk(kc)
pr(kc,jc,ic)=pr(kc,jc,ic)+dphhalo(kc,jc,ic)-be*( &
(dphhalo(kc,jc,ip) &
-2.0*dphhalo(kc,jc,ic) &
+dphhalo(kc,jc,im))*dzq+ &
(dphhalo(kc,jp,ic) &
-2.0*dphhalo(kc,jc,ic) &
+dphhalo(kc,jm,ic))*dyq+ &
(dphhalo(kp,jc,ic)*app &
+dphhalo(kc,jc,ic)*acc &
+dphhalo(km,jc,ic)*amm))
enddo
enddo
enddo
!$OMP END PARALLEL DO
return
end