-
Notifications
You must be signed in to change notification settings - Fork 39
/
TimeMarcher.F90
78 lines (60 loc) · 2.36 KB
/
TimeMarcher.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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! FILE: TimeMarcher.F90 !
! CONTAINS: subroutine TimeMarcher !
! !
! PURPOSE: Main time integrating routine, which calls !
! other subroutines for calculating the Navier-Stokes !
! equations and advancing velocity and temperature in !
! time !
! !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TimeMarcher
use param
use local_arrays
use mpih
use decomp_2d
implicit none
integer :: ns
integer :: j,k,i
beta=dt/ren*0.5d0
do ns=1,nsst
!RO Coefficients for time marching integration (alpha, gamma, rho)
al=alm(ns)
ga=gam(ns)
ro=rom(ns)
call ExplicitTermsVX
call ExplicitTermsVY
call ExplicitTermsVZ
call ExplicitTermsTemp
call ImplicitAndUpdateVX
call ImplicitAndUpdateVY
call ImplicitAndUpdateVZ
call ImplicitAndUpdateTemp
call update_halo(vy,lvlhalo)
call update_halo(vz,lvlhalo)
call CalcLocalDivergence
call SolvePressureCorrection
!EP this copy can be avoided by changing transpose_x_to_y_real and
!transpose_y_to_x_real so these routines can handle arrays with
!halo. This copy is a defacto array temporary. Using inferred size
!arrays in the transpose calls results in 5 more of these, and more
!memory usage. Time spent on this copy is 0.1% for 65^3 grid.
do i=xstart(3),xend(3)
do j=xstart(2),xend(2)
do k=1,nxm
dphhalo(k,j,i) = dph(k,j,i)
enddo
enddo
enddo
call update_halo(dphhalo,lvlhalo)
call CorrectVelocity
call CorrectPressure
call update_halo(vx,lvlhalo)
call update_halo(vy,lvlhalo)
call update_halo(vz,lvlhalo)
call update_halo(pr,lvlhalo)
call update_halo(temp,lvlhalo)
enddo
return
end