/[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.55 - (show annotations) (download)
Mon Jan 8 16:37:43 2001 UTC (23 years, 4 months ago) by heimbach
Branch: MAIN
Changes since 1.54: +64 -16 lines
Added or modified store directives for TAMC.
Updates are adopted from ecco_c32_e2.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.54 2000/11/13 16:32:57 heimbach Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6 C /==========================================================\
7 C | SUBROUTINE DYNAMICS |
8 C | o Controlling routine for the explicit part of the model |
9 C | dynamics. |
10 C |==========================================================|
11 C | This routine evaluates the "dynamics" terms for each |
12 C | block of ocean in turn. Because the blocks of ocean have |
13 C | overlap regions they are independent of one another. |
14 C | If terms involving lateral integrals are needed in this |
15 C | routine care will be needed. Similarly finite-difference |
16 C | operations with stencils wider than the overlap region |
17 C | require special consideration. |
18 C | Notes |
19 C | ===== |
20 C | C*P* comments indicating place holders for which code is |
21 C | presently being developed. |
22 C \==========================================================/
23 IMPLICIT NONE
24
25 C == Global variables ===
26 #include "SIZE.h"
27 #include "EEPARAMS.h"
28 #include "CG2D.h"
29 #include "PARAMS.h"
30 #include "DYNVARS.h"
31 #include "GRID.h"
32
33 #ifdef ALLOW_AUTODIFF_TAMC
34 # include "tamc.h"
35 # include "tamc_keys.h"
36 #endif /* ALLOW_AUTODIFF_TAMC */
37
38 #ifdef ALLOW_KPP
39 # include "KPP.h"
40 #endif
41
42 C == Routine arguments ==
43 C myTime - Current time in simulation
44 C myIter - Current iteration number in simulation
45 C myThid - Thread number for this instance of the routine.
46 _RL myTime
47 INTEGER myIter
48 INTEGER myThid
49
50 C == Local variables
51 C xA, yA - Per block temporaries holding face areas
52 C uTrans, vTrans, rTrans - Per block temporaries holding flow
53 C transport
54 C rVel o uTrans: Zonal transport
55 C o vTrans: Meridional transport
56 C o rTrans: Vertical transport
57 C o rVel: Vertical velocity at upper and
58 C lower cell faces.
59 C maskC,maskUp o maskC: land/water mask for tracer cells
60 C o maskUp: land/water mask for W points
61 C aTerm, xTerm, cTerm - Work arrays for holding separate terms in
62 C mTerm, pTerm, tendency equations.
63 C fZon, fMer, fVer[STUV] o aTerm: Advection term
64 C o xTerm: Mixing term
65 C o cTerm: Coriolis term
66 C o mTerm: Metric term
67 C o pTerm: Pressure term
68 C o fZon: Zonal flux term
69 C o fMer: Meridional flux term
70 C o fVer: Vertical flux term - note fVer
71 C is "pipelined" in the vertical
72 C so we need an fVer for each
73 C variable.
74 C rhoK, rhoKM1 - Density at current level, level above and level
75 C below.
76 C rhoKP1
77 C buoyK, buoyKM1 - Buoyancy at current level and level above.
78 C phiHyd - Hydrostatic part of the potential phiHydi.
79 C In z coords phiHydiHyd is the hydrostatic
80 C pressure anomaly
81 C In p coords phiHydiHyd is the geopotential
82 C surface height
83 C anomaly.
84 C etaSurfX, - Holds surface elevation gradient in X and Y.
85 C etaSurfY
86 C KappaRT, - Total diffusion in vertical for T and S.
87 C KappaRS (background + spatially varying, isopycnal term).
88 C iMin, iMax - Ranges and sub-block indices on which calculations
89 C jMin, jMax are applied.
90 C bi, bj
91 C k, kup, - Index for layer above and below. kup and kDown
92 C kDown, km1 are switched with layer to be the appropriate
93 C index into fVerTerm.
94 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100 _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
110 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
111 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
112 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
113 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
114 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115 _RL rhokp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 _RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 _RL rhotmp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
123 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
124 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
125 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
126 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
127 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
128 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
129
130 C This is currently also used by IVDC and Diagnostics
131 C #ifdef INCLUDE_CONVECT_CALL
132 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
133 C #endif
134
135 INTEGER iMin, iMax
136 INTEGER jMin, jMax
137 INTEGER bi, bj
138 INTEGER i, j
139 INTEGER k, km1, kup, kDown
140 LOGICAL BOTTOM_LAYER
141
142 #ifdef ALLOW_AUTODIFF_TAMC
143 INTEGER isbyte
144 PARAMETER( isbyte = 4 )
145
146 INTEGER act1, act2, act3, act4
147 INTEGER max1, max2, max3
148 INTEGER iikey, kkey
149 INTEGER maximpl
150 #endif /* ALLOW_AUTODIFF_TAMC */
151
152 C--- The algorithm...
153 C
154 C "Correction Step"
155 C =================
156 C Here we update the horizontal velocities with the surface
157 C pressure such that the resulting flow is either consistent
158 C with the free-surface evolution or the rigid-lid:
159 C U[n] = U* + dt x d/dx P
160 C V[n] = V* + dt x d/dy P
161 C
162 C "Calculation of Gs"
163 C ===================
164 C This is where all the accelerations and tendencies (ie.
165 C physics, parameterizations etc...) are calculated
166 C rVel = sum_r ( div. u[n] )
167 C rho = rho ( theta[n], salt[n] )
168 C b = b(rho, theta)
169 C K31 = K31 ( rho )
170 C Gu[n] = Gu( u[n], v[n], rVel, b, ... )
171 C Gv[n] = Gv( u[n], v[n], rVel, b, ... )
172 C Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
173 C Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
174 C
175 C "Time-stepping" or "Prediction"
176 C ================================
177 C The models variables are stepped forward with the appropriate
178 C time-stepping scheme (currently we use Adams-Bashforth II)
179 C - For momentum, the result is always *only* a "prediction"
180 C in that the flow may be divergent and will be "corrected"
181 C later with a surface pressure gradient.
182 C - Normally for tracers the result is the new field at time
183 C level [n+1} *BUT* in the case of implicit diffusion the result
184 C is also *only* a prediction.
185 C - We denote "predictors" with an asterisk (*).
186 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
187 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
188 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
189 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
190 C With implicit diffusion:
191 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
192 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
193 C (1 + dt * K * d_zz) theta[n] = theta*
194 C (1 + dt * K * d_zz) salt[n] = salt*
195 C---
196
197 #ifdef ALLOW_AUTODIFF_TAMC
198 C-- dummy statement to end declaration part
199 ikey = 1
200 #endif /* ALLOW_AUTODIFF_TAMC */
201
202
203 C-- Set up work arrays with valid (i.e. not NaN) values
204 C These inital values do not alter the numerical results. They
205 C just ensure that all memory references are to valid floating
206 C point numbers. This prevents spurious hardware signals due to
207 C uninitialised but inert locations.
208 DO j=1-OLy,sNy+OLy
209 DO i=1-OLx,sNx+OLx
210 xA(i,j) = 0. _d 0
211 yA(i,j) = 0. _d 0
212 uTrans(i,j) = 0. _d 0
213 vTrans(i,j) = 0. _d 0
214 aTerm(i,j) = 0. _d 0
215 xTerm(i,j) = 0. _d 0
216 cTerm(i,j) = 0. _d 0
217 mTerm(i,j) = 0. _d 0
218 pTerm(i,j) = 0. _d 0
219 fZon(i,j) = 0. _d 0
220 fMer(i,j) = 0. _d 0
221 DO k=1,Nr
222 phiHyd (i,j,k) = 0. _d 0
223 KappaRU(i,j,k) = 0. _d 0
224 KappaRV(i,j,k) = 0. _d 0
225 sigmaX(i,j,k) = 0. _d 0
226 sigmaY(i,j,k) = 0. _d 0
227 sigmaR(i,j,k) = 0. _d 0
228 ENDDO
229 rhoKM1 (i,j) = 0. _d 0
230 rhok (i,j) = 0. _d 0
231 rhoKP1 (i,j) = 0. _d 0
232 rhoTMP (i,j) = 0. _d 0
233 buoyKM1(i,j) = 0. _d 0
234 buoyK (i,j) = 0. _d 0
235 maskC (i,j) = 0. _d 0
236 ENDDO
237 ENDDO
238
239
240 #ifdef ALLOW_AUTODIFF_TAMC
241 C-- HPF directive to help TAMC
242 CHPF$ INDEPENDENT
243 #endif /* ALLOW_AUTODIFF_TAMC */
244
245 DO bj=myByLo(myThid),myByHi(myThid)
246
247 #ifdef ALLOW_AUTODIFF_TAMC
248 C-- HPF directive to help TAMC
249 CHPF$ INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
250 CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA
251 CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
252 CHPF$& )
253 #endif /* ALLOW_AUTODIFF_TAMC */
254
255 DO bi=myBxLo(myThid),myBxHi(myThid)
256
257 #ifdef ALLOW_AUTODIFF_TAMC
258 act1 = bi - myBxLo(myThid)
259 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
260
261 act2 = bj - myByLo(myThid)
262 max2 = myByHi(myThid) - myByLo(myThid) + 1
263
264 act3 = myThid - 1
265 max3 = nTx*nTy
266
267 act4 = ikey_dynamics - 1
268
269 ikey = (act1 + 1) + act2*max1
270 & + act3*max1*max2
271 & + act4*max1*max2*max3
272 #endif /* ALLOW_AUTODIFF_TAMC */
273
274 C-- Set up work arrays that need valid initial values
275 DO j=1-OLy,sNy+OLy
276 DO i=1-OLx,sNx+OLx
277 rTrans(i,j) = 0. _d 0
278 rVel (i,j,1) = 0. _d 0
279 rVel (i,j,2) = 0. _d 0
280 fVerT (i,j,1) = 0. _d 0
281 fVerT (i,j,2) = 0. _d 0
282 fVerS (i,j,1) = 0. _d 0
283 fVerS (i,j,2) = 0. _d 0
284 fVerU (i,j,1) = 0. _d 0
285 fVerU (i,j,2) = 0. _d 0
286 fVerV (i,j,1) = 0. _d 0
287 fVerV (i,j,2) = 0. _d 0
288 phiHyd(i,j,1) = 0. _d 0
289 ENDDO
290 ENDDO
291
292 DO k=1,Nr
293 DO j=1-OLy,sNy+OLy
294 DO i=1-OLx,sNx+OLx
295 #ifdef INCLUDE_CONVECT_CALL
296 ConvectCount(i,j,k) = 0.
297 #endif
298 KappaRT(i,j,k) = 0. _d 0
299 KappaRS(i,j,k) = 0. _d 0
300 ENDDO
301 ENDDO
302 ENDDO
303
304 iMin = 1-OLx+1
305 iMax = sNx+OLx
306 jMin = 1-OLy+1
307 jMax = sNy+OLy
308
309 k = 1
310 BOTTOM_LAYER = k .EQ. Nr
311
312 #ifdef DO_PIPELINED_CORRECTION_STEP
313 C-- Calculate gradient of surface pressure
314 CALL CALC_GRAD_ETA_SURF(
315 I bi,bj,iMin,iMax,jMin,jMax,
316 O etaSurfX,etaSurfY,
317 I myThid)
318 C-- Update fields in top level according to tendency terms
319 CALL CORRECTION_STEP(
320 I bi,bj,iMin,iMax,jMin,jMax,k,
321 I etaSurfX,etaSurfY,myTime,myThid)
322
323 #ifdef ALLOW_OBCS
324 IF (openBoundaries) THEN
325 #ifdef ALLOW_AUTODIFF_TAMC
326 CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
327 CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
328 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
329 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
330 #endif /* ALLOW_AUTODIFF_TAMC */
331 CALL APPLY_OBCS1( bi, bj, k, myThid )
332 END IF
333 #endif
334
335 IF ( .NOT. BOTTOM_LAYER ) THEN
336 C-- Update fields in layer below according to tendency terms
337 CALL CORRECTION_STEP(
338 I bi,bj,iMin,iMax,jMin,jMax,k+1,
339 I etaSurfX,etaSurfY,myTime,myThid)
340 #ifdef ALLOW_OBCS
341 IF (openBoundaries) THEN
342 #ifdef ALLOW_AUTODIFF_TAMC
343 CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
344 CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
345 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
346 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
347 #endif /* ALLOW_AUTODIFF_TAMC */
348 CALL APPLY_OBCS1( bi, bj, k+1, myThid )
349 END IF
350 #endif
351 ENDIF
352 #endif
353
354 C-- Density of 1st level (below W(1)) reference to level 1
355 #ifdef INCLUDE_FIND_RHO_CALL
356 #ifdef ALLOW_AUTODIFF_TAMC
357 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
358 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
359 #endif /* ALLOW_AUTODIFF_TAMC */
360 CALL FIND_RHO(
361 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
362 O rhoKm1,
363 I myThid )
364 #endif
365
366 IF (.NOT. BOTTOM_LAYER) THEN
367
368 C-- Check static stability with layer below
369 C-- and mix as needed.
370 #ifdef INCLUDE_FIND_RHO_CALL
371 #ifdef ALLOW_AUTODIFF_TAMC
372 CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj
373 CADJ & , key = ikey, byte = isbyte
374 CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj
375 CADJ & , key = ikey, byte = isbyte
376 #endif /* ALLOW_AUTODIFF_TAMC */
377 CALL FIND_RHO(
378 I bi, bj, iMin, iMax, jMin, jMax, k+1, k, eosType,
379 O rhoKp1,
380 I myThid )
381 #endif
382
383 #ifdef ALLOW_AUTODIFF_TAMC
384 CADJ STORE rhoKm1(:,:) = comlev1_bibj, key = ikey, byte = isbyte
385 CADJ STORE rhoKp1(:,:) = comlev1_bibj, key = ikey, byte = isbyte
386 #endif /* ALLOW_AUTODIFF_TAMC */
387
388 #ifdef INCLUDE_CONVECT_CALL
389
390 CALL CONVECT(
391 I bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,
392 U ConvectCount,
393 I myTime,myIter,myThid)
394
395 #ifdef ALLOW_AUTODIFF_TAMC
396 CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)
397 CADJ & = comlev1_bibj, key = ikey, byte = isbyte
398 CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)
399 CADJ & = comlev1_bibj, key = ikey, byte = isbyte
400 #endif /* ALLOW_AUTODIFF_TAMC */
401
402 #endif
403
404 C-- Implicit Vertical Diffusion for Convection
405 IF (ivdc_kappa.NE.0.) THEN
406 CALL CALC_IVDC(
407 I bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,
408 U ConvectCount, KappaRT, KappaRS,
409 I myTime,myIter,myThid)
410 ENDIF
411
412 C-- Recompute density after mixing
413 #ifdef INCLUDE_FIND_RHO_CALL
414 CALL FIND_RHO(
415 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
416 O rhoKm1,
417 I myThid )
418 #endif
419 ENDIF
420
421 C-- Calculate buoyancy
422 CALL CALC_BUOYANCY(
423 I bi,bj,iMin,iMax,jMin,jMax,k,rhoKm1,
424 O buoyKm1,
425 I myThid )
426
427 C-- Integrate hydrostatic balance for phiHyd with BC of
428 C-- phiHyd(z=0)=0
429 CALL CALC_PHI_HYD(
430 I bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyKm1,
431 U phiHyd,
432 I myThid )
433
434 #ifdef ALLOW_GMREDI
435 IF ( useGMRedi ) THEN
436 CALL GRAD_SIGMA(
437 I bi, bj, iMin, iMax, jMin, jMax, k,
438 I rhoKm1, rhoKm1, rhoKm1,
439 O sigmaX, sigmaY, sigmaR,
440 I myThid )
441 ELSE
442 DO j=1-OLy,sNy+OLy
443 DO i=1-OLx,sNx+OLx
444 sigmaX(i,j,k) = 0. _d 0
445 sigmaY(i,j,k) = 0. _d 0
446 sigmaR(i,j,k) = 0. _d 0
447 ENDDO
448 ENDDO
449 ENDIF
450 #endif
451
452 C-- Start of downward loop
453 DO k=2,Nr
454
455 #ifdef ALLOW_AUTODIFF_TAMC
456 kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1
457 #endif /* ALLOW_AUTODIFF_TAMC */
458
459 BOTTOM_LAYER = k .EQ. Nr
460
461 #ifdef DO_PIPELINED_CORRECTION_STEP
462 IF ( .NOT. BOTTOM_LAYER ) THEN
463 C-- Update fields in layer below according to tendency terms
464 CALL CORRECTION_STEP(
465 I bi,bj,iMin,iMax,jMin,jMax,k+1,
466 I etaSurfX,etaSurfY,myTime,myThid)
467 #ifdef ALLOW_OBCS
468 IF (openBoundaries) THEN
469 #ifdef ALLOW_AUTODIFF_TAMC
470 CADJ STORE uvel (:,:,k,bi,bj) = comlev1_bibj_k
471 CADJ & , key = kkey, byte = isbyte
472 CADJ STORE vvel (:,:,k,bi,bj) = comlev1_bibj_k
473 CADJ & , key = kkey, byte = isbyte
474 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k
475 CADJ & , key = kkey, byte = isbyte
476 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k
477 CADJ & , key = kkey, byte = isbyte
478 #endif /* ALLOW_AUTODIFF_TAMC */
479 CALL APPLY_OBCS1( bi, bj, k+1, myThid )
480 END IF
481 #endif
482 ENDIF
483 #endif /* DO_PIPELINED_CORRECTION_STEP */
484
485 C-- Density of k level (below W(k)) reference to k level
486 #ifdef INCLUDE_FIND_RHO_CALL
487 #ifdef ALLOW_AUTODIFF_TAMC
488 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k
489 CADJ & , key = kkey, byte = isbyte
490 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k
491 CADJ & , key = kkey, byte = isbyte
492 #endif /* ALLOW_AUTODIFF_TAMC */
493 CALL FIND_RHO(
494 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
495 O rhoK,
496 I myThid )
497
498 #ifdef ALLOW_AUTODIFF_TAMC
499 cph( storing not necessary
500 cphCADJ STORE rhoK(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
501 cph)
502 #endif /* ALLOW_AUTODIFF_TAMC */
503 #endif
504
505 IF (.NOT. BOTTOM_LAYER) THEN
506
507 C-- Check static stability with layer below and mix as needed.
508 C-- Density of k+1 level (below W(k+1)) reference to k level.
509 #ifdef INCLUDE_FIND_RHO_CALL
510 #ifdef ALLOW_AUTODIFF_TAMC
511 CADJ STORE theta(:,:,k+1,bi,bj) = comlev1_bibj_k
512 CADJ & , key = kkey, byte = isbyte
513 CADJ STORE salt (:,:,k+1,bi,bj) = comlev1_bibj_k
514 CADJ & , key = kkey, byte = isbyte
515 #endif /* ALLOW_AUTODIFF_TAMC */
516 CALL FIND_RHO(
517 I bi, bj, iMin, iMax, jMin, jMax, k+1, k, eosType,
518 O rhoKp1,
519 I myThid )
520 #ifdef ALLOW_AUTODIFF_TAMC
521 CADJ STORE rhoKp1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
522 #endif /* ALLOW_AUTODIFF_TAMC */
523 #endif
524
525 #ifdef INCLUDE_CONVECT_CALL
526 CALL CONVECT(
527 I bi,bj,iMin,iMax,jMin,jMax,k+1,rhoK,rhoKp1,
528 U ConvectCount,
529 I myTime,myIter,myThid)
530
531 #endif
532
533 C-- Implicit Vertical Diffusion for Convection
534 IF (ivdc_kappa.NE.0.) THEN
535 #ifdef ALLOW_AUTODIFF_TAMC
536 CADJ STORE rhoKm1(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
537 #endif /* ALLOW_AUTODIFF_TAMC */
538 CALL CALC_IVDC(
539 I bi,bj,iMin,iMax,jMin,jMax,k+1,rhoKm1,rhoKp1,
540 U ConvectCount, KappaRT, KappaRS,
541 I myTime,myIter,myThid)
542 END IF
543
544 C-- Recompute density after mixing
545 #ifdef INCLUDE_FIND_RHO_CALL
546 #ifdef ALLOW_AUTODIFF_TAMC
547 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k
548 CADJ & , key = kkey, byte = isbyte
549 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k
550 CADJ & , key = kkey, byte = isbyte
551 #endif /* ALLOW_AUTODIFF_TAMC */
552 CALL FIND_RHO(
553 I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,
554 O rhoK,
555 I myThid )
556 #endif
557
558 C-- IF (.NOT. BOTTOM_LAYER) ends here
559 ENDIF
560
561 C-- Calculate buoyancy
562 CALL CALC_BUOYANCY(
563 I bi,bj,iMin,iMax,jMin,jMax,k,rhoK,
564 O buoyK,
565 I myThid )
566
567 C-- Integrate hydrostatic balance for phiHyd with BC of
568 C-- phiHyd(z=0)=0
569 CALL CALC_PHI_HYD(
570 I bi,bj,iMin,iMax,jMin,jMax,k,buoyKm1,buoyK,
571 U phiHyd,
572 I myThid )
573
574 #ifdef INCLUDE_FIND_RHO_CALL
575 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
576
577 #ifdef ALLOW_AUTODIFF_TAMC
578 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k
579 CADJ & , key = kkey, byte = isbyte
580 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k
581 CADJ & , key = kkey, byte = isbyte
582 #endif /* ALLOW_AUTODIFF_TAMC */
583
584 CALL FIND_RHO(
585 I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,
586 O rhoTmp,
587 I myThid )
588 #endif
589
590
591 #ifdef ALLOW_GMREDI
592 IF ( useGMRedi ) THEN
593 CALL GRAD_SIGMA(
594 I bi, bj, iMin, iMax, jMin, jMax, k,
595 I rhoK, rhotmp, rhoK,
596 O sigmaX, sigmaY, sigmaR,
597 I myThid )
598 ELSE
599 DO j=1-OLy,sNy+OLy
600 DO i=1-OLx,sNx+OLx
601 sigmaX(i,j,k) = 0. _d 0
602 sigmaY(i,j,k) = 0. _d 0
603 sigmaR(i,j,k) = 0. _d 0
604 ENDDO
605 ENDDO
606 ENDIF
607 #endif
608
609 DO J=jMin,jMax
610 DO I=iMin,iMax
611 #ifdef INCLUDE_FIND_RHO_CALL
612 rhoKm1 (I,J) = rhoK(I,J)
613 #endif
614 buoyKm1(I,J) = buoyK(I,J)
615 ENDDO
616 ENDDO
617
618 C-- end of k loop
619 ENDDO
620
621 C Determines forcing terms based on external fields
622 C relaxation terms, etc.
623 CALL EXTERNAL_FORCING_SURF(
624 I bi, bj, iMin, iMax, jMin, jMax,
625 I myThid )
626
627 #ifdef ALLOW_AUTODIFF_TAMC
628
629 CADJ STORE surfacetendencyu(:,:,bi,bj)
630 CADJ & , surfacetendencyv(:,:,bi,bj)
631 CADJ & , surfacetendencys(:,:,bi,bj)
632 CADJ & , surfacetendencyt(:,:,bi,bj)
633 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
634
635 # ifdef ALLOW_GMREDI
636 CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte
637 CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte
638 CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte
639 # endif /* ALLOW_GMREDI */
640
641 #endif /* ALLOW_AUTODIFF_TAMC */
642
643 #ifdef ALLOW_GMREDI
644 IF (useGMRedi) THEN
645 DO k=1, Nr
646 CALL GMREDI_CALC_TENSOR(
647 I bi, bj, iMin, iMax, jMin, jMax, k,
648 I sigmaX, sigmaY, sigmaR,
649 I myThid )
650 ENDDO
651 #ifdef ALLOW_AUTODIFF_TAMC
652 ELSE
653 DO k=1, Nr
654 CALL GMREDI_CALC_TENSOR_DUMMY(
655 I bi, bj, iMin, iMax, jMin, jMax, k,
656 I sigmaX, sigmaY, sigmaR,
657 I myThid )
658 ENDDO
659 #endif /* ALLOW_AUTODIFF_TAMC */
660 ENDIF
661 #endif
662
663 #ifdef ALLOW_AUTODIFF_TAMC
664 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key=ikey, byte=isbyte
665 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key=ikey, byte=isbyte
666
667 C-- R.G. We need to define a new tape since Kw use mythid instead of bi,bj
668 CADJ STORE Kwx(:,:,:,myThid) = comlev1_bibj, key=ikey, byte=isbyte
669 CADJ STORE Kwy(:,:,:,myThid) = comlev1_bibj, key=ikey, byte=isbyte
670 CADJ STORE Kwz(:,:,:,myThid) = comlev1_bibj, key=ikey, byte=isbyte
671
672 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
673 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
674 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
675 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
676
677 C-- dummy initialization to break data flow because
678 C-- calc_div_ghat has a condition for initialization
679 DO J=jMin,jMax
680 DO I=iMin,iMax
681 cg2d_b(i,j,bi,bj) = 0.0
682 ENDDO
683 ENDDO
684 #endif /* ALLOW_AUTODIFF_TAMC */
685
686 #ifdef ALLOW_KPP
687 C-- Compute KPP mixing coefficients
688 IF (useKPP) THEN
689
690 CALL TIMER_START('KPP_CALC [DYNAMICS]', myThid)
691 CALL KPP_CALC(
692 I bi, bj, myTime, myThid )
693 CALL TIMER_STOP ('KPP_CALC [DYNAMICS]', myThid)
694
695 #ifdef ALLOW_AUTODIFF_TAMC
696 ELSE
697 DO j=1-OLy,sNy+OLy
698 DO i=1-OLx,sNx+OLx
699 KPPhbl (i,j,bi,bj) = 1.0
700 KPPfrac(i,j,bi,bj) = 0.0
701 DO k = 1,Nr
702 KPPghat (i,j,k,bi,bj) = 0.0
703 KPPviscAz (i,j,k,bi,bj) = viscAz
704 KPPdiffKzT(i,j,k,bi,bj) = diffKzT
705 KPPdiffKzS(i,j,k,bi,bj) = diffKzS
706 ENDDO
707 ENDDO
708 ENDDO
709 #endif /* ALLOW_AUTODIFF_TAMC */
710 ENDIF
711
712 #ifdef ALLOW_AUTODIFF_TAMC
713 CADJ STORE KPPghat (:,:,:,bi,bj)
714 CADJ & , KPPviscAz (:,:,:,bi,bj)
715 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
716 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
717 CADJ & , KPPfrac (:,: ,bi,bj)
718 CADJ & = comlev1_bibj, key=ikey, byte=isbyte
719 #endif /* ALLOW_AUTODIFF_TAMC */
720
721 #endif /* ALLOW_KPP */
722
723 C-- Start of upward loop
724 DO k = Nr, 1, -1
725
726 C-- km1 Points to level above k (=k-1)
727 C-- kup Cycles through 1,2 to point to layer above
728 C-- kDown Cycles through 2,1 to point to current layer
729
730 km1 =max(1,k-1)
731 kup =1+MOD(k+1,2)
732 kDown=1+MOD(k,2)
733
734 iMin = 1-OLx+2
735 iMax = sNx+OLx-1
736 jMin = 1-OLy+2
737 jMax = sNy+OLy-1
738
739 #ifdef ALLOW_AUTODIFF_TAMC
740 kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1
741
742 CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte
743 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
744 CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
745 CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
746 #endif /* ALLOW_AUTODIFF_TAMC */
747
748 C-- Get temporary terms used by tendency routines
749 CALL CALC_COMMON_FACTORS (
750 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
751 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
752 I myThid)
753
754 #ifdef ALLOW_OBCS
755 IF (openBoundaries) THEN
756 CALL APPLY_OBCS3( bi, bj, k, kup, rTrans, rVel, myThid )
757 ENDIF
758 #endif
759
760 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
761 C-- Calculate the total vertical diffusivity
762 CALL CALC_DIFFUSIVITY(
763 I bi,bj,iMin,iMax,jMin,jMax,k,
764 I maskC,maskUp,
765 O KappaRT,KappaRS,KappaRU,KappaRV,
766 I myThid)
767 #endif
768 C-- Calculate accelerations in the momentum equations
769 IF ( momStepping ) THEN
770 CALL CALC_MOM_RHS(
771 I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown,
772 I xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
773 I phiHyd,KappaRU,KappaRV,
774 U aTerm,xTerm,cTerm,mTerm,pTerm,
775 U fZon, fMer, fVerU, fVerV,
776 I myTime, myThid)
777 #ifdef ALLOW_AUTODIFF_TAMC
778 #ifdef INCLUDE_CD_CODE
779 ELSE
780 DO j=1-OLy,sNy+OLy
781 DO i=1-OLx,sNx+OLx
782 guCD(i,j,k,bi,bj) = 0.0
783 gvCD(i,j,k,bi,bj) = 0.0
784 END DO
785 END DO
786 #endif
787 #endif /* ALLOW_AUTODIFF_TAMC */
788 ENDIF
789 C-- Calculate active tracer tendencies
790 IF ( tempStepping ) THEN
791 CALL CALC_GT(
792 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
793 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
794 I KappaRT,
795 U aTerm,xTerm,fZon,fMer,fVerT,
796 I myTime, myThid)
797 ENDIF
798 IF ( saltStepping ) THEN
799 CALL CALC_GS(
800 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
801 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
802 I KappaRS,
803 U aTerm,xTerm,fZon,fMer,fVerS,
804 I myTime, myThid)
805 ENDIF
806 #ifdef ALLOW_OBCS
807 C-- Calculate future values on open boundaries
808 IF (openBoundaries) THEN
809 Caja CALL CYCLE_OBCS( k, bi, bj, myThid )
810 CALL SET_OBCS( k, bi, bj, myTime+deltaTclock, myThid )
811 ENDIF
812 #endif
813 C-- Prediction step (step forward all model variables)
814 CALL TIMESTEP(
815 I bi,bj,iMin,iMax,jMin,jMax,k,
816 I myIter, myThid)
817 #ifdef ALLOW_OBCS
818 C-- Apply open boundary conditions
819 IF (openBoundaries) THEN
820 #ifdef ALLOW_AUTODIFF_TAMC
821 CADJ STORE gunm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
822 CADJ STORE gvnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
823 CADJ STORE gwnm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
824 #endif /* ALLOW_AUTODIFF_TAMC */
825
826 CALL APPLY_OBCS2( bi, bj, k, myThid )
827 END IF
828 #endif
829 C-- Freeze water
830 IF (allowFreezing) THEN
831 #ifdef ALLOW_AUTODIFF_TAMC
832 CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
833 #endif /* ALLOW_AUTODIFF_TAMC */
834 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
835 END IF
836
837 #ifdef DIVG_IN_DYNAMICS
838 C-- Diagnose barotropic divergence of predicted fields
839 CALL CALC_DIV_GHAT(
840 I bi,bj,iMin,iMax,jMin,jMax,k,
841 I xA,yA,
842 I myThid)
843 #endif /* DIVG_IN_DYNAMICS */
844
845 C-- Cumulative diagnostic calculations (ie. time-averaging)
846 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
847 IF (taveFreq.GT.0.) THEN
848 CALL DO_TIME_AVERAGES(
849 I myTime, myIter, bi, bj, k, kup, kDown,
850 I rVel, ConvectCount,
851 I myThid )
852 ENDIF
853 #endif
854
855
856 C-- k loop
857 ENDDO
858
859 #ifdef ALLOW_AUTODIFF_TAMC
860 maximpl = 6
861 iikey = (ikey-1)*maximpl
862 #endif /* ALLOW_AUTODIFF_TAMC */
863
864 C-- Implicit diffusion
865 IF (implicitDiffusion) THEN
866
867 IF (tempStepping) THEN
868 #ifdef ALLOW_AUTODIFF_TAMC
869 idkey = iikey + 1
870 CADJ STORE gTNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
871 #endif /* ALLOW_AUTODIFF_TAMC */
872 CALL IMPLDIFF(
873 I bi, bj, iMin, iMax, jMin, jMax,
874 I deltaTtracer, KappaRT,recip_HFacC,
875 U gTNm1,
876 I myThid )
877 END IF
878
879 IF (saltStepping) THEN
880 #ifdef ALLOW_AUTODIFF_TAMC
881 idkey = iikey + 2
882 CADJ STORE gSNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
883 #endif /* ALLOW_AUTODIFF_TAMC */
884 CALL IMPLDIFF(
885 I bi, bj, iMin, iMax, jMin, jMax,
886 I deltaTtracer, KappaRS,recip_HFacC,
887 U gSNm1,
888 I myThid )
889 END IF
890
891 C-- implicitDiffusion
892 ENDIF
893
894 C-- Implicit viscosity
895 IF (implicitViscosity) THEN
896
897 IF (momStepping) THEN
898 #ifdef ALLOW_AUTODIFF_TAMC
899 idkey = iikey + 3
900 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
901 #endif /* ALLOW_AUTODIFF_TAMC */
902 CALL IMPLDIFF(
903 I bi, bj, iMin, iMax, jMin, jMax,
904 I deltaTmom, KappaRU,recip_HFacW,
905 U gUNm1,
906 I myThid )
907 #ifdef ALLOW_AUTODIFF_TAMC
908 idkey = iikey + 4
909 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
910 #endif /* ALLOW_AUTODIFF_TAMC */
911 CALL IMPLDIFF(
912 I bi, bj, iMin, iMax, jMin, jMax,
913 I deltaTmom, KappaRV,recip_HFacS,
914 U gVNm1,
915 I myThid )
916
917 #ifdef INCLUDE_CD_CODE
918
919 #ifdef ALLOW_AUTODIFF_TAMC
920 idkey = iikey + 5
921 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
922 #endif /* ALLOW_AUTODIFF_TAMC */
923 CALL IMPLDIFF(
924 I bi, bj, iMin, iMax, jMin, jMax,
925 I deltaTmom, KappaRU,recip_HFacW,
926 U vVelD,
927 I myThid )
928 #ifdef ALLOW_AUTODIFF_TAMC
929 idkey = iikey + 6
930 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
931 #endif /* ALLOW_AUTODIFF_TAMC */
932 CALL IMPLDIFF(
933 I bi, bj, iMin, iMax, jMin, jMax,
934 I deltaTmom, KappaRV,recip_HFacS,
935 U uVelD,
936 I myThid )
937
938 #endif
939
940 C-- momStepping
941 ENDIF
942
943 C-- implicitViscosity
944 ENDIF
945
946 ENDDO
947 ENDDO
948
949 RETURN
950 END

  ViewVC Help
Powered by ViewVC 1.1.22