/[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.28 - (show annotations) (download)
Sun Jan 26 21:06:11 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48b_post
Changes since 1.27: +27 -11 lines
r* coordinate added in #ifdef NONLIN_FRSURF block.
 (modification to pressure gradient not yet implemented)

1 C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.27 2001/09/26 18:09:16 cnh Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: TIMESTEP
8 C !INTERFACE:
9 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, K,
10 I phiHyd, phiSurfX, phiSurfY,
11 I myIter, myThid )
12 C !DESCRIPTION: \bv
13 C *==========================================================*
14 C | S/R TIMESTEP
15 C | o Step model fields forward in time
16 C *==========================================================*
17 C \ev
18
19 C !USES:
20 IMPLICIT NONE
21 C == Global variables ==
22 #include "SIZE.h"
23 #include "DYNVARS.h"
24 #include "EEPARAMS.h"
25 #include "PARAMS.h"
26 #include "GRID.h"
27 #include "SURFACE.h"
28
29 C !INPUT/OUTPUT PARAMETERS:
30 C == Routine Arguments ==
31 C phiHyd - Hydrostatic Potential (ocean: pressure/rho)
32 C (atmos: geopotentiel)
33 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
34 C phiSurfY or geopotentiel (atmos) in X and Y direction
35 INTEGER bi,bj,iMin,iMax,jMin,jMax
36 INTEGER K
37 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
38 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
39 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
40 INTEGER myIter, myThid
41
42 C !LOCAL VARIABLES:
43 C == Local variables ==
44 INTEGER i,j
45 _RL ab15,ab05
46 _RL phxFac,phyFac, psFac
47 _RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
48 _RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
49 CEOP
50
51 C Adams-Bashforth timestepping weights
52 IF (myIter .EQ. 0) THEN
53 ab15=1.0
54 ab05=0.0
55 ELSE
56 ab15=1.5+abeps
57 ab05=-0.5-abeps
58 ENDIF
59
60 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
61 C- Compute effective gU term (including Adams-Bashforth weights) :
62 DO j=jMin,jMax
63 DO i=iMin,iMax
64 gUtmp(i,j) = ab15*gU(i,j,k,bi,bj)
65 & + ab05*gUNm1(i,j,k,bi,bj)
66 #ifdef INCLUDE_CD_CODE
67 & + guCD(i,j,k,bi,bj)
68 #endif
69 ENDDO
70 ENDDO
71
72 #ifdef NONLIN_FRSURF
73 IF (.NOT. vectorInvariantMomentum
74 & .AND. nonlinFreeSurf.GT.1) THEN
75 IF (select_rStar.GT.0) THEN
76 DO j=jMin,jMax
77 DO i=iMin,iMax
78 gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj)
79 ENDDO
80 ENDDO
81 ELSE
82 DO j=jMin,jMax
83 DO i=iMin,iMax
84 IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
85 gUtmp(i,j) = gUtmp(i,j)
86 & *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
87 ENDIF
88 ENDDO
89 ENDDO
90 ENDIF
91 ENDIF
92 #endif
93
94 C Step forward zonal velocity (store in Gu)
95 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
96 DO j=jMin,jMax
97 DO i=iMin,iMax
98 gUNm1(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
99 & +deltaTmom*(
100 & gUtmp(i,j)
101 & - psFac*phiSurfX(i,j)
102 & )*_maskW(i,j,k,bi,bj)
103 ENDDO
104 ENDDO
105
106 IF (staggerTimeStep) THEN
107 C-- -grad Phi_Hyd has not been incorporated to gU and is added here:
108 phxFac = pfFacMom*deltaTmom
109 DO j=jMin,jMax
110 DO i=iMin,iMax
111 gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj)
112 & - _recip_dxC(i,j,bi,bj)
113 & *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac
114 & *_maskW(i,j,k,bi,bj)
115 ENDDO
116 ENDDO
117 ENDIF
118
119 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
120 C- Compute effective gV term (including Adams-Bashforth weights) :
121 DO j=jMin,jMax
122 DO i=iMin,iMax
123 gVtmp(i,j) = ab15*gV(i,j,k,bi,bj)
124 & + ab05*gVNm1(i,j,k,bi,bj)
125 #ifdef INCLUDE_CD_CODE
126 & + gvCD(i,j,k,bi,bj)
127 #endif
128 ENDDO
129 ENDDO
130
131 #ifdef NONLIN_FRSURF
132 IF (.NOT. vectorInvariantMomentum
133 & .AND. nonlinFreeSurf.GT.1) THEN
134 IF (select_rStar.GT.0) THEN
135 DO j=jMin,jMax
136 DO i=iMin,iMax
137 gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj)
138 ENDDO
139 ENDDO
140 ELSE
141 DO j=jMin,jMax
142 DO i=iMin,iMax
143 IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
144 gVtmp(i,j) = gVtmp(i,j)
145 & *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
146 ENDIF
147 ENDDO
148 ENDDO
149 ENDIF
150 ENDIF
151 #endif
152
153 C Step forward meridional velocity (store in Gv)
154 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
155 DO j=jMin,jMax
156 DO i=iMin,iMax
157 gVNm1(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
158 & +deltaTmom*(
159 & gVtmp(i,j)
160 & - psFac*phiSurfY(i,j)
161 & )*_maskS(i,j,k,bi,bj)
162 ENDDO
163 ENDDO
164
165 IF (staggerTimeStep) THEN
166 C-- -grad Phi_Hyd has not been incorporated to gV and is added here:
167 phyFac = pfFacMom*deltaTmom
168 DO j=jMin,jMax
169 DO i=iMin,iMax
170 gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj)
171 & - _recip_dyC(i,j,bi,bj)
172 & *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac
173 & *_maskS(i,j,k,bi,bj)
174 ENDDO
175 ENDDO
176 ENDIF
177
178 RETURN
179 END

  ViewVC Help
Powered by ViewVC 1.1.22