Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:41:45 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
108a00eab9 Ryan*0001 #include "LAYERS_OPTIONS.h"
                0002 
b0b5937c88 Jean*0003 CBOP 0
                0004 C     !ROUTINE: LAYERS_OUTPUT
                0005 
                0006 C     !INTERFACE:
108a00eab9 Ryan*0007       SUBROUTINE LAYERS_OUTPUT( myTime, myIter, myThid )
                0008 
                0009 C     !DESCRIPTION: \bv
                0010 C     *==========================================================*
                0011 C     | SUBROUTINE LAYERS_OUTPUT
                0012 C     | o general routine for LAYERS output
                0013 C     *==========================================================*
                0014 C     |   write time-average & snap-shot output
                0015 C     *==========================================================*
                0016 C     \ev
                0017 
                0018 C     !USES:
                0019       IMPLICIT NONE
                0020 
                0021 C     === Global variables ===
                0022 #include "SIZE.h"
                0023 #include "EEPARAMS.h"
                0024 #include "PARAMS.h"
                0025 #include "LAYERS_SIZE.h"
                0026 #include "LAYERS.h"
                0027 
                0028 C     !INPUT PARAMETERS:
                0029 C     == Routine arguments ==
                0030 C     myTime :: Current time of simulation ( s )
                0031 C     myIter :: Iteration number
                0032 C     myThid :: my Thread Id number
                0033       _RL     myTime
                0034       INTEGER myIter
                0035       INTEGER myThid
                0036 CEOP
                0037 
                0038 #ifdef ALLOW_LAYERS
1f346f78e1 Jean*0039 #ifdef ALLOW_TIMEAVE
108a00eab9 Ryan*0040 C     !LOCAL VARIABLES:
                0041 C     == Local variables ==
                0042       LOGICAL  DIFFERENT_MULTIPLE
                0043       EXTERNAL DIFFERENT_MULTIPLE
df5a9764ba Jean*0044       CHARACTER*(10) suff
406891c1a3 Gael*0045       INTEGER iLa
b0b5937c88 Jean*0046       INTEGER bi, bj
1d0846c989 Jean*0047 
                0048 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
108a00eab9 Ryan*0049 
b0b5937c88 Jean*0050       IF ( layers_taveFreq.GT.0. ) THEN
1d0846c989 Jean*0051 cgf layers_maxNum loop and dimension would be needed for
                0052 cgf the following and tave output to work beyond iLa=1
                0053 c      DO iLa=1,layers_maxNum
                0054        iLa=1
108a00eab9 Ryan*0055 
b52376b6e8 Davi*0056 c set arrays to zero if first timestep
b0b5937c88 Jean*0057        IF ( myIter.EQ.nIter0 ) THEN
                0058         DO bj = myByLo(myThid), myByHi(myThid)
                0059          DO bi = myBxLo(myThid), myBxHi(myThid)
1f346f78e1 Jean*0060           layers_TimeAve(bi,bj) = 0.
b52376b6e8 Davi*0061 #ifdef LAYERS_UFLUX
1f346f78e1 Jean*0062           CALL TIMEAVE_RESET(layers_UH_T,  Nlayers,bi,bj,myThid)
b52376b6e8 Davi*0063 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0064           CALL TIMEAVE_RESET(layers_Hw_T, Nlayers,bi,bj,myThid)
b52376b6e8 Davi*0065 #endif /* LAYERS_THICKNESS */
                0066 #endif /* LAYERS_UFLUX */
                0067 
                0068 #ifdef LAYERS_VFLUX
1f346f78e1 Jean*0069           CALL TIMEAVE_RESET(layers_VH_T, Nlayers,bi,bj,myThid)
b52376b6e8 Davi*0070 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0071           CALL TIMEAVE_RESET(layers_Hs_T, Nlayers,bi,bj,myThid)
b52376b6e8 Davi*0072 #endif /* LAYERS_THICKNESS */
                0073 #endif /* LAYERS_VFLUX */
1f346f78e1 Jean*0074 
                0075 #ifdef LAYERS_PRHO_REF
                0076           CALL TIMEAVE_RESET(prho_tave,Nr,bi,bj,myThid)
                0077 #endif /* LAYERS_PRHO_REF */
b52376b6e8 Davi*0078          ENDDO
                0079         ENDDO
                0080 
108a00eab9 Ryan*0081 C     Dump files and restart average computation if needed
b0b5937c88 Jean*0082        ELSEIF (
b52376b6e8 Davi*0083      &  DIFFERENT_MULTIPLE(layers_taveFreq,myTime,deltaTClock)
b0b5937c88 Jean*0084      &        ) THEN
108a00eab9 Ryan*0085 
                0086 C      Normalize by integrated time
b0b5937c88 Jean*0087         DO bj = myByLo(myThid), myByHi(myThid)
                0088          DO bi = myBxLo(myThid), myBxHi(myThid)
108a00eab9 Ryan*0089 
                0090 #ifdef LAYERS_UFLUX
1f346f78e1 Jean*0091           CALL TIMEAVE_NORMALIZE( layers_UH_T, layers_TimeAve,
b0b5937c88 Jean*0092      &                            Nlayers, bi, bj, myThid )
108a00eab9 Ryan*0093 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0094           CALL TIMEAVE_NORMALIZE( layers_Hw_T, layers_TimeAve,
b0b5937c88 Jean*0095      &                            Nlayers, bi, bj, myThid )
108a00eab9 Ryan*0096 #endif /* LAYERS_THICKNESS */
                0097 #endif /* LAYERS_UFLUX */
                0098 
                0099 #ifdef LAYERS_VFLUX
1f346f78e1 Jean*0100           CALL TIMEAVE_NORMALIZE( layers_VH_T, layers_TimeAve,
b0b5937c88 Jean*0101      &                            Nlayers, bi, bj, myThid )
108a00eab9 Ryan*0102 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0103           CALL TIMEAVE_NORMALIZE( layers_Hs_T, layers_TimeAve,
b0b5937c88 Jean*0104      &                            Nlayers, bi, bj, myThid )
108a00eab9 Ryan*0105 #endif /* LAYERS_THICKNESS */
                0106 #endif /* LAYERS_VFLUX */
                0107 
17ce8d85dd Davi*0108 #ifdef LAYERS_PRHO_REF
1d0846c989 Jean*0109           IF ( layers_num(1).EQ.3 )
1f346f78e1 Jean*0110      &    CALL TIMEAVE_NORMALIZE( prho_tave, layers_TimeAve,
17ce8d85dd Davi*0111      &                            Nr, bi, bj, myThid )
                0112 #endif /* LAYERS_PRHO_REF */
                0113 
b0b5937c88 Jean*0114          ENDDO
108a00eab9 Ryan*0115         ENDDO
                0116 
b0b5937c88 Jean*0117         IF ( layers_MDSIO ) THEN
df5a9764ba Jean*0118          IF ( rwSuffixType.EQ.0 ) THEN
                0119            WRITE(suff,'(I10.10)') myIter
                0120          ELSE
                0121            CALL RW_GET_SUFFIX( suff, myTime, myIter, myThid )
                0122          ENDIF
108a00eab9 Ryan*0123 #ifdef LAYERS_UFLUX
14bb2ae926 Ryan*0124          CALL WRITE_FLD_3D_RL( 'layers_UH-tave.', suff, Nlayers,
                0125      &                          layers_UH_T, myIter, myThid )
108a00eab9 Ryan*0126 #ifdef LAYERS_THICKNESS
14bb2ae926 Ryan*0127          CALL WRITE_FLD_3D_RL( 'layers_Hw-tave.', suff, Nlayers,
                0128      &                          layers_Hw_T, myIter, myThid )
108a00eab9 Ryan*0129 #endif /* LAYERS_THICKNESS */
                0130 #endif /* LAYERS_UFLUX */
                0131 #ifdef LAYERS_VFLUX
14bb2ae926 Ryan*0132          CALL WRITE_FLD_3D_RL( 'layers_VH-tave.', suff, Nlayers,
                0133      &                          layers_VH_T, myIter, myThid )
108a00eab9 Ryan*0134 #ifdef LAYERS_THICKNESS
14bb2ae926 Ryan*0135          CALL WRITE_FLD_3D_RL( 'layers_Hs-tave.', suff, Nlayers,
                0136      &                          layers_Hs_T, myIter, myThid )
108a00eab9 Ryan*0137 #endif /* LAYERS_THICKNESS */
                0138 #endif /* LAYERS_VFLUX */
17ce8d85dd Davi*0139 
                0140 #ifdef LAYERS_PRHO_REF
1d0846c989 Jean*0141          IF ( layers_num(1).EQ.3 )
                0142      &   CALL WRITE_FLD_3D_RL( 'layers_prho-tave.', suff, Nr,
17ce8d85dd Davi*0143      &                          prho_tave, myIter, myThid )
                0144 #endif /* LAYERS_PRHO_REF */
                0145 
b0b5937c88 Jean*0146         ENDIF
108a00eab9 Ryan*0147 
1d0846c989 Jean*0148 c#ifdef ALLOW_MNC
108a00eab9 Ryan*0149 C     Do MNC output.
1d0846c989 Jean*0150 c#endif
108a00eab9 Ryan*0151 
                0152 C      Reset averages to zero
b0b5937c88 Jean*0153         DO bj = myByLo(myThid), myByHi(myThid)
                0154          DO bi = myBxLo(myThid), myBxHi(myThid)
1f346f78e1 Jean*0155           layers_TimeAve(bi,bj) = 0.
108a00eab9 Ryan*0156 #ifdef LAYERS_UFLUX
1f346f78e1 Jean*0157           CALL TIMEAVE_RESET(layers_UH_T,  Nlayers,bi,bj,myThid)
108a00eab9 Ryan*0158 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0159           CALL TIMEAVE_RESET(layers_Hw_T, Nlayers,bi,bj,myThid)
108a00eab9 Ryan*0160 #endif /* LAYERS_THICKNESS */
                0161 #endif /* LAYERS_UFLUX */
                0162 
                0163 #ifdef LAYERS_VFLUX
1f346f78e1 Jean*0164           CALL TIMEAVE_RESET(layers_VH_T, Nlayers,bi,bj,myThid)
108a00eab9 Ryan*0165 #ifdef LAYERS_THICKNESS
1f346f78e1 Jean*0166           CALL TIMEAVE_RESET(layers_Hs_T, Nlayers,bi,bj,myThid)
108a00eab9 Ryan*0167 #endif /* LAYERS_THICKNESS */
                0168 #endif /* LAYERS_VFLUX */
17ce8d85dd Davi*0169 
                0170 #ifdef LAYERS_PRHO_REF
1f346f78e1 Jean*0171           IF ( layers_num(1).EQ.3 )
                0172      &    CALL TIMEAVE_RESET(prho_tave,Nr,bi,bj,myThid)
17ce8d85dd Davi*0173 #endif /* LAYERS_PRHO_REF */
108a00eab9 Ryan*0174          ENDDO
                0175         ENDDO
                0176 
b0b5937c88 Jean*0177 C--   end of bloc: if time is a multiple of layers_taveFreq
                0178        ENDIF
108a00eab9 Ryan*0179 
b0b5937c88 Jean*0180       ENDIF
108a00eab9 Ryan*0181 #endif /* ALLOW_TIMEAVE */
                0182 #endif /* ALLOW_LAYERS */
                0183 
                0184       RETURN
                0185       END