1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.16 1998/11/06 22:44:49 cnh Exp $ |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
C /==========================================================\ |
6 |
C | S/R TIMESTEP | |
7 |
C | o Step model fields forward in time | |
8 |
C \==========================================================/ |
9 |
SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, |
10 |
I K, |
11 |
I myIter, myThid ) |
12 |
implicit none |
13 |
! Common |
14 |
#include "SIZE.h" |
15 |
#include "DYNVARS.h" |
16 |
#include "EEPARAMS.h" |
17 |
#include "PARAMS.h" |
18 |
#include "GRID.h" |
19 |
C == Routine Arguments == |
20 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
21 |
INTEGER K |
22 |
INTEGER myIter, myThid |
23 |
C == Local variables == |
24 |
INTEGER i,j |
25 |
_RL ab15,ab05 |
26 |
|
27 |
C Adams-Bashforth timestepping weights |
28 |
Caja IF (myIter .EQ. 0) THEN |
29 |
Caja ab15=1.0 |
30 |
Caja ab05=0.0 |
31 |
Caja ELSE |
32 |
ab15=1.5+abeps |
33 |
ab05=-0.5-abeps |
34 |
Caja ENDIF |
35 |
|
36 |
C Step forward zonal velocity (store in Gu) |
37 |
DO j=jMin,jMax |
38 |
DO i=iMin,iMax |
39 |
gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj) |
40 |
& +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj) |
41 |
#ifdef INCLUDE_CD_CODE |
42 |
& +guCD(i,j,k,bi,bj) |
43 |
#endif |
44 |
& )*_maskW(i,j,k,bi,bj) |
45 |
ENDDO |
46 |
ENDDO |
47 |
|
48 |
C Step forward meridional velocity (store in Gv) |
49 |
DO j=jMin,jMax |
50 |
DO i=iMin,iMax |
51 |
gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj) |
52 |
& +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj) |
53 |
#ifdef INCLUDE_CD_CODE |
54 |
& +gvCD(i,j,k,bi,bj) |
55 |
#endif |
56 |
& )*_maskS(i,j,k,bi,bj) |
57 |
ENDDO |
58 |
ENDDO |
59 |
|
60 |
C Step forward temperature |
61 |
IF (tempStepping) THEN |
62 |
DO j=jMin,jMax |
63 |
DO i=iMin,iMax |
64 |
gTNm1(i,j,k,bi,bj)=theta(i,j,k,bi,bj) |
65 |
& +deltaTtracer*(ab15*gT(i,j,k,bi,bj)+ab05*gTNm1(i,j,k,bi,bj)) |
66 |
ENDDO |
67 |
ENDDO |
68 |
ENDIF |
69 |
|
70 |
C Step forward salt |
71 |
IF (saltStepping) THEN |
72 |
DO j=jMin,jMax |
73 |
DO i=iMin,iMax |
74 |
gSNm1(i,j,k,bi,bj)=salt(i,j,k,bi,bj) |
75 |
& +deltaTtracer*(ab15*gS(i,j,k,bi,bj)+ab05*gSNm1(i,j,k,bi,bj)) |
76 |
ENDDO |
77 |
ENDDO |
78 |
ENDIF |
79 |
|
80 |
RETURN |
81 |
END |