-
Notifications
You must be signed in to change notification settings - Fork 39
/
ExplicitTermsTemp.F90
100 lines (95 loc) · 2.59 KB
/
ExplicitTermsTemp.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: ExplicitTermsTemp.F90 !
! CONTAINS: subroutine ExplicitTermsTemp !
! !
! PURPOSE: Compute the non-linear terms associated to !
! the temperature. !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine ExplicitTermsTemp
use param
use local_arrays, only: vy,vx,temp,vz,hro
use decomp_2d, only: xstart,xend
implicit none
integer :: jc,kc,ic
integer :: km,kp,jm,jp,im,ip
real :: htx,hty,htz,udy,udz
real :: udzq,udyq
real :: dyyt,dzzt
udz=dz*0.25
udy=dy*0.25
udzq=dzq/pec
udyq=dyq/pec
!$OMP PARALLEL DO &
!$OMP DEFAULT(none) &
!$OMP SHARED(xstart,xend,vz,vy,vx,nxm) &
!$OMP SHARED(kmv,kpv,am3sk,ac3sk,ap3sk,udz) &
!$OMP SHARED(udy,udzq,udyq,udx3c,temp,hro) &
!$OMP PRIVATE(ic,jc,kc,im,ip,km,kp,jm,jp) &
!$OMP PRIVATE(htx,hty,htz,dyyt,dzzt)
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=2,nxm
km=kc-1
kp=kc+1
!
!
! rho vz term
!
!
! d rho q_z
! -----------
! d z
!
htz=((vz(km,jc,ip)+vz(kc,jc,ip))*(temp(kc,jc,ip)+temp(kc,jc,ic))- &
(vz(km,jc,ic)+vz(kc,jc,ic))*(temp(kc,jc,ic)+temp(kc,jc,im)) &
)*udz
!
!
! rho vy term
!
!
! d rho q_y
! -----------
! d y
!
hty=((vy(kc,jp,ic)+vy(km,jp,ic))*(temp(kc,jp,ic)+temp(kc,jc,ic))- &
(vy(kc,jc,ic)+vy(km,jc,ic))*(temp(kc,jc,ic)+temp(kc,jm,ic)) &
)*udy
!
! rho vx term
!
!
! d rho q_x
! -----------
! d x
!
htx=((vx(kp,jc,ic)+vx(kc,jc,ic))*(temp(kp,jc,ic)+temp(kc,jc,ic))- &
(vx(kc,jc,ic)+vx(km,jc,ic))*(temp(kc,jc,ic)+temp(km,jc,ic)) &
)*udx3c(kc)*0.25d0
!
!
! zz second derivatives of temp
!
dzzt=(temp(kc,jc,ip) &
-2.0*temp(kc,jc,ic) &
+temp(kc,jc,im))*udzq
!
! yy second derivatives of temp
!
dyyt=(temp(kc,jp,ic) &
-2.0*temp(kc,jc,ic) &
+temp(kc,jm,ic))*udyq
!
hro(kc,jc,ic) = -(htx+hty+htz)+dyyt+dzzt
enddo
enddo
enddo
!$OMP END PARALLEL DO
return
end