/[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.37 - (show annotations) (download)
Fri Apr 15 14:17:31 2005 UTC (19 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_post, checkpoint57i_post, checkpoint57n_post, checkpoint57l_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57h_done, checkpoint57j_post, checkpoint57o_post, checkpoint57k_post
Changes since 1.36: +52 -33 lines
call adams_bashforth2 (or adams_bashforth3) instead of local calculation

1 C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.36 2004/11/10 03:02:00 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: TIMESTEP
9 C !INTERFACE:
10 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, k,
11 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
12 I guDissip, gvDissip,
13 I myTime, myIter, myThid )
14 C !DESCRIPTION: \bv
15 C *==========================================================*
16 C | S/R TIMESTEP
17 C | o Step model fields forward in time
18 C *==========================================================*
19 C \ev
20
21 C !USES:
22 IMPLICIT NONE
23 C == Global variables ==
24 #include "SIZE.h"
25 #include "DYNVARS.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #include "GRID.h"
29 #include "SURFACE.h"
30
31 C !INPUT/OUTPUT PARAMETERS:
32 C == Routine Arguments ==
33 C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential
34 C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean)
35 C phiSurfY :: or geopotential (atmos) in X and Y direction
36 C guDissip :: dissipation tendency (all explicit terms), u component
37 C gvDissip :: dissipation tendency (all explicit terms), v component
38
39 INTEGER bi,bj,iMin,iMax,jMin,jMax
40 INTEGER k
41 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
42 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
43 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
44 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
45 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
46 _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
47 _RL myTime
48 INTEGER myIter, myThid
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 LOGICAL momForcing_In_AB
53 LOGICAL momDissip_In_AB
54 INTEGER i,j
55 _RL ab15,ab05
56 _RL phxFac,phyFac, psFac
57 _RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
58 _RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
59 #ifdef ALLOW_CD_CODE
60 _RL guCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61 _RL gvCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62 #endif
63 CEOP
64
65 C Adams-Bashforth timestepping weights
66 IF (myIter .EQ. 0) THEN
67 ab15=1.0
68 ab05=0.0
69 ELSE
70 ab15=1.5+abeps
71 ab05=-0.5-abeps
72 ENDIF
73
74 C-- explicit part of the surface potential gradient is added in this S/R
75 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
76
77 C-- including or excluding momentum forcing from Adams-Bashforth:
78 momForcing_In_AB = forcing_In_AB
79 momForcing_In_AB = .TRUE.
80 momDissip_In_AB = .TRUE.
81
82 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
83
84 C- Initialize local arrays (not really necessary but safer)
85 DO j=1-Oly,sNy+Oly
86 DO i=1-Olx,sNx+Olx
87 gUtmp(i,j) = 0. _d 0
88 gVtmp(i,j) = 0. _d 0
89 #ifdef ALLOW_CD_CODE
90 guCor(i,j) = 0. _d 0
91 gvCor(i,j) = 0. _d 0
92 #endif
93 ENDDO
94 ENDDO
95
96 C-- Stagger time step: grad Phi_Hyp will be added later
97 IF (staggerTimeStep) THEN
98 phxFac = pfFacMom
99 phyFac = pfFacMom
100 ELSE
101 C-- Synchronous time step: add grad Phi_Hyp to gU,gV before doing Adams-Bashforth
102 C note: already done in S/R mom_vecinv and mom_fluxform but would be better
103 C to add it to gU,gV here.
104 c DO j=jMin,jMax
105 c DO i=iMin,iMax
106 c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) - pfFacMom*dPhiHydX(i,j)
107 c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) - pfFacMom*dPhiHydY(i,j)
108 c ENDDO
109 c ENDDO
110 phxFac = 0.
111 phyFac = 0.
112 ENDIF
113
114 #ifdef ALLOW_MOM_VECINV
115 C-- Dissipation term inside the Adams-Bashforth:
116 C note: already in gU,gV if using fluxform
117 IF ( momViscosity .AND. momDissip_In_AB
118 & .AND. vectorInvariantMomentum ) THEN
119 DO j=jMin,jMax
120 DO i=iMin,iMax
121 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) + guDissip(i,j)
122 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) + gvDissip(i,j)
123 ENDDO
124 ENDDO
125 ENDIF
126 #endif
127
128 C-- Forcing term inside the Adams-Bashforth:
129 IF (momForcing .AND. momForcing_In_AB) THEN
130 CALL EXTERNAL_FORCING_U(
131 I iMin,iMax,jMin,jMax,bi,bj,k,
132 I myTime,myThid)
133 CALL EXTERNAL_FORCING_V(
134 I iMin,iMax,jMin,jMax,bi,bj,k,
135 I myTime,myThid)
136 ENDIF
137
138 IF (useCDscheme) THEN
139 C- for CD-scheme, store gU,Vtmp = gU,V^n + forcing
140 DO j=jMin,jMax
141 DO i=iMin,iMax
142 gUtmp(i,j) = gU(i,j,k,bi,bj)
143 gVtmp(i,j) = gV(i,j,k,bi,bj)
144 ENDDO
145 ENDDO
146 ENDIF
147
148 C- Compute effective gU,gV_[n+1/2] terms (including Adams-Bashforth weights)
149 C and save gU,gV_[n] into guNm1,gvNm1 for the next time step.
150 #ifdef ALLOW_ADAMSBASHFORTH_3
151 CALL ADAMS_BASHFORTH3(
152 I bi, bj, k,
153 U gU, guNm,
154 I myIter, myThid )
155 CALL ADAMS_BASHFORTH3(
156 I bi, bj, k,
157 U gV, gvNm,
158 I myIter, myThid )
159 #else /* ALLOW_ADAMSBASHFORTH_3 */
160 CALL ADAMS_BASHFORTH2(
161 I bi, bj, k,
162 U gU, guNm1,
163 I myIter, myThid )
164 CALL ADAMS_BASHFORTH2(
165 I bi, bj, k,
166 U gV, gvNm1,
167 I myIter, myThid )
168 #endif /* ALLOW_ADAMSBASHFORTH_3 */
169
170 C-- Forcing term outside the Adams-Bashforth:
171 C (not recommanded with CD-scheme ON)
172 IF (momForcing .AND. .NOT.momForcing_In_AB) THEN
173 IF (useCDscheme) THEN
174 DO j=jMin,jMax
175 DO i=iMin,iMax
176 gUtmp(i,j) = gUtmp(i,j) - gU(i,j,k,bi,bj)
177 gVtmp(i,j) = gVtmp(i,j) - gV(i,j,k,bi,bj)
178 ENDDO
179 ENDDO
180 ENDIF
181 CALL EXTERNAL_FORCING_U(
182 I iMin,iMax,jMin,jMax,bi,bj,k,
183 I myTime,myThid)
184 CALL EXTERNAL_FORCING_V(
185 I iMin,iMax,jMin,jMax,bi,bj,k,
186 I myTime,myThid)
187
188 C- for CD-scheme, compute gU,Vtmp = gU,V^n + forcing
189 IF (useCDscheme) THEN
190 DO j=jMin,jMax
191 DO i=iMin,iMax
192 gUtmp(i,j) = gUtmp(i,j) + gU(i,j,k,bi,bj)
193 gVtmp(i,j) = gVtmp(i,j) + gV(i,j,k,bi,bj)
194 ENDDO
195 ENDDO
196 ENDIF
197 ENDIF
198
199 #ifdef ALLOW_CD_CODE
200 IF (useCDscheme) THEN
201 C- Step forward D-grid velocity using C-grid gU,Vtmp = gU,V^n + forcing
202 C and return coriolis terms on C-grid (guCor,gvCor)
203 CALL CD_CODE_SCHEME(
204 I bi,bj,k, dPhiHydX,dPhiHydY, gUtmp,gVtmp,
205 O guCor,gvCor,
206 I myTime, myIter, myThid)
207 DO j=jMin,jMax
208 DO i=iMin,iMax
209 gUtmp(i,j) = gU(i,j,k,bi,bj)
210 & + guCor(i,j)
211 gVtmp(i,j) = gV(i,j,k,bi,bj)
212 & + gvCor(i,j)
213 ENDDO
214 ENDDO
215 ELSE
216 #endif /* ALLOW_CD_CODE */
217 DO j=jMin,jMax
218 DO i=iMin,iMax
219 gUtmp(i,j) = gU(i,j,k,bi,bj)
220 gVtmp(i,j) = gV(i,j,k,bi,bj)
221 ENDDO
222 ENDDO
223 #ifdef ALLOW_CD_CODE
224 ENDIF
225 #endif
226
227 #ifdef NONLIN_FRSURF
228 IF (.NOT. vectorInvariantMomentum
229 & .AND. nonlinFreeSurf.GT.1) THEN
230 IF (select_rStar.GT.0) THEN
231 DO j=jMin,jMax
232 DO i=iMin,iMax
233 gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj)
234 gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj)
235 ENDDO
236 ENDDO
237 ELSE
238 DO j=jMin,jMax
239 DO i=iMin,iMax
240 IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
241 gUtmp(i,j) = gUtmp(i,j)
242 & *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
243 ENDIF
244 IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
245 gVtmp(i,j) = gVtmp(i,j)
246 & *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
247 ENDIF
248 ENDDO
249 ENDDO
250 ENDIF
251 ENDIF
252 #endif /* NONLIN_FRSURF */
253
254 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
255
256 #ifdef ALLOW_MOM_VECINV
257 C-- Dissipation term outside the Adams-Bashforth:
258 C note: only implemented with vecinv formulation
259 IF ( momViscosity .AND. .NOT.momDissip_In_AB
260 & .AND. vectorInvariantMomentum ) THEN
261 DO j=jMin,jMax
262 DO i=iMin,iMax
263 gUtmp(i,j) = gUtmp(i,j) + guDissip(i,j)
264 gVtmp(i,j) = gVtmp(i,j) + gvDissip(i,j)
265 ENDDO
266 ENDDO
267 ENDIF
268 #endif
269
270 C Step forward zonal velocity (store in Gu)
271 DO j=jMin,jMax
272 DO i=iMin,iMax
273 gU(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
274 & +deltaTmom*(
275 & gUtmp(i,j)
276 & - psFac*phiSurfX(i,j)
277 & - phxFac*dPhiHydX(i,j)
278 & )*_maskW(i,j,k,bi,bj)
279 ENDDO
280 ENDDO
281
282 C Step forward meridional velocity (store in Gv)
283 DO j=jMin,jMax
284 DO i=iMin,iMax
285 gV(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
286 & +deltaTmom*(
287 & gVtmp(i,j)
288 & - psFac*phiSurfY(i,j)
289 & - phyFac*dPhiHydY(i,j)
290 & )*_maskS(i,j,k,bi,bj)
291 ENDDO
292 ENDDO
293
294 RETURN
295 END

  ViewVC Help
Powered by ViewVC 1.1.22