/[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.45 - (hide annotations) (download)
Mon Sep 22 22:12:05 2003 UTC (20 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.44: +8 -7 lines
initialize tendency to zero : proposed by Martin to fix the bug when
 setting temp/salt advection=.false.

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

  ViewVC Help
Powered by ViewVC 1.1.22