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

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

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


Revision 1.41 - (show annotations) (download)
Mon Jun 23 22:32:02 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51, checkpoint51b_pre, checkpoint50i_post, checkpoint51a_post
Changes since 1.40: +15 -21 lines
Preparing next differentiable checkpoint and sync
of MAIN vs. ecco-branch
(updating store after changes in checkpoint50b_post,
plus still messing around with init. sequence).

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.40 2003/05/13 17:42:00 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #ifdef ALLOW_AUTODIFF_TAMC
6 # ifdef ALLOW_GMREDI
7 # include "GMREDI_OPTIONS.h"
8 # endif
9 # ifdef ALLOW_KPP
10 # include "KPP_OPTIONS.h"
11 # endif
12 cswdice --- add ----
13 #ifdef ALLOW_THERM_SEAICE
14 #include "ICE.h"
15 #endif
16 cswdice ------
17 #endif /* ALLOW_AUTODIFF_TAMC */
18
19 CBOP
20 C !ROUTINE: THERMODYNAMICS
21 C !INTERFACE:
22 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
23 C !DESCRIPTION: \bv
24 C *==========================================================*
25 C | SUBROUTINE THERMODYNAMICS
26 C | o Controlling routine for the prognostic part of the
27 C | thermo-dynamics.
28 C *===========================================================
29 C | The algorithm...
30 C |
31 C | "Correction Step"
32 C | =================
33 C | Here we update the horizontal velocities with the surface
34 C | pressure such that the resulting flow is either consistent
35 C | with the free-surface evolution or the rigid-lid:
36 C | U[n] = U* + dt x d/dx P
37 C | V[n] = V* + dt x d/dy P
38 C |
39 C | "Calculation of Gs"
40 C | ===================
41 C | This is where all the accelerations and tendencies (ie.
42 C | physics, parameterizations etc...) are calculated
43 C | rho = rho ( theta[n], salt[n] )
44 C | b = b(rho, theta)
45 C | K31 = K31 ( rho )
46 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
47 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
48 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
49 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
50 C |
51 C | "Time-stepping" or "Prediction"
52 C | ================================
53 C | The models variables are stepped forward with the appropriate
54 C | time-stepping scheme (currently we use Adams-Bashforth II)
55 C | - For momentum, the result is always *only* a "prediction"
56 C | in that the flow may be divergent and will be "corrected"
57 C | later with a surface pressure gradient.
58 C | - Normally for tracers the result is the new field at time
59 C | level [n+1} *BUT* in the case of implicit diffusion the result
60 C | is also *only* a prediction.
61 C | - We denote "predictors" with an asterisk (*).
62 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
63 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
64 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66 C | With implicit diffusion:
67 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
68 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69 C | (1 + dt * K * d_zz) theta[n] = theta*
70 C | (1 + dt * K * d_zz) salt[n] = salt*
71 C |
72 C *==========================================================*
73 C \ev
74
75 C !USES:
76 IMPLICIT NONE
77 C == Global variables ===
78 #include "SIZE.h"
79 #include "EEPARAMS.h"
80 #include "PARAMS.h"
81 #include "DYNVARS.h"
82 #include "GRID.h"
83 #include "GAD.h"
84 #ifdef ALLOW_PASSIVE_TRACER
85 #include "TR1.h"
86 #endif
87 #ifdef ALLOW_AUTODIFF_TAMC
88 # include "tamc.h"
89 # include "tamc_keys.h"
90 # include "FFIELDS.h"
91 # include "EOS.h"
92 # ifdef ALLOW_KPP
93 # include "KPP.h"
94 # endif
95 # ifdef ALLOW_GMREDI
96 # include "GMREDI.h"
97 # endif
98 #endif /* ALLOW_AUTODIFF_TAMC */
99 #ifdef ALLOW_TIMEAVE
100 #include "TIMEAVE_STATV.h"
101 #endif
102
103 C !INPUT/OUTPUT PARAMETERS:
104 C == Routine arguments ==
105 C myTime - Current time in simulation
106 C myIter - Current iteration number in simulation
107 C myThid - Thread number for this instance of the routine.
108 _RL myTime
109 INTEGER myIter
110 INTEGER myThid
111
112 C !LOCAL VARIABLES:
113 C == Local variables
114 C xA, yA - Per block temporaries holding face areas
115 C uTrans, vTrans, rTrans - Per block temporaries holding flow
116 C transport
117 C o uTrans: Zonal transport
118 C o vTrans: Meridional transport
119 C o rTrans: Vertical transport
120 C maskUp o maskUp: land/water mask for W points
121 C fVer[STUV] o fVer: Vertical flux term - note fVer
122 C is "pipelined" in the vertical
123 C so we need an fVer for each
124 C variable.
125 C rhoK, rhoKM1 - Density at current level, and level above
126 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
127 C phiSurfY or geopotentiel (atmos) in X and Y direction
128 C KappaRT, - Total diffusion in vertical for T and S.
129 C KappaRS (background + spatially varying, isopycnal term).
130 C useVariableK = T when vertical diffusion is not constant
131 C iMin, iMax - Ranges and sub-block indices on which calculations
132 C jMin, jMax are applied.
133 C bi, bj
134 C k, kup, - Index for layer above and below. kup and kDown
135 C kDown, km1 are switched with layer to be the appropriate
136 C index into fVerTerm.
137 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
141 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
142 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
143 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
144 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
145 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
151 _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
152 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
153 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
154 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
155 C This is currently used by IVDC and Diagnostics
156 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
157 LOGICAL useVariableK
158 INTEGER iMin, iMax
159 INTEGER jMin, jMax
160 INTEGER bi, bj
161 INTEGER i, j
162 INTEGER k, km1, kup, kDown
163
164 CEOP
165
166 #ifndef DISABLE_DEBUGMODE
167 IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
168 #endif
169
170 #ifdef ALLOW_AUTODIFF_TAMC
171 C-- dummy statement to end declaration part
172 ikey = 1
173 itdkey = 1
174 #endif /* ALLOW_AUTODIFF_TAMC */
175
176 #ifdef ALLOW_AUTODIFF_TAMC
177 C-- HPF directive to help TAMC
178 CHPF$ INDEPENDENT
179 #endif /* ALLOW_AUTODIFF_TAMC */
180
181 DO bj=myByLo(myThid),myByHi(myThid)
182
183 #ifdef ALLOW_AUTODIFF_TAMC
184 C-- HPF directive to help TAMC
185 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
186 CHPF$& ,utrans,vtrans,xA,yA
187 CHPF$& ,KappaRT,KappaRS
188 CHPF$& )
189 #endif /* ALLOW_AUTODIFF_TAMC */
190
191 DO bi=myBxLo(myThid),myBxHi(myThid)
192
193 #ifdef ALLOW_AUTODIFF_TAMC
194 act1 = bi - myBxLo(myThid)
195 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
196 act2 = bj - myByLo(myThid)
197 max2 = myByHi(myThid) - myByLo(myThid) + 1
198 act3 = myThid - 1
199 max3 = nTx*nTy
200 act4 = ikey_dynamics - 1
201 itdkey = (act1 + 1) + act2*max1
202 & + act3*max1*max2
203 & + act4*max1*max2*max3
204 #endif /* ALLOW_AUTODIFF_TAMC */
205
206 C-- Set up work arrays with valid (i.e. not NaN) values
207 C These inital values do not alter the numerical results. They
208 C just ensure that all memory references are to valid floating
209 C point numbers. This prevents spurious hardware signals due to
210 C uninitialised but inert locations.
211
212 DO j=1-OLy,sNy+OLy
213 DO i=1-OLx,sNx+OLx
214 xA(i,j) = 0. _d 0
215 yA(i,j) = 0. _d 0
216 uTrans(i,j) = 0. _d 0
217 vTrans(i,j) = 0. _d 0
218 rhok (i,j) = 0. _d 0
219 rhoKM1 (i,j) = 0. _d 0
220 phiSurfX(i,j) = 0. _d 0
221 phiSurfY(i,j) = 0. _d 0
222 rTrans (i,j) = 0. _d 0
223 fVerT (i,j,1) = 0. _d 0
224 fVerT (i,j,2) = 0. _d 0
225 fVerS (i,j,1) = 0. _d 0
226 fVerS (i,j,2) = 0. _d 0
227 fVerTr1(i,j,1) = 0. _d 0
228 fVerTr1(i,j,2) = 0. _d 0
229 ENDDO
230 ENDDO
231
232 DO k=1,Nr
233 DO j=1-OLy,sNy+OLy
234 DO i=1-OLx,sNx+OLx
235 C This is currently also used by IVDC and Diagnostics
236 sigmaX(i,j,k) = 0. _d 0
237 sigmaY(i,j,k) = 0. _d 0
238 sigmaR(i,j,k) = 0. _d 0
239 ConvectCount(i,j,k) = 0.
240 KappaRT(i,j,k) = 0. _d 0
241 KappaRS(i,j,k) = 0. _d 0
242 #ifdef ALLOW_AUTODIFF_TAMC
243 cph all the following init. are necessary for TAF
244 cph although some of these are re-initialised later.
245 gT(i,j,k,bi,bj) = 0. _d 0
246 gS(i,j,k,bi,bj) = 0. _d 0
247 # ifdef ALLOW_PASSIVE_TRACER
248 gTr1(i,j,k,bi,bj) = 0. _d 0
249 # endif
250 # ifdef ALLOW_GMREDI
251 Kwx(i,j,k,bi,bj) = 0. _d 0
252 Kwy(i,j,k,bi,bj) = 0. _d 0
253 Kwz(i,j,k,bi,bj) = 0. _d 0
254 # ifdef GM_NON_UNITY_DIAGONAL
255 Kux(i,j,k,bi,bj) = 0. _d 0
256 Kvy(i,j,k,bi,bj) = 0. _d 0
257 # endif
258 # ifdef GM_EXTRA_DIAGONAL
259 Kuz(i,j,k,bi,bj) = 0. _d 0
260 Kvz(i,j,k,bi,bj) = 0. _d 0
261 # endif
262 # ifdef GM_BOLUS_ADVEC
263 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
264 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
265 # endif
266 # ifdef GM_VISBECK_VARIABLE_K
267 VisbeckK(i,j,bi,bj) = 0. _d 0
268 # endif
269 # endif /* ALLOW_GMREDI */
270 #endif /* ALLOW_AUTODIFF_TAMC */
271 ENDDO
272 ENDDO
273 ENDDO
274
275 iMin = 1-OLx
276 iMax = sNx+OLx
277 jMin = 1-OLy
278 jMax = sNy+OLy
279
280
281 #ifdef ALLOW_AUTODIFF_TAMC
282 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
283 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
284 CADJ STORE totphihyd
285 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
286 #ifdef ALLOW_KPP
287 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
288 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
289 #endif
290 #endif /* ALLOW_AUTODIFF_TAMC */
291
292 #ifndef DISABLE_DEBUGMODE
293 IF (debugMode) CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
294 #endif
295
296 C-- Start of diagnostic loop
297 DO k=Nr,1,-1
298
299 #ifdef ALLOW_AUTODIFF_TAMC
300 C? Patrick, is this formula correct now that we change the loop range?
301 C? Do we still need this?
302 cph kkey formula corrected.
303 cph Needed for rhok, rhokm1, in the case useGMREDI.
304 kkey = (itdkey-1)*Nr + k
305 #endif /* ALLOW_AUTODIFF_TAMC */
306
307 C-- Integrate continuity vertically for vertical velocity
308 c CALL INTEGRATE_FOR_W(
309 c I bi, bj, k, uVel, vVel,
310 c O wVel,
311 c I myThid )
312
313 #ifdef ALLOW_OBCS
314 #ifdef ALLOW_NONHYDROSTATIC
315 C-- Apply OBC to W if in N-H mode
316 c IF (useOBCS.AND.nonHydrostatic) THEN
317 c CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
318 c ENDIF
319 #endif /* ALLOW_NONHYDROSTATIC */
320 #endif /* ALLOW_OBCS */
321
322 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
323 C-- MOST of THERMODYNAMICS will be disabled
324 #ifndef SINGLE_LAYER_MODE
325
326 C-- Calculate gradients of potential density for isoneutral
327 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
328 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
329 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
330 #ifndef DISABLE_DEBUGMODE
331 IF (debugMode) CALL DEBUG_CALL('FIND_RHO',myThid)
332 #endif
333 #ifdef ALLOW_AUTODIFF_TAMC
334 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
335 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
336 #endif /* ALLOW_AUTODIFF_TAMC */
337 CALL FIND_RHO(
338 I bi, bj, iMin, iMax, jMin, jMax, k, k,
339 I theta, salt,
340 O rhoK,
341 I myThid )
342
343 IF (k.GT.1) THEN
344 #ifdef ALLOW_AUTODIFF_TAMC
345 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
346 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
347 #endif /* ALLOW_AUTODIFF_TAMC */
348 CALL FIND_RHO(
349 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
350 I theta, salt,
351 O rhoKm1,
352 I myThid )
353 ENDIF
354 #ifndef DISABLE_DEBUGMODE
355 IF (debugMode) CALL DEBUG_CALL('GRAD_SIGMA',myThid)
356 #endif
357 CALL GRAD_SIGMA(
358 I bi, bj, iMin, iMax, jMin, jMax, k,
359 I rhoK, rhoKm1, rhoK,
360 O sigmaX, sigmaY, sigmaR,
361 I myThid )
362 ENDIF
363
364 #ifdef ALLOW_AUTODIFF_TAMC
365 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
366 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
367 #endif /* ALLOW_AUTODIFF_TAMC */
368 C-- Implicit Vertical Diffusion for Convection
369 c ==> should use sigmaR !!!
370 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
371 #ifndef DISABLE_DEBUGMODE
372 IF (debugMode) CALL DEBUG_CALL('CALC_IVDC',myThid)
373 #endif
374 CALL CALC_IVDC(
375 I bi, bj, iMin, iMax, jMin, jMax, k,
376 I rhoKm1, rhoK,
377 U ConvectCount, KappaRT, KappaRS,
378 I myTime, myIter, myThid)
379 ENDIF
380
381 #endif /* SINGLE_LAYER_MODE */
382
383 C-- end of diagnostic k loop (Nr:1)
384 ENDDO
385
386 #ifdef ALLOW_AUTODIFF_TAMC
387 cph avoids recomputation of integrate_for_w
388 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
389 #endif /* ALLOW_AUTODIFF_TAMC */
390
391 #ifdef ALLOW_OBCS
392 C-- Calculate future values on open boundaries
393 IF (useOBCS) THEN
394 #ifndef DISABLE_DEBUGMODE
395 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC',myThid)
396 #endif
397 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
398 I uVel, vVel, wVel, theta, salt,
399 I myThid )
400 ENDIF
401 #endif /* ALLOW_OBCS */
402
403
404 c********************************************
405 cswdice --- add ---
406 #ifdef ALLOW_THERM_SEAICE
407 #ifndef DISABLE_DEBUGMODE
408 IF (debugMode) CALL DEBUG_CALL('ICE_FORCING',myThid)
409 #endif
410 C-- Determines forcing terms based on external fields
411 c-- including effects from ice
412 CALL ICE_FORCING(
413 I bi, bj, iMin, iMax, jMin, jMax,
414 I myThid )
415 #else
416
417 cswdice --- end add ---
418
419 C-- Determines forcing terms based on external fields
420 C relaxation terms, etc.
421 #ifndef DISABLE_DEBUGMODE
422 IF (debugMode) CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
423 #endif
424 CALL EXTERNAL_FORCING_SURF(
425 I bi, bj, iMin, iMax, jMin, jMax,
426 I myThid )
427 cswdice --- add ----
428 #endif
429 cswdice --- end add ---
430 c******************************************
431
432
433
434
435 #ifdef ALLOW_AUTODIFF_TAMC
436 cph needed for KPP
437 CADJ STORE surfacetendencyU(:,:,bi,bj)
438 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
439 CADJ STORE surfacetendencyV(:,:,bi,bj)
440 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
441 CADJ STORE surfacetendencyS(:,:,bi,bj)
442 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
443 CADJ STORE surfacetendencyT(:,:,bi,bj)
444 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
445 #endif /* ALLOW_AUTODIFF_TAMC */
446
447 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
448 C-- MOST of THERMODYNAMICS will be disabled
449 #ifndef SINGLE_LAYER_MODE
450
451 #ifdef ALLOW_GMREDI
452
453 #ifdef ALLOW_AUTODIFF_TAMC
454 cph storing here is needed only for one GMREDI_OPTIONS:
455 cph define GM_BOLUS_ADVEC
456 cph but I've avoided the #ifdef for now, in case more things change
457 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
458 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
459 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
460 #endif /* ALLOW_AUTODIFF_TAMC */
461
462 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
463 IF (useGMRedi) THEN
464 #ifndef DISABLE_DEBUGMODE
465 IF (debugMode) CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
466 #endif
467 CALL GMREDI_CALC_TENSOR(
468 I bi, bj, iMin, iMax, jMin, jMax,
469 I sigmaX, sigmaY, sigmaR,
470 I myThid )
471 #ifdef ALLOW_AUTODIFF_TAMC
472 ELSE
473 CALL GMREDI_CALC_TENSOR_DUMMY(
474 I bi, bj, iMin, iMax, jMin, jMax,
475 I sigmaX, sigmaY, sigmaR,
476 I myThid )
477 #endif /* ALLOW_AUTODIFF_TAMC */
478 ENDIF
479
480 #ifdef ALLOW_AUTODIFF_TAMC
481 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
482 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
483 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
484 #endif /* ALLOW_AUTODIFF_TAMC */
485
486 #endif /* ALLOW_GMREDI */
487
488 #ifdef ALLOW_KPP
489 C-- Compute KPP mixing coefficients
490 IF (useKPP) THEN
491 #ifndef DISABLE_DEBUGMODE
492 IF (debugMode) CALL DEBUG_CALL('KPP_CALC',myThid)
493 #endif
494 CALL KPP_CALC(
495 I bi, bj, myTime, myThid )
496 #ifdef ALLOW_AUTODIFF_TAMC
497 ELSE
498 CALL KPP_CALC_DUMMY(
499 I bi, bj, myTime, myThid )
500 #endif /* ALLOW_AUTODIFF_TAMC */
501 ENDIF
502
503 #ifdef ALLOW_AUTODIFF_TAMC
504 CADJ STORE KPPghat (:,:,:,bi,bj)
505 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
506 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
507 CADJ & , KPPfrac (:,: ,bi,bj)
508 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
509 #endif /* ALLOW_AUTODIFF_TAMC */
510
511 #endif /* ALLOW_KPP */
512
513 #ifdef ALLOW_AUTODIFF_TAMC
514 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
515 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
516 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
517 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
518 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
519 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
520 #ifdef ALLOW_PASSIVE_TRACER
521 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
522 #endif
523 #endif /* ALLOW_AUTODIFF_TAMC */
524
525 #ifdef ALLOW_AIM
526 C AIM - atmospheric intermediate model, physics package code.
527 IF ( useAIM ) THEN
528 #ifndef DISABLE_DEBUGMODE
529 IF (debugMode) CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
530 #endif
531 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
532 CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
533 CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
534 ENDIF
535 #endif /* ALLOW_AIM */
536
537 #ifndef DISABLE_MULTIDIM_ADVECTION
538 C-- Some advection schemes are better calculated using a multi-dimensional
539 C method in the absence of any other terms and, if used, is done here.
540 C
541 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
542 C The default is to use multi-dimensinal advection for non-linear advection
543 C schemes. However, for the sake of efficiency of the adjoint it is necessary
544 C to be able to exclude this scheme to avoid excessive storage and
545 C recomputation. It *is* differentiable, if you need it.
546 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
547 C disable this section of code.
548 IF (tempMultiDimAdvec) THEN
549 #ifndef DISABLE_DEBUGMODE
550 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
551 #endif
552 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
553 U theta,gT,
554 I myTime,myIter,myThid)
555 ENDIF
556 IF (saltMultiDimAdvec) THEN
557 #ifndef DISABLE_DEBUGMODE
558 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
559 #endif
560 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
561 U salt,gS,
562 I myTime,myIter,myThid)
563 ENDIF
564 C Since passive tracers are configurable separately from T,S we
565 C call the multi-dimensional method for PTRACERS regardless
566 C of whether multiDimAdvection is set or not.
567 #ifdef ALLOW_PTRACERS
568 IF ( usePTRACERS ) THEN
569 #ifndef DISABLE_DEBUGMODE
570 IF (debugMode) CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
571 #endif
572 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
573 ENDIF
574 #endif /* ALLOW_PTRACERS */
575 #endif /* DISABLE_MULTIDIM_ADVECTION */
576
577 #ifndef DISABLE_DEBUGMODE
578 IF (debugMode) CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
579 #endif
580
581 C-- Start of thermodynamics loop
582 DO k=Nr,1,-1
583 #ifdef ALLOW_AUTODIFF_TAMC
584 C? Patrick Is this formula correct?
585 cph Yes, but I rewrote it.
586 cph Also, the KappaR? need the index and subscript k!
587 kkey = (itdkey-1)*Nr + k
588 #endif /* ALLOW_AUTODIFF_TAMC */
589
590 C-- km1 Points to level above k (=k-1)
591 C-- kup Cycles through 1,2 to point to layer above
592 C-- kDown Cycles through 2,1 to point to current layer
593
594 km1 = MAX(1,k-1)
595 kup = 1+MOD(k+1,2)
596 kDown= 1+MOD(k,2)
597
598 iMin = 1-OLx
599 iMax = sNx+OLx
600 jMin = 1-OLy
601 jMax = sNy+OLy
602
603 C-- Get temporary terms used by tendency routines
604 CALL CALC_COMMON_FACTORS (
605 I bi,bj,iMin,iMax,jMin,jMax,k,
606 O xA,yA,uTrans,vTrans,rTrans,maskUp,
607 I myThid)
608
609 #ifdef ALLOW_GMREDI
610
611 C-- Residual transp = Bolus transp + Eulerian transp
612 IF (useGMRedi) THEN
613 CALL GMREDI_CALC_UVFLOW(
614 & uTrans, vTrans, bi, bj, k, myThid)
615 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
616 & rTrans, bi, bj, k, myThid)
617 ENDIF
618
619 #ifdef ALLOW_AUTODIFF_TAMC
620 #ifdef GM_BOLUS_ADVEC
621 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
622 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
623 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
624 #endif
625 #endif /* ALLOW_AUTODIFF_TAMC */
626
627 #endif /* ALLOW_GMREDI */
628
629 #ifdef ALLOW_AUTODIFF_TAMC
630 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
631 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
632 #endif /* ALLOW_AUTODIFF_TAMC */
633
634 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
635 C-- Calculate the total vertical diffusivity
636 CALL CALC_DIFFUSIVITY(
637 I bi,bj,iMin,iMax,jMin,jMax,k,
638 I maskUp,
639 O KappaRT,KappaRS,
640 I myThid)
641 #endif
642
643 iMin = 1-OLx+2
644 iMax = sNx+OLx-1
645 jMin = 1-OLy+2
646 jMax = sNy+OLy-1
647
648 C-- Calculate active tracer tendencies (gT,gS,...)
649 C and step forward storing result in gTnm1, gSnm1, etc.
650 IF ( tempStepping ) THEN
651 CALL CALC_GT(
652 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
653 I xA,yA,uTrans,vTrans,rTrans,maskUp,
654 I KappaRT,
655 U fVerT,
656 I myTime,myIter,myThid)
657 CALL TIMESTEP_TRACER(
658 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
659 I theta, gT,
660 I myIter, myThid)
661 ENDIF
662 cswdice ---- add ---
663 #ifdef ALLOW_THERM_SEAICE
664 if (k.eq.1) then
665 call ICE_FREEZE(bi, bj, iMin, iMax, jMin, jMax, myThid )
666 endif
667 #endif
668 cswdice -- end add ---
669 IF ( saltStepping ) THEN
670 CALL CALC_GS(
671 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
672 I xA,yA,uTrans,vTrans,rTrans,maskUp,
673 I KappaRS,
674 U fVerS,
675 I myTime,myIter,myThid)
676 CALL TIMESTEP_TRACER(
677 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
678 I salt, gS,
679 I myIter, myThid)
680 ENDIF
681 #ifdef ALLOW_PASSIVE_TRACER
682 IF ( tr1Stepping ) THEN
683 CALL CALC_GTR1(
684 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
685 I xA,yA,uTrans,vTrans,rTrans,maskUp,
686 I KappaRT,
687 U fVerTr1,
688 I myTime,myIter,myThid)
689 CALL TIMESTEP_TRACER(
690 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
691 I Tr1, gTr1,
692 I myIter,myThid)
693 ENDIF
694 #endif
695 #ifdef ALLOW_PTRACERS
696 IF ( usePTRACERS ) THEN
697 CALL PTRACERS_INTEGERATE(
698 I bi,bj,k,
699 I xA,yA,uTrans,vTrans,rTrans,maskUp,
700 X KappaRS,
701 I myIter,myTime,myThid)
702 ENDIF
703 #endif /* ALLOW_PTRACERS */
704
705 #ifdef ALLOW_OBCS
706 C-- Apply open boundary conditions
707 IF (useOBCS) THEN
708 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
709 END IF
710 #endif /* ALLOW_OBCS */
711
712 C-- Freeze water
713 IF ( allowFreezing .AND. .NOT. useSEAICE ) THEN
714 #ifdef ALLOW_AUTODIFF_TAMC
715 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
716 CADJ & , key = kkey, byte = isbyte
717 #endif /* ALLOW_AUTODIFF_TAMC */
718 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
719 END IF
720
721 C-- end of thermodynamic k loop (Nr:1)
722 ENDDO
723
724 cswdice -- add ---
725 #ifdef ALLOW_THERM_SEAICE
726 c timeaveraging for ice model values
727 CALL ICE_AVE(bi,bj,iMin,iMax,jMin,jMax,myThid )
728 #endif
729 cswdice --- end add ---
730
731
732
733
734 C-- Implicit diffusion
735 IF (implicitDiffusion) THEN
736
737 IF (tempStepping) THEN
738 #ifdef ALLOW_AUTODIFF_TAMC
739 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
740 #endif /* ALLOW_AUTODIFF_TAMC */
741 CALL IMPLDIFF(
742 I bi, bj, iMin, iMax, jMin, jMax,
743 I deltaTtracer, KappaRT, recip_HFacC,
744 U gT,
745 I myThid )
746 ENDIF
747
748 IF (saltStepping) THEN
749 #ifdef ALLOW_AUTODIFF_TAMC
750 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
751 #endif /* ALLOW_AUTODIFF_TAMC */
752 CALL IMPLDIFF(
753 I bi, bj, iMin, iMax, jMin, jMax,
754 I deltaTtracer, KappaRS, recip_HFacC,
755 U gS,
756 I myThid )
757 ENDIF
758
759 #ifdef ALLOW_PASSIVE_TRACER
760 IF (tr1Stepping) THEN
761 #ifdef ALLOW_AUTODIFF_TAMC
762 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
763 #endif /* ALLOW_AUTODIFF_TAMC */
764 CALL IMPLDIFF(
765 I bi, bj, iMin, iMax, jMin, jMax,
766 I deltaTtracer, KappaRT, recip_HFacC,
767 U gTr1,
768 I myThid )
769 ENDIF
770 #endif
771
772 #ifdef ALLOW_PTRACERS
773 C Vertical diffusion (implicit) for passive tracers
774 IF ( usePTRACERS ) THEN
775 CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
776 ENDIF
777 #endif /* ALLOW_PTRACERS */
778
779 #ifdef ALLOW_OBCS
780 C-- Apply open boundary conditions
781 IF (useOBCS) THEN
782 DO K=1,Nr
783 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
784 ENDDO
785 END IF
786 #endif /* ALLOW_OBCS */
787
788 C-- End If implicitDiffusion
789 ENDIF
790
791 #ifdef ALLOW_TIMEAVE
792 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
793 CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
794 I Nr, deltaTclock, bi, bj, myThid)
795 ENDIF
796 useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
797 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
798 IF (implicitDiffusion) THEN
799 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
800 I Nr, 3, deltaTclock, bi, bj, myThid)
801 ELSE
802 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
803 I Nr, 3, deltaTclock, bi, bj, myThid)
804 ENDIF
805 ENDIF
806 #endif /* ALLOW_TIMEAVE */
807
808 #endif /* SINGLE_LAYER_MODE */
809
810 C-- end bi,bj loops.
811 ENDDO
812 ENDDO
813
814 #ifndef DISABLE_DEBUGMODE
815 If (debugMode) THEN
816 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
817 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
818 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
819 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
820 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
821 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
822 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
823 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
824 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
825 #ifdef ALLOW_PTRACERS
826 IF ( usePTRACERS ) THEN
827 CALL PTRACERS_DEBUG(myThid)
828 ENDIF
829 #endif /* ALLOW_PTRACERS */
830 ENDIF
831 #endif
832
833 #ifndef DISABLE_DEBUGMODE
834 IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
835 #endif
836
837 RETURN
838 END

  ViewVC Help
Powered by ViewVC 1.1.22