1 |
C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.59 2014/08/14 16:51:10 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
#ifdef ALLOW_CD_CODE |
7 |
#include "CD_CODE_OPTIONS.h" |
8 |
#endif |
9 |
|
10 |
CBOP |
11 |
C !ROUTINE: TIMESTEP |
12 |
C !INTERFACE: |
13 |
SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, k, |
14 |
I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY, |
15 |
I guDissip, gvDissip, |
16 |
I myTime, myIter, myThid ) |
17 |
C !DESCRIPTION: \bv |
18 |
C *==========================================================* |
19 |
C | S/R TIMESTEP |
20 |
C | o Step model fields forward in time |
21 |
C *==========================================================* |
22 |
C \ev |
23 |
|
24 |
C !USES: |
25 |
IMPLICIT NONE |
26 |
C == Global variables == |
27 |
#include "SIZE.h" |
28 |
#include "EEPARAMS.h" |
29 |
#include "PARAMS.h" |
30 |
#include "GRID.h" |
31 |
#include "SURFACE.h" |
32 |
#include "RESTART.h" |
33 |
#include "DYNVARS.h" |
34 |
#ifdef ALLOW_NONHYDROSTATIC |
35 |
#include "NH_VARS.h" |
36 |
#endif |
37 |
|
38 |
C !INPUT/OUTPUT PARAMETERS: |
39 |
C == Routine Arguments == |
40 |
C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential |
41 |
C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean) |
42 |
C phiSurfY :: or geopotential (atmos) in X and Y direction |
43 |
C guDissip :: dissipation tendency (all explicit terms), u component |
44 |
C gvDissip :: dissipation tendency (all explicit terms), v component |
45 |
|
46 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
47 |
INTEGER k |
48 |
_RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
49 |
_RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
50 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
51 |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
_RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
53 |
_RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
54 |
_RL myTime |
55 |
INTEGER myIter, myThid |
56 |
|
57 |
C !LOCAL VARIABLES: |
58 |
C == Local variables == |
59 |
C guExt :: forcing tendency, u component |
60 |
C gvExt :: forcing tendency, v component |
61 |
C gu_AB :: tendency increment from Adams-Bashforth, u component |
62 |
C gv_AB :: tendency increment from Adams-Bashforth, v component |
63 |
INTEGER i,j |
64 |
_RL phxFac,phyFac, psFac |
65 |
_RL guExt(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
66 |
_RL gvExt(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
67 |
_RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
68 |
_RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
69 |
_RL gu_AB(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
70 |
_RL gv_AB(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
71 |
#ifdef ALLOW_NONHYDROSTATIC |
72 |
_RL nhFac |
73 |
#endif |
74 |
#ifdef ALLOW_CD_CODE |
75 |
_RL guCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
76 |
_RL gvCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
77 |
#endif |
78 |
CEOP |
79 |
|
80 |
C-- explicit part of the surface potential gradient is added in this S/R |
81 |
psFac = pfFacMom*(1. _d 0 - implicSurfPress) |
82 |
& *recip_deepFacC(k)*recip_rhoFacC(k) |
83 |
|
84 |
C-- factors for gradient (X & Y directions) of Hydrostatic Potential |
85 |
phxFac = pfFacMom |
86 |
phyFac = pfFacMom |
87 |
|
88 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
89 |
|
90 |
C- Initialize local arrays (not really necessary for all but safer) |
91 |
DO j=1-OLy,sNy+OLy |
92 |
DO i=1-OLx,sNx+OLx |
93 |
guExt(i,j) = 0. _d 0 |
94 |
gvExt(i,j) = 0. _d 0 |
95 |
gUtmp(i,j) = 0. _d 0 |
96 |
gVtmp(i,j) = 0. _d 0 |
97 |
#ifdef ALLOW_CD_CODE |
98 |
guCor(i,j) = 0. _d 0 |
99 |
gvCor(i,j) = 0. _d 0 |
100 |
#endif |
101 |
ENDDO |
102 |
ENDDO |
103 |
|
104 |
IF ( momForcing ) THEN |
105 |
C-- Collect forcing term in local array guExt,gvExt: |
106 |
CALL APPLY_FORCING_U( |
107 |
U guExt, |
108 |
I iMin,iMax,jMin,jMax, k, bi,bj, |
109 |
I myTime, myIter, myThid ) |
110 |
CALL APPLY_FORCING_V( |
111 |
U gvExt, |
112 |
I iMin,iMax,jMin,jMax, k, bi,bj, |
113 |
I myTime, myIter, myThid ) |
114 |
ENDIF |
115 |
|
116 |
IF ( .NOT.staggerTimeStep .AND. .NOT. implicitIntGravWave ) THEN |
117 |
C-- Synchronous time step: add grad Phi_Hyp to gU,gV before doing Adams-Bashforth |
118 |
DO j=jMin,jMax |
119 |
DO i=iMin,iMax |
120 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) - phxFac*dPhiHydX(i,j) |
121 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) - phyFac*dPhiHydY(i,j) |
122 |
ENDDO |
123 |
ENDDO |
124 |
phxFac = 0. |
125 |
phyFac = 0. |
126 |
c ELSE |
127 |
C-- Stagger time step: grad Phi_Hyp will be added later |
128 |
ENDIF |
129 |
|
130 |
C-- Dissipation term inside the Adams-Bashforth: |
131 |
IF ( momViscosity .AND. momDissip_In_AB) THEN |
132 |
DO j=jMin,jMax |
133 |
DO i=iMin,iMax |
134 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) + guDissip(i,j) |
135 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) + gvDissip(i,j) |
136 |
ENDDO |
137 |
ENDDO |
138 |
ENDIF |
139 |
|
140 |
C-- Forcing term inside the Adams-Bashforth: |
141 |
IF ( momForcing .AND. momForcingOutAB.NE.1 ) THEN |
142 |
DO j=jMin,jMax |
143 |
DO i=iMin,iMax |
144 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) + guExt(i,j) |
145 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) + gvExt(i,j) |
146 |
ENDDO |
147 |
ENDDO |
148 |
ENDIF |
149 |
|
150 |
#ifdef CD_CODE_NO_AB_MOMENTUM |
151 |
IF ( useCDscheme ) THEN |
152 |
C- CD-scheme, before doing AB, store gU,Vtmp = gU,V^n (+dissip. +forcing) |
153 |
DO j=jMin,jMax |
154 |
DO i=iMin,iMax |
155 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
156 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
157 |
ENDDO |
158 |
ENDDO |
159 |
ENDIF |
160 |
#endif /* CD_CODE_NO_AB_MOMENTUM */ |
161 |
|
162 |
C- Compute effective gU,gV_[n+1/2] terms (including Adams-Bashforth weights) |
163 |
C and save gU,gV_[n] into guNm1,gvNm1 for the next time step. |
164 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
165 |
CALL ADAMS_BASHFORTH3( |
166 |
I bi, bj, k, Nr, |
167 |
U gU(1-OLx,1-OLy,1,bi,bj), guNm, |
168 |
O gu_AB, |
169 |
I mom_StartAB, myIter, myThid ) |
170 |
CALL ADAMS_BASHFORTH3( |
171 |
I bi, bj, k, Nr, |
172 |
U gV(1-OLx,1-OLy,1,bi,bj), gvNm, |
173 |
O gv_AB, |
174 |
I mom_StartAB, myIter, myThid ) |
175 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
176 |
CALL ADAMS_BASHFORTH2( |
177 |
I bi, bj, k, Nr, |
178 |
U gU(1-OLx,1-OLy,1,bi,bj), |
179 |
U guNm1(1-OLx,1-OLy,1,bi,bj), |
180 |
O gu_AB, |
181 |
I mom_StartAB, myIter, myThid ) |
182 |
CALL ADAMS_BASHFORTH2( |
183 |
I bi, bj, k, Nr, |
184 |
U gV(1-OLx,1-OLy,1,bi,bj), |
185 |
U gvNm1(1-OLx,1-OLy,1,bi,bj), |
186 |
O gv_AB, |
187 |
I mom_StartAB, myIter, myThid ) |
188 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
189 |
#ifdef ALLOW_DIAGNOSTICS |
190 |
IF ( useDiagnostics ) THEN |
191 |
CALL DIAGNOSTICS_FILL(gu_AB,'AB_gU ',k,1,2,bi,bj,myThid) |
192 |
CALL DIAGNOSTICS_FILL(gv_AB,'AB_gV ',k,1,2,bi,bj,myThid) |
193 |
ENDIF |
194 |
#endif /* ALLOW_DIAGNOSTICS */ |
195 |
|
196 |
C- Make a local copy in gU,Vtmp of gU,V^n+1/2 (+dissip. +forcing) |
197 |
#ifdef CD_CODE_NO_AB_MOMENTUM |
198 |
IF ( .NOT.useCDscheme ) THEN |
199 |
#endif |
200 |
DO j=jMin,jMax |
201 |
DO i=iMin,iMax |
202 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
203 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
204 |
ENDDO |
205 |
ENDDO |
206 |
#ifdef CD_CODE_NO_AB_MOMENTUM |
207 |
ENDIF |
208 |
#endif |
209 |
|
210 |
C-- Forcing term outside the Adams-Bashforth: |
211 |
IF ( momForcing .AND. momForcingOutAB.EQ.1 ) THEN |
212 |
DO j=jMin,jMax |
213 |
DO i=iMin,iMax |
214 |
gUtmp(i,j) = gUtmp(i,j) + guExt(i,j) |
215 |
gVtmp(i,j) = gVtmp(i,j) + gvExt(i,j) |
216 |
ENDDO |
217 |
ENDDO |
218 |
ENDIF |
219 |
|
220 |
C-- Dissipation term outside the Adams-Bashforth: |
221 |
IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN |
222 |
DO j=jMin,jMax |
223 |
DO i=iMin,iMax |
224 |
gUtmp(i,j) = gUtmp(i,j) + guDissip(i,j) |
225 |
gVtmp(i,j) = gVtmp(i,j) + gvDissip(i,j) |
226 |
ENDDO |
227 |
ENDDO |
228 |
ENDIF |
229 |
|
230 |
#ifdef ALLOW_CD_CODE |
231 |
IF ( useCDscheme ) THEN |
232 |
C- Step forward D-grid velocity using C-grid gU,Vtmp = gU,V +dissip +forcing |
233 |
C and return coriolis terms on C-grid (guCor,gvCor) |
234 |
CALL CD_CODE_SCHEME( |
235 |
I bi,bj,k, dPhiHydX,dPhiHydY, gUtmp,gVtmp, |
236 |
O guCor,gvCor, |
237 |
I myTime, myIter, myThid) |
238 |
|
239 |
#ifdef CD_CODE_NO_AB_MOMENTUM |
240 |
IF ( momForcing .AND. momForcingOutAB.EQ.1 ) THEN |
241 |
DO j=jMin,jMax |
242 |
DO i=iMin,iMax |
243 |
gUtmp(i,j) = ( gU(i,j,k,bi,bj) + guExt(i,j) ) + guCor(i,j) |
244 |
gVtmp(i,j) = ( gV(i,j,k,bi,bj) + gvExt(i,j) ) + gvCor(i,j) |
245 |
ENDDO |
246 |
ENDDO |
247 |
ELSE |
248 |
DO j=jMin,jMax |
249 |
DO i=iMin,iMax |
250 |
gUtmp(i,j) = gU(i,j,k,bi,bj) + guCor(i,j) |
251 |
gVtmp(i,j) = gV(i,j,k,bi,bj) + gvCor(i,j) |
252 |
ENDDO |
253 |
ENDDO |
254 |
ENDIF |
255 |
IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN |
256 |
DO j=jMin,jMax |
257 |
DO i=iMin,iMax |
258 |
gUtmp(i,j) = gUtmp(i,j) + guDissip(i,j) |
259 |
gVtmp(i,j) = gVtmp(i,j) + gvDissip(i,j) |
260 |
ENDDO |
261 |
ENDDO |
262 |
ENDIF |
263 |
#else /* CD_CODE_NO_AB_MOMENTUM */ |
264 |
DO j=jMin,jMax |
265 |
DO i=iMin,iMax |
266 |
gUtmp(i,j) = gUtmp(i,j) + guCor(i,j) |
267 |
gVtmp(i,j) = gVtmp(i,j) + gvCor(i,j) |
268 |
ENDDO |
269 |
ENDDO |
270 |
#endif /* CD_CODE_NO_AB_MOMENTUM */ |
271 |
ENDIF |
272 |
#endif /* ALLOW_CD_CODE */ |
273 |
|
274 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
275 |
|
276 |
#ifdef NONLIN_FRSURF |
277 |
IF ( .NOT. vectorInvariantMomentum |
278 |
& .AND. nonlinFreeSurf.GT.1 ) THEN |
279 |
IF ( select_rStar.GT.0 ) THEN |
280 |
# ifndef DISABLE_RSTAR_CODE |
281 |
DO j=jMin,jMax |
282 |
DO i=iMin,iMax |
283 |
gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj) |
284 |
gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj) |
285 |
ENDDO |
286 |
ENDDO |
287 |
# endif /* DISABLE_RSTAR_CODE */ |
288 |
ELSEIF ( selectSigmaCoord.NE.0 ) THEN |
289 |
# ifndef DISABLE_SIGMA_CODE |
290 |
DO j=jMin,jMax |
291 |
DO i=iMin,iMax |
292 |
gUtmp(i,j) = gUtmp(i,j) |
293 |
& /( 1. _d 0 + dEtaWdt(i,j,bi,bj)*deltaTFreeSurf |
294 |
& *dBHybSigF(k)*recip_drF(k) |
295 |
& *recip_hFacW(i,j,k,bi,bj) |
296 |
& ) |
297 |
gVtmp(i,j) = gVtmp(i,j) |
298 |
& /( 1. _d 0 + dEtaSdt(i,j,bi,bj)*deltaTFreeSurf |
299 |
& *dBHybSigF(k)*recip_drF(k) |
300 |
& *recip_hFacS(i,j,k,bi,bj) |
301 |
& ) |
302 |
ENDDO |
303 |
ENDDO |
304 |
# endif /* DISABLE_SIGMA_CODE */ |
305 |
ELSE |
306 |
DO j=jMin,jMax |
307 |
DO i=iMin,iMax |
308 |
IF ( k.EQ.kSurfW(i,j,bi,bj) ) THEN |
309 |
gUtmp(i,j) = gUtmp(i,j) |
310 |
& *_hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj) |
311 |
ENDIF |
312 |
IF ( k.EQ.kSurfS(i,j,bi,bj) ) THEN |
313 |
gVtmp(i,j) = gVtmp(i,j) |
314 |
& *_hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj) |
315 |
ENDIF |
316 |
ENDDO |
317 |
ENDDO |
318 |
ENDIF |
319 |
ENDIF |
320 |
#endif /* NONLIN_FRSURF */ |
321 |
|
322 |
#ifdef ALLOW_NONHYDROSTATIC |
323 |
C-- explicit part of the NH pressure gradient is added here |
324 |
IF ( use3Dsolver .AND. implicitNHPress.NE.1. _d 0 ) THEN |
325 |
nhFac = pfFacMom*(1. _d 0 - implicitNHPress) |
326 |
& *recip_deepFacC(k)*recip_rhoFacC(k) |
327 |
IF ( exactConserv ) THEN |
328 |
DO j=jMin,jMax |
329 |
DO i=iMin,iMax |
330 |
gUtmp(i,j) = gUtmp(i,j) |
331 |
& - nhFac*_recip_dxC(i,j,bi,bj)* |
332 |
& ( (phi_nh(i,j,k,bi,bj)-phi_nh(i-1,j,k,bi,bj)) |
333 |
& -( dPhiNH(i,j,bi,bj) - dPhiNH(i-1,j,bi,bj) ) |
334 |
& ) |
335 |
gVtmp(i,j) = gVtmp(i,j) |
336 |
& - nhFac*_recip_dyC(i,j,bi,bj)* |
337 |
& ( (phi_nh(i,j,k,bi,bj)-phi_nh(i,j-1,k,bi,bj)) |
338 |
& -( dPhiNH(i,j,bi,bj) - dPhiNH(i,j-1,bi,bj) ) |
339 |
& ) |
340 |
ENDDO |
341 |
ENDDO |
342 |
ELSE |
343 |
DO j=jMin,jMax |
344 |
DO i=iMin,iMax |
345 |
gUtmp(i,j) = gUtmp(i,j) |
346 |
& - nhFac*_recip_dxC(i,j,bi,bj)* |
347 |
& (phi_nh(i,j,k,bi,bj)-phi_nh(i-1,j,k,bi,bj)) |
348 |
gVtmp(i,j) = gVtmp(i,j) |
349 |
& - nhFac*_recip_dyC(i,j,bi,bj)* |
350 |
& (phi_nh(i,j,k,bi,bj)-phi_nh(i,j-1,k,bi,bj)) |
351 |
ENDDO |
352 |
ENDDO |
353 |
ENDIF |
354 |
ENDIF |
355 |
#endif /* ALLOW_NONHYDROSTATIC */ |
356 |
|
357 |
C Step forward zonal velocity (store in Gu) |
358 |
DO j=jMin,jMax |
359 |
DO i=iMin,iMax |
360 |
gU(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) |
361 |
& +deltaTMom*( |
362 |
& gUtmp(i,j) |
363 |
& - psFac*phiSurfX(i,j) |
364 |
& - phxFac*dPhiHydX(i,j) |
365 |
& )*_maskW(i,j,k,bi,bj) |
366 |
ENDDO |
367 |
ENDDO |
368 |
|
369 |
C Step forward meridional velocity (store in Gv) |
370 |
DO j=jMin,jMax |
371 |
DO i=iMin,iMax |
372 |
gV(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) |
373 |
& +deltaTMom*( |
374 |
& gVtmp(i,j) |
375 |
& - psFac*phiSurfY(i,j) |
376 |
& - phyFac*dPhiHydY(i,j) |
377 |
& )*_maskS(i,j,k,bi,bj) |
378 |
ENDDO |
379 |
ENDDO |
380 |
|
381 |
#ifdef ALLOW_DIAGNOSTICS |
382 |
IF ( momViscosity .AND. useDiagnostics ) THEN |
383 |
CALL DIAGNOSTICS_FILL( guDissip,'Um_Diss ',k,1,2,bi,bj,myThid ) |
384 |
CALL DIAGNOSTICS_FILL( gvDissip,'Vm_Diss ',k,1,2,bi,bj,myThid ) |
385 |
ENDIF |
386 |
IF ( momForcing .AND. useDiagnostics ) THEN |
387 |
CALL DIAGNOSTICS_FILL( guExt,'Um_Ext ',k,1,2,bi,bj,myThid ) |
388 |
CALL DIAGNOSTICS_FILL( gvExt,'Vm_Ext ',k,1,2,bi,bj,myThid ) |
389 |
ENDIF |
390 |
#ifdef ALLOW_CD_CODE |
391 |
IF ( useCDscheme .AND. useDiagnostics ) THEN |
392 |
CALL DIAGNOSTICS_FILL( guCor,'Um_Cori ',k,1,2,bi,bj,myThid ) |
393 |
CALL DIAGNOSTICS_FILL( gvCor,'Vm_Cori ',k,1,2,bi,bj,myThid ) |
394 |
ENDIF |
395 |
#endif /* ALLOW_CD_CODE */ |
396 |
#endif /* ALLOW_DIAGNOSTICS */ |
397 |
|
398 |
RETURN |
399 |
END |