/[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.65 - (show annotations) (download)
Sun Jan 25 00:31:52 2004 UTC (20 years, 3 months ago) by dimitri
Branch: MAIN
CVS Tags: hrcube4, checkpoint52j_pre, hrcube_2, hrcube_3
Changes since 1.64: +3 -1 lines
o limit timeave output for hi-res integrations

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

  ViewVC Help
Powered by ViewVC 1.1.22