/[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.73 - (show annotations) (download)
Wed Jul 7 01:32:11 2004 UTC (19 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint54b_post
Changes since 1.72: +5 -1 lines
try to keep an accurate diagnostic of timeave_surf_flux

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

  ViewVC Help
Powered by ViewVC 1.1.22