/[MITgcm]/MITgcm/model/src/thermodynamics.F
ViewVC logotype

Annotation of /MITgcm/model/src/thermodynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.42 - (hide annotations) (download)
Fri Jun 27 01:51:10 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51b_post
Changes since 1.41: +28 -11 lines
o disentangled ALLOW_PTRACERS using new ALLOW_GCHEM

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

  ViewVC Help
Powered by ViewVC 1.1.22