/[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.96 - (show annotations) (download)
Thu Dec 8 15:44:34 2005 UTC (18 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57y_post
Changes since 1.95: +18 -1 lines
First step for a NLFS adjoint
o initially suppress rStar (new flag DISABLE_RSTAR_CODE)
o new init. routines for calc_r_star, calc_surf_dr
o still need to deal with ini_masks_etc
o testreport seemed happy

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

  ViewVC Help
Powered by ViewVC 1.1.22