/[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.52 - (show annotations) (download)
Fri Oct 10 22:56:08 2003 UTC (20 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51l_post, checkpoint51j_post, checkpoint51n_pre, checkpoint51l_pre, checkpoint51m_post
Branch point for: tg2-branch
Changes since 1.51: +7 -8 lines
adjusted some flow directives

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

  ViewVC Help
Powered by ViewVC 1.1.22