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

  ViewVC Help
Powered by ViewVC 1.1.22