/[MITgcm]/MITgcm/verification/matrix_example/code/thermodynamics.F
ViewVC logotype

Contents of /MITgcm/verification/matrix_example/code/thermodynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (show annotations) (download)
Wed Apr 20 18:30:58 2005 UTC (19 years ago) by spk
Branch: MAIN
CVS Tags: HEAD
Changes since 1.2: +1 -1 lines
FILE REMOVED
Removing files that are now incorporated in main code tree

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

  ViewVC Help
Powered by ViewVC 1.1.22