/[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.38 - (show annotations) (download)
Fri Feb 28 02:20:52 2003 UTC (21 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint50c_pre, checkpoint48i_post, checkpoint50, checkpoint50b_pre, checkpoint50a_post, checkpoint49, checkpoint50b_post
Changes since 1.37: +3 -7 lines
Changes to restore differentiability of code w.r.t. previous tag
(mostly adding new routines to make list and replacing
pressure by totPhiHyd).

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

  ViewVC Help
Powered by ViewVC 1.1.22