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

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

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


Revision 1.12 - (show annotations) (download)
Wed Jun 10 16:05:39 1998 UTC (25 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint7
Branch point for: checkpoint7-4degree-ref
Changes since 1.11: +16 -1 lines
Added code to bring "salt" up-to-date with "theta".
One caveat is that implicit diffusion of salt is done with the
diffusivity of theta. We'll sort this out later. In explicit
mode, diffKzS is used.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.11 1998/06/09 16:48:03 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 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 myThid
23 C == Local variables ==
24 INTEGER i,j
25 _RL ab15,ab05
26
27 C Adams-Bashforth timestepping weights
28 ab15=1.5+abeps
29 ab05=-0.5-abeps
30
31 C Step forward zonal velocity (store in Gu)
32 DO j=jMin,jMax
33 DO i=iMin,iMax
34 gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
35 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
36 #ifdef ALLOW_CD
37 & +guCD(i,j,k,bi,bj)
38 #endif
39 & )*_maskW(i,j,k,bi,bj)
40 ENDDO
41 ENDDO
42
43 C Step forward meridional velocity (store in Gv)
44 DO j=jMin,jMax
45 DO i=iMin,iMax
46 gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
47 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
48 #ifdef ALLOW_CD
49 & +gvCD(i,j,k,bi,bj)
50 #endif
51 & )*_maskS(i,j,k,bi,bj)
52 ENDDO
53 ENDDO
54
55 C Step forward temperature
56 IF (tempStepping) THEN
57 DO j=jMin,jMax
58 DO i=iMin,iMax
59 theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj)
60 & +deltaTtracer*(ab15*gT(i,j,k,bi,bj)+ab05*gTNm1(i,j,k,bi,bj))
61 gTNm1(i,j,k,bi,bj)=gT(i,j,k,bi,bj)
62 ENDDO
63 ENDDO
64 ENDIF
65
66 C Step forward salt
67 IF (saltStepping) THEN
68 DO j=jMin,jMax
69 DO i=iMin,iMax
70 salt(i,j,k,bi,bj)=salt(i,j,k,bi,bj)
71 & +deltaTtracer*(ab15*gS(i,j,k,bi,bj)+ab05*gSNm1(i,j,k,bi,bj))
72 gSNm1(i,j,k,bi,bj)=gS(i,j,k,bi,bj)
73 ENDDO
74 ENDDO
75 ENDIF
76
77 RETURN
78 END

  ViewVC Help
Powered by ViewVC 1.1.22