/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.75 - (show annotations) (download)
Fri Aug 3 19:06:11 2001 UTC (22 years, 9 months ago) by adcroft
Branch: MAIN
Changes since 1.74: +16 -467 lines
Split dynamics.F into dynamics.F and thermodynamics.F
 - idea is to make algorithm more transparent???
 - probably less efficient
 - has exchanges at end of thermodynamics.F (which are needed
   if using staggered time-stepping with the cube OR using AIM.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.74 2001/07/30 20:37:45 heimbach Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
7 C /==========================================================\
8 C | SUBROUTINE DYNAMICS |
9 C | o Controlling routine for the explicit part of the model |
10 C | dynamics. |
11 C |==========================================================|
12 C | This routine evaluates the "dynamics" terms for each |
13 C | block of ocean in turn. Because the blocks of ocean have |
14 C | overlap regions they are independent of one another. |
15 C | If terms involving lateral integrals are needed in this |
16 C | routine care will be needed. Similarly finite-difference |
17 C | operations with stencils wider than the overlap region |
18 C | require special consideration. |
19 C | Notes |
20 C | ===== |
21 C | C*P* comments indicating place holders for which code is |
22 C | presently being developed. |
23 C \==========================================================/
24 IMPLICIT NONE
25
26 C == Global variables ===
27 #include "SIZE.h"
28 #include "EEPARAMS.h"
29 #include "PARAMS.h"
30 #include "DYNVARS.h"
31 #include "GRID.h"
32 #ifdef ALLOW_PASSIVE_TRACER
33 #include "TR1.h"
34 #endif
35
36 #ifdef ALLOW_AUTODIFF_TAMC
37 # include "tamc.h"
38 # include "tamc_keys.h"
39 # include "FFIELDS.h"
40 # ifdef ALLOW_KPP
41 # include "KPP.h"
42 # endif
43 # ifdef ALLOW_GMREDI
44 # include "GMREDI.h"
45 # endif
46 #endif /* ALLOW_AUTODIFF_TAMC */
47
48 #ifdef ALLOW_TIMEAVE
49 #include "TIMEAVE_STATV.h"
50 #endif
51
52 C == Routine arguments ==
53 C myTime - Current time in simulation
54 C myIter - Current iteration number in simulation
55 C myThid - Thread number for this instance of the routine.
56 _RL myTime
57 INTEGER myIter
58 INTEGER myThid
59
60 C == Local variables
61 C xA, yA - Per block temporaries holding face areas
62 C uTrans, vTrans, rTrans - Per block temporaries holding flow
63 C transport
64 C o uTrans: Zonal transport
65 C o vTrans: Meridional transport
66 C o rTrans: Vertical transport
67 C maskUp o maskUp: land/water mask for W points
68 C fVer[STUV] o fVer: Vertical flux term - note fVer
69 C is "pipelined" in the vertical
70 C so we need an fVer for each
71 C variable.
72 C rhoK, rhoKM1 - Density at current level, and level above
73 C phiHyd - Hydrostatic part of the potential phiHydi.
74 C In z coords phiHydiHyd is the hydrostatic
75 C Potential (=pressure/rho0) anomaly
76 C In p coords phiHydiHyd is the geopotential
77 C surface height anomaly.
78 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
79 C phiSurfY or geopotentiel (atmos) in X and Y direction
80 C KappaRT, - Total diffusion in vertical for T and S.
81 C KappaRS (background + spatially varying, isopycnal term).
82 C iMin, iMax - Ranges and sub-block indices on which calculations
83 C jMin, jMax are applied.
84 C bi, bj
85 C k, kup, - Index for layer above and below. kup and kDown
86 C kDown, km1 are switched with layer to be the appropriate
87 C index into fVerTerm.
88 C tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf.
89 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
96 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
97 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
98 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
99 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
106 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
107 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
108 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
109 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
110 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
112 _RL tauAB
113
114 C This is currently used by IVDC and Diagnostics
115 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
116
117 INTEGER iMin, iMax
118 INTEGER jMin, jMax
119 INTEGER bi, bj
120 INTEGER i, j
121 INTEGER k, km1, kup, kDown
122
123 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
124 c CHARACTER*(MAX_LEN_MBUF) suff
125 c LOGICAL DIFFERENT_MULTIPLE
126 c EXTERNAL DIFFERENT_MULTIPLE
127 Cjmc(end)
128
129 C--- The algorithm...
130 C
131 C "Correction Step"
132 C =================
133 C Here we update the horizontal velocities with the surface
134 C pressure such that the resulting flow is either consistent
135 C with the free-surface evolution or the rigid-lid:
136 C U[n] = U* + dt x d/dx P
137 C V[n] = V* + dt x d/dy P
138 C
139 C "Calculation of Gs"
140 C ===================
141 C This is where all the accelerations and tendencies (ie.
142 C physics, parameterizations etc...) are calculated
143 C rho = rho ( theta[n], salt[n] )
144 C b = b(rho, theta)
145 C K31 = K31 ( rho )
146 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
147 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
148 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
149 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
150 C
151 C "Time-stepping" or "Prediction"
152 C ================================
153 C The models variables are stepped forward with the appropriate
154 C time-stepping scheme (currently we use Adams-Bashforth II)
155 C - For momentum, the result is always *only* a "prediction"
156 C in that the flow may be divergent and will be "corrected"
157 C later with a surface pressure gradient.
158 C - Normally for tracers the result is the new field at time
159 C level [n+1} *BUT* in the case of implicit diffusion the result
160 C is also *only* a prediction.
161 C - We denote "predictors" with an asterisk (*).
162 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
163 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
164 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
165 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
166 C With implicit diffusion:
167 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
168 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
169 C (1 + dt * K * d_zz) theta[n] = theta*
170 C (1 + dt * K * d_zz) salt[n] = salt*
171 C---
172
173 DO bj=myByLo(myThid),myByHi(myThid)
174 DO bi=myBxLo(myThid),myBxHi(myThid)
175 Ccs-
176
177 C-- Start computation of dynamics
178 iMin = 1-OLx+2
179 iMax = sNx+OLx-1
180 jMin = 1-OLy+2
181 jMax = sNy+OLy-1
182
183 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
184 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
185 IF (implicSurfPress.NE.1.) THEN
186 CALL CALC_GRAD_PHI_SURF(
187 I bi,bj,iMin,iMax,jMin,jMax,
188 I etaN,
189 O phiSurfX,phiSurfY,
190 I myThid )
191 ENDIF
192
193 C-- Start of dynamics loop
194 DO k=1,Nr
195
196 C-- km1 Points to level above k (=k-1)
197 C-- kup Cycles through 1,2 to point to layer above
198 C-- kDown Cycles through 2,1 to point to current layer
199
200 km1 = MAX(1,k-1)
201 kup = 1+MOD(k+1,2)
202 kDown= 1+MOD(k,2)
203
204 C-- Integrate hydrostatic balance for phiHyd with BC of
205 C phiHyd(z=0)=0
206 C distinguishe between Stagger and Non Stagger time stepping
207 IF (staggerTimeStep) THEN
208 CALL CALC_PHI_HYD(
209 I bi,bj,iMin,iMax,jMin,jMax,k,
210 I gTnm1, gSnm1,
211 U phiHyd,
212 I myThid )
213 ELSE
214 CALL CALC_PHI_HYD(
215 I bi,bj,iMin,iMax,jMin,jMax,k,
216 I theta, salt,
217 U phiHyd,
218 I myThid )
219 ENDIF
220
221 #ifdef ALLOW_AUTODIFF_TAMC
222 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
223 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
224 #endif /* ALLOW_AUTODIFF_TAMC */
225
226 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
227 C-- Calculate the total vertical diffusivity
228 CALL CALC_DIFFUSIVITY(
229 I bi,bj,iMin,iMax,jMin,jMax,k,
230 I maskUp,
231 O KappaRT,KappaRS,KappaRU,KappaRV,
232 I myThid)
233 #endif
234
235 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
236 C and step forward storing the result in gUnm1, gVnm1, etc...
237 IF ( momStepping ) THEN
238 CALL CALC_MOM_RHS(
239 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
240 I phiHyd,KappaRU,KappaRV,
241 U fVerU, fVerV,
242 I myTime, myThid)
243 CALL TIMESTEP(
244 I bi,bj,iMin,iMax,jMin,jMax,k,
245 I phiHyd, phiSurfX, phiSurfY,
246 I myIter, myThid)
247
248 #ifdef ALLOW_OBCS
249 C-- Apply open boundary conditions
250 IF (useOBCS) THEN
251 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
252 END IF
253 #endif /* ALLOW_OBCS */
254
255 #ifdef ALLOW_AUTODIFF_TAMC
256 #ifdef INCLUDE_CD_CODE
257 ELSE
258 DO j=1-OLy,sNy+OLy
259 DO i=1-OLx,sNx+OLx
260 guCD(i,j,k,bi,bj) = 0.0
261 gvCD(i,j,k,bi,bj) = 0.0
262 END DO
263 END DO
264 #endif /* INCLUDE_CD_CODE */
265 #endif /* ALLOW_AUTODIFF_TAMC */
266 ENDIF
267
268
269 C-- end of dynamics k loop (1:Nr)
270 ENDDO
271
272
273
274 C-- Implicit viscosity
275 IF (implicitViscosity.AND.momStepping) THEN
276 #ifdef ALLOW_AUTODIFF_TAMC
277 idkey = iikey + 3
278 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
279 #endif /* ALLOW_AUTODIFF_TAMC */
280 CALL IMPLDIFF(
281 I bi, bj, iMin, iMax, jMin, jMax,
282 I deltaTmom, KappaRU,recip_HFacW,
283 U gUNm1,
284 I myThid )
285 #ifdef ALLOW_AUTODIFF_TAMC
286 idkey = iikey + 4
287 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
288 #endif /* ALLOW_AUTODIFF_TAMC */
289 CALL IMPLDIFF(
290 I bi, bj, iMin, iMax, jMin, jMax,
291 I deltaTmom, KappaRV,recip_HFacS,
292 U gVNm1,
293 I myThid )
294
295 #ifdef ALLOW_OBCS
296 C-- Apply open boundary conditions
297 IF (useOBCS) THEN
298 DO K=1,Nr
299 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
300 ENDDO
301 END IF
302 #endif /* ALLOW_OBCS */
303
304 #ifdef INCLUDE_CD_CODE
305 #ifdef ALLOW_AUTODIFF_TAMC
306 idkey = iikey + 5
307 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
308 #endif /* ALLOW_AUTODIFF_TAMC */
309 CALL IMPLDIFF(
310 I bi, bj, iMin, iMax, jMin, jMax,
311 I deltaTmom, KappaRU,recip_HFacW,
312 U vVelD,
313 I myThid )
314 #ifdef ALLOW_AUTODIFF_TAMC
315 idkey = iikey + 6
316 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
317 #endif /* ALLOW_AUTODIFF_TAMC */
318 CALL IMPLDIFF(
319 I bi, bj, iMin, iMax, jMin, jMax,
320 I deltaTmom, KappaRV,recip_HFacS,
321 U uVelD,
322 I myThid )
323 #endif /* INCLUDE_CD_CODE */
324 C-- End If implicitViscosity.AND.momStepping
325 ENDIF
326
327 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
328 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
329 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
330 c WRITE(suff,'(I10.10)') myIter+1
331 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
332 c ENDIF
333 Cjmc(end)
334
335 #ifdef ALLOW_TIMEAVE
336 IF (taveFreq.GT.0.) THEN
337 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
338 I deltaTclock, bi, bj, myThid)
339 IF (ivdc_kappa.NE.0.) THEN
340 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
341 I deltaTclock, bi, bj, myThid)
342 ENDIF
343 ENDIF
344 #endif /* ALLOW_TIMEAVE */
345
346 ENDDO
347 ENDDO
348
349 #ifndef EXCLUDE_DEBUGMODE
350 If (debugMode) THEN
351 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
352 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
353 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
354 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
355 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
356 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
357 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
358 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
359 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
360 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
361 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
362 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
363 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
364 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
365 ENDIF
366 #endif
367
368 RETURN
369 END

  ViewVC Help
Powered by ViewVC 1.1.22