/[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.75 - (show annotations) (download)
Wed Sep 1 16:25:06 2004 UTC (19 years, 9 months ago) by stephd
Branch: MAIN
Changes since 1.74: +19 -1 lines
o adding offline package

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

  ViewVC Help
Powered by ViewVC 1.1.22