/[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.63 - (show annotations) (download)
Sat Jan 3 01:01:34 2004 UTC (20 years, 5 months ago) by jmc
Branch: MAIN
Changes since 1.62: +55 -30 lines
add calls for implicit vertical direction (advection & diffusion)
    but keep impldiff for implicit diffusion & viscosity only.

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.62 2003/12/30 02:33:42 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 #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 C !LOCAL VARIABLES:
117 C == Local variables
118 C xA, yA - Per block temporaries holding face areas
119 C uTrans, vTrans, rTrans - Per block temporaries holding flow
120 C transport
121 C o uTrans: Zonal transport
122 C o vTrans: Meridional transport
123 C o rTrans: Vertical transport
124 C maskUp o maskUp: land/water mask for W points
125 C fVer[STUV] o fVer: Vertical flux term - note fVer
126 C is "pipelined" in the vertical
127 C so we need an fVer for each
128 C variable.
129 C rhoK, rhoKM1 - Density at current level, and level above
130 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
131 C phiSurfY or geopotentiel (atmos) in X and Y direction
132 C KappaRT, - Total diffusion in vertical for T and S.
133 C KappaRS (background + spatially varying, isopycnal term).
134 C useVariableK = T when vertical diffusion is not constant
135 C iMin, iMax - Ranges and sub-block indices on which calculations
136 C jMin, jMax are applied.
137 C bi, bj
138 C k, kup, - Index for layer above and below. kup and kDown
139 C kDown, km1 are switched with layer to be the appropriate
140 C index into fVerTerm.
141 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
142 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
143 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
144 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
146 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
148 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149 #ifdef ALLOW_PASSIVE_TRACER
150 _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
151 #endif
152 #ifdef ALLOW_PTRACERS
153 _RL fVerP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
154 #endif
155 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
158 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
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 C This is currently used by IVDC and Diagnostics
165 _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
166 LOGICAL useVariableK
167 INTEGER iMin, iMax
168 INTEGER jMin, jMax
169 INTEGER bi, bj
170 INTEGER i, j
171 INTEGER k, km1, kup, kDown
172 INTEGER iTracer, ip
173
174 CEOP
175
176 #ifdef ALLOW_DEBUG
177 IF ( debugLevel .GE. debLevB )
178 & CALL DEBUG_ENTER('THERMODYNAMICS',myThid)
179 #endif
180
181 #ifdef ALLOW_AUTODIFF_TAMC
182 C-- dummy statement to end declaration part
183 ikey = 1
184 itdkey = 1
185 #endif /* ALLOW_AUTODIFF_TAMC */
186
187 #ifdef ALLOW_AUTODIFF_TAMC
188 C-- HPF directive to help TAMC
189 CHPF$ INDEPENDENT
190 #endif /* ALLOW_AUTODIFF_TAMC */
191
192 DO bj=myByLo(myThid),myByHi(myThid)
193
194 #ifdef ALLOW_AUTODIFF_TAMC
195 C-- HPF directive to help TAMC
196 CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS
197 CHPF$& ,utrans,vtrans,xA,yA
198 CHPF$& ,KappaRT,KappaRS
199 CHPF$& )
200 #endif /* ALLOW_AUTODIFF_TAMC */
201
202 DO bi=myBxLo(myThid),myBxHi(myThid)
203
204 #ifdef ALLOW_AUTODIFF_TAMC
205 act1 = bi - myBxLo(myThid)
206 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
207 act2 = bj - myByLo(myThid)
208 max2 = myByHi(myThid) - myByLo(myThid) + 1
209 act3 = myThid - 1
210 max3 = nTx*nTy
211 act4 = ikey_dynamics - 1
212 itdkey = (act1 + 1) + act2*max1
213 & + act3*max1*max2
214 & + act4*max1*max2*max3
215 #endif /* ALLOW_AUTODIFF_TAMC */
216
217 C-- Set up work arrays with valid (i.e. not NaN) values
218 C These inital values do not alter the numerical results. They
219 C just ensure that all memory references are to valid floating
220 C point numbers. This prevents spurious hardware signals due to
221 C uninitialised but inert locations.
222
223 DO j=1-OLy,sNy+OLy
224 DO i=1-OLx,sNx+OLx
225 xA(i,j) = 0. _d 0
226 yA(i,j) = 0. _d 0
227 uTrans(i,j) = 0. _d 0
228 vTrans(i,j) = 0. _d 0
229 rhok (i,j) = 0. _d 0
230 rhoKM1 (i,j) = 0. _d 0
231 phiSurfX(i,j) = 0. _d 0
232 phiSurfY(i,j) = 0. _d 0
233 rTrans (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 #ifdef ALLOW_PASSIVE_TRACER
239 fVerTr1(i,j,1) = 0. _d 0
240 fVerTr1(i,j,2) = 0. _d 0
241 #endif
242 #ifdef ALLOW_PTRACERS
243 DO ip=1,PTRACERS_num
244 fVerP (i,j,1,ip) = 0. _d 0
245 fVerP (i,j,2,ip) = 0. _d 0
246 ENDDO
247 #endif
248 ENDDO
249 ENDDO
250
251 DO k=1,Nr
252 DO j=1-OLy,sNy+OLy
253 DO i=1-OLx,sNx+OLx
254 C This is currently also used by IVDC and Diagnostics
255 sigmaX(i,j,k) = 0. _d 0
256 sigmaY(i,j,k) = 0. _d 0
257 sigmaR(i,j,k) = 0. _d 0
258 ConvectCount(i,j,k) = 0.
259 KappaRT(i,j,k) = 0. _d 0
260 KappaRS(i,j,k) = 0. _d 0
261 C- tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
262 gT(i,j,k,bi,bj) = 0. _d 0
263 gS(i,j,k,bi,bj) = 0. _d 0
264 # ifdef ALLOW_PASSIVE_TRACER
265 ceh3 needs an IF ( use PASSIVE_TRACER) THEN
266 gTr1(i,j,k,bi,bj) = 0. _d 0
267 # endif
268 # ifdef ALLOW_PTRACERS
269 ceh3 this should have an IF ( usePTRACERS ) THEN
270 DO iTracer=1,PTRACERS_numInUse
271 gPTr(i,j,k,bi,bj,itracer) = 0. _d 0
272 ENDDO
273 # endif
274 #ifdef ALLOW_AUTODIFF_TAMC
275 cph all the following init. are necessary for TAF
276 cph although some of these are re-initialised later.
277 # ifdef ALLOW_GMREDI
278 Kwx(i,j,k,bi,bj) = 0. _d 0
279 Kwy(i,j,k,bi,bj) = 0. _d 0
280 Kwz(i,j,k,bi,bj) = 0. _d 0
281 # ifdef GM_NON_UNITY_DIAGONAL
282 Kux(i,j,k,bi,bj) = 0. _d 0
283 Kvy(i,j,k,bi,bj) = 0. _d 0
284 # endif
285 # ifdef GM_EXTRA_DIAGONAL
286 Kuz(i,j,k,bi,bj) = 0. _d 0
287 Kvz(i,j,k,bi,bj) = 0. _d 0
288 # endif
289 # ifdef GM_BOLUS_ADVEC
290 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
291 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
292 # endif
293 # ifdef GM_VISBECK_VARIABLE_K
294 VisbeckK(i,j,bi,bj) = 0. _d 0
295 # endif
296 # endif /* ALLOW_GMREDI */
297 #endif /* ALLOW_AUTODIFF_TAMC */
298 ENDDO
299 ENDDO
300 ENDDO
301
302 iMin = 1-OLx
303 iMax = sNx+OLx
304 jMin = 1-OLy
305 jMax = sNy+OLy
306
307 #ifdef ALLOW_AUTODIFF_TAMC
308 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
309 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
310 CADJ STORE totphihyd
311 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
312 #ifdef ALLOW_KPP
313 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
314 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
315 #endif
316 #endif /* ALLOW_AUTODIFF_TAMC */
317
318 #ifdef ALLOW_DEBUG
319 IF ( debugLevel .GE. debLevB )
320 & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
321 #endif
322
323 C-- Start of diagnostic loop
324 DO k=Nr,1,-1
325
326 #ifdef ALLOW_AUTODIFF_TAMC
327 C? Patrick, is this formula correct now that we change the loop range?
328 C? Do we still need this?
329 cph kkey formula corrected.
330 cph Needed for rhok, rhokm1, in the case useGMREDI.
331 kkey = (itdkey-1)*Nr + k
332 #endif /* ALLOW_AUTODIFF_TAMC */
333
334 C-- Integrate continuity vertically for vertical velocity
335 c CALL INTEGRATE_FOR_W(
336 c I bi, bj, k, uVel, vVel,
337 c O wVel,
338 c I myThid )
339
340 #ifdef ALLOW_OBCS
341 #ifdef ALLOW_NONHYDROSTATIC
342 C-- Apply OBC to W if in N-H mode
343 c IF (useOBCS.AND.nonHydrostatic) THEN
344 c CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
345 c ENDIF
346 #endif /* ALLOW_NONHYDROSTATIC */
347 #endif /* ALLOW_OBCS */
348
349 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
350 C-- MOST of THERMODYNAMICS will be disabled
351 #ifndef SINGLE_LAYER_MODE
352
353 C-- Calculate gradients of potential density for isoneutral
354 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
355 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
356 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
357 #ifdef ALLOW_DEBUG
358 IF ( debugLevel .GE. debLevB )
359 & CALL DEBUG_CALL('FIND_RHO',myThid)
360 #endif
361 #ifdef ALLOW_AUTODIFF_TAMC
362 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
363 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
364 #endif /* ALLOW_AUTODIFF_TAMC */
365 CALL FIND_RHO(
366 I bi, bj, iMin, iMax, jMin, jMax, k, k,
367 I theta, salt,
368 O rhoK,
369 I myThid )
370
371 IF (k.GT.1) THEN
372 #ifdef ALLOW_AUTODIFF_TAMC
373 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
374 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
375 #endif /* ALLOW_AUTODIFF_TAMC */
376 CALL FIND_RHO(
377 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
378 I theta, salt,
379 O rhoKm1,
380 I myThid )
381 ENDIF
382 #ifdef ALLOW_DEBUG
383 IF ( debugLevel .GE. debLevB )
384 & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
385 #endif
386 CALL GRAD_SIGMA(
387 I bi, bj, iMin, iMax, jMin, jMax, k,
388 I rhoK, rhoKm1, rhoK,
389 O sigmaX, sigmaY, sigmaR,
390 I myThid )
391 ENDIF
392
393 #ifdef ALLOW_AUTODIFF_TAMC
394 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
395 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
396 #endif /* ALLOW_AUTODIFF_TAMC */
397 C-- Implicit Vertical Diffusion for Convection
398 c ==> should use sigmaR !!!
399 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
400 #ifdef ALLOW_DEBUG
401 IF ( debugLevel .GE. debLevB )
402 & CALL DEBUG_CALL('CALC_IVDC',myThid)
403 #endif
404 CALL CALC_IVDC(
405 I bi, bj, iMin, iMax, jMin, jMax, k,
406 I rhoKm1, rhoK,
407 U ConvectCount, KappaRT, KappaRS,
408 I myTime, myIter, myThid)
409 ENDIF
410
411 #endif /* SINGLE_LAYER_MODE */
412
413 C-- end of diagnostic k loop (Nr:1)
414 ENDDO
415
416 #ifdef ALLOW_AUTODIFF_TAMC
417 cph avoids recomputation of integrate_for_w
418 CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
419 #endif /* ALLOW_AUTODIFF_TAMC */
420
421 #ifdef ALLOW_OBCS
422 C-- Calculate future values on open boundaries
423 IF (useOBCS) THEN
424 #ifdef ALLOW_DEBUG
425 IF ( debugLevel .GE. debLevB )
426 & CALL DEBUG_CALL('OBCS_CALC',myThid)
427 #endif
428 CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
429 I uVel, vVel, wVel, theta, salt,
430 I myThid )
431 ENDIF
432 #endif /* ALLOW_OBCS */
433
434 IF ( buoyancyRelation(1:7) .EQ. 'OCEANIC' ) THEN
435 C-- Determines forcing terms based on external fields
436 C relaxation terms, etc.
437 #ifdef ALLOW_DEBUG
438 IF ( debugLevel .GE. debLevB )
439 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
440 #endif
441 CALL EXTERNAL_FORCING_SURF(
442 I bi, bj, iMin, iMax, jMin, jMax,
443 I myTime, myIter, myThid )
444
445 #ifdef ALLOW_AUTODIFF_TAMC
446 cph needed for KPP
447 CADJ STORE surfacetendencyU(:,:,bi,bj)
448 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
449 CADJ STORE surfacetendencyV(:,:,bi,bj)
450 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
451 CADJ STORE surfacetendencyS(:,:,bi,bj)
452 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
453 CADJ STORE surfacetendencyT(:,:,bi,bj)
454 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
455 # ifdef ALLOW_SEAICE
456 CADJ STORE surfacetendencyTice(:,:,bi,bj)
457 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
458 # endif
459 #endif /* ALLOW_AUTODIFF_TAMC */
460 ENDIF
461
462 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
463 C-- MOST of THERMODYNAMICS will be disabled
464 #ifndef SINGLE_LAYER_MODE
465
466 #ifdef ALLOW_GMREDI
467
468 #ifdef ALLOW_AUTODIFF_TAMC
469 cph storing here is needed only for one GMREDI_OPTIONS:
470 cph define GM_BOLUS_ADVEC
471 cph but I've avoided the #ifdef for now, in case more things change
472 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
473 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
474 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
475 #endif /* ALLOW_AUTODIFF_TAMC */
476
477 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
478 IF (useGMRedi) THEN
479 #ifdef ALLOW_DEBUG
480 IF ( debugLevel .GE. debLevB )
481 & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
482 #endif
483 CALL GMREDI_CALC_TENSOR(
484 I bi, bj, iMin, iMax, jMin, jMax,
485 I sigmaX, sigmaY, sigmaR,
486 I myThid )
487 #ifdef ALLOW_AUTODIFF_TAMC
488 ELSE
489 CALL GMREDI_CALC_TENSOR_DUMMY(
490 I bi, bj, iMin, iMax, jMin, jMax,
491 I sigmaX, sigmaY, sigmaR,
492 I myThid )
493 #endif /* ALLOW_AUTODIFF_TAMC */
494 ENDIF
495
496 #ifdef ALLOW_AUTODIFF_TAMC
497 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
498 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
499 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
500 #endif /* ALLOW_AUTODIFF_TAMC */
501
502 #endif /* ALLOW_GMREDI */
503
504 #ifdef ALLOW_KPP
505 C-- Compute KPP mixing coefficients
506 IF (useKPP) THEN
507 #ifdef ALLOW_DEBUG
508 IF ( debugLevel .GE. debLevB )
509 & CALL DEBUG_CALL('KPP_CALC',myThid)
510 #endif
511 CALL KPP_CALC(
512 I bi, bj, myTime, myThid )
513 #ifdef ALLOW_AUTODIFF_TAMC
514 ELSE
515 CALL KPP_CALC_DUMMY(
516 I bi, bj, myTime, myThid )
517 #endif /* ALLOW_AUTODIFF_TAMC */
518 ENDIF
519
520 #ifdef ALLOW_AUTODIFF_TAMC
521 CADJ STORE KPPghat (:,:,:,bi,bj)
522 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
523 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
524 CADJ & , KPPfrac (:,: ,bi,bj)
525 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
526 #endif /* ALLOW_AUTODIFF_TAMC */
527
528 #endif /* ALLOW_KPP */
529
530 #ifdef ALLOW_AUTODIFF_TAMC
531 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
532 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
533 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
534 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
535 #ifdef ALLOW_PASSIVE_TRACER
536 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
537 #endif
538 #ifdef ALLOW_PTRACERS
539 cph-- moved to forward_step to avoid key computation
540 cphCADJ STORE ptracer(:,:,:,bi,bj,itracer) = comlev1_bibj,
541 cphCADJ & key=itdkey, byte=isbyte
542 #endif
543 #endif /* ALLOW_AUTODIFF_TAMC */
544
545 #ifdef ALLOW_AIM
546 C AIM - atmospheric intermediate model, physics package code.
547 IF ( useAIM ) THEN
548 #ifdef ALLOW_DEBUG
549 IF ( debugLevel .GE. debLevB )
550 & CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
551 #endif
552 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
553 CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
554 CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
555 ENDIF
556 #endif /* ALLOW_AIM */
557
558 #ifndef DISABLE_MULTIDIM_ADVECTION
559 C-- Some advection schemes are better calculated using a multi-dimensional
560 C method in the absence of any other terms and, if used, is done here.
561 C
562 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
563 C The default is to use multi-dimensinal advection for non-linear advection
564 C schemes. However, for the sake of efficiency of the adjoint it is necessary
565 C to be able to exclude this scheme to avoid excessive storage and
566 C recomputation. It *is* differentiable, if you need it.
567 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
568 C disable this section of code.
569 IF (tempMultiDimAdvec) THEN
570 #ifdef ALLOW_DEBUG
571 IF ( debugLevel .GE. debLevB )
572 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
573 #endif
574 CALL GAD_ADVECTION(
575 I tempImplVertAdv,tempAdvScheme,GAD_TEMPERATURE,
576 I uVel, vVel, wVel, theta,
577 O gT,
578 I bi,bj,myTime,myIter,myThid)
579 ENDIF
580 IF (saltMultiDimAdvec) THEN
581 #ifdef ALLOW_DEBUG
582 IF ( debugLevel .GE. debLevB )
583 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
584 #endif
585 CALL GAD_ADVECTION(
586 I saltImplVertAdv,saltAdvScheme,GAD_SALINITY,
587 I uVel, vVel, wVel, salt,
588 O gS,
589 I bi,bj,myTime,myIter,myThid)
590 ENDIF
591 C Since passive tracers are configurable separately from T,S we
592 C call the multi-dimensional method for PTRACERS regardless
593 C of whether multiDimAdvection is set or not.
594 #ifdef ALLOW_PTRACERS
595 IF ( usePTRACERS ) THEN
596 #ifdef ALLOW_DEBUG
597 IF ( debugLevel .GE. debLevB )
598 & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
599 #endif
600 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
601 ENDIF
602 #endif /* ALLOW_PTRACERS */
603 #endif /* DISABLE_MULTIDIM_ADVECTION */
604
605 #ifdef ALLOW_DEBUG
606 IF ( debugLevel .GE. debLevB )
607 & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
608 #endif
609
610 C-- Start of thermodynamics loop
611 DO k=Nr,1,-1
612 #ifdef ALLOW_AUTODIFF_TAMC
613 C? Patrick Is this formula correct?
614 cph Yes, but I rewrote it.
615 cph Also, the KappaR? need the index and subscript k!
616 kkey = (itdkey-1)*Nr + k
617 #endif /* ALLOW_AUTODIFF_TAMC */
618
619 C-- km1 Points to level above k (=k-1)
620 C-- kup Cycles through 1,2 to point to layer above
621 C-- kDown Cycles through 2,1 to point to current layer
622
623 km1 = MAX(1,k-1)
624 kup = 1+MOD(k+1,2)
625 kDown= 1+MOD(k,2)
626
627 iMin = 1-OLx
628 iMax = sNx+OLx
629 jMin = 1-OLy
630 jMax = sNy+OLy
631
632 C-- Get temporary terms used by tendency routines
633 CALL CALC_COMMON_FACTORS (
634 I bi,bj,iMin,iMax,jMin,jMax,k,
635 O xA,yA,uTrans,vTrans,rTrans,maskUp,
636 I myThid)
637
638 #ifdef ALLOW_GMREDI
639
640 C-- Residual transp = Bolus transp + Eulerian transp
641 IF (useGMRedi) THEN
642 CALL GMREDI_CALC_UVFLOW(
643 & uTrans, vTrans, bi, bj, k, myThid)
644 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
645 & rTrans, bi, bj, k, myThid)
646 ENDIF
647
648 #ifdef ALLOW_AUTODIFF_TAMC
649 #ifdef GM_BOLUS_ADVEC
650 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
651 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
652 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
653 #endif
654 #endif /* ALLOW_AUTODIFF_TAMC */
655
656 #endif /* ALLOW_GMREDI */
657
658 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
659 C-- Calculate the total vertical diffusivity
660 CALL CALC_DIFFUSIVITY(
661 I bi,bj,iMin,iMax,jMin,jMax,k,
662 I maskUp,
663 O KappaRT,KappaRS,
664 I myThid)
665 # ifdef ALLOW_AUTODIFF_TAMC
666 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
667 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
668 # endif /* ALLOW_AUTODIFF_TAMC */
669 #endif
670
671 iMin = 1-OLx+2
672 iMax = sNx+OLx-1
673 jMin = 1-OLy+2
674 jMax = sNy+OLy-1
675
676 C-- Calculate active tracer tendencies (gT,gS,...)
677 C and step forward storing result in gTnm1, gSnm1, etc.
678 IF ( tempStepping ) THEN
679 CALL CALC_GT(
680 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
681 I xA,yA,uTrans,vTrans,rTrans,maskUp,
682 I KappaRT,
683 U fVerT,
684 I myTime,myIter,myThid)
685 CALL TIMESTEP_TRACER(
686 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
687 I theta, gT,
688 I myIter, myThid)
689 ENDIF
690
691 IF ( saltStepping ) THEN
692 CALL CALC_GS(
693 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
694 I xA,yA,uTrans,vTrans,rTrans,maskUp,
695 I KappaRS,
696 U fVerS,
697 I myTime,myIter,myThid)
698 CALL TIMESTEP_TRACER(
699 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
700 I salt, gS,
701 I myIter, myThid)
702 ENDIF
703 #ifdef ALLOW_PASSIVE_TRACER
704 ceh3 needs an IF ( usePASSIVE_TRACER ) THEN
705 IF ( tr1Stepping ) THEN
706 CALL CALC_GTR1(
707 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
708 I xA,yA,uTrans,vTrans,rTrans,maskUp,
709 I KappaRT,
710 U fVerTr1,
711 I myTime,myIter,myThid)
712 CALL TIMESTEP_TRACER(
713 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
714 I Tr1, gTr1,
715 I myIter,myThid)
716 ENDIF
717 #endif
718 #ifdef ALLOW_PTRACERS
719 IF ( usePTRACERS ) THEN
720 CALL PTRACERS_INTEGRATE(
721 I bi,bj,k,
722 I xA,yA,uTrans,vTrans,rTrans,maskUp,
723 X fVerP, KappaRS,
724 I myIter,myTime,myThid)
725 ENDIF
726 #endif /* ALLOW_PTRACERS */
727
728 #ifdef ALLOW_OBCS
729 C-- Apply open boundary conditions
730 IF (useOBCS) THEN
731 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
732 END IF
733 #endif /* ALLOW_OBCS */
734
735 C-- Freeze water
736 C this bit of code is left here for backward compatibility.
737 C freezing at surface level has been moved to FORWARD_STEP
738 IF ( useOldFreezing .AND. .NOT. useSEAICE
739 & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
740 #ifdef ALLOW_AUTODIFF_TAMC
741 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
742 CADJ & , key = kkey, byte = isbyte
743 #endif /* ALLOW_AUTODIFF_TAMC */
744 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
745 ENDIF
746
747 C-- end of thermodynamic k loop (Nr:1)
748 ENDDO
749
750
751 C-- Implicit vertical advection & diffusion
752 #ifdef INCLUDE_IMPLVERTADV_CODE
753 IF ( tempImplVertAdv ) THEN
754 CALL GAD_IMPLICIT_R(
755 I tempImplVertAdv, tempAdvScheme, GAD_TEMPERATURE,
756 I kappaRT, wVel, theta,
757 U gT,
758 I bi, bj, myTime, myIter, myThid )
759 ELSEIF ( tempStepping .AND. implicitDiffusion ) THEN
760 #else /* INCLUDE_IMPLVERTADV_CODE */
761 IF ( tempStepping .AND. implicitDiffusion ) THEN
762 #endif /* INCLUDE_IMPLVERTADV_CODE */
763 #ifdef ALLOW_AUTODIFF_TAMC
764 CADJ STORE KappaRT(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
765 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
766 #endif /* ALLOW_AUTODIFF_TAMC */
767 CALL IMPLDIFF(
768 I bi, bj, iMin, iMax, jMin, jMax,
769 I deltaTtracer, KappaRT, recip_HFacC,
770 U gT,
771 I myThid )
772 ENDIF
773
774 #ifdef INCLUDE_IMPLVERTADV_CODE
775 IF ( saltImplVertAdv ) THEN
776 CALL GAD_IMPLICIT_R(
777 I saltImplVertAdv, saltAdvScheme, GAD_SALINITY,
778 I kappaRS, wVel, salt,
779 U gS,
780 I bi, bj, myTime, myIter, myThid )
781 ELSEIF ( saltStepping .AND. implicitDiffusion ) THEN
782 #else /* INCLUDE_IMPLVERTADV_CODE */
783 IF ( saltStepping .AND. implicitDiffusion ) THEN
784 #endif /* INCLUDE_IMPLVERTADV_CODE */
785 #ifdef ALLOW_AUTODIFF_TAMC
786 CADJ STORE KappaRS(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
787 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
788 #endif /* ALLOW_AUTODIFF_TAMC */
789 CALL IMPLDIFF(
790 I bi, bj, iMin, iMax, jMin, jMax,
791 I deltaTtracer, KappaRS, recip_HFacC,
792 U gS,
793 I myThid )
794 ENDIF
795
796 #ifdef ALLOW_PASSIVE_TRACER
797 IF ( tr1Stepping .AND. implicitDiffusion ) THEN
798 #ifdef ALLOW_AUTODIFF_TAMC
799 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
800 #endif /* ALLOW_AUTODIFF_TAMC */
801 CALL IMPLDIFF(
802 I bi, bj, iMin, iMax, jMin, jMax,
803 I deltaTtracer, KappaRT, recip_HFacC,
804 U gTr1,
805 I myThid )
806 ENDIF
807 #endif
808
809 #ifdef ALLOW_PTRACERS
810 c #ifdef INCLUDE_IMPLVERTADV_CODE
811 c IF ( usePTRACERS .AND. ptracerImplVertAdv ) THEN
812 c ELSEIF ( usePTRACERS .AND. implicitDiffusion ) THEN
813 c #else
814 IF ( usePTRACERS .AND. implicitDiffusion ) THEN
815 C-- Vertical diffusion (implicit) for passive tracers
816 CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
817 ENDIF
818 #endif /* ALLOW_PTRACERS */
819
820 #ifdef ALLOW_OBCS
821 C-- Apply open boundary conditions
822 IF ( ( implicitDiffusion
823 & .OR. tempImplVertAdv
824 & .OR. saltImplVertAdv
825 & ) .AND. useOBCS ) THEN
826 DO K=1,Nr
827 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
828 ENDDO
829 ENDIF
830 #endif /* ALLOW_OBCS */
831
832 #ifdef ALLOW_TIMEAVE
833 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
834 CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
835 I Nr, deltaTclock, bi, bj, myThid)
836 ENDIF
837 useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
838 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
839 IF (implicitDiffusion) THEN
840 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
841 I Nr, 3, deltaTclock, bi, bj, myThid)
842 ELSE
843 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
844 I Nr, 3, deltaTclock, bi, bj, myThid)
845 ENDIF
846 ENDIF
847 #endif /* ALLOW_TIMEAVE */
848
849 #endif /* SINGLE_LAYER_MODE */
850
851 C-- end bi,bj loops.
852 ENDDO
853 ENDDO
854
855 #ifdef ALLOW_DEBUG
856 If (debugMode) THEN
857 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
858 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
859 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
860 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
861 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
862 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
863 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
864 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
865 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
866 #ifdef ALLOW_PTRACERS
867 IF ( usePTRACERS ) THEN
868 CALL PTRACERS_DEBUG(myThid)
869 ENDIF
870 #endif /* ALLOW_PTRACERS */
871 ENDIF
872 #endif
873
874 #ifdef ALLOW_DEBUG
875 IF ( debugLevel .GE. debLevB )
876 & CALL DEBUG_LEAVE('THERMODYNAMICS',myThid)
877 #endif
878
879 RETURN
880 END

  ViewVC Help
Powered by ViewVC 1.1.22