/[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.84 - (show annotations) (download)
Sat Dec 4 00:12:14 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
Changes since 1.83: +11 -5 lines
depth convergence accelerator: replace deltaTtracer by dTtracerLev(k)

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

  ViewVC Help
Powered by ViewVC 1.1.22