/[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.77 - (show annotations) (download)
Mon Aug 13 23:26:56 2001 UTC (22 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.76: +16 -16 lines
Moved call to calc_viscosity before k=1,Nr loop.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.76 2001/08/13 18:05:26 heimbach Exp $
2 C $Name: checkpoint40pre6 $
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 fVer[STUV] o fVer: Vertical flux term - note fVer
62 C is "pipelined" in the vertical
63 C so we need an fVer for each
64 C variable.
65 C rhoK, rhoKM1 - Density at current level, and level above
66 C phiHyd - Hydrostatic part of the potential phiHydi.
67 C In z coords phiHydiHyd is the hydrostatic
68 C Potential (=pressure/rho0) anomaly
69 C In p coords phiHydiHyd is the geopotential
70 C surface height anomaly.
71 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
72 C phiSurfY or geopotentiel (atmos) in X and Y direction
73 C iMin, iMax - Ranges and sub-block indices on which calculations
74 C jMin, jMax are applied.
75 C bi, bj
76 C k, kup, - Index for layer above and below. kup and kDown
77 C kDown, km1 are switched with layer to be the appropriate
78 C index into fVerTerm.
79 C tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf.
80 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
81 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
82 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
83 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
88 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
89 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
92 _RL tauAB
93
94 C This is currently used by IVDC and Diagnostics
95 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96
97 INTEGER iMin, iMax
98 INTEGER jMin, jMax
99 INTEGER bi, bj
100 INTEGER i, j
101 INTEGER k, km1, kp1, kup, kDown
102
103 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
104 c CHARACTER*(MAX_LEN_MBUF) suff
105 c LOGICAL DIFFERENT_MULTIPLE
106 c EXTERNAL DIFFERENT_MULTIPLE
107 Cjmc(end)
108
109 C--- The algorithm...
110 C
111 C "Correction Step"
112 C =================
113 C Here we update the horizontal velocities with the surface
114 C pressure such that the resulting flow is either consistent
115 C with the free-surface evolution or the rigid-lid:
116 C U[n] = U* + dt x d/dx P
117 C V[n] = V* + dt x d/dy P
118 C
119 C "Calculation of Gs"
120 C ===================
121 C This is where all the accelerations and tendencies (ie.
122 C physics, parameterizations etc...) are calculated
123 C rho = rho ( theta[n], salt[n] )
124 C b = b(rho, theta)
125 C K31 = K31 ( rho )
126 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
127 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
128 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
129 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
130 C
131 C "Time-stepping" or "Prediction"
132 C ================================
133 C The models variables are stepped forward with the appropriate
134 C time-stepping scheme (currently we use Adams-Bashforth II)
135 C - For momentum, the result is always *only* a "prediction"
136 C in that the flow may be divergent and will be "corrected"
137 C later with a surface pressure gradient.
138 C - Normally for tracers the result is the new field at time
139 C level [n+1} *BUT* in the case of implicit diffusion the result
140 C is also *only* a prediction.
141 C - We denote "predictors" with an asterisk (*).
142 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
143 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
144 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
145 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
146 C With implicit diffusion:
147 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
148 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
149 C (1 + dt * K * d_zz) theta[n] = theta*
150 C (1 + dt * K * d_zz) salt[n] = salt*
151 C---
152
153 C-- Set up work arrays with valid (i.e. not NaN) values
154 C These inital values do not alter the numerical results. They
155 C just ensure that all memory references are to valid floating
156 C point numbers. This prevents spurious hardware signals due to
157 C uninitialised but inert locations.
158 DO j=1-OLy,sNy+OLy
159 DO i=1-OLx,sNx+OLx
160 DO k=1,Nr
161 phiHyd(i,j,k) = 0. _d 0
162 cph KappaRU(i,j,k) = 0. _d 0
163 cph KappaRV(i,j,k) = 0. _d 0
164 sigmaX(i,j,k) = 0. _d 0
165 sigmaY(i,j,k) = 0. _d 0
166 sigmaR(i,j,k) = 0. _d 0
167 ENDDO
168 rhoKM1 (i,j) = 0. _d 0
169 rhok (i,j) = 0. _d 0
170 phiSurfX(i,j) = 0. _d 0
171 phiSurfY(i,j) = 0. _d 0
172 ENDDO
173 ENDDO
174
175 #ifdef ALLOW_AUTODIFF_TAMC
176 C-- HPF directive to help TAMC
177 CHPF$ INDEPENDENT
178 #endif /* ALLOW_AUTODIFF_TAMC */
179
180 DO bj=myByLo(myThid),myByHi(myThid)
181
182 #ifdef ALLOW_AUTODIFF_TAMC
183 C-- HPF directive to help TAMC
184 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
185 CHPF$& ,phiHyd
186 CHPF$& ,KappaRU,KappaRV
187 CHPF$& )
188 #endif /* ALLOW_AUTODIFF_TAMC */
189
190 DO bi=myBxLo(myThid),myBxHi(myThid)
191
192 #ifdef ALLOW_AUTODIFF_TAMC
193 act1 = bi - myBxLo(myThid)
194 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
195
196 act2 = bj - myByLo(myThid)
197 max2 = myByHi(myThid) - myByLo(myThid) + 1
198
199 act3 = myThid - 1
200 max3 = nTx*nTy
201
202 act4 = ikey_dynamics - 1
203
204 ikey = (act1 + 1) + act2*max1
205 & + act3*max1*max2
206 & + act4*max1*max2*max3
207 #endif /* ALLOW_AUTODIFF_TAMC */
208
209 C-- Set up work arrays that need valid initial values
210 DO j=1-OLy,sNy+OLy
211 DO i=1-OLx,sNx+OLx
212 fVerU (i,j,1) = 0. _d 0
213 fVerU (i,j,2) = 0. _d 0
214 fVerV (i,j,1) = 0. _d 0
215 fVerV (i,j,2) = 0. _d 0
216 ENDDO
217 ENDDO
218
219 C-- Start computation of dynamics
220 iMin = 1-OLx+2
221 iMax = sNx+OLx-1
222 jMin = 1-OLy+2
223 jMax = sNy+OLy-1
224
225 #ifdef ALLOW_AUTODIFF_TAMC
226 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
227 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
228 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
229 #endif /* ALLOW_AUTODIFF_TAMC */
230
231 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
232 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
233 IF (implicSurfPress.NE.1.) THEN
234 CALL CALC_GRAD_PHI_SURF(
235 I bi,bj,iMin,iMax,jMin,jMax,
236 I etaN,
237 O phiSurfX,phiSurfY,
238 I myThid )
239 ENDIF
240
241 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
242 C-- Calculate the total vertical diffusivity
243 DO k=1,Nr
244 CALL CALC_VISCOSITY(
245 I bi,bj,iMin,iMax,jMin,jMax,k,
246 O KappaRU,KappaRV,
247 I myThid)
248 ENDDO
249 #endif
250
251 C-- Start of dynamics loop
252 DO k=1,Nr
253
254 C-- km1 Points to level above k (=k-1)
255 C-- kup Cycles through 1,2 to point to layer above
256 C-- kDown Cycles through 2,1 to point to current layer
257
258 km1 = MAX(1,k-1)
259 kp1 = MIN(k+1,Nr)
260 kup = 1+MOD(k+1,2)
261 kDown= 1+MOD(k,2)
262
263 #ifdef ALLOW_AUTODIFF_TAMC
264 kkey = (ikey-1)*Nr + k
265 #endif /* ALLOW_AUTODIFF_TAMC */
266
267 C-- Integrate hydrostatic balance for phiHyd with BC of
268 C phiHyd(z=0)=0
269 C distinguishe between Stagger and Non Stagger time stepping
270 IF (staggerTimeStep) THEN
271 CALL CALC_PHI_HYD(
272 I bi,bj,iMin,iMax,jMin,jMax,k,
273 I gTnm1, gSnm1,
274 U phiHyd,
275 I myThid )
276 ELSE
277 CALL CALC_PHI_HYD(
278 I bi,bj,iMin,iMax,jMin,jMax,k,
279 I theta, salt,
280 U phiHyd,
281 I myThid )
282 ENDIF
283
284 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
285 C and step forward storing the result in gUnm1, gVnm1, etc...
286 IF ( momStepping ) THEN
287 CALL CALC_MOM_RHS(
288 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
289 I phiHyd,KappaRU,KappaRV,
290 U fVerU, fVerV,
291 I myTime, myThid)
292 CALL TIMESTEP(
293 I bi,bj,iMin,iMax,jMin,jMax,k,
294 I phiHyd, phiSurfX, phiSurfY,
295 I myIter, myThid)
296
297 #ifdef ALLOW_OBCS
298 C-- Apply open boundary conditions
299 IF (useOBCS) THEN
300 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
301 END IF
302 #endif /* ALLOW_OBCS */
303
304 #ifdef ALLOW_AUTODIFF_TAMC
305 #ifdef INCLUDE_CD_CODE
306 ELSE
307 DO j=1-OLy,sNy+OLy
308 DO i=1-OLx,sNx+OLx
309 guCD(i,j,k,bi,bj) = 0.0
310 gvCD(i,j,k,bi,bj) = 0.0
311 END DO
312 END DO
313 #endif /* INCLUDE_CD_CODE */
314 #endif /* ALLOW_AUTODIFF_TAMC */
315 ENDIF
316
317
318 C-- end of dynamics k loop (1:Nr)
319 ENDDO
320
321
322
323 C-- Implicit viscosity
324 IF (implicitViscosity.AND.momStepping) THEN
325 #ifdef ALLOW_AUTODIFF_TAMC
326 idkey = iikey + 3
327 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
328 #endif /* ALLOW_AUTODIFF_TAMC */
329 CALL IMPLDIFF(
330 I bi, bj, iMin, iMax, jMin, jMax,
331 I deltaTmom, KappaRU,recip_HFacW,
332 U gUNm1,
333 I myThid )
334 #ifdef ALLOW_AUTODIFF_TAMC
335 idkey = iikey + 4
336 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
337 #endif /* ALLOW_AUTODIFF_TAMC */
338 CALL IMPLDIFF(
339 I bi, bj, iMin, iMax, jMin, jMax,
340 I deltaTmom, KappaRV,recip_HFacS,
341 U gVNm1,
342 I myThid )
343
344 #ifdef ALLOW_OBCS
345 C-- Apply open boundary conditions
346 IF (useOBCS) THEN
347 DO K=1,Nr
348 CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
349 ENDDO
350 END IF
351 #endif /* ALLOW_OBCS */
352
353 #ifdef INCLUDE_CD_CODE
354 #ifdef ALLOW_AUTODIFF_TAMC
355 idkey = iikey + 5
356 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
357 #endif /* ALLOW_AUTODIFF_TAMC */
358 CALL IMPLDIFF(
359 I bi, bj, iMin, iMax, jMin, jMax,
360 I deltaTmom, KappaRU,recip_HFacW,
361 U vVelD,
362 I myThid )
363 #ifdef ALLOW_AUTODIFF_TAMC
364 idkey = iikey + 6
365 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
366 #endif /* ALLOW_AUTODIFF_TAMC */
367 CALL IMPLDIFF(
368 I bi, bj, iMin, iMax, jMin, jMax,
369 I deltaTmom, KappaRV,recip_HFacS,
370 U uVelD,
371 I myThid )
372 #endif /* INCLUDE_CD_CODE */
373 C-- End If implicitViscosity.AND.momStepping
374 ENDIF
375
376 Cjmc : add for phiHyd output <- but not working if multi tile per CPU
377 c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
378 c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
379 c WRITE(suff,'(I10.10)') myIter+1
380 c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
381 c ENDIF
382 Cjmc(end)
383
384 #ifdef ALLOW_TIMEAVE
385 IF (taveFreq.GT.0.) THEN
386 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
387 I deltaTclock, bi, bj, myThid)
388 IF (ivdc_kappa.NE.0.) THEN
389 CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
390 I deltaTclock, bi, bj, myThid)
391 ENDIF
392 ENDIF
393 #endif /* ALLOW_TIMEAVE */
394
395 ENDDO
396 ENDDO
397
398 #ifndef EXCLUDE_DEBUGMODE
399 If (debugMode) THEN
400 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
401 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
402 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
403 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
404 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
405 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
406 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
407 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
408 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
409 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
410 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
411 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
412 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
413 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
414 ENDIF
415 #endif
416
417 RETURN
418 END

  ViewVC Help
Powered by ViewVC 1.1.22