Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:36:25 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
bcc58d6972 Mich*0001 #include "PACKAGES_CONFIG.h"
                0002 #include "CPP_OPTIONS.h"
                0003 
                0004 #ifdef ALLOW_GMREDI
                0005 # include "GMREDI_OPTIONS.h"
                0006 #endif
                0007 
434645f3d4 Jean*0008       SUBROUTINE CALC_EDDY_STRESS( bi, bj, myThid )
bcc58d6972 Mich*0009 
                0010 C     !DESCRIPTION: \bv
                0011 C     *==========================================================*
                0012 C     | S/R CALC_EDDY_STRESS
                0013 C     | o Calculates the eddy stress when running a residual
                0014 C     |   ocean model
                0015 C     *==========================================================*
                0016 C     | Calculates the eddy stress.  Later this will be added to
                0017 C     | gU the same as external sources (e.g. wind stress, bottom
                0018 C     | friction, etc.
                0019 C     *==========================================================*
                0020 C     \ev
                0021 
                0022       IMPLICIT NONE
                0023 
                0024 #include "SIZE.h"
                0025 #include "EEPARAMS.h"
                0026 #include "PARAMS.h"
                0027 #include "GRID.h"
                0028 #include "FFIELDS.h"
                0029 #ifdef ALLOW_GMREDI
                0030 # include "GMREDI.h"
                0031 #endif
                0032 
                0033 C     !INPUT/OUTPUT PARAMETERS:
                0034 C     == Routine arguments ==
                0035       INTEGER bi,bj
                0036       INTEGER myThid
                0037 #ifdef ALLOW_EDDYPSI
                0038 
                0039 C     !LOCAL VARIABLES:
                0040 C     == Local variables ==
                0041 C     Loop counters
                0042       INTEGER i,j,k
                0043 C     Interpolated stream function and coriolis
                0044       _RL psix, psiy, coriU, coriV
                0045 
                0046 C     Calculate the eddy stress from the eddy induced streamfunction
                0047 #ifdef ALLOW_GMREDI
                0048       IF ( GM_InMomAsStress ) THEN
                0049 #endif
                0050         DO k=1,Nr
434645f3d4 Jean*0051 
                0052          DO j=1-OLy,sNy+OLy-1
                0053           DO i=1-OLx+1,sNx+OLx
bcc58d6972 Mich*0054 #ifdef ALLOW_GMREDI
                0055            psiy = op25*(GM_PsiY(i,  j  ,k,bi,bj)
                0056      &                 +GM_PsiY(i,  j+1,k,bi,bj)
                0057      &                 +GM_PsiY(i-1,j  ,k,bi,bj)
                0058      &                 +GM_PsiY(i-1,j+1,k,bi,bj))
                0059 #else
                0060            psiy = op25*(eddyPsiY(i,  j  ,k,bi,bj)
                0061      &                 +eddyPsiY(i,  j+1,k,bi,bj)
                0062      &                 +eddyPsiY(i-1,j  ,k,bi,bj)
                0063      &                 +eddyPsiY(i-1,j+1,k,bi,bj))
                0064 #endif
                0065            coriU = op5*(fcori(i-1,j,bi,bj)
                0066      &                 +fCori(i  ,j,bi,bj))
c68a19be70 Mich*0067            tauxEddy(i,j,k,bi,bj) =  rhoConst*coriU*psiy
bcc58d6972 Mich*0068           ENDDO
                0069          ENDDO
                0070 
434645f3d4 Jean*0071          DO j=1-OLy+1,sNy+OLy
                0072           DO i=1-OLx,sNx+OLx-1
bcc58d6972 Mich*0073 #ifdef ALLOW_GMREDI
                0074            psix = op25*(GM_PsiX(i,  j  ,k,bi,bj)
                0075      &                 +GM_PsiX(i+1,j  ,k,bi,bj)
                0076      &                 +GM_PsiX(i  ,j-1,k,bi,bj)
                0077      &                 +GM_PsiX(i+1,j-1,k,bi,bj))
                0078 #else
                0079            psix = op25*(eddyPsiX(i,  j  ,k,bi,bj)
                0080      &                 +eddyPsiX(i+1,j  ,k,bi,bj)
                0081      &                 +eddyPsiX(i  ,j-1,k,bi,bj)
                0082      &                 +eddyPsiX(i+1,j-1,k,bi,bj))
                0083 #endif
                0084            coriV = op5*(fcori(i,j-1,bi,bj)
                0085      &                 +fCori(i,j  ,bi,bj))
                0086            tauyEddy(i,j,k,bi,bj) = -rhoConst*coriV*psix
                0087           ENDDO
                0088          ENDDO
                0089 
434645f3d4 Jean*0090 C-    end k loop
bcc58d6972 Mich*0091         ENDDO
                0092 
                0093 #ifdef ALLOW_DIAGNOSTICS
434645f3d4 Jean*0094         IF ( useDiagnostics ) THEN
                0095           CALL DIAGNOSTICS_FILL( tauxEddy, 'TAUXEDDY',
                0096      &                           0, Nr, 1, bi, bj, myThid )
                0097           CALL DIAGNOSTICS_FILL( tauyEddy, 'TAUYEDDY',
                0098      &                           0, Nr, 1, bi, bj, myThid )
                0099         ENDIF
bcc58d6972 Mich*0100 #endif
                0101 #ifdef ALLOW_GMREDI
                0102       ENDIF
                0103 #endif
                0104 #endif /* ALLOW_EDDYPSI */
                0105       RETURN
                0106       END