/[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.40 - (show annotations) (download)
Tue May 13 17:42:00 2003 UTC (21 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint50f_post, checkpoint50f_pre, checkpoint50g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50e_post
Changes since 1.39: +53 -1 lines
Extended pkg/debug and instrumented main code to help track down fatal
errors.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/thermodynamics.F,v 1.39 2003/05/01 22:30:33 jmc 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 C-- Set up work arrays with valid (i.e. not NaN) values
177 C These inital values do not alter the numerical results. They
178 C just ensure that all memory references are to valid floating
179 C point numbers. This prevents spurious hardware signals due to
180 C uninitialised but inert locations.
181 DO j=1-OLy,sNy+OLy
182 DO i=1-OLx,sNx+OLx
183 xA(i,j) = 0. _d 0
184 yA(i,j) = 0. _d 0
185 uTrans(i,j) = 0. _d 0
186 vTrans(i,j) = 0. _d 0
187 rhok (i,j) = 0. _d 0
188 phiSurfX(i,j) = 0. _d 0
189 phiSurfY(i,j) = 0. _d 0
190 ENDDO
191 ENDDO
192
193
194 #ifdef ALLOW_AUTODIFF_TAMC
195 C-- HPF directive to help TAMC
196 CHPF$ INDEPENDENT
197 #endif /* ALLOW_AUTODIFF_TAMC */
198
199 DO bj=myByLo(myThid),myByHi(myThid)
200
201 #ifdef ALLOW_AUTODIFF_TAMC
202 C-- HPF directive to help TAMC
203 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
204 CHPF$& ,utrans,vtrans,xA,yA
205 CHPF$& ,KappaRT,KappaRS
206 CHPF$& )
207 #endif /* ALLOW_AUTODIFF_TAMC */
208
209 DO bi=myBxLo(myThid),myBxHi(myThid)
210
211 #ifdef ALLOW_AUTODIFF_TAMC
212 act1 = bi - myBxLo(myThid)
213 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
214 act2 = bj - myByLo(myThid)
215 max2 = myByHi(myThid) - myByLo(myThid) + 1
216 act3 = myThid - 1
217 max3 = nTx*nTy
218 act4 = ikey_dynamics - 1
219 itdkey = (act1 + 1) + act2*max1
220 & + act3*max1*max2
221 & + act4*max1*max2*max3
222 #endif /* ALLOW_AUTODIFF_TAMC */
223
224 C-- Set up work arrays that need valid initial values
225 DO j=1-OLy,sNy+OLy
226 DO i=1-OLx,sNx+OLx
227 rTrans (i,j) = 0. _d 0
228 fVerT (i,j,1) = 0. _d 0
229 fVerT (i,j,2) = 0. _d 0
230 fVerS (i,j,1) = 0. _d 0
231 fVerS (i,j,2) = 0. _d 0
232 fVerTr1(i,j,1) = 0. _d 0
233 fVerTr1(i,j,2) = 0. _d 0
234 rhoKM1 (i,j) = 0. _d 0
235 ENDDO
236 ENDDO
237
238 DO k=1,Nr
239 DO j=1-OLy,sNy+OLy
240 DO i=1-OLx,sNx+OLx
241 C This is currently also used by IVDC and Diagnostics
242 sigmaX(i,j,k) = 0. _d 0
243 sigmaY(i,j,k) = 0. _d 0
244 sigmaR(i,j,k) = 0. _d 0
245 ConvectCount(i,j,k) = 0.
246 KappaRT(i,j,k) = 0. _d 0
247 KappaRS(i,j,k) = 0. _d 0
248 #ifdef ALLOW_AUTODIFF_TAMC
249 cph all the following init. are necessary for TAF
250 cph although some of these are re-initialised later.
251 gT(i,j,k,bi,bj) = 0. _d 0
252 gS(i,j,k,bi,bj) = 0. _d 0
253 # ifdef ALLOW_PASSIVE_TRACER
254 gTr1(i,j,k,bi,bj) = 0. _d 0
255 # endif
256 # ifdef ALLOW_GMREDI
257 Kwx(i,j,k,bi,bj) = 0. _d 0
258 Kwy(i,j,k,bi,bj) = 0. _d 0
259 Kwz(i,j,k,bi,bj) = 0. _d 0
260 # ifdef GM_NON_UNITY_DIAGONAL
261 Kux(i,j,k,bi,bj) = 0. _d 0
262 Kvy(i,j,k,bi,bj) = 0. _d 0
263 # endif
264 # ifdef GM_EXTRA_DIAGONAL
265 Kuz(i,j,k,bi,bj) = 0. _d 0
266 Kvz(i,j,k,bi,bj) = 0. _d 0
267 # endif
268 # ifdef GM_BOLUS_ADVEC
269 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
270 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
271 # endif
272 # ifdef GM_VISBECK_VARIABLE_K
273 VisbeckK(i,j,bi,bj) = 0. _d 0
274 # endif
275 # endif /* ALLOW_GMREDI */
276 #endif /* ALLOW_AUTODIFF_TAMC */
277 ENDDO
278 ENDDO
279 ENDDO
280
281 iMin = 1-OLx
282 iMax = sNx+OLx
283 jMin = 1-OLy
284 jMax = sNy+OLy
285
286
287 #ifdef ALLOW_AUTODIFF_TAMC
288 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
289 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
290 CADJ STORE totphihyd
291 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
292 #ifdef ALLOW_KPP
293 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
294 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
295 #endif
296 #endif /* ALLOW_AUTODIFF_TAMC */
297
298 #ifndef DISABLE_DEBUGMODE
299 IF (debugMode) CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
300 #endif
301
302 C-- Start of diagnostic loop
303 DO k=Nr,1,-1
304
305 #ifdef ALLOW_AUTODIFF_TAMC
306 C? Patrick, is this formula correct now that we change the loop range?
307 C? Do we still need this?
308 cph kkey formula corrected.
309 cph Needed for rhok, rhokm1, in the case useGMREDI.
310 kkey = (itdkey-1)*Nr + k
311 #endif /* ALLOW_AUTODIFF_TAMC */
312
313 C-- Integrate continuity vertically for vertical velocity
314 c CALL INTEGRATE_FOR_W(
315 c I bi, bj, k, uVel, vVel,
316 c O wVel,
317 c I myThid )
318
319 #ifdef ALLOW_OBCS
320 #ifdef ALLOW_NONHYDROSTATIC
321 C-- Apply OBC to W if in N-H mode
322 c IF (useOBCS.AND.nonHydrostatic) THEN
323 c CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
324 c ENDIF
325 #endif /* ALLOW_NONHYDROSTATIC */
326 #endif /* ALLOW_OBCS */
327
328 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
329 C-- MOST of THERMODYNAMICS will be disabled
330 #ifndef SINGLE_LAYER_MODE
331
332 C-- Calculate gradients of potential density for isoneutral
333 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
334 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
335 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
336 #ifndef DISABLE_DEBUGMODE
337 IF (debugMode) CALL DEBUG_CALL('FIND_RHO',myThid)
338 #endif
339 #ifdef ALLOW_AUTODIFF_TAMC
340 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
341 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
342 #endif /* ALLOW_AUTODIFF_TAMC */
343 CALL FIND_RHO(
344 I bi, bj, iMin, iMax, jMin, jMax, k, k,
345 I theta, salt,
346 O rhoK,
347 I myThid )
348
349 IF (k.GT.1) THEN
350 #ifdef ALLOW_AUTODIFF_TAMC
351 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
352 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
353 #endif /* ALLOW_AUTODIFF_TAMC */
354 CALL FIND_RHO(
355 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
356 I theta, salt,
357 O rhoKm1,
358 I myThid )
359 ENDIF
360 #ifndef DISABLE_DEBUGMODE
361 IF (debugMode) CALL DEBUG_CALL('GRAD_SIGMA',myThid)
362 #endif
363 CALL GRAD_SIGMA(
364 I bi, bj, iMin, iMax, jMin, jMax, k,
365 I rhoK, rhoKm1, rhoK,
366 O sigmaX, sigmaY, sigmaR,
367 I myThid )
368 ENDIF
369
370 #ifdef ALLOW_AUTODIFF_TAMC
371 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
372 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
373 #endif /* ALLOW_AUTODIFF_TAMC */
374 C-- Implicit Vertical Diffusion for Convection
375 c ==> should use sigmaR !!!
376 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
377 #ifndef DISABLE_DEBUGMODE
378 IF (debugMode) CALL DEBUG_CALL('CALC_IVDC',myThid)
379 #endif
380 CALL CALC_IVDC(
381 I bi, bj, iMin, iMax, jMin, jMax, k,
382 I rhoKm1, rhoK,
383 U ConvectCount, KappaRT, KappaRS,
384 I myTime, myIter, myThid)
385 ENDIF
386
387 #endif /* SINGLE_LAYER_MODE */
388
389 C-- end of diagnostic k loop (Nr:1)
390 ENDDO
391
392 #ifdef ALLOW_AUTODIFF_TAMC
393 cph avoids recomputation of integrate_for_w
394 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
395 #endif /* ALLOW_AUTODIFF_TAMC */
396
397 #ifdef ALLOW_OBCS
398 C-- Calculate future values on open boundaries
399 IF (useOBCS) THEN
400 #ifndef DISABLE_DEBUGMODE
401 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC',myThid)
402 #endif
403 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
404 I uVel, vVel, wVel, theta, salt,
405 I myThid )
406 ENDIF
407 #endif /* ALLOW_OBCS */
408
409
410 c********************************************
411 cswdice --- add ---
412 #ifdef ALLOW_THERM_SEAICE
413 #ifndef DISABLE_DEBUGMODE
414 IF (debugMode) CALL DEBUG_CALL('ICE_FORCING',myThid)
415 #endif
416 C-- Determines forcing terms based on external fields
417 c-- including effects from ice
418 CALL ICE_FORCING(
419 I bi, bj, iMin, iMax, jMin, jMax,
420 I myThid )
421 #else
422
423 cswdice --- end add ---
424
425 C-- Determines forcing terms based on external fields
426 C relaxation terms, etc.
427 #ifndef DISABLE_DEBUGMODE
428 IF (debugMode) CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
429 #endif
430 CALL EXTERNAL_FORCING_SURF(
431 I bi, bj, iMin, iMax, jMin, jMax,
432 I myThid )
433 cswdice --- add ----
434 #endif
435 cswdice --- end add ---
436 c******************************************
437
438
439
440
441 #ifdef ALLOW_AUTODIFF_TAMC
442 cph needed for KPP
443 CADJ STORE surfacetendencyU(:,:,bi,bj)
444 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
445 CADJ STORE surfacetendencyV(:,:,bi,bj)
446 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
447 CADJ STORE surfacetendencyS(:,:,bi,bj)
448 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
449 CADJ STORE surfacetendencyT(:,:,bi,bj)
450 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
451 #endif /* ALLOW_AUTODIFF_TAMC */
452
453 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
454 C-- MOST of THERMODYNAMICS will be disabled
455 #ifndef SINGLE_LAYER_MODE
456
457 #ifdef ALLOW_GMREDI
458
459 #ifdef ALLOW_AUTODIFF_TAMC
460 cph storing here is needed only for one GMREDI_OPTIONS:
461 cph define GM_BOLUS_ADVEC
462 cph but I've avoided the #ifdef for now, in case more things change
463 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
464 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
465 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
466 #endif /* ALLOW_AUTODIFF_TAMC */
467
468 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
469 IF (useGMRedi) THEN
470 #ifndef DISABLE_DEBUGMODE
471 IF (debugMode) CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
472 #endif
473 CALL GMREDI_CALC_TENSOR(
474 I bi, bj, iMin, iMax, jMin, jMax,
475 I sigmaX, sigmaY, sigmaR,
476 I myThid )
477 #ifdef ALLOW_AUTODIFF_TAMC
478 ELSE
479 CALL GMREDI_CALC_TENSOR_DUMMY(
480 I bi, bj, iMin, iMax, jMin, jMax,
481 I sigmaX, sigmaY, sigmaR,
482 I myThid )
483 #endif /* ALLOW_AUTODIFF_TAMC */
484 ENDIF
485
486 #ifdef ALLOW_AUTODIFF_TAMC
487 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
488 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
489 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
490 #endif /* ALLOW_AUTODIFF_TAMC */
491
492 #endif /* ALLOW_GMREDI */
493
494 #ifdef ALLOW_KPP
495 C-- Compute KPP mixing coefficients
496 IF (useKPP) THEN
497 #ifndef DISABLE_DEBUGMODE
498 IF (debugMode) CALL DEBUG_CALL('KPP_CALC',myThid)
499 #endif
500 CALL KPP_CALC(
501 I bi, bj, myTime, myThid )
502 #ifdef ALLOW_AUTODIFF_TAMC
503 ELSE
504 CALL KPP_CALC_DUMMY(
505 I bi, bj, myTime, myThid )
506 #endif /* ALLOW_AUTODIFF_TAMC */
507 ENDIF
508
509 #ifdef ALLOW_AUTODIFF_TAMC
510 CADJ STORE KPPghat (:,:,:,bi,bj)
511 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
512 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
513 CADJ & , KPPfrac (:,: ,bi,bj)
514 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
515 #endif /* ALLOW_AUTODIFF_TAMC */
516
517 #endif /* ALLOW_KPP */
518
519 #ifdef ALLOW_AUTODIFF_TAMC
520 CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
521 CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
522 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
523 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
524 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
525 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
526 #ifdef ALLOW_PASSIVE_TRACER
527 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
528 #endif
529 #endif /* ALLOW_AUTODIFF_TAMC */
530
531 #ifdef ALLOW_AIM
532 C AIM - atmospheric intermediate model, physics package code.
533 IF ( useAIM ) THEN
534 #ifndef DISABLE_DEBUGMODE
535 IF (debugMode) CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
536 #endif
537 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
538 CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
539 CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
540 ENDIF
541 #endif /* ALLOW_AIM */
542
543 #ifndef DISABLE_MULTIDIM_ADVECTION
544 C-- Some advection schemes are better calculated using a multi-dimensional
545 C method in the absence of any other terms and, if used, is done here.
546 C
547 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
548 C The default is to use multi-dimensinal advection for non-linear advection
549 C schemes. However, for the sake of efficiency of the adjoint it is necessary
550 C to be able to exclude this scheme to avoid excessive storage and
551 C recomputation. It *is* differentiable, if you need it.
552 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
553 C disable this section of code.
554 IF (tempMultiDimAdvec) THEN
555 #ifndef DISABLE_DEBUGMODE
556 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
557 #endif
558 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
559 U theta,gT,
560 I myTime,myIter,myThid)
561 ENDIF
562 IF (saltMultiDimAdvec) THEN
563 #ifndef DISABLE_DEBUGMODE
564 IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
565 #endif
566 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
567 U salt,gS,
568 I myTime,myIter,myThid)
569 ENDIF
570 C Since passive tracers are configurable separately from T,S we
571 C call the multi-dimensional method for PTRACERS regardless
572 C of whether multiDimAdvection is set or not.
573 #ifdef ALLOW_PTRACERS
574 IF ( usePTRACERS ) THEN
575 #ifndef DISABLE_DEBUGMODE
576 IF (debugMode) CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
577 #endif
578 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
579 ENDIF
580 #endif /* ALLOW_PTRACERS */
581 #endif /* DISABLE_MULTIDIM_ADVECTION */
582
583 #ifndef DISABLE_DEBUGMODE
584 IF (debugMode) CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
585 #endif
586
587 C-- Start of thermodynamics loop
588 DO k=Nr,1,-1
589 #ifdef ALLOW_AUTODIFF_TAMC
590 C? Patrick Is this formula correct?
591 cph Yes, but I rewrote it.
592 cph Also, the KappaR? need the index and subscript k!
593 kkey = (itdkey-1)*Nr + k
594 #endif /* ALLOW_AUTODIFF_TAMC */
595
596 C-- km1 Points to level above k (=k-1)
597 C-- kup Cycles through 1,2 to point to layer above
598 C-- kDown Cycles through 2,1 to point to current layer
599
600 km1 = MAX(1,k-1)
601 kup = 1+MOD(k+1,2)
602 kDown= 1+MOD(k,2)
603
604 iMin = 1-OLx
605 iMax = sNx+OLx
606 jMin = 1-OLy
607 jMax = sNy+OLy
608
609 C-- Get temporary terms used by tendency routines
610 CALL CALC_COMMON_FACTORS (
611 I bi,bj,iMin,iMax,jMin,jMax,k,
612 O xA,yA,uTrans,vTrans,rTrans,maskUp,
613 I myThid)
614
615 #ifdef ALLOW_GMREDI
616
617 C-- Residual transp = Bolus transp + Eulerian transp
618 IF (useGMRedi) THEN
619 CALL GMREDI_CALC_UVFLOW(
620 & uTrans, vTrans, bi, bj, k, myThid)
621 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
622 & rTrans, bi, bj, k, myThid)
623 ENDIF
624
625 #ifdef ALLOW_AUTODIFF_TAMC
626 #ifdef GM_BOLUS_ADVEC
627 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
628 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
629 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
630 #endif
631 #endif /* ALLOW_AUTODIFF_TAMC */
632
633 #endif /* ALLOW_GMREDI */
634
635 #ifdef ALLOW_AUTODIFF_TAMC
636 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
637 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
638 #endif /* ALLOW_AUTODIFF_TAMC */
639
640 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
641 C-- Calculate the total vertical diffusivity
642 CALL CALC_DIFFUSIVITY(
643 I bi,bj,iMin,iMax,jMin,jMax,k,
644 I maskUp,
645 O KappaRT,KappaRS,
646 I myThid)
647 #endif
648
649 iMin = 1-OLx+2
650 iMax = sNx+OLx-1
651 jMin = 1-OLy+2
652 jMax = sNy+OLy-1
653
654 C-- Calculate active tracer tendencies (gT,gS,...)
655 C and step forward storing result in gTnm1, gSnm1, etc.
656 IF ( tempStepping ) THEN
657 CALL CALC_GT(
658 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
659 I xA,yA,uTrans,vTrans,rTrans,maskUp,
660 I KappaRT,
661 U fVerT,
662 I myTime,myIter,myThid)
663 CALL TIMESTEP_TRACER(
664 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
665 I theta, gT,
666 I myIter, myThid)
667 ENDIF
668 cswdice ---- add ---
669 #ifdef ALLOW_THERM_SEAICE
670 if (k.eq.1) then
671 call ICE_FREEZE(bi, bj, iMin, iMax, jMin, jMax, myThid )
672 endif
673 #endif
674 cswdice -- end add ---
675 IF ( saltStepping ) THEN
676 CALL CALC_GS(
677 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
678 I xA,yA,uTrans,vTrans,rTrans,maskUp,
679 I KappaRS,
680 U fVerS,
681 I myTime,myIter,myThid)
682 CALL TIMESTEP_TRACER(
683 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
684 I salt, gS,
685 I myIter, myThid)
686 ENDIF
687 #ifdef ALLOW_PASSIVE_TRACER
688 IF ( tr1Stepping ) THEN
689 CALL CALC_GTR1(
690 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
691 I xA,yA,uTrans,vTrans,rTrans,maskUp,
692 I KappaRT,
693 U fVerTr1,
694 I myTime,myIter,myThid)
695 CALL TIMESTEP_TRACER(
696 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
697 I Tr1, gTr1,
698 I myIter,myThid)
699 ENDIF
700 #endif
701 #ifdef ALLOW_PTRACERS
702 IF ( usePTRACERS ) THEN
703 CALL PTRACERS_INTEGERATE(
704 I bi,bj,k,
705 I xA,yA,uTrans,vTrans,rTrans,maskUp,
706 X KappaRS,
707 I myIter,myTime,myThid)
708 ENDIF
709 #endif /* ALLOW_PTRACERS */
710
711 #ifdef ALLOW_OBCS
712 C-- Apply open boundary conditions
713 IF (useOBCS) THEN
714 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
715 END IF
716 #endif /* ALLOW_OBCS */
717
718 C-- Freeze water
719 IF ( allowFreezing .AND. .NOT. useSEAICE ) THEN
720 #ifdef ALLOW_AUTODIFF_TAMC
721 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
722 CADJ & , key = kkey, byte = isbyte
723 #endif /* ALLOW_AUTODIFF_TAMC */
724 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
725 END IF
726
727 C-- end of thermodynamic k loop (Nr:1)
728 ENDDO
729
730 cswdice -- add ---
731 #ifdef ALLOW_THERM_SEAICE
732 c timeaveraging for ice model values
733 CALL ICE_AVE(bi,bj,iMin,iMax,jMin,jMax,myThid )
734 #endif
735 cswdice --- end add ---
736
737
738
739
740 C-- Implicit diffusion
741 IF (implicitDiffusion) THEN
742
743 IF (tempStepping) THEN
744 #ifdef ALLOW_AUTODIFF_TAMC
745 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
746 #endif /* ALLOW_AUTODIFF_TAMC */
747 CALL IMPLDIFF(
748 I bi, bj, iMin, iMax, jMin, jMax,
749 I deltaTtracer, KappaRT, recip_HFacC,
750 U gT,
751 I myThid )
752 ENDIF
753
754 IF (saltStepping) THEN
755 #ifdef ALLOW_AUTODIFF_TAMC
756 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
757 #endif /* ALLOW_AUTODIFF_TAMC */
758 CALL IMPLDIFF(
759 I bi, bj, iMin, iMax, jMin, jMax,
760 I deltaTtracer, KappaRS, recip_HFacC,
761 U gS,
762 I myThid )
763 ENDIF
764
765 #ifdef ALLOW_PASSIVE_TRACER
766 IF (tr1Stepping) THEN
767 #ifdef ALLOW_AUTODIFF_TAMC
768 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
769 #endif /* ALLOW_AUTODIFF_TAMC */
770 CALL IMPLDIFF(
771 I bi, bj, iMin, iMax, jMin, jMax,
772 I deltaTtracer, KappaRT, recip_HFacC,
773 U gTr1,
774 I myThid )
775 ENDIF
776 #endif
777
778 #ifdef ALLOW_PTRACERS
779 C Vertical diffusion (implicit) for passive tracers
780 IF ( usePTRACERS ) THEN
781 CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
782 ENDIF
783 #endif /* ALLOW_PTRACERS */
784
785 #ifdef ALLOW_OBCS
786 C-- Apply open boundary conditions
787 IF (useOBCS) THEN
788 DO K=1,Nr
789 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
790 ENDDO
791 END IF
792 #endif /* ALLOW_OBCS */
793
794 C-- End If implicitDiffusion
795 ENDIF
796
797 #ifdef ALLOW_TIMEAVE
798 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
799 CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
800 I Nr, deltaTclock, bi, bj, myThid)
801 ENDIF
802 useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
803 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
804 IF (implicitDiffusion) THEN
805 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
806 I Nr, 3, deltaTclock, bi, bj, myThid)
807 ELSE
808 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
809 I Nr, 3, deltaTclock, bi, bj, myThid)
810 ENDIF
811 ENDIF
812 #endif /* ALLOW_TIMEAVE */
813
814 #endif /* SINGLE_LAYER_MODE */
815
816 C-- end bi,bj loops.
817 ENDDO
818 ENDDO
819
820 #ifndef DISABLE_DEBUGMODE
821 If (debugMode) THEN
822 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
823 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
824 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
825 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
826 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
827 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
828 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
829 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
830 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
831 #ifdef ALLOW_PTRACERS
832 IF ( usePTRACERS ) THEN
833 CALL PTRACERS_DEBUG(myThid)
834 ENDIF
835 #endif /* ALLOW_PTRACERS */
836 ENDIF
837 #endif
838
839 #ifndef DISABLE_DEBUGMODE
840 IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
841 #endif
842
843 RETURN
844 END

  ViewVC Help
Powered by ViewVC 1.1.22