/[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.82 - (show annotations) (download)
Tue Nov 2 20:41:30 2004 UTC (19 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint56b_post, checkpoint56, checkpoint55j_post, checkpoint56a_post, checkpoint56c_post
Changes since 1.81: +1 -4 lines
removed unused arrays

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.81 2004/10/28 00:30:57 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_AUTODIFF_TAMC
8 # ifdef ALLOW_GMREDI
9 # include "GMREDI_OPTIONS.h"
10 # endif
11 # ifdef ALLOW_KPP
12 # include "KPP_OPTIONS.h"
13 # endif
14 #endif /* ALLOW_AUTODIFF_TAMC */
15
16 CBOP
17 C !ROUTINE: THERMODYNAMICS
18 C !INTERFACE:
19 SUBROUTINE THERMODYNAMICS(myTime, myIter, myThid)
20 C !DESCRIPTION: \bv
21 C *==========================================================*
22 C | SUBROUTINE THERMODYNAMICS
23 C | o Controlling routine for the prognostic part of the
24 C | thermo-dynamics.
25 C *===========================================================
26 C | The algorithm...
27 C |
28 C | "Correction Step"
29 C | =================
30 C | Here we update the horizontal velocities with the surface
31 C | pressure such that the resulting flow is either consistent
32 C | with the free-surface evolution or the rigid-lid:
33 C | U[n] = U* + dt x d/dx P
34 C | V[n] = V* + dt x d/dy P
35 C |
36 C | "Calculation of Gs"
37 C | ===================
38 C | This is where all the accelerations and tendencies (ie.
39 C | physics, parameterizations etc...) are calculated
40 C | rho = rho ( theta[n], salt[n] )
41 C | b = b(rho, theta)
42 C | K31 = K31 ( rho )
43 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
44 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
45 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
46 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
47 C |
48 C | "Time-stepping" or "Prediction"
49 C | ================================
50 C | The models variables are stepped forward with the appropriate
51 C | time-stepping scheme (currently we use Adams-Bashforth II)
52 C | - For momentum, the result is always *only* a "prediction"
53 C | in that the flow may be divergent and will be "corrected"
54 C | later with a surface pressure gradient.
55 C | - Normally for tracers the result is the new field at time
56 C | level [n+1} *BUT* in the case of implicit diffusion the result
57 C | is also *only* a prediction.
58 C | - We denote "predictors" with an asterisk (*).
59 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
60 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
61 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63 C | With implicit diffusion:
64 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66 C | (1 + dt * K * d_zz) theta[n] = theta*
67 C | (1 + dt * K * d_zz) salt[n] = salt*
68 C |
69 C *==========================================================*
70 C \ev
71
72 C !USES:
73 IMPLICIT NONE
74 C == Global variables ===
75 #include "SIZE.h"
76 #include "EEPARAMS.h"
77 #include "PARAMS.h"
78 #include "DYNVARS.h"
79 #include "GRID.h"
80 #include "GAD.h"
81 #ifdef ALLOW_OFFLINE
82 #include "OFFLINE.h"
83 #endif
84 #ifdef ALLOW_PTRACERS
85 #include "PTRACERS_SIZE.h"
86 #include "PTRACERS.h"
87 #endif
88 #ifdef ALLOW_TIMEAVE
89 #include "TIMEAVE_STATV.h"
90 #endif
91
92 #ifdef ALLOW_AUTODIFF_TAMC
93 # include "tamc.h"
94 # include "tamc_keys.h"
95 # include "FFIELDS.h"
96 # include "EOS.h"
97 # ifdef ALLOW_KPP
98 # include "KPP.h"
99 # endif
100 # ifdef ALLOW_GMREDI
101 # include "GMREDI.h"
102 # endif
103 # ifdef ALLOW_EBM
104 # include "EBM.h"
105 # endif
106 #endif /* ALLOW_AUTODIFF_TAMC */
107
108
109 C !INPUT/OUTPUT PARAMETERS:
110 C == Routine arguments ==
111 C myTime - Current time in simulation
112 C myIter - Current iteration number in simulation
113 C myThid - Thread number for this instance of the routine.
114 _RL myTime
115 INTEGER myIter
116 INTEGER myThid
117
118 C !LOCAL VARIABLES:
119 C == Local variables
120 C xA, yA - Per block temporaries holding face areas
121 C uTrans, vTrans, rTrans - Per block temporaries holding flow
122 C transport
123 C o uTrans: Zonal transport
124 C o vTrans: Meridional transport
125 C o rTrans: Vertical transport
126 C rTransKp1 o vertical volume transp. at interface k+1
127 C maskUp o maskUp: land/water mask for W points
128 C fVer[STUV] o fVer: Vertical flux term - note fVer
129 C is "pipelined" in the vertical
130 C so we need an fVer for each
131 C variable.
132 C kappaRT, - Total diffusion in vertical at level k, for T and S
133 C kappaRS (background + spatially varying, isopycnal term).
134 C kappaRTr - Total diffusion in vertical at level k,
135 C for each passive Tracer
136 C kappaRk - Total diffusion in vertical, all levels, 1 tracer
137 C useVariableK = T when vertical diffusion is not constant
138 C iMin, iMax - Ranges and sub-block indices on which calculations
139 C jMin, jMax are applied.
140 C bi, bj
141 C k, kup, - Index for layer above and below. kup and kDown
142 C kDown, km1 are switched with layer to be the appropriate
143 C index into fVerTerm.
144 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
152 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
153 _RL kappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
154 _RL kappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
155 #ifdef ALLOW_PTRACERS
156 _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
157 _RL kappaRTr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,PTRACERS_num)
158 #endif
159 _RL kappaRk (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
160 _RL kp1Msk
161 LOGICAL useVariableK
162 INTEGER iMin, iMax
163 INTEGER jMin, jMax
164 INTEGER bi, bj
165 INTEGER i, j
166 INTEGER k, km1, kup, kDown
167 INTEGER iTracer, ip
168
169 CEOP
170
171 #ifdef ALLOW_DEBUG
172 IF ( debugLevel .GE. debLevB )
173 & CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
174 #endif
175
176 #ifdef ALLOW_AUTODIFF_TAMC
177 C-- dummy statement to end declaration part
178 ikey = 1
179 itdkey = 1
180 #endif /* ALLOW_AUTODIFF_TAMC */
181
182 #ifdef ALLOW_AUTODIFF_TAMC
183 C-- HPF directive to help TAMC
184 CHPF$ INDEPENDENT
185 #endif /* ALLOW_AUTODIFF_TAMC */
186
187 DO bj=myByLo(myThid),myByHi(myThid)
188
189 #ifdef ALLOW_AUTODIFF_TAMC
190 C-- HPF directive to help TAMC
191 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
192 CHPF$& ,utrans,vtrans,xA,yA
193 CHPF$& ,kappaRT,kappaRS
194 CHPF$& )
195 #endif /* ALLOW_AUTODIFF_TAMC */
196
197 DO bi=myBxLo(myThid),myBxHi(myThid)
198
199 #ifdef ALLOW_AUTODIFF_TAMC
200 act1 = bi - myBxLo(myThid)
201 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
202 act2 = bj - myByLo(myThid)
203 max2 = myByHi(myThid) - myByLo(myThid) + 1
204 act3 = myThid - 1
205 max3 = nTx*nTy
206 act4 = ikey_dynamics - 1
207 itdkey = (act1 + 1) + act2*max1
208 & + act3*max1*max2
209 & + act4*max1*max2*max3
210 #endif /* ALLOW_AUTODIFF_TAMC */
211
212 C-- Set up work arrays with valid (i.e. not NaN) values
213 C These inital values do not alter the numerical results. They
214 C just ensure that all memory references are to valid floating
215 C point numbers. This prevents spurious hardware signals due to
216 C uninitialised but inert locations.
217
218 DO j=1-OLy,sNy+OLy
219 DO i=1-OLx,sNx+OLx
220 xA(i,j) = 0. _d 0
221 yA(i,j) = 0. _d 0
222 uTrans(i,j) = 0. _d 0
223 vTrans(i,j) = 0. _d 0
224 rTrans (i,j) = 0. _d 0
225 rTransKp1(i,j) = 0. _d 0
226 fVerT (i,j,1) = 0. _d 0
227 fVerT (i,j,2) = 0. _d 0
228 fVerS (i,j,1) = 0. _d 0
229 fVerS (i,j,2) = 0. _d 0
230 kappaRT(i,j) = 0. _d 0
231 kappaRS(i,j) = 0. _d 0
232 #ifdef ALLOW_PTRACERS
233 DO ip=1,PTRACERS_num
234 fVerP (i,j,1,ip) = 0. _d 0
235 fVerP (i,j,2,ip) = 0. _d 0
236 kappaRTr(i,j,ip) = 0. _d 0
237 ENDDO
238 #endif
239 ENDDO
240 ENDDO
241
242 DO k=1,Nr
243 DO j=1-OLy,sNy+OLy
244 DO i=1-OLx,sNx+OLx
245 C This is currently also used by IVDC and Diagnostics
246 kappaRk(i,j,k) = 0. _d 0
247 C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
248 gT(i,j,k,bi,bj) = 0. _d 0
249 gS(i,j,k,bi,bj) = 0. _d 0
250 # ifdef ALLOW_PTRACERS
251 DO iTracer=1,PTRACERS_numInUse
252 gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
253 ENDDO
254 # endif
255 ENDDO
256 ENDDO
257 ENDDO
258
259 c iMin = 1-OLx
260 c iMax = sNx+OLx
261 c jMin = 1-OLy
262 c jMax = sNy+OLy
263
264 #ifdef ALLOW_AUTODIFF_TAMC
265 cph avoids recomputation of integrate_for_w
266 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
267 #endif /* ALLOW_AUTODIFF_TAMC */
268
269 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
270 C-- MOST of THERMODYNAMICS will be disabled
271 #ifndef SINGLE_LAYER_MODE
272
273 #ifdef ALLOW_AUTODIFF_TAMC
274 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
275 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
276 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
277 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
278 #ifdef ALLOW_PTRACERS
279 cph-- moved to forward_step to avoid key computation
280 cphCADJ STORE ptracer(:,:,:,bi,bj,itracer) = comlev1_bibj,
281 cphCADJ & key=itdkey, byte=isbyte
282 #endif
283 #endif /* ALLOW_AUTODIFF_TAMC */
284
285 #ifndef DISABLE_MULTIDIM_ADVECTION
286 C-- Some advection schemes are better calculated using a multi-dimensional
287 C method in the absence of any other terms and, if used, is done here.
288 C
289 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
290 C The default is to use multi-dimensinal advection for non-linear advection
291 C schemes. However, for the sake of efficiency of the adjoint it is necessary
292 C to be able to exclude this scheme to avoid excessive storage and
293 C recomputation. It *is* differentiable, if you need it.
294 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
295 C disable this section of code.
296 #ifndef ALLOW_OFFLINE
297 IF (tempMultiDimAdvec) THEN
298 #ifdef ALLOW_DEBUG
299 IF ( debugLevel .GE. debLevB )
300 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
301 #endif
302 CALL GAD_ADVECTION(
303 I tempImplVertAdv, tempAdvScheme, tempVertAdvScheme,
304 I GAD_TEMPERATURE,
305 I uVel, vVel, wVel, theta,
306 O gT,
307 I bi,bj,myTime,myIter,myThid)
308 ENDIF
309 #endif
310 #ifndef ALLOW_OFFLINE
311 IF (saltMultiDimAdvec) THEN
312 #ifdef ALLOW_DEBUG
313 IF ( debugLevel .GE. debLevB )
314 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
315 #endif
316 CALL GAD_ADVECTION(
317 I saltImplVertAdv, saltAdvScheme, saltVertAdvScheme,
318 I GAD_SALINITY,
319 I uVel, vVel, wVel, salt,
320 O gS,
321 I bi,bj,myTime,myIter,myThid)
322 ENDIF
323 #endif
324 C Since passive tracers are configurable separately from T,S we
325 C call the multi-dimensional method for PTRACERS regardless
326 C of whether multiDimAdvection is set or not.
327 #ifdef ALLOW_PTRACERS
328 IF ( usePTRACERS ) THEN
329 #ifdef ALLOW_DEBUG
330 IF ( debugLevel .GE. debLevB )
331 & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
332 #endif
333 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
334 ENDIF
335 #endif /* ALLOW_PTRACERS */
336 #endif /* DISABLE_MULTIDIM_ADVECTION */
337
338 #ifdef ALLOW_DEBUG
339 IF ( debugLevel .GE. debLevB )
340 & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
341 #endif
342
343 C-- Start of thermodynamics loop
344 DO k=Nr,1,-1
345 #ifdef ALLOW_AUTODIFF_TAMC
346 C? Patrick Is this formula correct?
347 cph Yes, but I rewrote it.
348 cph Also, the kappaR? need the index and subscript k!
349 kkey = (itdkey-1)*Nr + k
350 #endif /* ALLOW_AUTODIFF_TAMC */
351
352 C-- km1 Points to level above k (=k-1)
353 C-- kup Cycles through 1,2 to point to layer above
354 C-- kDown Cycles through 2,1 to point to current layer
355
356 km1 = MAX(1,k-1)
357 kup = 1+MOD(k+1,2)
358 kDown= 1+MOD(k,2)
359
360 iMin = 1-OLx
361 iMax = sNx+OLx
362 jMin = 1-OLy
363 jMax = sNy+OLy
364
365 kp1Msk=1.
366 IF (k.EQ.Nr) kp1Msk=0.
367 DO j=1-Oly,sNy+Oly
368 DO i=1-Olx,sNx+Olx
369 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
370 ENDDO
371 ENDDO
372 #ifdef ALLOW_AUTODIFF_TAMC
373 CADJ STORE rTransKp1(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
374 #endif
375
376 C-- Get temporary terms used by tendency routines
377 CALL CALC_COMMON_FACTORS (
378 I bi,bj,iMin,iMax,jMin,jMax,k,
379 O xA,yA,uTrans,vTrans,rTrans,maskUp,
380 I myThid)
381
382 IF (k.EQ.1) THEN
383 C- Surface interface :
384 DO j=1-Oly,sNy+Oly
385 DO i=1-Olx,sNx+Olx
386 rTrans(i,j) = 0.
387 ENDDO
388 ENDDO
389 ELSE
390 C- Interior interface :
391 DO j=1-Oly,sNy+Oly
392 DO i=1-Olx,sNx+Olx
393 rTrans(i,j) = rTrans(i,j)*maskC(i,j,k-1,bi,bj)
394 ENDDO
395 ENDDO
396 ENDIF
397
398 #ifdef ALLOW_GMREDI
399
400 C-- Residual transp = Bolus transp + Eulerian transp
401 IF (useGMRedi) THEN
402 CALL GMREDI_CALC_UVFLOW(
403 & uTrans, vTrans, bi, bj, k, myThid)
404 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
405 & rTrans, bi, bj, k, myThid)
406 ENDIF
407
408 #ifdef ALLOW_AUTODIFF_TAMC
409 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
410 #ifdef GM_BOLUS_ADVEC
411 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
412 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
413 #endif
414 #endif /* ALLOW_AUTODIFF_TAMC */
415
416 #endif /* ALLOW_GMREDI */
417
418 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
419 C-- Calculate the total vertical diffusivity
420 IF ( .NOT.implicitDiffusion ) THEN
421 CALL CALC_DIFFUSIVITY(
422 I bi,bj,iMin,iMax,jMin,jMax,k,
423 I maskUp,
424 O kappaRT,kappaRS,
425 I myThid)
426 ENDIF
427 # ifdef ALLOW_AUTODIFF_TAMC
428 CADJ STORE kappaRT(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
429 CADJ STORE kappaRS(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
430 # endif /* ALLOW_AUTODIFF_TAMC */
431 #endif
432
433 iMin = 1-OLx+2
434 iMax = sNx+OLx-1
435 jMin = 1-OLy+2
436 jMax = sNy+OLy-1
437
438 C-- Calculate active tracer tendencies (gT,gS,...)
439 C and step forward storing result in gTnm1, gSnm1, etc.
440 #ifndef ALLOW_OFFLINE
441 IF ( tempStepping ) THEN
442 CALL CALC_GT(
443 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
444 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
445 I kappaRT,
446 U fVerT,
447 I myTime,myIter,myThid)
448 CALL TIMESTEP_TRACER(
449 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
450 I theta, gT,
451 I myIter, myThid)
452 ENDIF
453 #endif
454
455 #ifndef ALLOW_OFFLINE
456 IF ( saltStepping ) THEN
457 CALL CALC_GS(
458 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
459 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
460 I kappaRS,
461 U fVerS,
462 I myTime,myIter,myThid)
463 CALL TIMESTEP_TRACER(
464 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
465 I salt, gS,
466 I myIter, myThid)
467 ENDIF
468 #endif
469 #ifdef ALLOW_PTRACERS
470 IF ( usePTRACERS ) THEN
471 IF ( .NOT.implicitDiffusion ) THEN
472 CALL PTRACERS_CALC_DIFF(
473 I bi,bj,iMin,iMax,jMin,jMax,k,
474 I maskUp,
475 O kappaRTr,
476 I myThid)
477 ENDIF
478 CALL PTRACERS_INTEGRATE(
479 I bi,bj,k,
480 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
481 U fVerP,
482 I kappaRTr,
483 I myIter,myTime,myThid)
484 ENDIF
485 #endif /* ALLOW_PTRACERS */
486
487 #ifdef ALLOW_OBCS
488 C-- Apply open boundary conditions
489 IF (useOBCS) THEN
490 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
491 END IF
492 #endif /* ALLOW_OBCS */
493
494 C-- Freeze water
495 C this bit of code is left here for backward compatibility.
496 C freezing at surface level has been moved to FORWARD_STEP
497 #ifndef ALLOW_OFFLINE
498 IF ( useOldFreezing .AND. .NOT. useSEAICE
499 & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
500 #ifdef ALLOW_AUTODIFF_TAMC
501 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
502 CADJ & , key = kkey, byte = isbyte
503 #endif /* ALLOW_AUTODIFF_TAMC */
504 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
505 ENDIF
506 #endif
507
508 C-- end of thermodynamic k loop (Nr:1)
509 ENDDO
510
511 iMin = 1
512 iMax = sNx
513 jMin = 1
514 jMax = sNy
515
516 C-- Implicit vertical advection & diffusion
517 #ifndef ALLOW_OFFLINE
518 IF ( tempStepping .AND. implicitDiffusion ) THEN
519 CALL CALC_3D_DIFFUSIVITY(
520 I bi,bj,iMin,iMax,jMin,jMax,
521 I GAD_TEMPERATURE, useGMredi, useKPP,
522 O kappaRk,
523 I myThid)
524 ENDIF
525 #ifdef INCLUDE_IMPLVERTADV_CODE
526 IF ( tempImplVertAdv ) THEN
527 CALL GAD_IMPLICIT_R(
528 I tempImplVertAdv, tempAdvScheme, GAD_TEMPERATURE,
529 I kappaRk, wVel, theta,
530 U gT,
531 I bi, bj, myTime, myIter, myThid )
532 ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN
533 #else /* INCLUDE_IMPLVERTADV_CODE */
534 IF ( tempStepping .AND. implicitDiffusion ) THEN
535 #endif /* INCLUDE_IMPLVERTADV_CODE */
536 #ifdef ALLOW_AUTODIFF_TAMC
537 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
538 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
539 #endif /* ALLOW_AUTODIFF_TAMC */
540 CALL IMPLDIFF(
541 I bi, bj, iMin, iMax, jMin, jMax,
542 I deltaTtracer, kappaRk, recip_HFacC,
543 U gT,
544 I myThid )
545 ENDIF
546 #endif /* ndef ALLOW_OFFLINE */
547
548 #ifdef ALLOW_TIMEAVE
549 #ifndef HRCUBE
550 useVariableK = useKPP .OR. usePP81 .OR. useMY82 .OR. useGGL90
551 & .OR. useGMredi .OR. ivdc_kappa.NE.0.
552 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
553 IF (implicitDiffusion) THEN
554 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRk,
555 I Nr, 3, deltaTclock, bi, bj, myThid)
556 c ELSE
557 c CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
558 c I Nr, 3, deltaTclock, bi, bj, myThid)
559 ENDIF
560 ENDIF
561 #endif /* ndef HRCUBE */
562 #endif /* ALLOW_TIMEAVE */
563
564 #ifndef ALLOW_OFFLINE
565 IF ( saltStepping .AND. implicitDiffusion ) THEN
566 CALL CALC_3D_DIFFUSIVITY(
567 I bi,bj,iMin,iMax,jMin,jMax,
568 I GAD_SALINITY, useGMredi, useKPP,
569 O kappaRk,
570 I myThid)
571 ENDIF
572
573 #ifdef INCLUDE_IMPLVERTADV_CODE
574 IF ( saltImplVertAdv ) THEN
575 CALL GAD_IMPLICIT_R(
576 I saltImplVertAdv, saltAdvScheme, GAD_SALINITY,
577 I kappaRk, wVel, salt,
578 U gS,
579 I bi, bj, myTime, myIter, myThid )
580 ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN
581 #else /* INCLUDE_IMPLVERTADV_CODE */
582 IF ( saltStepping .AND. implicitDiffusion ) THEN
583 #endif /* INCLUDE_IMPLVERTADV_CODE */
584 #ifdef ALLOW_AUTODIFF_TAMC
585 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
586 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
587 #endif /* ALLOW_AUTODIFF_TAMC */
588 CALL IMPLDIFF(
589 I bi, bj, iMin, iMax, jMin, jMax,
590 I deltaTtracer, kappaRk, recip_HFacC,
591 U gS,
592 I myThid )
593 ENDIF
594 #endif
595
596 #ifdef ALLOW_PTRACERS
597 IF ( usePTRACERS .AND. implicitDiffusion ) THEN
598 C-- Vertical diffusion (implicit) for passive tracers
599 CALL PTRACERS_IMPLDIFF( bi,bj,kappaRk,myThid )
600 ENDIF
601 #endif /* ALLOW_PTRACERS */
602
603 #ifdef ALLOW_OBCS
604 C-- Apply open boundary conditions
605 IF ( ( implicitDiffusion
606 & .OR. tempImplVertAdv
607 & .OR. saltImplVertAdv
608 & ) .AND. useOBCS ) THEN
609 DO K=1,Nr
610 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
611 ENDDO
612 ENDIF
613 #endif /* ALLOW_OBCS */
614
615 #ifdef ALLOW_TIMEAVE
616 IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
617 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
618 ENDIF
619 #ifndef HRCUBE
620 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
621 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
622 I Nr, deltaTclock, bi, bj, myThid)
623 ENDIF
624 #endif /* ndef HRCUBE */
625 #endif /* ALLOW_TIMEAVE */
626
627 #endif /* SINGLE_LAYER_MODE */
628
629 C-- end bi,bj loops.
630 ENDDO
631 ENDDO
632
633 #ifdef ALLOW_DEBUG
634 If (debugMode) THEN
635 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
636 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
637 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
638 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
639 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
640 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
641 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
642 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
643 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
644 #ifdef ALLOW_PTRACERS
645 IF ( usePTRACERS ) THEN
646 CALL PTRACERS_DEBUG(myThid)
647 ENDIF
648 #endif /* ALLOW_PTRACERS */
649 ENDIF
650 #endif
651
652 #ifdef ALLOW_DEBUG
653 IF ( debugLevel .GE. debLevB )
654 & CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
655 #endif
656
657 RETURN
658 END

  ViewVC Help
Powered by ViewVC 1.1.22