/[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.106 - (show annotations) (download)
Thu Jun 22 21:31:09 2006 UTC (18 years ago) by heimbach
Branch: MAIN
Changes since 1.105: +13 -4 lines
Remove a "hidden" recomputation in bottom_ctrl_5x5/
But does not reconcile G r a d accuracy w.r.t. reference.

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

  ViewVC Help
Powered by ViewVC 1.1.22