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

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

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

revision 1.18 by adcroft, Tue Mar 5 14:15:34 2002 UTC revision 1.42 by heimbach, Fri Jun 27 01:51:10 2003 UTC
# Line 2  C $Header$ Line 2  C $Header$
2  C $Name$  C $Name$
3    
4  #include "CPP_OPTIONS.h"  #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    # ifdef ALLOW_PTRACERS
13    #  include "PTRACERS_OPTIONS.h"
14    # endif
15    #endif /* ALLOW_AUTODIFF_TAMC */
16    
17  CBOP  CBOP
18  C     !ROUTINE: THERMODYNAMICS  C     !ROUTINE: THERMODYNAMICS
# Line 71  C     == Global variables === Line 82  C     == Global variables ===
82  #ifdef ALLOW_PASSIVE_TRACER  #ifdef ALLOW_PASSIVE_TRACER
83  #include "TR1.h"  #include "TR1.h"
84  #endif  #endif
85    #ifdef ALLOW_TIMEAVE
86    #include "TIMEAVE_STATV.h"
87    #endif
88    
89  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
90  # include "tamc.h"  # include "tamc.h"
91  # include "tamc_keys.h"  # include "tamc_keys.h"
92  # include "FFIELDS.h"  # include "FFIELDS.h"
93    # include "EOS.h"
94  # ifdef ALLOW_KPP  # ifdef ALLOW_KPP
95  #  include "KPP.h"  #  include "KPP.h"
96  # endif  # endif
97  # ifdef ALLOW_GMREDI  # ifdef ALLOW_GMREDI
98  #  include "GMREDI.h"  #  include "GMREDI.h"
99  # endif  # endif
100    # 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  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
 #ifdef ALLOW_TIMEAVE  
 #include "TIMEAVE_STATV.h"  
 #endif  
109    
110  C     !INPUT/OUTPUT PARAMETERS:  C     !INPUT/OUTPUT PARAMETERS:
111  C     == Routine arguments ==  C     == Routine arguments ==
# Line 109  C                                      i Line 130  C                                      i
130  C                                      so we need an fVer for each  C                                      so we need an fVer for each
131  C                                      variable.  C                                      variable.
132  C     rhoK, rhoKM1   - Density at current level, and level above  C     rhoK, rhoKM1   - Density at current level, and level above
 C     phiHyd         - Hydrostatic part of the potential phiHydi.  
 C                      In z coords phiHydiHyd is the hydrostatic  
 C                      Potential (=pressure/rho0) anomaly  
 C                      In p coords phiHydiHyd is the geopotential  
 C                      surface height anomaly.  
133  C     phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)  C     phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
134  C     phiSurfY             or geopotentiel (atmos) in X and Y direction  C     phiSurfY             or geopotentiel (atmos) in X and Y direction
135  C     KappaRT,       - Total diffusion in vertical for T and S.  C     KappaRT,       - Total diffusion in vertical for T and S.
136  C     KappaRS          (background + spatially varying, isopycnal term).  C     KappaRS          (background + spatially varying, isopycnal term).
137    C     useVariableK   = T when vertical diffusion is not constant
138  C     iMin, iMax     - Ranges and sub-block indices on which calculations  C     iMin, iMax     - Ranges and sub-block indices on which calculations
139  C     jMin, jMax       are applied.  C     jMin, jMax       are applied.
140  C     bi, bj  C     bi, bj
# Line 133  C                      index into fVerTe Line 150  C                      index into fVerTe
150        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
151        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerS   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
152        _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerTr1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
       _RL phiHyd  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
153        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
# Line 145  C                      index into fVerTe Line 161  C                      index into fVerTe
161        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
162  C     This is currently used by IVDC and Diagnostics  C     This is currently used by IVDC and Diagnostics
163        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
164          LOGICAL useVariableK
165        INTEGER iMin, iMax        INTEGER iMin, iMax
166        INTEGER jMin, jMax        INTEGER jMin, jMax
167        INTEGER bi, bj        INTEGER bi, bj
168        INTEGER i, j        INTEGER i, j
169        INTEGER k, km1, kup, kDown        INTEGER k, km1, kup, kDown
170          INTEGER iTracer
171    
172  CEOP  CEOP
173    
174    #ifndef DISABLE_DEBUGMODE
175             IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
176    #endif
177    
178  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
179  C--   dummy statement to end declaration part  C--   dummy statement to end declaration part
180        ikey = 1        ikey = 1
181          itdkey = 1
182  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
183    
 C--   Set up work arrays with valid (i.e. not NaN) values  
 C     These inital values do not alter the numerical results. They  
 C     just ensure that all memory references are to valid floating  
 C     point numbers. This prevents spurious hardware signals due to  
 C     uninitialised but inert locations.  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         xA(i,j)      = 0. _d 0  
         yA(i,j)      = 0. _d 0  
         uTrans(i,j)  = 0. _d 0  
         vTrans(i,j)  = 0. _d 0  
         DO k=1,Nr  
          phiHyd(i,j,k)  = 0. _d 0  
          sigmaX(i,j,k) = 0. _d 0  
          sigmaY(i,j,k) = 0. _d 0  
          sigmaR(i,j,k) = 0. _d 0  
         ENDDO  
         rhoKM1 (i,j) = 0. _d 0  
         rhok   (i,j) = 0. _d 0  
         phiSurfX(i,j) = 0. _d 0  
         phiSurfY(i,j) = 0. _d 0  
        ENDDO  
       ENDDO  
   
   
184  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
185  C--   HPF directive to help TAMC  C--   HPF directive to help TAMC
186  CHPF$ INDEPENDENT  CHPF$ INDEPENDENT
# Line 193  CHPF$ INDEPENDENT Line 191  CHPF$ INDEPENDENT
191  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
192  C--    HPF directive to help TAMC  C--    HPF directive to help TAMC
193  CHPF$  INDEPENDENT, NEW (rTrans,fVerT,fVerS  CHPF$  INDEPENDENT, NEW (rTrans,fVerT,fVerS
194  CHPF$&                  ,phiHyd,utrans,vtrans,xA,yA  CHPF$&                  ,utrans,vtrans,xA,yA
195  CHPF$&                  ,KappaRT,KappaRS  CHPF$&                  ,KappaRT,KappaRS
196  CHPF$&                  )  CHPF$&                  )
197  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
# Line 208  CHPF$&                  ) Line 206  CHPF$&                  )
206            act3 = myThid - 1            act3 = myThid - 1
207            max3 = nTx*nTy            max3 = nTx*nTy
208            act4 = ikey_dynamics - 1            act4 = ikey_dynamics - 1
209            ikey = (act1 + 1) + act2*max1            itdkey = (act1 + 1) + act2*max1
210       &                      + act3*max1*max2       &                      + act3*max1*max2
211       &                      + act4*max1*max2*max3       &                      + act4*max1*max2*max3
212  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
213    
214  C--     Set up work arrays that need valid initial values  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          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
221           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
222              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            rTrans (i,j)   = 0. _d 0            rTrans (i,j)   = 0. _d 0
231            fVerT  (i,j,1) = 0. _d 0            fVerT  (i,j,1) = 0. _d 0
232            fVerT  (i,j,2) = 0. _d 0            fVerT  (i,j,2) = 0. _d 0
# Line 230  C--     Set up work arrays that need val Line 241  C--     Set up work arrays that need val
241           DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
242            DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
243  C This is currently also used by IVDC and Diagnostics  C This is currently also used by IVDC and Diagnostics
244               sigmaX(i,j,k) = 0. _d 0
245               sigmaY(i,j,k) = 0. _d 0
246               sigmaR(i,j,k) = 0. _d 0
247             ConvectCount(i,j,k) = 0.             ConvectCount(i,j,k) = 0.
248             KappaRT(i,j,k) = 0. _d 0             KappaRT(i,j,k)    = 0. _d 0
249             KappaRS(i,j,k) = 0. _d 0             KappaRS(i,j,k)    = 0. _d 0
250  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
251             gT(i,j,k,bi,bj) = 0. _d 0  cph all the following init. are necessary for TAF
252             gS(i,j,k,bi,bj) = 0. _d 0  cph although some of these are re-initialised later.
253  #ifdef ALLOW_PASSIVE_TRACER             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             gTr1(i,j,k,bi,bj) = 0. _d 0             gTr1(i,j,k,bi,bj) = 0. _d 0
257  #endif  # endif
258  #endif  # 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    # 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    #  endif
279    #  ifdef GM_VISBECK_VARIABLE_K
280               VisbeckK(i,j,bi,bj)   = 0. _d 0
281    #  endif
282    # endif /* ALLOW_GMREDI */
283    #endif /* ALLOW_AUTODIFF_TAMC */
284            ENDDO            ENDDO
285           ENDDO           ENDDO
286          ENDDO          ENDDO
287    
288          iMin = 1-OLx+1          iMin = 1-OLx
289          iMax = sNx+OLx          iMax = sNx+OLx
290          jMin = 1-OLy+1          jMin = 1-OLy
291          jMax = sNy+OLy          jMax = sNy+OLy
292    
   
293  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
294  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
295  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
296    CADJ STORE totphihyd
297    CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
298  #ifdef ALLOW_KPP  #ifdef ALLOW_KPP
299  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
300  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
301  #endif  #endif
302  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
303    
304    #ifndef DISABLE_DEBUGMODE
305            IF (debugMode) CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
306    #endif
307    
308  C--     Start of diagnostic loop  C--     Start of diagnostic loop
309          DO k=Nr,1,-1          DO k=Nr,1,-1
310    
# Line 267  C? Patrick, is this formula correct now Line 313  C? Patrick, is this formula correct now
313  C? Do we still need this?  C? Do we still need this?
314  cph kkey formula corrected.  cph kkey formula corrected.
315  cph Needed for rhok, rhokm1, in the case useGMREDI.  cph Needed for rhok, rhokm1, in the case useGMREDI.
316           kkey = (ikey-1)*Nr + k           kkey = (itdkey-1)*Nr + k
 CADJ STORE rhokm1(:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
 CADJ STORE rhok  (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
317  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
318    
319  C--       Integrate continuity vertically for vertical velocity  C--       Integrate continuity vertically for vertical velocity
320            CALL INTEGRATE_FOR_W(  c         CALL INTEGRATE_FOR_W(
321       I                         bi, bj, k, uVel, vVel,  c    I                         bi, bj, k, uVel, vVel,
322       O                         wVel,  c    O                         wVel,
323       I                         myThid )  c    I                         myThid )
324    
325  #ifdef    ALLOW_OBCS  #ifdef    ALLOW_OBCS
326  #ifdef    ALLOW_NONHYDROSTATIC  #ifdef    ALLOW_NONHYDROSTATIC
327  C--       Apply OBC to W if in N-H mode  C--       Apply OBC to W if in N-H mode
328            IF (useOBCS.AND.nonHydrostatic) THEN  c         IF (useOBCS.AND.nonHydrostatic) THEN
329              CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )  c           CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid )
330            ENDIF  c         ENDIF
331  #endif    /* ALLOW_NONHYDROSTATIC */  #endif    /* ALLOW_NONHYDROSTATIC */
332  #endif    /* ALLOW_OBCS */  #endif    /* ALLOW_OBCS */
333    
334    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  C--       Calculate gradients of potential density for isoneutral  C--       Calculate gradients of potential density for isoneutral
339  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)  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  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            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN
342    #ifndef DISABLE_DEBUGMODE
343                IF (debugMode) CALL DEBUG_CALL('FIND_RHO',myThid)
344    #endif
345  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
346  CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte  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  CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
348  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
349              CALL FIND_RHO(              CALL FIND_RHO(
350       I        bi, bj, iMin, iMax, jMin, jMax, k, k, eosType,       I        bi, bj, iMin, iMax, jMin, jMax, k, k,
351       I        theta, salt,       I        theta, salt,
352       O        rhoK,       O        rhoK,
353       I        myThid )       I        myThid )
354    
355              IF (k.GT.1) THEN              IF (k.GT.1) THEN
356  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
357  CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte  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  CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
359  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
360               CALL FIND_RHO(               CALL FIND_RHO(
361       I        bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType,       I        bi, bj, iMin, iMax, jMin, jMax, k-1, k,
362       I        theta, salt,       I        theta, salt,
363       O        rhoKm1,       O        rhoKm1,
364       I        myThid )       I        myThid )
365              ENDIF              ENDIF
366    #ifndef DISABLE_DEBUGMODE
367                IF (debugMode) CALL DEBUG_CALL('GRAD_SIGMA',myThid)
368    #endif
369              CALL GRAD_SIGMA(              CALL GRAD_SIGMA(
370       I             bi, bj, iMin, iMax, jMin, jMax, k,       I             bi, bj, iMin, iMax, jMin, jMax, k,
371       I             rhoK, rhoKm1, rhoK,       I             rhoK, rhoKm1, rhoK,
# Line 318  CADJ STORE salt (:,:,k-1,bi,bj) = comlev Line 373  CADJ STORE salt (:,:,k-1,bi,bj) = comlev
373       I             myThid )       I             myThid )
374            ENDIF            ENDIF
375    
376    #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  C--       Implicit Vertical Diffusion for Convection  C--       Implicit Vertical Diffusion for Convection
381  c ==> should use sigmaR !!!  c ==> should use sigmaR !!!
382            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
383    #ifndef DISABLE_DEBUGMODE
384                IF (debugMode) CALL DEBUG_CALL('CALC_IVDC',myThid)
385    #endif
386              CALL CALC_IVDC(              CALL CALC_IVDC(
387       I        bi, bj, iMin, iMax, jMin, jMax, k,       I        bi, bj, iMin, iMax, jMin, jMax, k,
388       I        rhoKm1, rhoK,       I        rhoKm1, rhoK,
# Line 328  c ==> should use sigmaR !!! Line 390  c ==> should use sigmaR !!!
390       I        myTime, myIter, myThid)       I        myTime, myIter, myThid)
391            ENDIF            ENDIF
392    
393    #endif /* SINGLE_LAYER_MODE */
394    
395  C--     end of diagnostic k loop (Nr:1)  C--     end of diagnostic k loop (Nr:1)
396          ENDDO          ENDDO
397    
398  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
399  cph avoids recomputation of integrate_for_w  cph avoids recomputation of integrate_for_w
400  CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
401  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
402    
403  #ifdef  ALLOW_OBCS  #ifdef  ALLOW_OBCS
404  C--     Calculate future values on open boundaries  C--     Calculate future values on open boundaries
405          IF (useOBCS) THEN          IF (useOBCS) THEN
406    #ifndef DISABLE_DEBUGMODE
407              IF (debugMode) CALL DEBUG_CALL('OBCS_CALC',myThid)
408    #endif
409            CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,            CALL OBCS_CALC( bi, bj, myTime+deltaT, myIter+1,
410       I            uVel, vVel, wVel, theta, salt,       I            uVel, vVel, wVel, theta, salt,
411       I            myThid )       I            myThid )
412          ENDIF          ENDIF
413  #endif  /* ALLOW_OBCS */  #endif  /* ALLOW_OBCS */
414    
415    
416    c********************************************
417    cswdice --- add ---
418    #ifdef ALLOW_THERM_SEAICE
419    #ifndef DISABLE_DEBUGMODE
420            IF (debugMode) CALL DEBUG_CALL('ICE_FORCING',myThid)
421    #endif
422    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  C--     Determines forcing terms based on external fields  C--     Determines forcing terms based on external fields
432  C       relaxation terms, etc.  C       relaxation terms, etc.
433          CALL EXTERNAL_FORCING_SURF(  #ifndef DISABLE_DEBUGMODE
434            IF (debugMode) CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
435    #endif
436            CALL EXTERNAL_FORCING_SURF(
437       I             bi, bj, iMin, iMax, jMin, jMax,       I             bi, bj, iMin, iMax, jMin, jMax,
438       I             myThid )       I             myThid )
439    cswdice --- add ----
440    #endif
441    cswdice --- end add ---
442    c******************************************
443    
444    
445    
446    
447  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
448  cph needed for KPP  cph needed for KPP
449  CADJ STORE surfacetendencyU(:,:,bi,bj)  CADJ STORE surfacetendencyU(:,:,bi,bj)
450  CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
451  CADJ STORE surfacetendencyV(:,:,bi,bj)  CADJ STORE surfacetendencyV(:,:,bi,bj)
452  CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
453  CADJ STORE surfacetendencyS(:,:,bi,bj)  CADJ STORE surfacetendencyS(:,:,bi,bj)
454  CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
455  CADJ STORE surfacetendencyT(:,:,bi,bj)  CADJ STORE surfacetendencyT(:,:,bi,bj)
456  CADJ &     = comlev1_bibj, key=ikey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
457  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
458    
459    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  #ifdef  ALLOW_GMREDI  #ifdef  ALLOW_GMREDI
464    
465  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
466  CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte  cph storing here is needed only for one GMREDI_OPTIONS:
467  CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte  cph define GM_BOLUS_ADVEC
468  CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte  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  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
473    
474  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
475          IF (useGMRedi) THEN          IF (useGMRedi) THEN
476    #ifndef DISABLE_DEBUGMODE
477              IF (debugMode) CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
478    #endif
479            CALL GMREDI_CALC_TENSOR(            CALL GMREDI_CALC_TENSOR(
480       I             bi, bj, iMin, iMax, jMin, jMax,       I             bi, bj, iMin, iMax, jMin, jMax,
481       I             sigmaX, sigmaY, sigmaR,       I             sigmaX, sigmaY, sigmaR,
# Line 385  C--     Calculate iso-neutral slopes for Line 490  C--     Calculate iso-neutral slopes for
490          ENDIF          ENDIF
491    
492  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
493  CADJ STORE Kwx(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE Kwx(:,:,:,bi,bj)   = comlev1_bibj, key=itdkey, byte=isbyte
494  CADJ STORE Kwy(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE Kwy(:,:,:,bi,bj)   = comlev1_bibj, key=itdkey, byte=isbyte
495  CADJ STORE Kwz(:,:,:,bi,bj)   = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE Kwz(:,:,:,bi,bj)   = comlev1_bibj, key=itdkey, byte=isbyte
496  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
497    
498  #endif  /* ALLOW_GMREDI */  #endif  /* ALLOW_GMREDI */
# Line 395  CADJ STORE Kwz(:,:,:,bi,bj)   = comlev1_ Line 500  CADJ STORE Kwz(:,:,:,bi,bj)   = comlev1_
500  #ifdef  ALLOW_KPP  #ifdef  ALLOW_KPP
501  C--     Compute KPP mixing coefficients  C--     Compute KPP mixing coefficients
502          IF (useKPP) THEN          IF (useKPP) THEN
503    #ifndef DISABLE_DEBUGMODE
504              IF (debugMode) CALL DEBUG_CALL('KPP_CALC',myThid)
505    #endif
506            CALL KPP_CALC(            CALL KPP_CALC(
507       I                  bi, bj, myTime, myThid )       I                  bi, bj, myTime, myThid )
508  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 409  CADJ STORE KPPghat   (:,:,:,bi,bj) Line 517  CADJ STORE KPPghat   (:,:,:,bi,bj)
517  CADJ &   , KPPdiffKzT(:,:,:,bi,bj)  CADJ &   , KPPdiffKzT(:,:,:,bi,bj)
518  CADJ &   , KPPdiffKzS(:,:,:,bi,bj)  CADJ &   , KPPdiffKzS(:,:,:,bi,bj)
519  CADJ &   , KPPfrac   (:,:  ,bi,bj)  CADJ &   , KPPfrac   (:,:  ,bi,bj)
520  CADJ &                 = comlev1_bibj, key=ikey, byte=isbyte  CADJ &                 = comlev1_bibj, key=itdkey, byte=isbyte
521  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
522    
523  #endif  /* ALLOW_KPP */  #endif  /* ALLOW_KPP */
524    
525  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
526  CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE KappaRT(:,:,:)     = comlev1_bibj, key=itdkey, byte=isbyte
527  CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE KappaRS(:,:,:)     = comlev1_bibj, key=itdkey, byte=isbyte
528  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
529  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
530  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
531  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
532  #ifdef ALLOW_PASSIVE_TRACER  #ifdef ALLOW_PASSIVE_TRACER
533  CADJ STORE tr1  (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte  CADJ STORE tr1  (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
534    #endif
535    #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  #endif
540  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
541    
542  #ifdef ALLOW_AIM  #ifdef ALLOW_AIM
543  C       AIM - atmospheric intermediate model, physics package code.  C       AIM - atmospheric intermediate model, physics package code.
 C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics  
544          IF ( useAIM ) THEN          IF ( useAIM ) THEN
545           CALL TIMER_START('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)  #ifndef DISABLE_DEBUGMODE
546           CALL AIM_DO_ATMOS_PHYSICS( phiHyd, bi, bj, myTime, myThid )            IF (debugMode) CALL DEBUG_CALL('AIM_DO_PHYSICS',myThid)
547           CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS      [DYNAMICS]', myThid)  #endif
548             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          ENDIF          ENDIF
552  #endif /* ALLOW_AIM */  #endif /* ALLOW_AIM */
553    
 #ifdef ALLOW_TIMEAVE  
         IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN  
           CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,  
      I                               deltaTclock, bi, bj, myThid)  
         ENDIF  
 #endif /* ALLOW_TIMEAVE */  
   
554  #ifndef DISABLE_MULTIDIM_ADVECTION  #ifndef DISABLE_MULTIDIM_ADVECTION
555  C--     Some advection schemes are better calculated using a multi-dimensional  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.  C       method in the absence of any other terms and, if used, is done here.
# Line 454  C to be able to exclude this scheme to a Line 562  C to be able to exclude this scheme to a
562  C recomputation. It *is* differentiable, if you need it.  C recomputation. It *is* differentiable, if you need it.
563  C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to  C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
564  C disable this section of code.  C disable this section of code.
565          IF (multiDimAdvection) THEN          IF (tempMultiDimAdvec) THEN
566           IF (tempStepping .AND.  #ifndef DISABLE_DEBUGMODE
567       &       tempAdvScheme.NE.ENUM_CENTERED_2ND .AND.            IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
568       &       tempAdvScheme.NE.ENUM_UPWIND_3RD .AND.  #endif
      &       tempAdvScheme.NE.ENUM_CENTERED_4TH ) THEN  
569            CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,            CALL GAD_ADVECTION(bi,bj,tempAdvScheme,GAD_TEMPERATURE,
570       U                      theta,gT,       U                      theta,gT,
571       I                      myTime,myIter,myThid)       I                      myTime,myIter,myThid)
572           ENDIF          ENDIF
573           IF (saltStepping .AND.          IF (saltMultiDimAdvec) THEN
574       &       saltAdvScheme.NE.ENUM_CENTERED_2ND .AND.  #ifndef DISABLE_DEBUGMODE
575       &       saltAdvScheme.NE.ENUM_UPWIND_3RD .AND.            IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
576       &       saltAdvScheme.NE.ENUM_CENTERED_4TH ) THEN  #endif
577            CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,            CALL GAD_ADVECTION(bi,bj,saltAdvScheme,GAD_SALINITY,
578       U                      salt,gS,       U                      salt,gS,
579       I                      myTime,myIter,myThid)       I                      myTime,myIter,myThid)
          ENDIF  
580          ENDIF          ENDIF
581  C Since passive tracers are configurable separately from T,S we  C Since passive tracers are configurable separately from T,S we
582  C call the multi-dimensional method for PTRACERS regardless  C call the multi-dimensional method for PTRACERS regardless
583  C of whether multiDimAdvection is set or not.  C of whether multiDimAdvection is set or not.
584  #ifdef ALLOW_PTRACERS  #ifdef ALLOW_PTRACERS
585          IF ( usePTRACERS ) THEN          IF ( usePTRACERS ) THEN
586    #ifndef DISABLE_DEBUGMODE
587              IF (debugMode) CALL DEBUG_CALL('PTRACERS_ADVECTION',myThid)
588    #endif
589           CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )           CALL PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
590          ENDIF          ENDIF
591  #endif /* ALLOW_PTRACERS */  #endif /* ALLOW_PTRACERS */
592  #endif /* DISABLE_MULTIDIM_ADVECTION */  #endif /* DISABLE_MULTIDIM_ADVECTION */
593    
594    #ifndef DISABLE_DEBUGMODE
595           IF (debugMode) CALL DEBUG_MSG('ENTERING DOWNWARD K LOOP',myThid)
596    #endif
597    
598  C--     Start of thermodynamics loop  C--     Start of thermodynamics loop
599          DO k=Nr,1,-1          DO k=Nr,1,-1
600  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
601  C? Patrick Is this formula correct?  C? Patrick Is this formula correct?
602  cph Yes, but I rewrote it.  cph Yes, but I rewrote it.
603  cph Also, the KappaR? need the index and subscript k!  cph Also, the KappaR? need the index and subscript k!
604           kkey = (ikey-1)*Nr + k           kkey = (itdkey-1)*Nr + k
605  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
606    
607  C--       km1    Points to level above k (=k-1)  C--       km1    Points to level above k (=k-1)
# Line 510  C--       Get temporary terms used by te Line 623  C--       Get temporary terms used by te
623       O         xA,yA,uTrans,vTrans,rTrans,maskUp,       O         xA,yA,uTrans,vTrans,rTrans,maskUp,
624       I         myThid)       I         myThid)
625    
626    #ifdef ALLOW_GMREDI
627    
628    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    
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    #endif /* ALLOW_GMREDI */
645    
646  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
647  CADJ STORE KappaRT(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte  CADJ STORE KappaRT(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte
648  CADJ STORE KappaRS(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte  CADJ STORE KappaRS(:,:,k)    = comlev1_bibj_k, key=kkey, byte=isbyte
# Line 543  C        and step forward storing result Line 676  C        and step forward storing result
676       I         theta, gT,       I         theta, gT,
677       I         myIter, myThid)       I         myIter, myThid)
678           ENDIF           ENDIF
679    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           IF ( saltStepping ) THEN           IF ( saltStepping ) THEN
687             CALL CALC_GS(             CALL CALC_GS(
688       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown,
# Line 571  C        and step forward storing result Line 711  C        and step forward storing result
711  #endif  #endif
712  #ifdef ALLOW_PTRACERS  #ifdef ALLOW_PTRACERS
713           IF ( usePTRACERS ) THEN           IF ( usePTRACERS ) THEN
714             CALL PTRACERS_INTEGERATE(             CALL PTRACERS_INTEGRATE(
715       I         bi,bj,k,       I         bi,bj,k,
716       I         xA,yA,uTrans,vTrans,rTrans,maskUp,       I         xA,yA,uTrans,vTrans,rTrans,maskUp,
717       X         KappaRS,       X         KappaRS,
# Line 587  C--      Apply open boundary conditions Line 727  C--      Apply open boundary conditions
727  #endif   /* ALLOW_OBCS */  #endif   /* ALLOW_OBCS */
728    
729  C--      Freeze water  C--      Freeze water
730           IF (allowFreezing) THEN           IF ( allowFreezing .AND. .NOT. useSEAICE ) THEN
731  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
732  CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k  CADJ STORE gT(:,:,k,bi,bj) = comlev1_bibj_k
733  CADJ &   , key = kkey, byte = isbyte  CADJ &   , key = kkey, byte = isbyte
# Line 598  CADJ &   , key = kkey, byte = isbyte Line 738  CADJ &   , key = kkey, byte = isbyte
738  C--     end of thermodynamic k loop (Nr:1)  C--     end of thermodynamic k loop (Nr:1)
739          ENDDO          ENDDO
740    
741    cswdice -- add ---
742    #ifdef ALLOW_THERM_SEAICE
743    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    
 #ifdef ALLOW_AUTODIFF_TAMC  
 C? Patrick? What about this one?  
 cph Keys iikey and idkey dont seem to be needed  
 cph since storing occurs on different tape for each  
 cph impldiff call anyways.  
 cph Thus, common block comlev1_impl isnt needed either.  
 cph Storing below needed in the case useGMREDI.  
         iikey = (ikey-1)*maximpl  
 #endif /* ALLOW_AUTODIFF_TAMC */  
750    
751  C--     Implicit diffusion  C--     Implicit diffusion
752          IF (implicitDiffusion) THEN          IF (implicitDiffusion) THEN
753    
754           IF (tempStepping) THEN           IF (tempStepping) THEN
755  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
756              idkey = iikey + 1  CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
 CADJ STORE gT(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  
757  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
758              CALL IMPLDIFF(              CALL IMPLDIFF(
759       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
# Line 626  CADJ STORE gT(:,:,:,bi,bj) = comlev1_bib Line 764  CADJ STORE gT(:,:,:,bi,bj) = comlev1_bib
764    
765           IF (saltStepping) THEN           IF (saltStepping) THEN
766  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
767           idkey = iikey + 2  CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
 CADJ STORE gS(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  
768  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
769              CALL IMPLDIFF(              CALL IMPLDIFF(
770       I         bi, bj, iMin, iMax, jMin, jMax,       I         bi, bj, iMin, iMax, jMin, jMax,
# Line 639  CADJ STORE gS(:,:,:,bi,bj) = comlev1_bib Line 776  CADJ STORE gS(:,:,:,bi,bj) = comlev1_bib
776  #ifdef ALLOW_PASSIVE_TRACER  #ifdef ALLOW_PASSIVE_TRACER
777           IF (tr1Stepping) THEN           IF (tr1Stepping) THEN
778  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
779  CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte  CADJ STORE gTr1(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
780  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
781            CALL IMPLDIFF(            CALL IMPLDIFF(
782       I      bi, bj, iMin, iMax, jMin, jMax,       I      bi, bj, iMin, iMax, jMin, jMax,
# Line 668  C--      Apply open boundary conditions Line 805  C--      Apply open boundary conditions
805  C--     End If implicitDiffusion  C--     End If implicitDiffusion
806          ENDIF          ENDIF
807    
808  Ccs-  #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    #endif /* SINGLE_LAYER_MODE */
826    
827    C--   end bi,bj loops.
828         ENDDO         ENDDO
829        ENDDO        ENDDO
830    
 #ifdef ALLOW_AIM  
       IF ( useAIM ) THEN  
        CALL AIM_AIM2DYN_EXCHANGES( myTime, myThid )  
       ENDIF  
        _EXCH_XYZ_R8(gT,myThid)  
        _EXCH_XYZ_R8(gS,myThid)  
 #else  
       IF (staggerTimeStep.AND.useCubedSphereExchange) THEN  
        _EXCH_XYZ_R8(gT,myThid)  
        _EXCH_XYZ_R8(gS,myThid)  
       ENDIF  
 #endif /* ALLOW_AIM */  
   
831  #ifndef DISABLE_DEBUGMODE  #ifndef DISABLE_DEBUGMODE
832        If (debugMode) THEN        If (debugMode) THEN
833         CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)         CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (THERMODYNAMICS)',myThid)
# Line 704  Ccs- Line 847  Ccs-
847        ENDIF        ENDIF
848  #endif  #endif
849    
850    #ifndef DISABLE_DEBUGMODE
851             IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
852    #endif
853    
854        RETURN        RETURN
855        END        END

Legend:
Removed from v.1.18  
changed lines
  Added in v.1.42

  ViewVC Help
Powered by ViewVC 1.1.22