Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-09 06:09:16 UTC

view on githubraw file Latest commit 84f265d4 on 2018-03-08 18:49:29 UTC
108a00eab9 Ryan*0001 #include "LAYERS_OPTIONS.h"
910fe2b443 Jean*0002 #ifdef ALLOW_GMREDI
a60f4dd950 Davi*0003 #include "GMREDI_OPTIONS.h"
910fe2b443 Jean*0004 #endif
108a00eab9 Ryan*0005 
4e07ab4926 Jean*0006 CBOP 0
                0007 C !ROUTINE: LAYERS_CALC
108a00eab9 Ryan*0008 
4e07ab4926 Jean*0009 C !INTERFACE:
108a00eab9 Ryan*0010       SUBROUTINE LAYERS_CALC(
62b1f3e57f Jean*0011      I                        myTime, myIter, myThid )
108a00eab9 Ryan*0012 
4e07ab4926 Jean*0013 C !DESCRIPTION:
108a00eab9 Ryan*0014 C ===================================================================
                0015 C     Calculate the transport in isopycnal layers.
aa668e01f1 Gael*0016 C     This was the meat of the LAYERS package, which
                0017 C     has been moved to S/R LAYERS_FLUXCALC.F
108a00eab9 Ryan*0018 C ===================================================================
                0019 
4e07ab4926 Jean*0020 C !USES:
108a00eab9 Ryan*0021       IMPLICIT NONE
                0022 #include "SIZE.h"
                0023 #include "EEPARAMS.h"
                0024 #include "PARAMS.h"
4e07ab4926 Jean*0025 #include "GRID.h"
                0026 #include "DYNVARS.h"
108a00eab9 Ryan*0027 #include "LAYERS_SIZE.h"
                0028 #include "LAYERS.h"
910fe2b443 Jean*0029 #ifdef ALLOW_GMREDI
                0030 # include "GMREDI.h"
a60f4dd950 Davi*0031 #endif
108a00eab9 Ryan*0032 
4e07ab4926 Jean*0033 C !INPUT PARAMETERS:
                0034 C     myTime :: Current time in simulation
                0035 C     myIter :: Current iteration number
                0036 C     myThid :: my Thread Id number
                0037       _RL     myTime
                0038       INTEGER myIter
                0039       INTEGER myThid
                0040 CEOP
108a00eab9 Ryan*0041 
                0042 #ifdef ALLOW_LAYERS
1f346f78e1 Jean*0043 C !FUNCTIONS:
                0044       LOGICAL  DIFFERENT_MULTIPLE
                0045       EXTERNAL DIFFERENT_MULTIPLE
108a00eab9 Ryan*0046 
4e07ab4926 Jean*0047 C !LOCAL VARIABLES:
62b1f3e57f Jean*0048 C --  3D Layers fields. The vertical dimension in these fields is Nlayers,
                0049 C     i.e. the isopycnal coordinate.
                0050 C      layers_UH  :: U integrated over layer (m^2/s)
                0051 C      layers_VH  :: V integrated over layer (m^2/s)
                0052 C      layers_Hw  :: Layer thickness at the U point (m)
                0053 C      layers_Hs  :: Layer thickness at the V point (m)
                0054 C      layers_PIw :: 1 if layer exists, 0 otherwise
                0055 C      layers_PIs :: 1 if layer exists, 0 otherwise
                0056 C      layers_U   :: mean zonal velocity in layer (only if layer exists) (m/s)
                0057 C      layers_V   :: mean meridional velocity in layer (only if layer exists) (m/s)
                0058 #ifdef LAYERS_UFLUX
                0059       _RL layers_UH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0060 # ifdef LAYERS_THICKNESS
                0061       _RL layers_Hw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0062       _RL layers_PIw(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0063       _RL layers_U  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0064 # endif /* LAYERS_THICKNESS */
                0065 #endif /* LAYERS_UFLUX */
                0066 #ifdef LAYERS_VFLUX
                0067       _RL layers_VH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0068 # ifdef LAYERS_THICKNESS
                0069       _RL layers_Hs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0070       _RL layers_PIs(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0071       _RL layers_V  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
                0072 # endif /* LAYERS_THICKNESS */
                0073 #endif /* LAYERS_VFLUX */
84f265d4b3 dfer*0074 #if (defined LAYERS_PRHO_REF) || (defined LAYERS_MSE)
62b1f3e57f Jean*0075       _RL prho(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
84f265d4b3 dfer*0076 #endif
                0077 #ifdef LAYERS_PRHO_REF
62b1f3e57f Jean*0078       _RL rhoShift
84f265d4b3 dfer*0079 #endif
                0080 #ifdef LAYERS_MSE
                0081       _RL conv_theta2T
62b1f3e57f Jean*0082 #endif
                0083 C --  other local variables:
4e07ab4926 Jean*0084 C     bi, bj   :: tile indices
108a00eab9 Ryan*0085 C     i,j      :: horizontal indices
62b1f3e57f Jean*0086 C     iLa      :: layer-type index
108a00eab9 Ryan*0087 C     k        :: vertical index for model grid
d21fbd4209 Jean*0088       INTEGER bi, bj, iLa
df5a9764ba Jean*0089       CHARACTER*(10) sufx
                0090       CHARACTER*(13) suff
84f265d4b3 dfer*0091       INTEGER i, j, k
50d8304171 Ryan*0092 
d21fbd4209 Jean*0093 #ifdef LAYERS_THERMODYNAMICS
                0094       INTEGER iTracer
                0095 #endif
1f346f78e1 Jean*0096 #ifdef ALLOW_DIAGNOSTICS
                0097       CHARACTER*8    diagName
                0098 #endif
                0099 c#ifdef ALLOW_MNC
                0100 c      CHARACTER*(1) pf
                0101 c#endif
108a00eab9 Ryan*0102 
05eacaf171 Jean*0103 #ifndef LAYERS_UFLUX
1f346f78e1 Jean*0104       _RL layers_UH(1)
05eacaf171 Jean*0105 #endif
                0106 #ifndef LAYERS_VFLUX
1f346f78e1 Jean*0107       _RL layers_VH(1)
05eacaf171 Jean*0108 #endif
                0109 #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_UFLUX)
1f346f78e1 Jean*0110       _RL layers_Hw(1), layers_PIw(1), layers_U(1)
05eacaf171 Jean*0111 #endif
                0112 #if !(defined LAYERS_THICKNESS) || !(defined LAYERS_VFLUX)
1f346f78e1 Jean*0113       _RL layers_Hs(1), layers_PIs(1), layers_V(1)
05eacaf171 Jean*0114 #endif
                0115 
4e07ab4926 Jean*0116 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0117 
1f346f78e1 Jean*0118       IF ( myIter.EQ.nIter0 ) RETURN
                0119 
cf336ab6c5 Ryan*0120 #ifdef LAYERS_THERMODYNAMICS
                0121       CALL LAYERS_CALC_RHS(myThid)
                0122 #endif
                0123 
406891c1a3 Gael*0124       DO iLa=1,layers_maxNum
77bf6311ca Jean*0125 
                0126        IF ( layers_num(iLa) .EQ. 1 ) THEN
                0127         CALL LAYERS_FLUXCALC( uVel,vVel,theta,iLa,
1f346f78e1 Jean*0128      &              layers_UH, layers_VH,
                0129      &              layers_Hw, layers_Hs,
                0130      &              layers_PIw,layers_PIs,
                0131      &              layers_U,  layers_V,
76c15ea089 Jean*0132      &              myThid )
cf336ab6c5 Ryan*0133 #ifdef LAYERS_THERMODYNAMICS
                0134         CALL LAYERS_DIAPYCNAL( theta,iLa,
50d8304171 Ryan*0135      &              layers_TtendSurf,
                0136      &              layers_TtendDiffh, layers_TtendDiffr,
                0137      &              layers_TtendAdvh, layers_TtendAdvr,
                0138      &              layers_Ttendtot,
                0139      &              layers_StendSurf,
cf336ab6c5 Ryan*0140      &              layers_StendDiffh, layers_StendDiffr,
50d8304171 Ryan*0141      &              layers_StendAdvh, layers_StendAdvr,
                0142      &              layers_Stendtot,
cf336ab6c5 Ryan*0143      &              layers_Hc, layers_PIc,
                0144      &              myThid)
                0145 #endif
77bf6311ca Jean*0146        ELSEIF ( layers_num(iLa) .EQ. 2 ) THEN
                0147         CALL LAYERS_FLUXCALC( uVel,vVel,salt,iLa,
1f346f78e1 Jean*0148      &              layers_UH, layers_VH,
                0149      &              layers_Hw, layers_Hs,
                0150      &              layers_PIw,layers_PIs,
                0151      &              layers_U,  layers_V,
77bf6311ca Jean*0152      &              myThid )
cf336ab6c5 Ryan*0153 #ifdef LAYERS_THERMODYNAMICS
                0154         CALL LAYERS_DIAPYCNAL( salt,iLa,
50d8304171 Ryan*0155      &              layers_TtendSurf,
                0156      &              layers_TtendDiffh, layers_TtendDiffr,
                0157      &              layers_TtendAdvh, layers_TtendAdvr,
                0158      &              layers_Ttendtot,
                0159      &              layers_StendSurf,
cf336ab6c5 Ryan*0160      &              layers_StendDiffh, layers_StendDiffr,
50d8304171 Ryan*0161      &              layers_StendAdvh, layers_StendAdvr,
                0162      &              layers_Stendtot,
cf336ab6c5 Ryan*0163      &              layers_Hc, layers_PIc,
                0164      &              myThid)
                0165 #endif
77bf6311ca Jean*0166        ELSEIF ( layers_num(iLa) .EQ. 3 ) THEN
17ce8d85dd Davi*0167 #ifdef LAYERS_PRHO_REF
9643385778 Jean*0168 C     For layers_num(iLa) = 3, calculate the potential density (minus 1000)
                0169 C     referenced to the model level given by layers_krho.
                0170         rhoShift = rhoConst - 1000. _d 0
77bf6311ca Jean*0171         DO bj=myByLo(myThid),myByHi(myThid)
                0172          DO bi=myBxLo(myThid),myBxHi(myThid)
62b1f3e57f Jean*0173 C     Initialise layers variable prho:
                0174           DO k=1,Nr
                0175            DO j=1-OLy,sNy+OLy
                0176             DO i=1-OLx,sNx+OLx
                0177              prho(i,j,k,bi,bj) = 0. _d 0
                0178             ENDDO
                0179            ENDDO
                0180           ENDDO
77bf6311ca Jean*0181           DO k = 1,Nr
                0182            CALL FIND_RHO_2D( 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
                0183      &                       layers_krho(iLa),
                0184      &                       theta(1-OLx,1-OLy,k,bi,bj),
                0185      &                       salt(1-OLx,1-OLy,k,bi,bj),
1f346f78e1 Jean*0186      &                       prho(1-OLx,1-OLy,k,bi,bj),
77bf6311ca Jean*0187      &                       k, bi, bj, myThid )
50d8304171 Ryan*0188 #ifdef LAYERS_THERMODYNAMICS
                0189 C -- it might be more memory efficient not to store alpha and beta
                0190 C    but to multiply the fluxes in place here
                0191            CALL FIND_ALPHA( bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
                0192      &                      k, layers_krho(iLa),
                0193      &                      layers_alpha(1-OLx,1-OLy,k,bi,bj), myThid )
                0194            CALL FIND_BETA(  bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
                0195      &                      k, layers_krho(iLa),
62b1f3e57f Jean*0196      &                      layers_beta(1-OLx,1-OLy,k,bi,bj), myThid )
                0197 #endif /* LAYERS_THERMODYNAMICS */
77bf6311ca Jean*0198            DO j = 1-OLy,sNy+OLy
                0199             DO i = 1-OLx,sNx+OLx
9643385778 Jean*0200              prho(i,j,k,bi,bj) = prho(i,j,k,bi,bj) + rhoShift
77bf6311ca Jean*0201             ENDDO
                0202            ENDDO
                0203           ENDDO
de7e561253 Jean*0204          ENDDO
                0205         ENDDO
1f346f78e1 Jean*0206         CALL LAYERS_FLUXCALC( uVel,vVel, prho, iLa,
                0207      &              layers_UH, layers_VH,
                0208      &              layers_Hw, layers_Hs,
                0209      &              layers_PIw,layers_PIs,
                0210      &              layers_U,  layers_V,
77bf6311ca Jean*0211      &              myThid )
cf336ab6c5 Ryan*0212 #ifdef LAYERS_THERMODYNAMICS
                0213         CALL LAYERS_DIAPYCNAL( prho,iLa,
50d8304171 Ryan*0214      &              layers_TtendSurf,
                0215      &              layers_TtendDiffh, layers_TtendDiffr,
                0216      &              layers_TtendAdvh, layers_TtendAdvr,
                0217      &              layers_Ttendtot,
                0218      &              layers_StendSurf,
cf336ab6c5 Ryan*0219      &              layers_StendDiffh, layers_StendDiffr,
50d8304171 Ryan*0220      &              layers_StendAdvh, layers_StendAdvr,
                0221      &              layers_Stendtot,
cf336ab6c5 Ryan*0222      &              layers_Hc, layers_PIc,
                0223      &              myThid)
17ce8d85dd Davi*0224 #endif
cf336ab6c5 Ryan*0225 #endif /* LAYERS_PRHO_REF */
84f265d4b3 dfer*0226        ELSEIF ( layers_num(iLa) .EQ. 4 ) THEN
                0227 #ifdef LAYERS_MSE
                0228 C     For layers_num(iLa) = 4, calculate the absolute temperature referenced to 1000 mb
                0229         DO bj=myByLo(myThid),myByHi(myThid)
                0230          DO bi=myBxLo(myThid),myBxHi(myThid)
                0231           DO k = 1,Nr
                0232            conv_theta2T = (rC(k)/atm_Po)**atm_kappa
                0233            DO j = 1-OLy,sNy+OLy
                0234             DO i = 1-OLx,sNx+OLx
                0235              prho(i,j,k,bi,bj) = atm_Cp*conv_theta2T*theta(i,j,k,bi,bj)
                0236      &                         +  totPhiHyd(i,j,k,bi,bj) +  phiRef(2*k)
                0237      &                         +  2501. _d 0 * salt(i,j,k,bi,bj)
                0238             ENDDO
                0239            ENDDO
                0240           ENDDO
                0241          ENDDO
                0242         ENDDO
                0243         CALL LAYERS_FLUXCALC( uVel,vVel, prho, iLa,
                0244      &              layers_UH, layers_VH,
                0245      &              layers_Hw, layers_Hs,
                0246      &              layers_PIw,layers_PIs,
                0247      &              layers_U,  layers_V,
                0248      &              myThid )
                0249 #endif /* LAYERS_MSE */
77bf6311ca Jean*0250        ENDIF
108a00eab9 Ryan*0251 
1f346f78e1 Jean*0252 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0253 C--   Direct Snap-shot output
                0254        IF ( DIFFERENT_MULTIPLE(layers_diagFreq,myTime,deltaTClock)
                0255      &    .AND. layers_num(iLa).NE.0 ) THEN
                0256 
                0257         IF ( layers_MDSIO ) THEN
df5a9764ba Jean*0258           IF ( rwSuffixType.EQ.0 ) THEN
                0259             WRITE(suff,'(I2.2,A1,I10.10)') iLa, '.', myIter
                0260           ELSE
                0261             CALL RW_GET_SUFFIX( sufx, myTime, myIter, myThid )
                0262             WRITE(suff,'(I2.2,A1,A)') iLa, '.', sufx
                0263           ENDIF
1f346f78e1 Jean*0264 #ifdef LAYERS_UFLUX
                0265           CALL WRITE_FLD_3D_RL( 'layers_UH.', suff, Nlayers,
                0266      &                           layers_UH, myIter, myThid )
                0267 #ifdef LAYERS_THICKNESS
                0268           CALL WRITE_FLD_3D_RL( 'layers_Hw.', suff, Nlayers,
                0269      &                           layers_Hw, myIter, myThid )
                0270 #endif /* LAYERS_THICKNESS */
                0271 #endif /* LAYERS_UFLUX */
                0272 #ifdef LAYERS_VFLUX
                0273           CALL WRITE_FLD_3D_RL( 'layers_VH.', suff, Nlayers,
                0274      &                           layers_VH, myIter, myThid )
                0275 #ifdef LAYERS_THICKNESS
ee16a2cae4 Ryan*0276           CALL WRITE_FLD_3D_RL( 'layers_Hs.', suff, Nlayers,
3719fa93a7 Jean*0277      &                           layers_Hs, myIter, myThid )
1f346f78e1 Jean*0278 #endif /* LAYERS_THICKNESS */
                0279 #endif /* LAYERS_VFLUX */
                0280 #ifdef LAYERS_PRHO_REF
5bfb3f2224 Jean*0281           IF ( layers_num(iLa).EQ.3 ) THEN
1f346f78e1 Jean*0282            CALL WRITE_FLD_3D_RL( 'layers_prho.', suff, Nr,
                0283      &                           prho, myIter, myThid )
                0284           ENDIF
                0285 #endif /* LAYERS_PRHO_REF */
50d8304171 Ryan*0286 
62b1f3e57f Jean*0287 #ifdef LAYERS_THERMODYNAMICS
                0288           CALL WRITE_FLD_3D_RL( 'layers_Ttottend.', suff, 2*Nr,
50d8304171 Ryan*0289      &       layers_tottend, myIter, myThid )
                0290 #ifdef SHORTWAVE_HEATING
62b1f3e57f Jean*0291           CALL WRITE_FLD_3D_RL( 'layers_sw.', suff, Nr,
50d8304171 Ryan*0292      &       layers_sw, myIter, myThid )
                0293 #endif /* LAYERS_SHORTWAVE */
62b1f3e57f Jean*0294           CALL WRITE_FLD_3D_RL( 'layers_surfflux.', suff, 2,
ee16a2cae4 Ryan*0295      &                           layers_surfflux, myIter, myThid )
62b1f3e57f Jean*0296           CALL WRITE_FLD_3D_RL( 'layers_dfx.', suff, 2*Nr,
ee16a2cae4 Ryan*0297      &                           layers_dfx, myIter, myThid )
50d8304171 Ryan*0298           CALL WRITE_FLD_3D_RL( 'layers_dfy.', suff, 2*Nr,
ee16a2cae4 Ryan*0299      &                           layers_dfy, myIter, myThid )
50d8304171 Ryan*0300           CALL WRITE_FLD_3D_RL( 'layers_dfr.', suff, 2*Nr,
ee16a2cae4 Ryan*0301      &                           layers_dfr, myIter, myThid )
50d8304171 Ryan*0302           CALL WRITE_FLD_3D_RL( 'layers_afx.', suff, 2*Nr,
                0303      &                           layers_afx, myIter, myThid )
                0304           CALL WRITE_FLD_3D_RL( 'layers_afy.', suff, 2*Nr,
                0305      &                           layers_afy, myIter, myThid )
                0306           CALL WRITE_FLD_3D_RL( 'layers_afr.', suff, 2*Nr,
                0307      &                           layers_afr, myIter, myThid )
                0308 #endif /* LAYERS_THERMODYNAMICS */
1f346f78e1 Jean*0309         ENDIF
                0310 
                0311 c#ifdef ALLOW_MNC
                0312 c#ifdef LAYERS_MNC
                0313 c      IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
                0314 c        pf(1:1) = 'D'
                0315 c      ELSE
                0316 c        pf(1:1) = 'R'
                0317 c      ENDIF
                0318 c        IF ( layers_MNC) THEN
                0319 C           Do MNC output...  But how?
                0320 c        ENDIF
                0321 c#endif /* LAYERS_MNC */
                0322 c#endif /* ALLOW_MNC */
                0323 
                0324        ENDIF
                0325 
                0326 #ifdef ALLOW_DIAGNOSTICS
                0327 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0328 C--   Fill-in diagnostics
                0329        IF ( useDiagnostics .AND. layers_num(iLa).NE.0 ) THEN
                0330 
5bfb3f2224 Jean*0331 #ifdef LAYERS_PRHO_REF
                0332          IF ( layers_num(iLa).EQ.3 ) THEN
55f64449a3 Ryan*0333           WRITE(diagName,'(A4,I1,A3)') 'LaTr',iLa,layers_name(iLa)
5bfb3f2224 Jean*0334           CALL DIAGNOSTICS_FILL( prho,
                0335      &                           diagName, 0, Nr, 0, 1, 1, myThid )
                0336          ENDIF
                0337 #endif /* LAYERS_PRHO_REF */
                0338 
1f346f78e1 Jean*0339 #ifdef LAYERS_UFLUX
                0340          WRITE(diagName,'(A4,I1,A3)') 'LaUH',iLa,layers_name(iLa)
                0341          CALL DIAGNOSTICS_FILL( layers_UH,
                0342      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0343 # ifdef LAYERS_THICKNESS
                0344          WRITE(diagName,'(A4,I1,A3)') 'LaHw',iLa,layers_name(iLa)
                0345          CALL DIAGNOSTICS_FILL( layers_Hw,
                0346      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0347          WRITE(diagName,'(A4,I1,A3)') 'LaPw',iLa,layers_name(iLa)
                0348          CALL DIAGNOSTICS_FILL( layers_PIw,
                0349      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0350          WRITE(diagName,'(A4,I1,A3)') 'LaUa',iLa,layers_name(iLa)
                0351          CALL DIAGNOSTICS_FILL( layers_U,
                0352      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0353 # endif
                0354 #endif /* LAYERS_UFLUX */
                0355 
                0356 #ifdef LAYERS_VFLUX
                0357          WRITE(diagName,'(A4,I1,A3)') 'LaVH',iLa,layers_name(iLa)
                0358          CALL DIAGNOSTICS_FILL( layers_VH,
                0359      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0360 # ifdef LAYERS_THICKNESS
                0361          WRITE(diagName,'(A4,I1,A3)') 'LaHs',iLa,layers_name(iLa)
                0362          CALL DIAGNOSTICS_FILL( layers_Hs,
                0363      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0364          WRITE(diagName,'(A4,I1,A3)') 'LaPs',iLa,layers_name(iLa)
                0365          CALL DIAGNOSTICS_FILL( layers_PIs,
                0366      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0367          WRITE(diagName,'(A4,I1,A3)') 'LaVa',iLa,layers_name(iLa)
                0368          CALL DIAGNOSTICS_FILL( layers_V,
                0369      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0370 # endif
                0371 #endif /* LAYERS_VFLUX */
                0372 
cf336ab6c5 Ryan*0373 #ifdef LAYERS_THERMODYNAMICS
                0374          WRITE(diagName,'(A4,I1,A3)') 'LaHc',iLa,layers_name(iLa)
                0375          CALL DIAGNOSTICS_FILL( layers_Hc,
                0376      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0377          WRITE(diagName,'(A4,I1,A3)') 'LaPc',iLa,layers_name(iLa)
                0378          CALL DIAGNOSTICS_FILL( layers_PIc,
                0379      &                          diagName,0,Nlayers, 0, 1, 1, myThid )
                0380          WRITE(diagName,'(A4,I1,A3)') 'LaTs',iLa,layers_name(iLa)
                0381          CALL DIAGNOSTICS_FILL( layers_TtendSurf,
50d8304171 Ryan*0382      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
cf336ab6c5 Ryan*0383          WRITE(diagName,'(A4,I1,A3)') 'LaTh',iLa,layers_name(iLa)
                0384          CALL DIAGNOSTICS_FILL( layers_TtendDiffh,
50d8304171 Ryan*0385      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
55f64449a3 Ryan*0386          WRITE(diagName,'(A4,I1,A3)') 'LaTz',iLa,layers_name(iLa)
cf336ab6c5 Ryan*0387          CALL DIAGNOSTICS_FILL( layers_TtendDiffr,
50d8304171 Ryan*0388      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
                0389          WRITE(diagName,'(A4,I1,A3)') 'LTha',iLa,layers_name(iLa)
                0390          CALL DIAGNOSTICS_FILL( layers_TtendAdvh,
                0391      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
7139ed4ee5 Ryan*0392          WRITE(diagName,'(A4,I1,A3)') 'LTza',iLa,layers_name(iLa)
50d8304171 Ryan*0393          CALL DIAGNOSTICS_FILL( layers_TtendAdvr,
                0394      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
                0395          WRITE(diagName,'(A4,I1,A3)') 'LTto',iLa,layers_name(iLa)
                0396          CALL DIAGNOSTICS_FILL( layers_Ttendtot,
                0397      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
cf336ab6c5 Ryan*0398          WRITE(diagName,'(A4,I1,A3)') 'LaSs',iLa,layers_name(iLa)
                0399          CALL DIAGNOSTICS_FILL( layers_StendSurf,
50d8304171 Ryan*0400      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
cf336ab6c5 Ryan*0401          WRITE(diagName,'(A4,I1,A3)') 'LaSh',iLa,layers_name(iLa)
                0402          CALL DIAGNOSTICS_FILL( layers_StendDiffh,
50d8304171 Ryan*0403      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
7139ed4ee5 Ryan*0404          WRITE(diagName,'(A4,I1,A3)') 'LaSz',iLa,layers_name(iLa)
cf336ab6c5 Ryan*0405          CALL DIAGNOSTICS_FILL( layers_StendDiffr,
50d8304171 Ryan*0406      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
                0407          WRITE(diagName,'(A4,I1,A3)') 'LSha',iLa,layers_name(iLa)
                0408          CALL DIAGNOSTICS_FILL( layers_StendAdvh,
                0409      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
7139ed4ee5 Ryan*0410          WRITE(diagName,'(A4,I1,A3)') 'LSza',iLa,layers_name(iLa)
50d8304171 Ryan*0411          CALL DIAGNOSTICS_FILL( layers_StendAdvr,
                0412      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
62b1f3e57f Jean*0413          WRITE(diagName,'(A4,I1,A3)') 'LSto',iLa,layers_name(iLa)
50d8304171 Ryan*0414          CALL DIAGNOSTICS_FILL( layers_Stendtot,
                0415      &                          diagName,0,Nlayers-1, 0, 1, 1, myThid )
cf336ab6c5 Ryan*0416 #endif /* LAYERS_THERMODYNAMICS */
                0417 
1f346f78e1 Jean*0418        ENDIF
                0419 #endif /* ALLOW_DIAGNOSTICS */
                0420 
108a00eab9 Ryan*0421 #ifdef ALLOW_TIMEAVE
1f346f78e1 Jean*0422 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
108a00eab9 Ryan*0423 C--   Time-average
406891c1a3 Gael*0424 cgf layers_maxNum loop and dimension would be needed for
                0425 cgf the following and tave output to work beyond iLa.EQ.1
77bf6311ca Jean*0426        IF ( layers_taveFreq.GT.0. .AND. iLa.EQ.1 ) THEN
aa668e01f1 Gael*0427 C --- The tile loops
77bf6311ca Jean*0428         DO bj=myByLo(myThid),myByHi(myThid)
                0429          DO bi=myBxLo(myThid),myBxHi(myThid)
108a00eab9 Ryan*0430 
                0431 #ifdef LAYERS_UFLUX
77bf6311ca Jean*0432           CALL TIMEAVE_CUMULATE( layers_UH_T, layers_UH, Nlayers,
df5a9764ba Jean*0433      &                           deltaTClock, bi, bj, myThid )
108a00eab9 Ryan*0434 #ifdef LAYERS_THICKNESS
77bf6311ca Jean*0435           CALL TIMEAVE_CUMULATE( layers_Hw_T, layers_Hw, Nlayers,
df5a9764ba Jean*0436      &                           deltaTClock, bi, bj, myThid )
108a00eab9 Ryan*0437 #endif /* LAYERS_THICKNESS */
                0438 #endif /* LAYERS_UFLUX */
                0439 #ifdef LAYERS_VFLUX
77bf6311ca Jean*0440           CALL TIMEAVE_CUMULATE( layers_VH_T, layers_VH, Nlayers,
df5a9764ba Jean*0441      &                           deltaTClock, bi, bj, myThid )
108a00eab9 Ryan*0442 #ifdef LAYERS_THICKNESS
77bf6311ca Jean*0443           CALL TIMEAVE_CUMULATE( layers_Hs_T, layers_Hs, Nlayers,
df5a9764ba Jean*0444      &                           deltaTClock, bi, bj, myThid )
108a00eab9 Ryan*0445 #endif /* LAYERS_THICKNESS */
                0446 #endif /* LAYERS_VFLUX */
                0447 
17ce8d85dd Davi*0448 #ifdef LAYERS_PRHO_REF
77bf6311ca Jean*0449           IF ( layers_num(iLa) .EQ. 3 )
                0450      &    CALL TIMEAVE_CUMULATE( prho_tave, prho, Nr,
df5a9764ba Jean*0451      &                           deltaTClock, bi, bj, myThid )
17ce8d85dd Davi*0452 #endif /* LAYERS_PRHO_REF */
                0453 
df5a9764ba Jean*0454           layers_TimeAve(bi,bj)=layers_TimeAve(bi,bj)+deltaTClock
108a00eab9 Ryan*0455 
                0456 C --- End bi,bj loop
77bf6311ca Jean*0457          ENDDO
                0458         ENDDO
                0459        ENDIF
aa668e01f1 Gael*0460 #endif /* ALLOW_TIMEAVE */
108a00eab9 Ryan*0461 
406891c1a3 Gael*0462       ENDDO !DO iLa=1,layers_maxNum
77bf6311ca Jean*0463 
a3d5c99c9c Ryan*0464 #ifdef LAYERS_THERMODYNAMICS
                0465 C--   Reset temporary flux arrays to zero
                0466       DO bj = myByLo(myThid), myByHi(myThid)
                0467        DO bi = myBxLo(myThid), myBxHi(myThid)
                0468         DO iTracer = 1,2
                0469          DO j=1-OLy,sNy+OLy
                0470           DO i=1-OLx,sNx+OLx
                0471            layers_surfflux(i,j,1,iTracer,bi,bj) = 0. _d 0
                0472           ENDDO
                0473          ENDDO
                0474          DO k=1,Nr
                0475           DO j=1-OLy,sNy+OLy
                0476            DO i=1-OLx,sNx+OLx
                0477             layers_dfx     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0478             layers_dfy     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0479             layers_dfr     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0480             layers_afx     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0481             layers_afy     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0482             layers_afr     (i,j,k,iTracer,bi,bj) = 0. _d 0
                0483             layers_tottend (i,j,k,iTracer,bi,bj) = 0. _d 0
                0484 #ifdef SHORTWAVE_HEATING
                0485             layers_sw       (i,j,k,1      ,bi,bj) = 0. _d 0
                0486 #endif /* SHORTWAVE_HEATING */
                0487            ENDDO
                0488           ENDDO
                0489          ENDDO
                0490         ENDDO
                0491        ENDDO
                0492       ENDDO
                0493 #endif /* LAYERS_THERMODYNAMICS */
                0494 
108a00eab9 Ryan*0495 #endif /* ALLOW_LAYERS */
                0496 
                0497       RETURN
                0498       END