/[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.18 - (show annotations) (download)
Mon Sep 11 20:53:25 2000 UTC (23 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: branch-atmos-merge-start, checkpoint33, checkpoint32, checkpoint31, checkpoint34, branch-atmos-merge-phase1
Branch point for: branch-atmos-merge
Changes since 1.17: +5 -2 lines
Cosmetic.

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

  ViewVC Help
Powered by ViewVC 1.1.22