/[MITgcm]/MITgcm/model/src/timestep.F
ViewVC logotype

Annotation of /MITgcm/model/src/timestep.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.17 - (hide annotations) (download)
Mon Aug 30 18:25:33 1999 UTC (24 years, 9 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint28, checkpoint29, checkpoint25, checkpoint27, checkpoint26, checkpoint30
Changes since 1.16: +10 -5 lines
Added myIter as an argument to timestep() to allow proper
timestepping at myIter=0

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

  ViewVC Help
Powered by ViewVC 1.1.22