/[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.61 - (show annotations) (download)
Sun Nov 23 01:28:05 2003 UTC (20 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint52d_pre, checkpoint52b_post, checkpoint52c_post, checkpoint52d_post, branch-netcdf
Branch point for: netcdf-sm0
Changes since 1.60: +2 -38 lines
use the new thermodynamic Sea-Ice pkg: thSIce

1 C $Header: /u/gcmpack/MITgcm/model/src/thermodynamics.F,v 1.60 2003/11/20 03:46:07 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('FORWARD_STEP',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 C-- Determines forcing terms based on external fields
435 C relaxation terms, etc.
436 #ifdef ALLOW_DEBUG
437 IF ( debugLevel .GE. debLevB )
438 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
439 #endif
440 CALL EXTERNAL_FORCING_SURF(
441 I bi, bj, iMin, iMax, jMin, jMax,
442 I myTime, myIter, myThid )
443
444 #ifdef ALLOW_AUTODIFF_TAMC
445 cph needed for KPP
446 CADJ STORE surfacetendencyU(:,:,bi,bj)
447 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
448 CADJ STORE surfacetendencyV(:,:,bi,bj)
449 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
450 CADJ STORE surfacetendencyS(:,:,bi,bj)
451 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
452 CADJ STORE surfacetendencyT(:,:,bi,bj)
453 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
454 # ifdef ALLOW_SEAICE
455 CADJ STORE surfacetendencyTice(:,:,bi,bj)
456 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
457 # endif
458 #endif /* ALLOW_AUTODIFF_TAMC */
459
460 C-- Attention: by defining "SINGLE_LAYER_MODE" in CPP_OPTIONS.h
461 C-- MOST of THERMODYNAMICS will be disabled
462 #ifndef SINGLE_LAYER_MODE
463
464 #ifdef ALLOW_GMREDI
465
466 #ifdef ALLOW_AUTODIFF_TAMC
467 cph storing here is needed only for one GMREDI_OPTIONS:
468 cph define GM_BOLUS_ADVEC
469 cph but I've avoided the #ifdef for now, in case more things change
470 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
471 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
472 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
473 #endif /* ALLOW_AUTODIFF_TAMC */
474
475 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
476 IF (useGMRedi) THEN
477 #ifdef ALLOW_DEBUG
478 IF ( debugLevel .GE. debLevB )
479 & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
480 #endif
481 CALL GMREDI_CALC_TENSOR(
482 I bi, bj, iMin, iMax, jMin, jMax,
483 I sigmaX, sigmaY, sigmaR,
484 I myThid )
485 #ifdef ALLOW_AUTODIFF_TAMC
486 ELSE
487 CALL GMREDI_CALC_TENSOR_DUMMY(
488 I bi, bj, iMin, iMax, jMin, jMax,
489 I sigmaX, sigmaY, sigmaR,
490 I myThid )
491 #endif /* ALLOW_AUTODIFF_TAMC */
492 ENDIF
493
494 #ifdef ALLOW_AUTODIFF_TAMC
495 CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
496 CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
497 CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
498 #endif /* ALLOW_AUTODIFF_TAMC */
499
500 #endif /* ALLOW_GMREDI */
501
502 #ifdef ALLOW_KPP
503 C-- Compute KPP mixing coefficients
504 IF (useKPP) THEN
505 #ifdef ALLOW_DEBUG
506 IF ( debugLevel .GE. debLevB )
507 & CALL DEBUG_CALL('KPP_CALC',myThid)
508 #endif
509 CALL KPP_CALC(
510 I bi, bj, myTime, myThid )
511 #ifdef ALLOW_AUTODIFF_TAMC
512 ELSE
513 CALL KPP_CALC_DUMMY(
514 I bi, bj, myTime, myThid )
515 #endif /* ALLOW_AUTODIFF_TAMC */
516 ENDIF
517
518 #ifdef ALLOW_AUTODIFF_TAMC
519 CADJ STORE KPPghat (:,:,:,bi,bj)
520 CADJ & , KPPdiffKzT(:,:,:,bi,bj)
521 CADJ & , KPPdiffKzS(:,:,:,bi,bj)
522 CADJ & , KPPfrac (:,: ,bi,bj)
523 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
524 #endif /* ALLOW_AUTODIFF_TAMC */
525
526 #endif /* ALLOW_KPP */
527
528 #ifdef ALLOW_AUTODIFF_TAMC
529 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
530 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
531 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
532 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
533 #ifdef ALLOW_PASSIVE_TRACER
534 CADJ STORE tr1 (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
535 #endif
536 #ifdef ALLOW_PTRACERS
537 cph-- moved to forward_step to avoid key computation
538 cphCADJ STORE ptracer(:,:,:,bi,bj,itracer) = comlev1_bibj,
539 cphCADJ & key=itdkey, byte=isbyte
540 #endif
541 #endif /* ALLOW_AUTODIFF_TAMC */
542
543 #ifdef ALLOW_AIM
544 C AIM - atmospheric intermediate model, physics package code.
545 IF ( useAIM ) THEN
546 #ifdef ALLOW_DEBUG
547 IF ( debugLevel .GE. debLevB )
548 & CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
549 #endif
550 CALL TIMER_START('AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
551 CALL AIM_DO_PHYSICS( bi, bj, myTime, myIter, myThid )
552 CALL TIMER_STOP( 'AIM_DO_PHYSICS [THERMODYNAMICS]', myThid)
553 ENDIF
554 #endif /* ALLOW_AIM */
555
556 #ifndef DISABLE_MULTIDIM_ADVECTION
557 C-- Some advection schemes are better calculated using a multi-dimensional
558 C method in the absence of any other terms and, if used, is done here.
559 C
560 C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
561 C The default is to use multi-dimensinal advection for non-linear advection
562 C schemes. However, for the sake of efficiency of the adjoint it is necessary
563 C to be able to exclude this scheme to avoid excessive storage and
564 C recomputation. It *is* differentiable, if you need it.
565 C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
566 C disable this section of code.
567 IF (tempMultiDimAdvec) THEN
568 #ifdef ALLOW_DEBUG
569 IF ( debugLevel .GE. debLevB )
570 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
571 #endif
572 CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
573 U theta,gT,
574 I myTime,myIter,myThid)
575 ENDIF
576 IF (saltMultiDimAdvec) THEN
577 #ifdef ALLOW_DEBUG
578 IF ( debugLevel .GE. debLevB )
579 & CALL DEBUG_CALL('GAD_ADVECTION',myThid)
580 #endif
581 CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
582 U salt,gS,
583 I myTime,myIter,myThid)
584 ENDIF
585 C Since passive tracers are configurable separately from T,S we
586 C call the multi-dimensional method for PTRACERS regardless
587 C of whether multiDimAdvection is set or not.
588 #ifdef ALLOW_PTRACERS
589 IF ( usePTRACERS ) THEN
590 #ifdef ALLOW_DEBUG
591 IF ( debugLevel .GE. debLevB )
592 & CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
593 #endif
594 CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
595 ENDIF
596 #endif /* ALLOW_PTRACERS */
597 #endif /* DISABLE_MULTIDIM_ADVECTION */
598
599 #ifdef ALLOW_DEBUG
600 IF ( debugLevel .GE. debLevB )
601 & CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
602 #endif
603
604 C-- Start of thermodynamics loop
605 DO k=Nr,1,-1
606 #ifdef ALLOW_AUTODIFF_TAMC
607 C? Patrick Is this formula correct?
608 cph Yes, but I rewrote it.
609 cph Also, the KappaR? need the index and subscript k!
610 kkey = (itdkey-1)*Nr + k
611 #endif /* ALLOW_AUTODIFF_TAMC */
612
613 C-- km1 Points to level above k (=k-1)
614 C-- kup Cycles through 1,2 to point to layer above
615 C-- kDown Cycles through 2,1 to point to current layer
616
617 km1 = MAX(1,k-1)
618 kup = 1+MOD(k+1,2)
619 kDown= 1+MOD(k,2)
620
621 iMin = 1-OLx
622 iMax = sNx+OLx
623 jMin = 1-OLy
624 jMax = sNy+OLy
625
626 C-- Get temporary terms used by tendency routines
627 CALL CALC_COMMON_FACTORS (
628 I bi,bj,iMin,iMax,jMin,jMax,k,
629 O xA,yA,uTrans,vTrans,rTrans,maskUp,
630 I myThid)
631
632 #ifdef ALLOW_GMREDI
633
634 C-- Residual transp = Bolus transp + Eulerian transp
635 IF (useGMRedi) THEN
636 CALL GMREDI_CALC_UVFLOW(
637 & uTrans, vTrans, bi, bj, k, myThid)
638 IF (K.GE.2) CALL GMREDI_CALC_WFLOW(
639 & rTrans, bi, bj, k, myThid)
640 ENDIF
641
642 #ifdef ALLOW_AUTODIFF_TAMC
643 #ifdef GM_BOLUS_ADVEC
644 CADJ STORE uTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
645 CADJ STORE vTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
646 CADJ STORE rTrans(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
647 #endif
648 #endif /* ALLOW_AUTODIFF_TAMC */
649
650 #endif /* ALLOW_GMREDI */
651
652 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
653 C-- Calculate the total vertical diffusivity
654 CALL CALC_DIFFUSIVITY(
655 I bi,bj,iMin,iMax,jMin,jMax,k,
656 I maskUp,
657 O KappaRT,KappaRS,
658 I myThid)
659 # ifdef ALLOW_AUTODIFF_TAMC
660 CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
661 CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte
662 # endif /* ALLOW_AUTODIFF_TAMC */
663 #endif
664
665 iMin = 1-OLx+2
666 iMax = sNx+OLx-1
667 jMin = 1-OLy+2
668 jMax = sNy+OLy-1
669
670 C-- Calculate active tracer tendencies (gT,gS,...)
671 C and step forward storing result in gTnm1, gSnm1, etc.
672 IF ( tempStepping ) THEN
673 CALL CALC_GT(
674 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
675 I xA,yA,uTrans,vTrans,rTrans,maskUp,
676 I KappaRT,
677 U fVerT,
678 I myTime,myIter,myThid)
679 CALL TIMESTEP_TRACER(
680 I bi,bj,iMin,iMax,jMin,jMax,k,tempAdvScheme,
681 I theta, gT,
682 I myIter, myThid)
683 ENDIF
684
685 IF ( saltStepping ) THEN
686 CALL CALC_GS(
687 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
688 I xA,yA,uTrans,vTrans,rTrans,maskUp,
689 I KappaRS,
690 U fVerS,
691 I myTime,myIter,myThid)
692 CALL TIMESTEP_TRACER(
693 I bi,bj,iMin,iMax,jMin,jMax,k,saltAdvScheme,
694 I salt, gS,
695 I myIter, myThid)
696 ENDIF
697 #ifdef ALLOW_PASSIVE_TRACER
698 ceh3 needs an IF ( usePASSIVE_TRACER ) THEN
699 IF ( tr1Stepping ) THEN
700 CALL CALC_GTR1(
701 I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
702 I xA,yA,uTrans,vTrans,rTrans,maskUp,
703 I KappaRT,
704 U fVerTr1,
705 I myTime,myIter,myThid)
706 CALL TIMESTEP_TRACER(
707 I bi,bj,iMin,iMax,jMin,jMax,k,tracerAdvScheme,
708 I Tr1, gTr1,
709 I myIter,myThid)
710 ENDIF
711 #endif
712 #ifdef ALLOW_PTRACERS
713 IF ( usePTRACERS ) THEN
714 CALL PTRACERS_INTEGRATE(
715 I bi,bj,k,
716 I xA,yA,uTrans,vTrans,rTrans,maskUp,
717 X fVerP, KappaRS,
718 I myIter,myTime,myThid)
719 ENDIF
720 #endif /* ALLOW_PTRACERS */
721
722 #ifdef ALLOW_OBCS
723 C-- Apply open boundary conditions
724 IF (useOBCS) THEN
725 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
726 END IF
727 #endif /* ALLOW_OBCS */
728
729 C-- Freeze water
730 C this bit of code is left here for backward compatibility.
731 C freezing at surface level has been moved to FORWARD_STEP
732 IF ( useOldFreezing .AND. .NOT. useSEAICE
733 & .AND. .NOT.(useThSIce.AND.k.EQ.1) ) THEN
734 #ifdef ALLOW_AUTODIFF_TAMC
735 CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
736 CADJ & , key = kkey, byte = isbyte
737 #endif /* ALLOW_AUTODIFF_TAMC */
738 CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid )
739 ENDIF
740
741 C-- end of thermodynamic k loop (Nr:1)
742 ENDDO
743
744
745 C-- Implicit diffusion
746 IF (implicitDiffusion) THEN
747
748 IF (tempStepping) THEN
749 #ifdef ALLOW_AUTODIFF_TAMC
750 CADJ STORE KappaRT(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
751 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
752 #endif /* ALLOW_AUTODIFF_TAMC */
753 CALL IMPLDIFF(
754 I bi, bj, iMin, iMax, jMin, jMax,
755 I deltaTtracer, KappaRT, recip_HFacC,
756 U gT,
757 I myThid )
758 ENDIF
759
760 IF (saltStepping) THEN
761 #ifdef ALLOW_AUTODIFF_TAMC
762 CADJ STORE KappaRS(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
763 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
764 #endif /* ALLOW_AUTODIFF_TAMC */
765 CALL IMPLDIFF(
766 I bi, bj, iMin, iMax, jMin, jMax,
767 I deltaTtracer, KappaRS, recip_HFacC,
768 U gS,
769 I myThid )
770 ENDIF
771
772 #ifdef ALLOW_PASSIVE_TRACER
773 IF (tr1Stepping) THEN
774 #ifdef ALLOW_AUTODIFF_TAMC
775 CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
776 #endif /* ALLOW_AUTODIFF_TAMC */
777 CALL IMPLDIFF(
778 I bi, bj, iMin, iMax, jMin, jMax,
779 I deltaTtracer, KappaRT, recip_HFacC,
780 U gTr1,
781 I myThid )
782 ENDIF
783 #endif
784
785 #ifdef ALLOW_PTRACERS
786 C Vertical diffusion (implicit) for passive tracers
787 IF ( usePTRACERS ) THEN
788 CALL PTRACERS_IMPLDIFF( bi,bj,KappaRS,myThid )
789 ENDIF
790 #endif /* ALLOW_PTRACERS */
791
792 #ifdef ALLOW_OBCS
793 C-- Apply open boundary conditions
794 IF (useOBCS) THEN
795 DO K=1,Nr
796 CALL OBCS_APPLY_TS( bi, bj, k, gT, gS, myThid )
797 ENDDO
798 END IF
799 #endif /* ALLOW_OBCS */
800
801 C-- End If implicitDiffusion
802 ENDIF
803
804 #ifdef ALLOW_TIMEAVE
805 ceh3 needs an IF ( useTIMEAVE ) THEN
806 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
807 CALL TIMEAVE_CUMUL_1T(ConvectCountTave, ConvectCount,
808 I Nr, deltaTclock, bi, bj, myThid)
809 ENDIF
810 useVariableK = useKPP .OR. useGMredi .OR. ivdc_kappa.NE.0.
811 IF (taveFreq.GT.0. .AND. useVariableK ) THEN
812 IF (implicitDiffusion) THEN
813 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, gT, kappaRT,
814 I Nr, 3, deltaTclock, bi, bj, myThid)
815 ELSE
816 CALL TIMEAVE_CUMUL_DIF_1T(TdiffRtave, theta, kappaRT,
817 I Nr, 3, deltaTclock, bi, bj, myThid)
818 ENDIF
819 ENDIF
820 #endif /* ALLOW_TIMEAVE */
821
822 #endif /* SINGLE_LAYER_MODE */
823
824 C-- end bi,bj loops.
825 ENDDO
826 ENDDO
827
828 #ifdef ALLOW_DEBUG
829 If (debugMode) THEN
830 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
831 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (THERMODYNAMICS)',myThid)
832 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (THERMODYNAMICS)',myThid)
833 CALL DEBUG_STATS_RL(Nr,theta,'Theta (THERMODYNAMICS)',myThid)
834 CALL DEBUG_STATS_RL(Nr,salt,'Salt (THERMODYNAMICS)',myThid)
835 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (THERMODYNAMICS)',myThid)
836 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (THERMODYNAMICS)',myThid)
837 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (THERMODYNAMICS)',myThid)
838 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (THERMODYNAMICS)',myThid)
839 #ifdef ALLOW_PTRACERS
840 IF ( usePTRACERS ) THEN
841 CALL PTRACERS_DEBUG(myThid)
842 ENDIF
843 #endif /* ALLOW_PTRACERS */
844 ENDIF
845 #endif
846
847 #ifdef ALLOW_DEBUG
848 IF ( debugLevel .GE. debLevB )
849 & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
850 #endif
851
852 RETURN
853 END

  ViewVC Help
Powered by ViewVC 1.1.22