1 |
C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.46 2006/06/07 01:55:13 heimbach 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 "RESTART.h" |
30 |
#include "SURFACE.h" |
31 |
|
32 |
C !INPUT/OUTPUT PARAMETERS: |
33 |
C == Routine Arguments == |
34 |
C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential |
35 |
C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean) |
36 |
C phiSurfY :: or geopotential (atmos) in X and Y direction |
37 |
C guDissip :: dissipation tendency (all explicit terms), u component |
38 |
C gvDissip :: dissipation tendency (all explicit terms), v component |
39 |
|
40 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
41 |
INTEGER k |
42 |
_RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
43 |
_RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
44 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
45 |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
46 |
_RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
47 |
_RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
48 |
_RL myTime |
49 |
INTEGER myIter, myThid |
50 |
|
51 |
C !LOCAL VARIABLES: |
52 |
C == Local variables == |
53 |
INTEGER i,j |
54 |
_RL phxFac,phyFac, psFac |
55 |
_RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
56 |
_RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
57 |
#ifdef ALLOW_DIAGNOSTICS |
58 |
C Allow diagnosis of external forcing |
59 |
LOGICAL externForcDiagIsOn |
60 |
LOGICAL DIAGNOSTICS_IS_ON |
61 |
EXTERNAL DIAGNOSTICS_IS_ON |
62 |
_RL gUext(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
63 |
_RL gVext(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
#endif |
65 |
#ifdef ALLOW_CD_CODE |
66 |
_RL guCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
67 |
_RL gvCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
68 |
#endif |
69 |
CEOP |
70 |
|
71 |
C-- explicit part of the surface potential gradient is added in this S/R |
72 |
psFac = pfFacMom*(1. _d 0 - implicSurfPress) |
73 |
|
74 |
C-- factors for gradient (X & Y directions) of Hydrostatic Potential |
75 |
phxFac = pfFacMom |
76 |
phyFac = pfFacMom |
77 |
|
78 |
#ifdef ALLOW_DIAGNOSTICS |
79 |
externForcDiagIsOn = useDiagnostics .AND. momForcing |
80 |
IF ( externForcDiagIsOn ) THEN |
81 |
externForcDiagIsOn = DIAGNOSTICS_IS_ON('Um_Ext ',myThid) |
82 |
& .OR. DIAGNOSTICS_IS_ON('Vm_Ext ',myThid) |
83 |
ENDIF |
84 |
#endif /* ALLOW_DIAGNOSTICS */ |
85 |
|
86 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
87 |
|
88 |
C- Initialize local arrays (not really necessary but safer) |
89 |
DO j=1-Oly,sNy+Oly |
90 |
DO i=1-Olx,sNx+Olx |
91 |
gUtmp(i,j) = 0. _d 0 |
92 |
gVtmp(i,j) = 0. _d 0 |
93 |
#ifdef ALLOW_CD_CODE |
94 |
guCor(i,j) = 0. _d 0 |
95 |
gvCor(i,j) = 0. _d 0 |
96 |
#endif |
97 |
ENDDO |
98 |
ENDDO |
99 |
|
100 |
IF ( .NOT.staggerTimeStep .AND. .NOT. implicitIntGravWave ) THEN |
101 |
C-- Synchronous time step: add grad Phi_Hyp to gU,gV before doing Adams-Bashforth |
102 |
DO j=jMin,jMax |
103 |
DO i=iMin,iMax |
104 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) - phxFac*dPhiHydX(i,j) |
105 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) - phyFac*dPhiHydY(i,j) |
106 |
ENDDO |
107 |
ENDDO |
108 |
phxFac = 0. |
109 |
phyFac = 0. |
110 |
c ELSE |
111 |
C-- Stagger time step: grad Phi_Hyp will be added later |
112 |
ENDIF |
113 |
|
114 |
C-- Dissipation term inside the Adams-Bashforth: |
115 |
IF ( momViscosity .AND. momDissip_In_AB) THEN |
116 |
DO j=jMin,jMax |
117 |
DO i=iMin,iMax |
118 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) + guDissip(i,j) |
119 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) + gvDissip(i,j) |
120 |
ENDDO |
121 |
ENDDO |
122 |
ENDIF |
123 |
|
124 |
C-- Forcing term inside the Adams-Bashforth: |
125 |
IF ( momForcing .AND. momForcingOutAB.NE.1 ) THEN |
126 |
#ifdef ALLOW_DIAGNOSTICS |
127 |
IF ( useDiagnostics .AND. externForcDiagIsOn ) THEN |
128 |
DO j=1,sNy+1 |
129 |
DO i=1,sNx+1 |
130 |
gUext(i,j) = gU(i,j,k,bi,bj) |
131 |
gVext(i,j) = gV(i,j,k,bi,bj) |
132 |
ENDDO |
133 |
ENDDO |
134 |
ENDIF |
135 |
#endif /* ALLOW_DIAGNOSTICS */ |
136 |
|
137 |
CALL EXTERNAL_FORCING_U( |
138 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
139 |
I myTime,myThid) |
140 |
CALL EXTERNAL_FORCING_V( |
141 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
142 |
I myTime,myThid) |
143 |
|
144 |
#ifdef ALLOW_DIAGNOSTICS |
145 |
IF ( useDiagnostics .AND. externForcDiagIsOn ) THEN |
146 |
DO j=1,sNy+1 |
147 |
DO i=1,sNx+1 |
148 |
gUext(i,j) = gU(i,j,k,bi,bj)-gUext(i,j) |
149 |
gVext(i,j) = gV(i,j,k,bi,bj)-gVext(i,j) |
150 |
ENDDO |
151 |
ENDDO |
152 |
ENDIF |
153 |
#endif /* ALLOW_DIAGNOSTICS */ |
154 |
ENDIF |
155 |
|
156 |
IF (useCDscheme) THEN |
157 |
C- for CD-scheme, store gU,Vtmp = gU,V^n + dissip. + forcing |
158 |
IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN |
159 |
DO j=jMin,jMax |
160 |
DO i=iMin,iMax |
161 |
gUtmp(i,j) = gU(i,j,k,bi,bj) + guDissip(i,j) |
162 |
gVtmp(i,j) = gV(i,j,k,bi,bj) + gvDissip(i,j) |
163 |
ENDDO |
164 |
ENDDO |
165 |
ELSE |
166 |
DO j=jMin,jMax |
167 |
DO i=iMin,iMax |
168 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
169 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
170 |
ENDDO |
171 |
ENDDO |
172 |
ENDIF |
173 |
ENDIF |
174 |
|
175 |
C- Compute effective gU,gV_[n+1/2] terms (including Adams-Bashforth weights) |
176 |
C and save gU,gV_[n] into guNm1,gvNm1 for the next time step. |
177 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
178 |
CALL ADAMS_BASHFORTH3( |
179 |
I bi, bj, k, |
180 |
U gU, guNm, |
181 |
I mom_StartAB, myIter, myThid ) |
182 |
CALL ADAMS_BASHFORTH3( |
183 |
I bi, bj, k, |
184 |
U gV, gvNm, |
185 |
I mom_StartAB, myIter, myThid ) |
186 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
187 |
CALL ADAMS_BASHFORTH2( |
188 |
I bi, bj, k, |
189 |
U gU, guNm1, |
190 |
I mom_StartAB, myIter, myThid ) |
191 |
CALL ADAMS_BASHFORTH2( |
192 |
I bi, bj, k, |
193 |
U gV, gvNm1, |
194 |
I mom_StartAB, myIter, myThid ) |
195 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
196 |
|
197 |
C-- Forcing term outside the Adams-Bashforth: |
198 |
C (not recommended with CD-scheme ON) |
199 |
IF ( momForcing .AND. momForcingOutAB.EQ.1 ) THEN |
200 |
IF (useCDscheme) THEN |
201 |
DO j=jMin,jMax |
202 |
DO i=iMin,iMax |
203 |
gUtmp(i,j) = gUtmp(i,j) - gU(i,j,k,bi,bj) |
204 |
gVtmp(i,j) = gVtmp(i,j) - gV(i,j,k,bi,bj) |
205 |
ENDDO |
206 |
ENDDO |
207 |
ENDIF |
208 |
#ifdef ALLOW_DIAGNOSTICS |
209 |
IF ( useDiagnostics .AND. externForcDiagIsOn ) THEN |
210 |
DO j=1,sNy+1 |
211 |
DO i=1,sNx+1 |
212 |
gUext(i,j) = gU(i,j,k,bi,bj) |
213 |
gVext(i,j) = gV(i,j,k,bi,bj) |
214 |
ENDDO |
215 |
ENDDO |
216 |
ENDIF |
217 |
#endif /* ALLOW_DIAGNOSTICS */ |
218 |
|
219 |
CALL EXTERNAL_FORCING_U( |
220 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
221 |
I myTime,myThid) |
222 |
CALL EXTERNAL_FORCING_V( |
223 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
224 |
I myTime,myThid) |
225 |
|
226 |
#ifdef ALLOW_DIAGNOSTICS |
227 |
IF ( useDiagnostics .AND. externForcDiagIsOn ) THEN |
228 |
DO j=1,sNy+1 |
229 |
DO i=1,sNx+1 |
230 |
gUext(i,j) = gU(i,j,k,bi,bj)-gUext(i,j) |
231 |
gVext(i,j) = gV(i,j,k,bi,bj)-gVext(i,j) |
232 |
ENDDO |
233 |
ENDDO |
234 |
ENDIF |
235 |
#endif /* ALLOW_DIAGNOSTICS */ |
236 |
|
237 |
C- for CD-scheme, compute gU,Vtmp = gU,V^n + forcing |
238 |
IF (useCDscheme) THEN |
239 |
DO j=jMin,jMax |
240 |
DO i=iMin,iMax |
241 |
gUtmp(i,j) = gUtmp(i,j) + gU(i,j,k,bi,bj) |
242 |
gVtmp(i,j) = gVtmp(i,j) + gV(i,j,k,bi,bj) |
243 |
ENDDO |
244 |
ENDDO |
245 |
ENDIF |
246 |
ENDIF |
247 |
|
248 |
#ifdef ALLOW_CD_CODE |
249 |
IF (useCDscheme) THEN |
250 |
C- Step forward D-grid velocity using C-grid gU,Vtmp = gU,V^n + forcing |
251 |
C and return coriolis terms on C-grid (guCor,gvCor) |
252 |
CALL CD_CODE_SCHEME( |
253 |
I bi,bj,k, dPhiHydX,dPhiHydY, gUtmp,gVtmp, |
254 |
O guCor,gvCor, |
255 |
I myTime, myIter, myThid) |
256 |
DO j=jMin,jMax |
257 |
DO i=iMin,iMax |
258 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
259 |
& + guCor(i,j) |
260 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
261 |
& + gvCor(i,j) |
262 |
ENDDO |
263 |
ENDDO |
264 |
ELSE |
265 |
#endif /* ALLOW_CD_CODE */ |
266 |
DO j=jMin,jMax |
267 |
DO i=iMin,iMax |
268 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
269 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
270 |
ENDDO |
271 |
ENDDO |
272 |
#ifdef ALLOW_CD_CODE |
273 |
ENDIF |
274 |
#endif |
275 |
|
276 |
#ifdef NONLIN_FRSURF |
277 |
IF (.NOT. vectorInvariantMomentum |
278 |
& .AND. nonlinFreeSurf.GT.1) THEN |
279 |
IF (select_rStar.GT.0) THEN |
280 |
DO j=jMin,jMax |
281 |
DO i=iMin,iMax |
282 |
gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj) |
283 |
gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj) |
284 |
ENDDO |
285 |
ENDDO |
286 |
ELSE |
287 |
DO j=jMin,jMax |
288 |
DO i=iMin,iMax |
289 |
IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN |
290 |
gUtmp(i,j) = gUtmp(i,j) |
291 |
& *_hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj) |
292 |
ENDIF |
293 |
IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN |
294 |
gVtmp(i,j) = gVtmp(i,j) |
295 |
& *_hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj) |
296 |
ENDIF |
297 |
ENDDO |
298 |
ENDDO |
299 |
ENDIF |
300 |
ENDIF |
301 |
#endif /* NONLIN_FRSURF */ |
302 |
|
303 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
304 |
|
305 |
C-- Dissipation term outside the Adams-Bashforth: |
306 |
IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN |
307 |
DO j=jMin,jMax |
308 |
DO i=iMin,iMax |
309 |
gUtmp(i,j) = gUtmp(i,j) + guDissip(i,j) |
310 |
gVtmp(i,j) = gVtmp(i,j) + gvDissip(i,j) |
311 |
ENDDO |
312 |
ENDDO |
313 |
ENDIF |
314 |
|
315 |
C Step forward zonal velocity (store in Gu) |
316 |
DO j=jMin,jMax |
317 |
DO i=iMin,iMax |
318 |
gU(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) |
319 |
& +deltaTmom*( |
320 |
& gUtmp(i,j) |
321 |
& - psFac*phiSurfX(i,j) |
322 |
& - phxFac*dPhiHydX(i,j) |
323 |
& )*_maskW(i,j,k,bi,bj) |
324 |
ENDDO |
325 |
ENDDO |
326 |
|
327 |
C Step forward meridional velocity (store in Gv) |
328 |
DO j=jMin,jMax |
329 |
DO i=iMin,iMax |
330 |
gV(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) |
331 |
& +deltaTmom*( |
332 |
& gVtmp(i,j) |
333 |
& - psFac*phiSurfY(i,j) |
334 |
& - phyFac*dPhiHydY(i,j) |
335 |
& )*_maskS(i,j,k,bi,bj) |
336 |
ENDDO |
337 |
ENDDO |
338 |
|
339 |
#ifdef ALLOW_DIAGNOSTICS |
340 |
IF ( useDiagnostics .AND. externForcDiagIsOn ) THEN |
341 |
CALL DIAGNOSTICS_FILL(gUext,'Um_Ext ',k,1,2,bi,bj,myThid) |
342 |
CALL DIAGNOSTICS_FILL(gVext,'Vm_Ext ',k,1,2,bi,bj,myThid) |
343 |
ENDIF |
344 |
#ifdef ALLOW_CD_CODE |
345 |
IF ( useCDscheme .AND. useDiagnostics ) THEN |
346 |
CALL DIAGNOSTICS_FILL(guCor,'Um_Cori ',k,1,2,bi,bj,myThid) |
347 |
CALL DIAGNOSTICS_FILL(gvCor,'Vm_Cori ',k,1,2,bi,bj,myThid) |
348 |
ENDIF |
349 |
#endif /* ALLOW_CD_CODE */ |
350 |
#endif /* ALLOW_DIAGNOSTICS */ |
351 |
|
352 |
RETURN |
353 |
END |