Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:37:06 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
45e6cba2ac Jean*0001 #include "PACKAGES_CONFIG.h"
                0002 #include "CPP_OPTIONS.h"
                0003 #ifdef ALLOW_GMREDI
                0004 # include "GMREDI_OPTIONS.h"
                0005 #endif
434645f3d4 Jean*0006 
45e6cba2ac Jean*0007 C--  File taueddy_tendency_apply.F: Routines to apply TAUEDDY tendencies
                0008 C--   Contents
                0009 C--   o TAUEDDY_TENDENCY_APPLY_U
                0010 C--   o TAUEDDY_TENDENCY_APPLY_V
                0011 
                0012 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0013 CBOP
                0014 C     !ROUTINE: TAUEDDY_TENDENCY_APPLY_U
                0015 C     !INTERFACE:
                0016       SUBROUTINE TAUEDDY_TENDENCY_APPLY_U(
                0017      U                     gU_arr,
                0018      I                     iMin,iMax,jMin,jMax, k, bi, bj,
                0019      I                     myTime, myIter, myThid )
                0020 C     !DESCRIPTION: \bv
                0021 C     *==========================================================*
                0022 C     | S/R TAUEDDY_TENDENCY_APPLY_U
                0023 C     | o Contains problem specific forcing for zonal velocity.
                0024 C     *==========================================================*
                0025 C     | Adds terms to gU for forcing by external sources
                0026 C     | e.g. wind stress, bottom friction etc..................
                0027 C     *==========================================================*
                0028 C     \ev
                0029 
                0030 C     !USES:
                0031       IMPLICIT NONE
                0032 C     == Global data ==
                0033 #include "SIZE.h"
                0034 #include "EEPARAMS.h"
                0035 #include "PARAMS.h"
                0036 #include "GRID.h"
                0037 #include "FFIELDS.h"
                0038 #ifdef ALLOW_GMREDI
                0039 # include "GMREDI.h"
                0040 #endif
                0041 
                0042 C     !INPUT/OUTPUT PARAMETERS:
                0043 C     gU_arr    :: the tendency array
                0044 C     iMin,iMax :: Working range of x-index for applying forcing.
                0045 C     jMin,jMax :: Working range of y-index for applying forcing.
                0046 C     k         :: Current vertical level index
                0047 C     bi,bj     :: Current tile indices
                0048 C     myTime    :: Current time in simulation
                0049 C     myIter    :: Current iteration number
                0050 C     myThid    :: my Thread Id number
                0051       _RL     gU_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0052       INTEGER iMin, iMax, jMin, jMax
                0053       INTEGER k, bi, bj
                0054       _RL     myTime
                0055       INTEGER myIter
                0056       INTEGER myThid
                0057 CEOP
                0058 
                0059 #ifdef ALLOW_EDDYPSI
                0060 C     !LOCAL VARIABLES:
                0061 C     i, j      :: Loop counters
                0062       INTEGER i, j
                0063       INTEGER kp1
                0064       _RL maskm1, maskp1
                0065 
                0066 C     Add zonal eddy momentum impulse into the layer
                0067 #ifdef ALLOW_GMREDI
                0068       IF ( GM_InMomAsStress ) THEN
                0069 #endif
                0070       kp1 = MIN(k+1,Nr)
                0071       maskp1 = 1.
                0072       maskm1 = 1.
                0073       IF (k.EQ.Nr) maskp1 = 0.
                0074       IF (k.EQ.1)  maskm1 = 0.
                0075       DO j=jMin,jMax
                0076        DO i=iMin,iMax
                0077         gU_arr(i,j) = gU_arr(i,j)
                0078      &  +foFacMom*recip_rhoConst*
                0079      &  ( maskm1*_maskW(i,j, k ,bi,bj)*tauxEddy(i,j, k ,bi,bj)
                0080      &  - maskp1*_maskW(i,j,kp1,bi,bj)*tauxEddy(i,j,kp1,bi,bj) )
                0081      &  *recip_drF(k)*_recip_hFacW(i,j,k,bi,bj)
                0082        ENDDO
                0083       ENDDO
                0084 #ifdef ALLOW_GMREDI
                0085       ENDIF
                0086 #endif
                0087 
434645f3d4 Jean*0088 #endif /* ALLOW_EDDYPSI */
45e6cba2ac Jean*0089 
                0090       RETURN
                0091       END
                0092 
                0093 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0094 CBOP
                0095 C     !ROUTINE: TAUEDDY_TENDENCY_APPLY_V
                0096 C     !INTERFACE:
                0097       SUBROUTINE TAUEDDY_TENDENCY_APPLY_V(
                0098      U                     gV_arr,
                0099      I                     iMin,iMax,jMin,jMax, k, bi, bj,
                0100      I                     myTime, myIter, myThid )
                0101 C     !DESCRIPTION: \bv
                0102 C     *==========================================================*
                0103 C     | S/R TAUEDDY_TENDENCY_APPLY_V
                0104 C     | o Contains problem specific forcing for merid velocity.
                0105 C     *==========================================================*
                0106 C     | Adds terms to gV for forcing by external sources
                0107 C     | e.g. wind stress, bottom friction etc..................
                0108 C     *==========================================================*
                0109 C     \ev
                0110 
                0111 C     !USES:
                0112       IMPLICIT NONE
                0113 C     == Global data ==
                0114 #include "SIZE.h"
                0115 #include "EEPARAMS.h"
                0116 #include "PARAMS.h"
                0117 #include "GRID.h"
                0118 #include "FFIELDS.h"
                0119 #ifdef ALLOW_GMREDI
                0120 #include "GMREDI.h"
                0121 #endif
                0122 
                0123 C     !INPUT/OUTPUT PARAMETERS:
                0124 C     gV_arr    :: the tendency array
                0125 C     iMin,iMax :: Working range of x-index for applying forcing.
                0126 C     jMin,jMax :: Working range of y-index for applying forcing.
                0127 C     k         :: Current vertical level index
                0128 C     bi,bj     :: Current tile indices
                0129 C     myTime    :: Current time in simulation
                0130 C     myIter    :: Current iteration number
                0131 C     myThid    :: my Thread Id number
                0132       _RL     gV_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0133       INTEGER iMin, iMax, jMin, jMax
                0134       INTEGER k, bi, bj
                0135       _RL     myTime
                0136       INTEGER myIter
                0137       INTEGER myThid
                0138 CEOP
                0139 
                0140 #ifdef ALLOW_EDDYPSI
                0141 C     !LOCAL VARIABLES:
                0142 C     i, j      :: Loop counters
                0143       INTEGER i, j
                0144       INTEGER kp1
                0145       _RL maskm1, maskp1
                0146 
                0147 C     Add meridional eddy momentum impulse into the layer
                0148 #ifdef ALLOW_GMREDI
                0149       IF ( GM_InMomAsStress ) THEN
                0150 #endif
                0151       kp1 = MIN(k+1,Nr)
                0152       maskp1 = 1.
                0153       maskm1 = 1.
                0154       IF (k.EQ.Nr) maskp1 = 0.
                0155       IF (k.EQ.1)  maskm1 = 0.
                0156       DO j=jMin,jMax
                0157        DO i=iMin,iMax
                0158         gV_arr(i,j) = gV_arr(i,j)
                0159      &  +foFacMom*recip_rhoConst*
                0160      &  ( maskm1*_maskS(i,j, k ,bi,bj)*tauyEddy(i,j, k ,bi,bj)
                0161      &  - maskp1*_maskS(i,j,kp1,bi,bj)*tauyEddy(i,j,kp1,bi,bj) )
                0162      &  *recip_drF(k)*_recip_hFacS(i,j,k,bi,bj)
                0163        ENDDO
                0164       ENDDO
                0165 #ifdef ALLOW_GMREDI
                0166       ENDIF
                0167 #endif
                0168 
434645f3d4 Jean*0169 #endif /* ALLOW_EDDYPSI */
45e6cba2ac Jean*0170 
                0171       RETURN
                0172       END