Back to home page

MITgcm

 
 

    


File indexing completed on 2024-10-18 05:11:27 UTC

view on githubraw file Latest commit 5bb179dd on 2024-10-17 18:00:27 UTC
5bb179ddc2 Mart*0001 #include "SEAICE_OPTIONS.h"
                0002 #ifdef ALLOW_OBCS
                0003 # include "OBCS_OPTIONS.h"
                0004 #endif
                0005 #ifdef ALLOW_AUTODIFF
                0006 # include "AUTODIFF_OPTIONS.h"
                0007 #endif
                0008 #undef SEAICE_ALTERNATIVE_MODULUS_U_AND_V
                0009 
                0010 CBOP
                0011 C     !ROUTINE: SEAICE_SIDEDRAG
                0012 C     !INTERFACE:
                0013       SUBROUTINE SEAICE_SIDEDRAG_STRESS(
                0014      I     uIceArg, vIceArg,
                0015      I     coastRoughUarg, coastRoughVarg, AREAarg,
                0016      O     sideDragUarg,sideDragVarg,
                0017      I     iStep, myTime, myIter, myThid )
                0018 
                0019 C     !DESCRIPTION: \bv
                0020 C     *==========================================================*
                0021 C     | SUBROUTINE SEAICE_BOTTOMDRAG_COEFFS
                0022 C     | o Compute the side drag stress for ice-side drag,
                0023 C     |    as a parameterization for island-drag fastice
                0024 C     |    (Liu et al 2022)
                0025 C     *==========================================================*
                0026 C     | written by Yuqing Liu, Jul 2020
                0027 C     *==========================================================*
                0028 C     \ev
                0029 
                0030 C     !USES:
                0031       IMPLICIT NONE
                0032 
                0033 C     === Global variables ===
                0034 #include "SIZE.h"
                0035 #include "EEPARAMS.h"
                0036 #include "PARAMS.h"
                0037 #include "GRID.h"
                0038 #include "DYNVARS.h"
                0039 #include "SEAICE_SIZE.h"
                0040 #include "SEAICE_PARAMS.h"
                0041 #include "SEAICE.h"
                0042 
                0043 C     !INPUT/OUTPUT PARAMETERS:
                0044 C     === Routine arguments ===
                0045 C     u/vIceArg      :: local copies of the current ice velocity
                0046 C     coastRoughUarg :: local copy of normalized coast length
                0047 C     coastRoughVarg :: local copy of normalized coast length
                0048 C     SideDraXarg    :: side drag in x direction
                0049 C     SideDraYarg    :: side drag in y direction
                0050 C     iStep     :: current sub-time step iterate
                0051 C     myTime    :: Simulation time
                0052 C     myIter    :: Simulation timestep number
                0053 C     myThid    :: my Thread Id. number
                0054       _RL uIceArg       (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0055       _RL vIceArg       (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0056       _RL AREAarg       (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0057       _RL sideDragUarg  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0058       _RL sideDragVarg  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0059       _RL coastRoughUarg(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0060       _RL coastRoughVarg(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0061 
                0062       INTEGER iStep
                0063       _RL     myTime
                0064       INTEGER myIter
                0065       INTEGER myThid
                0066 
                0067 #ifdef SEAICE_ALLOW_SIDEDRAG
                0068 C     === local variables ===
                0069 C     i,j,bi,bj,ksrf :: loop indices
                0070       INTEGER i,j,bi,bj
                0071       _RL    tmp, tmpx, tmpy
                0072       _RL    uSpd(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0073       _RL     u_0
                0074 
                0075       u_0 = 0.0005 _d 0
                0076       DO bj=myByLo(myThid),myByHi(myThid)
                0077        DO bi=myBxLo(myThid),myBxHi(myThid)
                0078 C     initialize fields
                0079         DO j=1-OLy,sNy+OLy
                0080          DO i=1-OLx,sNx+OLx
                0081           sideDragUarg(i,j,bi,bj) = 0. _d 0
                0082           sideDragVarg(i,j,bi,bj) = 0. _d 0
                0083           uSpd(i,j)               = 0. _d 0
                0084          ENDDO
                0085         ENDDO
                0086 
                0087         DO j=1-OLy,sNy+OLy-1
                0088          DO i=1-OLx,sNx+OLx-1
                0089           uSpd(i,j) = 0.5 _d 0 *
                0090      &          SQRT((uIceArg(i,  j,bi,bj)+uIceArg(i+1,j,  bi,bj))**2
                0091      &              +(vIceArg(i,  j,bi,bj)+vIceArg(i,  j+1,bi,bj))**2)
                0092          ENDDO
                0093         ENDDO
                0094 C     calculate the sidedrag coeff for landfast ice parameterisations
                0095 C     (Liu et al 2022 use SEAICEsideDrag 2e-4)
                0096         IF (SEAICEsideDrag .GT. 0. _d 0 ) THEN
                0097          DO j=1-OLy+1,sNy+OLy-1
                0098           DO i=1-OLx+1,sNx+OLx-1
                0099            IF ( AREAarg(i,j,bi,bj) .GT. 0.01 _d 0 ) THEN
                0100 # ifndef SEAICE_ALTERNATIVE_MODULUS_U_AND_V
                0101             tmpx = 0.5 _d 0 * ( uSpd(i,j) + uSpd(i-1,j) )
                0102             tmpy = 0.5 _d 0 * ( uSpd(i,j) + uSpd(i,j-1) )
                0103 # else /* ifdef SEAICE_ALTERNATIVE_MODULUS_U_AND_V */
                0104 C     A different way to compute |u|, |v|, deprecated
                0105             tmpx = sqrt(uIceArg(i,j,bi,bj)**2+ (0.25 _d 0*
                0106      &           (vIceArg(i,  j,  bi,bj)+vIceArg(i-1,j,bi,bj)+
                0107      &            vIceArg(i-1,j+1,bi,bj)+vIceArg(i,  j+1,bi,bj)))**2)
                0108             tmpy = sqrt(vIceArg(i,j,bi,bj)**2+ (0.25 _d 0*
                0109      &           (uIceArg(i,  j,  bi,bj)+uIceArg(i,  j-1,bi,bj)+
                0110      &            uIceArg(i+1,j-1,bi,bj)+uIceArg(i+1,j,bi,bj)))**2)
                0111 #endif
                0112 
                0113             sideDragUarg(i,j,bi,bj) =
                0114      &           SEAICEsideDrag/(tmpx + u_0) * coastRoughUarg(i,j,bi,bj)
                0115      &                                       *    seaiceMassU(i,j,bi,bj)
                0116             sideDragVarg(i,j,bi,bj) =
                0117      &           SEAICEsideDrag/(tmpy + u_0) * coastRoughVarg(i,j,bi,bj)
                0118      &                                       *    seaiceMassV(i,j,bi,bj)
                0119            ENDIF
                0120           ENDDO
                0121          ENDDO
                0122         ELSEIF ( SEAICEsideDrag .LT. 0. _d 0 ) THEN
                0123 C     this is the quadratic drag law version, deprecated
                0124          DO j=1-OLy+1,sNy+OLy-1
                0125           DO i=1-OLx+1,sNx+OLx-1
                0126            IF ( AREAarg(i,j,bi,bj) .GT. 0.01 _d 0 ) THEN
                0127 # ifndef SEAICE_ALTERNATIVE_MODULUS_U_AND_V
                0128             tmpx = 0.5 _d 0 * ( uSpd(i,j) + uSpd(i-1,j) )
                0129             tmpy = 0.5 _d 0 * ( uSpd(i,j) + uSpd(i,j-1) )
                0130 # else /* ifdef SEAICE_ALTERNATIVE_MODULUS_U_AND_V */
                0131 C     A different way to compute |u|, |v|, deprecated
                0132             tmpx = sqrt(uIceArg(i,j,bi,bj)**2+ (0.25 _d 0*
                0133      &           (vIceArg(i,  j,  bi,bj)+vIceArg(i-1,j,bi,bj)+
                0134      &            vIceArg(i-1,j+1,bi,bj)+vIceArg(i,  j+1,bi,bj)))**2)
                0135             tmpy = sqrt(vIceArg(i,j,bi,bj)**2+ (0.25 _d 0*
                0136      &           (uIceArg(i,  j,  bi,bj)+uIceArg(i,  j-1,bi,bj)+
                0137      &            uIceArg(i+1,j-1,bi,bj)+uIceArg(i+1,j,bi,bj)))**2)
                0138 # endif
                0139 
                0140             sideDragUarg(i,j,bi,bj) =
                0141      &           - SEAICEsideDrag * tmpx * coastRoughUarg(i,j,bi,bj)
                0142      &                                   *    seaiceMassU(i,j,bi,bj)
                0143             sideDragVarg(i,j,bi,bj) =
                0144      &           - SEAICEsideDrag * tmpy * coastRoughVarg(i,j,bi,bj)
                0145      &                                   *    seaiceMassV(i,j,bi,bj)
                0146            ENDIF
                0147           ENDDO
                0148          ENDDO
                0149         ENDIF
                0150        ENDDO
                0151       ENDDO
                0152 
                0153 #endif /* SEAICE_ALLOW_SIDEDRAG */
                0154 
                0155       RETURN
                0156       END