Back to home page

darwin3

 
 

    


File indexing completed on 2025-11-15 13:24:08 UTC

view on githubraw file Latest commit b7411f1a on 2025-11-06 19:05:26 UTC
cec2469d72 Alis*0001 #include "MOM_VECINV_OPTIONS.h"
a340904e5a Ou W*0002 #include "MOM_COMMON_OPTIONS.h"
88b07ffa37 Jean*0003 #ifdef ALLOW_AUTODIFF
                0004 # include "AUTODIFF_OPTIONS.h"
                0005 #endif
9293d3c672 Hajo*0006 #ifdef ALLOW_GGL90
                0007 # include "GGL90_OPTIONS.h"
f0501c53d1 Jean*0008 #endif
aea29c8517 Alis*0009 
039fe61d35 Jean*0010       SUBROUTINE MOM_VECINV(
07e17fa184 Jean*0011      I        bi,bj,k,iMin,iMax,jMin,jMax,
1833b564cb Jean*0012      I        kappaRU, kappaRV,
07e17fa184 Jean*0013      I        fVerUkm, fVerVkm,
                0014      O        fVerUkp, fVerVkp,
4667e4066d Jean*0015      O        guDiss, gvDiss,
07e17fa184 Jean*0016      I        myTime, myIter, myThid )
                0017 C     *==========================================================*
66d5e1666c Alis*0018 C     | S/R MOM_VECINV                                           |
aea29c8517 Alis*0019 C     | o Form the right hand-side of the momentum equation.     |
07e17fa184 Jean*0020 C     *==========================================================*
aea29c8517 Alis*0021 C     | Terms are evaluated one layer at a time working from     |
                0022 C     | the bottom to the top. The vertically integrated         |
                0023 C     | barotropic flow tendency term is evluated by summing the |
                0024 C     | tendencies.                                              |
                0025 C     | Notes:                                                   |
                0026 C     | We have not sorted out an entirely satisfactory formula  |
                0027 C     | for the diffusion equation bc with lopping. The present  |
                0028 C     | form produces a diffusive flux that does not scale with  |
                0029 C     | open-area. Need to do something to solidfy this and to   |
                0030 C     | deal "properly" with thin walls.                         |
07e17fa184 Jean*0031 C     *==========================================================*
aea29c8517 Alis*0032       IMPLICIT NONE
                0033 
                0034 C     == Global variables ==
                0035 #include "SIZE.h"
                0036 #include "EEPARAMS.h"
                0037 #include "PARAMS.h"
                0038 #include "GRID.h"
01f860d49e Jean*0039 #include "SURFACE.h"
f0501c53d1 Jean*0040 #include "DYNVARS.h"
79074ef66b Jean*0041 #include "FFIELDS.h"
a340904e5a Ou W*0042 #include "MOM_VISC.h"
9293d3c672 Hajo*0043 #if ( defined ALLOW_GGL90 && defined ALLOW_GGL90_LANGMUIR )
                0044 # include "GGL90.h"
                0045 #endif
f0501c53d1 Jean*0046 #ifdef ALLOW_MNC
                0047 # include "MNC_PARAMS.h"
cab1a69b8a Jean*0048 #endif
cd3c16aeda Patr*0049 #ifdef ALLOW_AUTODIFF_TAMC
                0050 # include "tamc.h"
                0051 #endif
aea29c8517 Alis*0052 
                0053 C     == Routine arguments ==
07e17fa184 Jean*0054 C     bi,bj   :: current tile indices
                0055 C     k       :: current vertical level
                0056 C     iMin,iMax,jMin,jMax :: loop ranges
                0057 C     fVerU   :: Flux of momentum in the vertical direction, out of the upper
f1a4cec01a Jean*0058 C     fVerV   :: face of a cell k ( flux into the cell above ).
07e17fa184 Jean*0059 C     fVerUkm :: vertical viscous flux of U, interface above (k-1/2)
                0060 C     fVerVkm :: vertical viscous flux of V, interface above (k-1/2)
                0061 C     fVerUkp :: vertical viscous flux of U, interface below (k+1/2)
                0062 C     fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
                0063 
                0064 C     guDiss  :: dissipation tendency (all explicit terms), u component
                0065 C     gvDiss  :: dissipation tendency (all explicit terms), v component
                0066 C     myTime  :: current time
                0067 C     myIter  :: current time-step number
                0068 C     myThid  :: my Thread Id number
                0069       INTEGER bi,bj,k
                0070       INTEGER iMin,iMax,jMin,jMax
1833b564cb Jean*0071       _RL kappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
                0072       _RL kappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
07e17fa184 Jean*0073       _RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0074       _RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0075       _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0076       _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
4667e4066d Jean*0077       _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0078       _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
f7d48db11c Jean*0079       _RL     myTime
3a279374db Alis*0080       INTEGER myIter
                0081       INTEGER myThid
aea29c8517 Alis*0082 
99bc607d7a Ed H*0083 #ifdef ALLOW_MOM_VECINV
3a279374db Alis*0084 C     == Functions ==
94a46dfe0d Jean*0085       LOGICAL  DIFFERENT_MULTIPLE
                0086       EXTERNAL DIFFERENT_MULTIPLE
a340904e5a Ou W*0087 #ifdef ALLOW_DIAGNOSTICS
                0088       LOGICAL  DIAGNOSTICS_IS_ON
                0089       EXTERNAL DIAGNOSTICS_IS_ON
                0090 #endif
3a279374db Alis*0091 
aea29c8517 Alis*0092 C     == Local variables ==
ed2dbbe83d Jean*0093 C     strainBC :: same as strain but account for no-slip BC
                0094 C     vort3BC  :: same as vort3  but account for no-slip BC
a340904e5a Ou W*0095 C     i,j      :: Loop counters
                0096       INTEGER i,j
aea29c8517 Alis*0097       _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
df1ec59c8a Jean*0098       _RL      vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0099       _RL      uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0100       _RL      vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0101       _RS hFacZ   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
01f860d49e Jean*0102       _RS h0FacZ  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
df1ec59c8a Jean*0103       _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0104       _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0105       _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
9293d3c672 Hajo*0106 #if ( defined ALLOW_GGL90 && defined ALLOW_GGL90_LANGMUIR )
                0107       _RL uRes    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0108       _RL vRes    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0109 #endif
df1ec59c8a Jean*0110       _RL del2u   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0111       _RL del2v   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0112       _RL dStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0113       _RL zStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0114       _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0115       _RL strain  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
ed2dbbe83d Jean*0116       _RL strainBC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
f59d76b0dd Ed D*0117       _RL stretching(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
06244a5e4f Jean*0118 #ifdef ALLOW_LEITH_QG
f59d76b0dd Ed D*0119       _RL Nsquare (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
06244a5e4f Jean*0120 #endif
79074ef66b Jean*0121       _RL cDrag   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
df1ec59c8a Jean*0122       _RL KE      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0123       _RL omega3  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0124       _RL vort3   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
ed2dbbe83d Jean*0125       _RL vort3BC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
df1ec59c8a Jean*0126       _RL hDiv    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0127       _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0128       _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0129       _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0130       _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
07e17fa184 Jean*0131 C     xxxFac :: On-off tracer parameters used for switching terms off.
aea29c8517 Alis*0132       _RL  ArDudrFac
                0133       _RL  ArDvdrFac
df1ec59c8a Jean*0134       _RL  sideMaskFac
aea29c8517 Alis*0135       LOGICAL bottomDragTerms
f7d48db11c Jean*0136       LOGICAL writeDiag
a340904e5a Ou W*0137 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0138       LOGICAL botDragDiagIsOn
                0139       LOGICAL shelfDragDiagIsOn
                0140 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
cd3c16aeda Patr*0141 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0142 C     kkey :: tape key (depends on level and tile indices)
                0143       INTEGER kkey
cd3c16aeda Patr*0144 #endif
cb356b4c5f Ed H*0145 #ifdef ALLOW_MNC
                0146       INTEGER offsets(9)
b22b541fe9 Ed H*0147       CHARACTER*(1) pf
cb356b4c5f Ed H*0148 #endif
                0149 
88b07ffa37 Jean*0150 #ifdef ALLOW_AUTODIFF
7d3b758ca2 Patr*0151 C--   only the kDown part of fverU/V is set in this subroutine
                0152 C--   the kUp is still required
f1a4cec01a Jean*0153 C--   In the case of mom_fluxform kUp is set as well
7d3b758ca2 Patr*0154 C--   (at least in part)
07e17fa184 Jean*0155       fVerUkm(1,1) = fVerUkm(1,1)
                0156       fVerVkm(1,1) = fVerVkm(1,1)
7d3b758ca2 Patr*0157 #endif
                0158 
cd3c16aeda Patr*0159 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0160       kkey = bi + (bj-1)*nSx + (ikey_dynamics-1)*nSx*nSy
                0161       kkey = k  + (kkey-1)*Nr
cd3c16aeda Patr*0162 #endif /* ALLOW_AUTODIFF_TAMC */
                0163 
94a46dfe0d Jean*0164       writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
aea29c8517 Alis*0165 
ef1e764710 Ed H*0166 #ifdef ALLOW_MNC
                0167       IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
b22b541fe9 Ed H*0168         IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
                0169           pf(1:1) = 'D'
                0170         ELSE
                0171           pf(1:1) = 'R'
                0172         ENDIF
cb356b4c5f Ed H*0173         IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
                0174           CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
987ff12cb6 Ed H*0175           CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
cb356b4c5f Ed H*0176           CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
987ff12cb6 Ed H*0177           CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
cb356b4c5f Ed H*0178         ENDIF
                0179         DO i = 1,9
                0180           offsets(i) = 0
                0181         ENDDO
                0182         offsets(3) = k
9c98e82bbb Jean*0183 c       write(*,*) 'offsets = ',(offsets(i),i=1,9)
ef1e764710 Ed H*0184       ENDIF
                0185 #endif /*  ALLOW_MNC  */
                0186 
9c98e82bbb Jean*0187 C--   Initialise intermediate terms
                0188       DO j=1-OLy,sNy+OLy
                0189        DO i=1-OLx,sNx+OLx
4667e4066d Jean*0190         vF(i,j)    = 0.
                0191         vrF(i,j)   = 0.
aea29c8517 Alis*0192         uCf(i,j)   = 0.
                0193         vCf(i,j)   = 0.
                0194         del2u(i,j) = 0.
                0195         del2v(i,j) = 0.
                0196         dStar(i,j) = 0.
                0197         zStar(i,j) = 0.
4667e4066d Jean*0198         guDiss(i,j)= 0.
                0199         gvDiss(i,j)= 0.
9293d3c672 Hajo*0200 #if ( defined ALLOW_GGL90 && defined ALLOW_GGL90_LANGMUIR )
                0201 c       uRes(i,j)  = 0.
                0202 c       vRes(i,j)  = 0.
                0203 #endif
aea29c8517 Alis*0204         vort3(i,j) = 0.
4667e4066d Jean*0205         omega3(i,j)= 0.
df1ec59c8a Jean*0206         KE(i,j)    = 0.
9c98e82bbb Jean*0207 C-    need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
                0208         hDiv(i,j)  = 0.
34e137f064 Jean*0209 c       viscAh_Z(i,j) = 0.
                0210 c       viscAh_D(i,j) = 0.
                0211 c       viscA4_Z(i,j) = 0.
                0212 c       viscA4_D(i,j) = 0.
629d440dd4 Patr*0213         strain(i,j)  = 0. _d 0
ed2dbbe83d Jean*0214         strainBC(i,j)= 0. _d 0
f59d76b0dd Ed D*0215         stretching(i,j) = 0. _d 0
06244a5e4f Jean*0216 #ifdef ALLOW_LEITH_QG
f59d76b0dd Ed D*0217         Nsquare(i,j) = 0. _d 0
06244a5e4f Jean*0218 #endif
629d440dd4 Patr*0219         tension(i,j) = 0. _d 0
88b07ffa37 Jean*0220 #ifdef ALLOW_AUTODIFF
cdc9f269ae Patr*0221         hFacZ(i,j)   = 0. _d 0
629d440dd4 Patr*0222 #endif
aea29c8517 Alis*0223        ENDDO
                0224       ENDDO
                0225 
                0226 C--   Term by term tracer parmeters
                0227 C     o U momentum equation
                0228       ArDudrFac    = vfFacMom*1.
                0229 C     o V momentum equation
                0230       ArDvdrFac    = vfFacMom*1.
                0231 
df1ec59c8a Jean*0232 C note: using standard stencil (no mask) results in under-estimating
                0233 C       vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
                0234       IF ( no_slip_sides ) THEN
                0235         sideMaskFac = sideDragFactor
                0236       ELSE
                0237         sideMaskFac = 0. _d 0
                0238       ENDIF
                0239 
99731c717f Jean*0240       IF ( selectImplicitDrag.EQ.0 .AND.
                0241      &      (  no_slip_bottom
932b38363b Jean*0242      &    .OR. selectBotDragQuadr.GE.0
99731c717f Jean*0243      &    .OR. bottomDragLinear.NE.0. ) ) THEN
aea29c8517 Alis*0244        bottomDragTerms=.TRUE.
                0245       ELSE
                0246        bottomDragTerms=.FALSE.
                0247       ENDIF
                0248 
a340904e5a Ou W*0249 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0250       botDragDiagIsOn = .FALSE.
                0251       shelfDragDiagIsOn = .FALSE.
                0252       IF ( useDiagnostics ) THEN
                0253        IF ( bottomDragTerms ) botDragDiagIsOn =
                0254      &            DIAGNOSTICS_IS_ON( 'UBotDrag', myThid )
                0255      &       .OR. DIAGNOSTICS_IS_ON( 'VBotDrag', myThid )
                0256        IF ( useShelfIce ) shelfDragDiagIsOn =
                0257      &            DIAGNOSTICS_IS_ON( 'UShIDrag', myThid )
                0258      &       .OR. DIAGNOSTICS_IS_ON( 'VShIDrag', myThid )
                0259       ENDIF
                0260 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
                0261 
aea29c8517 Alis*0262 C--   Calculate open water fraction at vorticity points
                0263       CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
                0264 
                0265 C     Make local copies of horizontal flow field
                0266       DO j=1-OLy,sNy+OLy
                0267        DO i=1-OLx,sNx+OLx
                0268         uFld(i,j) = uVel(i,j,k,bi,bj)
                0269         vFld(i,j) = vVel(i,j,k,bi,bj)
                0270        ENDDO
                0271       ENDDO
                0272 
d09cd10c60 Gael*0273 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0274 CADJ STORE hFacZ(:,:)   = comlev1_bibj_k, key = kkey, byte = isbyte
                0275 CADJ STORE r_hFacZ(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0276 CADJ STORE fverukm(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0277 CADJ STORE fvervkm(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0278 #endif
                0279 
cab1a69b8a Jean*0280 C note (jmc) : Dissipation and Vort3 advection do not necesary
                0281 C              use the same maskZ (and hFacZ)  => needs 2 call(s)
                0282 c     CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
                0283 
b8082fc644 Jean*0284       CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
aea29c8517 Alis*0285 
7c7b0b4a46 Alis*0286       CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
aea29c8517 Alis*0287 
ed2dbbe83d Jean*0288 C-    mask vort3 and account for no-slip / free-slip BC in vort3BC:
                0289       DO j=1-OLy,sNy+OLy
                0290        DO i=1-OLx,sNx+OLx
                0291          vort3BC(i,j) = vort3(i,j)
                0292          IF ( hFacZ(i,j).EQ.zeroRS ) THEN
                0293            vort3BC(i,j) = sideMaskFac*vort3BC(i,j)
                0294            vort3(i,j)   = 0.
                0295          ENDIF
                0296        ENDDO
                0297       ENDDO
                0298 
d09cd10c60 Gael*0299 #ifdef ALLOW_AUTODIFF_TAMC
7448700841 Mart*0300 CADJ STORE KE(:,:)      = comlev1_bibj_k, key = kkey, byte = isbyte
edb6656069 Mart*0301 CADJ STORE vort3(:,:)   = comlev1_bibj_k, key = kkey, byte = isbyte
                0302 CADJ STORE vort3bc(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0303 #endif
                0304 
aea29c8517 Alis*0305       IF (momViscosity) THEN
039fe61d35 Jean*0306 C--    For viscous term, compute horizontal divergence, tension & strain
df1ec59c8a Jean*0307 C      and mask relative vorticity (free-slip case):
                0308 
01f860d49e Jean*0309        DO j=1-OLy,sNy+OLy
                0310         DO i=1-OLx,sNx+OLx
                0311           h0FacZ(i,j) = hFacZ(i,j)
                0312         ENDDO
                0313        ENDDO
                0314 #ifdef NONLIN_FRSURF
                0315        IF ( no_slip_sides .AND. nonlinFreeSurf.GT.0 ) THEN
                0316         DO j=2-OLy,sNy+OLy
                0317          DO i=2-OLx,sNx+OLx
                0318           h0FacZ(i,j) = MIN(
                0319      &       MIN( h0FacW(i,j,k,bi,bj), h0FacW(i,j-1,k,bi,bj) ),
                0320      &       MIN( h0FacS(i,j,k,bi,bj), h0FacS(i-1,j,k,bi,bj) ) )
                0321          ENDDO
                0322         ENDDO
                0323        ENDIF
7448700841 Mart*0324 # ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0325 CADJ STORE h0FacZ(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
7448700841 Mart*0326 # endif
                0327 #endif /* NONLIN_FRSURF */
d09cd10c60 Gael*0328 
df1ec59c8a Jean*0329        CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
                0330 
ed2dbbe83d Jean*0331        IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
                0332         CALL MOM_CALC_TENSION( bi,bj,k,uFld,vFld,tension,myThid )
                0333         CALL MOM_CALC_STRAIN( bi,bj,k,uFld,vFld,hFacZ,strain,myThid )
                0334 C-    mask strain and account for no-slip / free-slip BC in strainBC:
                0335         DO j=1-OLy,sNy+OLy
                0336          DO i=1-OLx,sNx+OLx
                0337            strainBC(i,j) = strain(i,j)
                0338            IF ( hFacZ(i,j).EQ.zeroRS ) THEN
                0339              strainBC(i,j) = sideMaskFac*strainBC(i,j)
                0340              strain(i,j)   = 0.
                0341            ENDIF
                0342          ENDDO
df1ec59c8a Jean*0343         ENDDO
f59d76b0dd Ed D*0344 #ifdef ALLOW_LEITH_QG
a125ab7be7 Jean*0345         IF ( viscC2LeithQG.NE.zeroRL ) THEN
f59d76b0dd Ed D*0346           CALL MOM_VISC_QGL_STRETCH(bi,bj,k,
                0347      O                            stretching, Nsquare,
                0348      I                            myTime, myIter, myThid )
                0349           CALL MOM_VISC_QGL_LIMIT(bi,bj,k,
                0350      O                          stretching,
                0351      I                          Nsquare, uFld,vFld,vort3,
                0352      I                          myTime, myIter, myThid )
                0353         ENDIF
                0354 #endif /* ALLOW_LEITH_QG */
ed2dbbe83d Jean*0355        ENDIF
df1ec59c8a Jean*0356 
d09cd10c60 Gael*0357 #ifdef ALLOW_AUTODIFF_TAMC
7448700841 Mart*0358 CADJ STORE hDiv(:,:)     = comlev1_bibj_k, key = kkey, byte = isbyte
edb6656069 Mart*0359 CADJ STORE tension(:,:)  = comlev1_bibj_k, key = kkey, byte = isbyte
                0360 CADJ STORE strain(:,:)   = comlev1_bibj_k, key = kkey, byte = isbyte
7448700841 Mart*0361 CADJ STORE strainBC(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0362 #endif
                0363 
34e137f064 Jean*0364 C--    Calculate Lateral Viscosities
                0365        DO j=1-OLy,sNy+OLy
                0366         DO i=1-OLx,sNx+OLx
                0367          viscAh_D(i,j) = viscAhD
                0368          viscAh_Z(i,j) = viscAhZ
                0369          viscA4_D(i,j) = viscA4D
                0370          viscA4_Z(i,j) = viscA4Z
                0371         ENDDO
                0372        ENDDO
                0373        IF ( useVariableVisc ) THEN
ed2dbbe83d Jean*0374 C-     uses vort3BC & strainBC which account for no-slip / free-slip BC
79074ef66b Jean*0375 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0376 CADJ STORE KE(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
79074ef66b Jean*0377 #endif
34e137f064 Jean*0378          CALL MOM_CALC_VISC( bi, bj, k,
                0379      O            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
f59d76b0dd Ed D*0380      I            hDiv, vort3BC, tension, strainBC, stretching,
                0381      I            KE, hfacZ, myThid )
d09cd10c60 Gael*0382 #ifdef ALLOW_AUTODIFF_TAMC
4240547d2d Mart*0383 # ifndef AUTODIFF_ALLOW_VISCFACADJ
                0384 C     These store directive must not be here, if you want to recompute
                0385 C     these viscosity coefficients with a modified viscFacAdj = viscFacInAd
                0386 C     because the store directives intentionally prevent the recomputation.
edb6656069 Mart*0387 CADJ STORE viscAh_Z(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0388 CADJ STORE viscAh_D(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0389 CADJ STORE viscA4_Z(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0390 CADJ STORE viscA4_D(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
4240547d2d Mart*0391 # endif /* AUTODIFF_ALLOW_VISCFACADJ */
                0392 #endif /* ALLOW_AUTODIFF_TAMC */
                0393        ENDIF
d09cd10c60 Gael*0394 
                0395 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0396 CADJ STORE hDiv(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0397 #endif
                0398 
aea29c8517 Alis*0399 C      Calculate del^2 u and del^2 v for bi-harmonic term
f0501c53d1 Jean*0400        IF (useBiharmonicVisc) THEN
3a279374db Alis*0401          CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
                0402      O                      del2u,del2v,
ed2dbbe83d Jean*0403      I                      myThid)
88e5e58941 Jean*0404          CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
                0405          CALL MOM_CALC_RELVORT3(bi,bj,k,
                0406      &                          del2u,del2v,hFacZ,zStar,myThid)
df1ec59c8a Jean*0407        ENDIF
                0408 
d09cd10c60 Gael*0409 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0410 CADJ STORE del2u(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0411 CADJ STORE del2v(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0412 CADJ STORE dStar(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0413 CADJ STORE zStar(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0414 #endif
                0415 
df1ec59c8a Jean*0416 C---   Calculate dissipation terms for U and V equations
b0c3bd7ab0 Bayl*0417 
ed2dbbe83d Jean*0418 C-     in terms of tension and strain
b0c3bd7ab0 Bayl*0419        IF (useStrainTensionVisc) THEN
ed2dbbe83d Jean*0420 C      use masked strain as if free-slip since side-drag is computed separately
f0501c53d1 Jean*0421          CALL MOM_HDISSIP( bi, bj, k,
ed2dbbe83d Jean*0422      I            tension, strain, hFacZ,
f0501c53d1 Jean*0423      I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
                0424      I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
                0425      O            guDiss, gvDiss,
                0426      I            myThid )
b0c3bd7ab0 Bayl*0427        ELSE
ed2dbbe83d Jean*0428 C-     in terms of vorticity and divergence
f0501c53d1 Jean*0429          CALL MOM_VI_HDISSIP( bi, bj, k,
ed2dbbe83d Jean*0430      I            hDiv, vort3, dStar, zStar, hFacZ,
f0501c53d1 Jean*0431      I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
                0432      I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
                0433      O            guDiss, gvDiss,
ed2dbbe83d Jean*0434      I            myThid )
07cc642809 Alis*0435        ENDIF
cab1a69b8a Jean*0436 
df1ec59c8a Jean*0437 C---  Other dissipation terms in Zonal momentum equation
aea29c8517 Alis*0438 
                0439 C--   Vertical flux (fVer is at upper face of "u" cell)
                0440 C     Eddy component of vertical flux (interior component only) -> vrF
79074ef66b Jean*0441        IF ( .NOT.implicitViscosity ) THEN
                0442         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,kappaRU,vrF,myThid)
aea29c8517 Alis*0443 C     Combine fluxes
79074ef66b Jean*0444         DO j=jMin,jMax
                0445          DO i=iMin,iMax
                0446           fVerUkp(i,j) = ArDudrFac*vrF(i,j)
                0447          ENDDO
4667e4066d Jean*0448         ENDDO
d09cd10c60 Gael*0449 
                0450 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0451 CADJ STORE fVerUkp(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0452 #endif
                0453 
4667e4066d Jean*0454 C--   Tendency is minus divergence of the fluxes
f1a4cec01a Jean*0455 C     vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
79074ef66b Jean*0456         DO j=jMin,jMax
                0457          DO i=iMin,iMax
                0458           guDiss(i,j) = guDiss(i,j)
                0459      &    -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
                0460      &    *recip_rAw(i,j,bi,bj)
                0461      &    *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
                0462      &    *recip_deepFac2C(k)*recip_rhoFacC(k)
                0463          ENDDO
4667e4066d Jean*0464         ENDDO
79074ef66b Jean*0465        ENDIF
aea29c8517 Alis*0466 
039fe61d35 Jean*0467 C-- No-slip and drag BCs appear as body forces in cell abutting topography
79074ef66b Jean*0468        IF ( no_slip_sides ) THEN
aea29c8517 Alis*0469 C-     No-slip BCs impose a drag at walls...
79074ef66b Jean*0470         CALL MOM_U_SIDEDRAG( bi, bj, k,
                0471      I           uFld, del2u, h0FacZ,
                0472      I           viscAh_Z, viscA4_Z,
                0473      I           useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
                0474      O           vF,
                0475      I           myThid )
                0476         DO j=jMin,jMax
                0477          DO i=iMin,iMax
                0478           guDiss(i,j) = guDiss(i,j)+vF(i,j)
                0479          ENDDO
aea29c8517 Alis*0480         ENDDO
79074ef66b Jean*0481        ENDIF
34e137f064 Jean*0482 
aea29c8517 Alis*0483 C-    No-slip BCs impose a drag at bottom
79074ef66b Jean*0484        IF ( bottomDragTerms ) THEN
                0485 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0486 CADJ STORE KE(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
79074ef66b Jean*0487 #endif
a125ab7be7 Jean*0488         CALL MOM_U_BOTDRAG_COEFF( bi, bj, k, .TRUE.,
79074ef66b Jean*0489      I                  uFld, vFld, kappaRU, KE,
                0490      O                  cDrag,
                0491      I                  myIter, myThid )
7448700841 Mart*0492 #ifdef ALLOW_AUTODIFF_TAMC
                0493 CADJ STORE cDrag(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0494 #endif
79074ef66b Jean*0495         DO j=jMin,jMax
                0496          DO i=iMin,iMax
a340904e5a Ou W*0497            vF(i,j) = -cDrag(i,j)*uFld(i,j)
                0498      &             *_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
                0499            gUdiss(i,j) = gUdiss(i,j) + vF(i,j)
79074ef66b Jean*0500          ENDDO
aea29c8517 Alis*0501         ENDDO
79074ef66b Jean*0502         IF ( useDiagnostics ) THEN
                0503          DO j=jMin,jMax
                0504           DO i=iMin,iMax
                0505            botDragU(i,j,bi,bj) = botDragU(i,j,bi,bj)
                0506      &                         - cDrag(i,j)*uFld(i,j)*rUnit2mass
                0507           ENDDO
                0508          ENDDO
                0509         ENDIF
a340904e5a Ou W*0510 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0511         IF ( botDragDiagIsOn ) THEN
                0512          CALL DIAGNOSTICS_FILL( vF, 'UBotDrag', k,1,2,bi,bj, myThid )
                0513         ENDIF
                0514 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
79074ef66b Jean*0515        ENDIF
dd49642a1d Mart*0516 #ifdef ALLOW_SHELFICE
e2d988bd46 Jean*0517        IF ( useShelfIce .AND. selectImplicitDrag.EQ.0 ) THEN
79074ef66b Jean*0518 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0519 CADJ STORE KE(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
79074ef66b Jean*0520 #endif
e2d988bd46 Jean*0521         CALL SHELFICE_U_DRAG_COEFF( bi, bj, k, .TRUE.,
                0522      I                  uFld, vFld, kappaRU, KE,
                0523      O                  cDrag,
                0524      I                  myIter, myThid )
7448700841 Mart*0525 #ifdef ALLOW_AUTODIFF_TAMC
                0526 CADJ STORE cDrag(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0527 #endif
79074ef66b Jean*0528         DO j=jMin,jMax
                0529          DO i=iMin,iMax
e2d988bd46 Jean*0530            gUdiss(i,j) = gUdiss(i,j)
                0531      &                 - cDrag(i,j)*uFld(i,j)
                0532      &                 *_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
79074ef66b Jean*0533          ENDDO
dd49642a1d Mart*0534         ENDDO
a340904e5a Ou W*0535 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0536         IF ( shelfDragDiagIsOn ) THEN
                0537          DO j=jMin,jMax
                0538           DO i=iMin,iMax
                0539            vrF(i,j) = -cDrag(i,j)*uFld(i,j)
                0540      &              *_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
                0541           ENDDO
                0542          ENDDO
                0543          CALL DIAGNOSTICS_FILL( vrF, 'UShIDrag', k,1,2,bi,bj, myThid )
                0544         ENDIF
                0545 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
79074ef66b Jean*0546        ENDIF
dd49642a1d Mart*0547 #endif /* ALLOW_SHELFICE */
                0548 
df1ec59c8a Jean*0549 C---  Other dissipation terms in Meridional momentum equation
aea29c8517 Alis*0550 
                0551 C--   Vertical flux (fVer is at upper face of "v" cell)
                0552 C     Eddy component of vertical flux (interior component only) -> vrF
79074ef66b Jean*0553        IF ( .NOT.implicitViscosity ) THEN
                0554         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,kappaRV,vrF,myThid)
aea29c8517 Alis*0555 C     Combine fluxes -> fVerV
79074ef66b Jean*0556         DO j=jMin,jMax
                0557          DO i=iMin,iMax
                0558           fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
                0559          ENDDO
4667e4066d Jean*0560         ENDDO
d09cd10c60 Gael*0561 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0562 CADJ STORE fVerVkp(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0563 #endif
4667e4066d Jean*0564 C--   Tendency is minus divergence of the fluxes
f1a4cec01a Jean*0565 C     vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
79074ef66b Jean*0566         DO j=jMin,jMax
                0567          DO i=iMin,iMax
                0568           gvDiss(i,j) = gvDiss(i,j)
                0569      &    -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
                0570      &    *recip_rAs(i,j,bi,bj)
                0571      &    *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
                0572      &    *recip_deepFac2C(k)*recip_rhoFacC(k)
                0573          ENDDO
4667e4066d Jean*0574         ENDDO
79074ef66b Jean*0575        ENDIF
aea29c8517 Alis*0576 
039fe61d35 Jean*0577 C-- No-slip and drag BCs appear as body forces in cell abutting topography
79074ef66b Jean*0578        IF ( no_slip_sides ) THEN
aea29c8517 Alis*0579 C-     No-slip BCs impose a drag at walls...
79074ef66b Jean*0580         CALL MOM_V_SIDEDRAG( bi, bj, k,
                0581      I           vFld, del2v, h0FacZ,
                0582      I           viscAh_Z, viscA4_Z,
                0583      I           useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
                0584      O           vF,
                0585      I           myThid )
                0586         DO j=jMin,jMax
                0587          DO i=iMin,iMax
                0588           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
                0589          ENDDO
aea29c8517 Alis*0590         ENDDO
79074ef66b Jean*0591        ENDIF
34e137f064 Jean*0592 
aea29c8517 Alis*0593 C-    No-slip BCs impose a drag at bottom
79074ef66b Jean*0594        IF ( bottomDragTerms ) THEN
                0595 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0596 CADJ STORE KE(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
79074ef66b Jean*0597 #endif
a125ab7be7 Jean*0598         CALL MOM_V_BOTDRAG_COEFF( bi, bj, k, .TRUE.,
79074ef66b Jean*0599      I                  uFld, vFld, kappaRV, KE,
                0600      O                  cDrag,
                0601      I                  myIter, myThid )
7448700841 Mart*0602 #ifdef ALLOW_AUTODIFF_TAMC
                0603 CADJ STORE cDrag(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0604 #endif
79074ef66b Jean*0605         DO j=jMin,jMax
                0606          DO i=iMin,iMax
a340904e5a Ou W*0607            vF(i,j) = -cDrag(i,j)*vFld(i,j)
                0608      &             *_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
                0609            gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
79074ef66b Jean*0610          ENDDO
aea29c8517 Alis*0611         ENDDO
79074ef66b Jean*0612         IF ( useDiagnostics ) THEN
                0613          DO j=jMin,jMax
                0614           DO i=iMin,iMax
                0615            botDragV(i,j,bi,bj) = botDragV(i,j,bi,bj)
                0616      &                         - cDrag(i,j)*vFld(i,j)*rUnit2mass
                0617           ENDDO
                0618          ENDDO
                0619         ENDIF
a340904e5a Ou W*0620 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0621         IF ( botDragDiagIsOn ) THEN
                0622          CALL DIAGNOSTICS_FILL( vF, 'VBotDrag', k,1,2,bi,bj, myThid )
                0623         ENDIF
                0624 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
79074ef66b Jean*0625        ENDIF
dd49642a1d Mart*0626 #ifdef ALLOW_SHELFICE
e2d988bd46 Jean*0627        IF ( useShelfIce .AND. selectImplicitDrag.EQ.0 ) THEN
79074ef66b Jean*0628 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0629 CADJ STORE KE(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
79074ef66b Jean*0630 #endif
e2d988bd46 Jean*0631         CALL SHELFICE_V_DRAG_COEFF( bi, bj, k, .TRUE.,
                0632      I                  uFld, vFld, kappaRV, KE,
                0633      O                  cDrag,
                0634      I                  myIter, myThid )
7448700841 Mart*0635 #ifdef ALLOW_AUTODIFF_TAMC
                0636 CADJ STORE cDrag(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0637 #endif
79074ef66b Jean*0638         DO j=jMin,jMax
                0639          DO i=iMin,iMax
e2d988bd46 Jean*0640            gvDiss(i,j) = gvDiss(i,j)
                0641      &                 - cDrag(i,j)*vFld(i,j)
                0642      &                 *_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
79074ef66b Jean*0643          ENDDO
932b38363b Jean*0644         ENDDO
a340904e5a Ou W*0645 #if ( defined ALLOW_MOM_TEND_EXTRA_DIAGS && defined ALLOW_DIAGNOSTICS )
                0646         IF ( shelfDragDiagIsOn ) THEN
                0647          DO j=jMin,jMax
                0648           DO i=iMin,iMax
                0649            vrF(i,j) = -cDrag(i,j)*vFld(i,j)
                0650      &              *_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
                0651           ENDDO
                0652          ENDDO
                0653          CALL DIAGNOSTICS_FILL( vrF, 'VShIDrag', k,1,2,bi,bj, myThid )
                0654         ENDIF
                0655 #endif /* ALLOW_MOM_TEND_EXTRA_DIAGS && ALLOW_DIAGNOSTICS */
79074ef66b Jean*0656        ENDIF
dd49642a1d Mart*0657 #endif /* ALLOW_SHELFICE */
                0658 
34e137f064 Jean*0659 C--   if (momViscosity) end of block.
                0660       ENDIF
                0661 
                0662 C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:
                0663 c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
aea29c8517 Alis*0664 
df1ec59c8a Jean*0665 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0666 
                0667 C---  Prepare for Advection & Coriolis terms:
ed2dbbe83d Jean*0668 C-    calculate absolute vorticity
df1ec59c8a Jean*0669       IF (useAbsVorticity)
                0670      &  CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
aea29c8517 Alis*0671 
d09cd10c60 Gael*0672 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0673 CADJ STORE omega3(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0674 #endif
                0675 
aea29c8517 Alis*0676 C--   Horizontal Coriolis terms
a75214ff78 Jean*0677 c     IF (useCoriolis .AND. .NOT.useCDscheme
                0678 c    &    .AND. .NOT. useAbsVorticity) THEN
                0679 C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
711329b07b Jean*0680       IF ( useCoriolis .AND.
a75214ff78 Jean*0681      &     .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
                0682      &   ) THEN
                0683        IF (useAbsVorticity) THEN
7c3c2cec96 Jean*0684         CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0685      &                         vFld,omega3,hFacZ,r_hFacZ,
a75214ff78 Jean*0686      &                         uCf,myThid)
7c3c2cec96 Jean*0687         CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0688      &                         uFld,omega3,hFacZ,r_hFacZ,
a75214ff78 Jean*0689      &                         vCf,myThid)
9293d3c672 Hajo*0690 #if ( defined ALLOW_GGL90 && defined ALLOW_GGL90_LANGMUIR )
                0691        ELSEIF ( useLANGMUIR ) THEN
                0692         CALL GGL90_ADD_STOKESDRIFT(
                0693      O                 uRes, vRes,
                0694      I                 uFld, vFld, k, bi, bj, myThid )
                0695         CALL MOM_VI_CORIOLIS(bi,bj,k,uRes,vRes,hFacZ,r_hFacZ,
                0696      &                       uCf,vCf,myThid)
                0697 #endif
a75214ff78 Jean*0698        ELSE
                0699         CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
                0700      &                       uCf,vCf,myThid)
                0701        ENDIF
481f592dad Jean*0702        DO j=jMin,jMax
                0703         DO i=iMin,iMax
1aff67ca82 Jean*0704          gU(i,j,k,bi,bj) = uCf(i,j)
                0705          gV(i,j,k,bi,bj) = vCf(i,j)
481f592dad Jean*0706         ENDDO
aea29c8517 Alis*0707        ENDDO
f7d48db11c Jean*0708        IF ( writeDiag ) THEN
ef1e764710 Ed H*0709          IF (snapshot_mdsio) THEN
                0710            CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
                0711            CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
                0712          ENDIF
                0713 #ifdef ALLOW_MNC
                0714          IF (useMNC .AND. snapshot_mnc) THEN
b22b541fe9 Ed H*0715            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
cb356b4c5f Ed H*0716      &          offsets, myThid)
b22b541fe9 Ed H*0717            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
cb356b4c5f Ed H*0718      &          offsets, myThid)
ef1e764710 Ed H*0719          ENDIF
                0720 #endif /*  ALLOW_MNC  */
f7d48db11c Jean*0721        ENDIF
711329b07b Jean*0722 #ifdef ALLOW_DIAGNOSTICS
                0723        IF ( useDiagnostics ) THEN
                0724          CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
                0725          CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
                0726        ENDIF
                0727 #endif /* ALLOW_DIAGNOSTICS */
4667e4066d Jean*0728       ELSE
                0729        DO j=jMin,jMax
                0730         DO i=iMin,iMax
1aff67ca82 Jean*0731          gU(i,j,k,bi,bj) = 0. _d 0
                0732          gV(i,j,k,bi,bj) = 0. _d 0
4667e4066d Jean*0733         ENDDO
                0734        ENDDO
481f592dad Jean*0735       ENDIF
                0736 
d09cd10c60 Gael*0737 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0738 CADJ STORE ucf(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0739 CADJ STORE vcf(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0740 #endif
                0741 
481f592dad Jean*0742       IF (momAdvection) THEN
3add9696e1 Jean*0743 C--   Horizontal advection of relative (or absolute) vorticity
7fe6343684 Jean*0744        IF ( (highOrderVorticity.OR.upwindVorticity)
                0745      &     .AND.useAbsVorticity ) THEN
79074ef66b Jean*0746         CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
3370e71df9 Mart*0747      &                         highOrderVorticity,upwindVorticity,
                0748      &                         vFld,omega3,r_hFacZ,
d4c99033f5 Jean*0749      &                         uCf,myThid)
7fe6343684 Jean*0750        ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
79074ef66b Jean*0751         CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
3370e71df9 Mart*0752      &                         highOrderVorticity, upwindVorticity,
                0753      &                         vFld,vort3, r_hFacZ,
3add9696e1 Jean*0754      &                         uCf,myThid)
7fe6343684 Jean*0755        ELSEIF ( useAbsVorticity ) THEN
3370e71df9 Mart*0756         CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0757      &                         vFld,omega3,hFacZ,r_hFacZ,
5d5fff3a81 Alis*0758      &                         uCf,myThid)
                0759        ELSE
3370e71df9 Mart*0760         CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0761      &                         vFld,vort3, hFacZ,r_hFacZ,
5d5fff3a81 Alis*0762      &                         uCf,myThid)
                0763        ENDIF
481f592dad Jean*0764        DO j=jMin,jMax
                0765         DO i=iMin,iMax
                0766          gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
                0767         ENDDO
aea29c8517 Alis*0768        ENDDO
7fe6343684 Jean*0769        IF ( (highOrderVorticity.OR.upwindVorticity)
                0770      &     .AND.useAbsVorticity ) THEN
79074ef66b Jean*0771         CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
3370e71df9 Mart*0772      &                         highOrderVorticity, upwindVorticity,
                0773      &                         uFld,omega3,r_hFacZ,
d4c99033f5 Jean*0774      &                         vCf,myThid)
7fe6343684 Jean*0775        ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
79074ef66b Jean*0776         CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
3370e71df9 Mart*0777      &                         highOrderVorticity, upwindVorticity,
                0778      &                         uFld,vort3, r_hFacZ,
3add9696e1 Jean*0779      &                         vCf,myThid)
7fe6343684 Jean*0780        ELSEIF ( useAbsVorticity ) THEN
3370e71df9 Mart*0781         CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0782      &                         uFld,omega3,hFacZ,r_hFacZ,
5d5fff3a81 Alis*0783      &                         vCf,myThid)
                0784        ELSE
3370e71df9 Mart*0785         CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
                0786      &                         uFld,vort3, hFacZ,r_hFacZ,
5d5fff3a81 Alis*0787      &                         vCf,myThid)
                0788        ENDIF
481f592dad Jean*0789        DO j=jMin,jMax
                0790         DO i=iMin,iMax
                0791          gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
                0792         ENDDO
aea29c8517 Alis*0793        ENDDO
                0794 
d09cd10c60 Gael*0795 #ifdef ALLOW_AUTODIFF_TAMC
edb6656069 Mart*0796 CADJ STORE ucf(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
                0797 CADJ STORE vcf(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte
d09cd10c60 Gael*0798 #endif
                0799 
f7d48db11c Jean*0800        IF ( writeDiag ) THEN
ef1e764710 Ed H*0801          IF (snapshot_mdsio) THEN
                0802            CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
                0803            CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
                0804          ENDIF
                0805 #ifdef ALLOW_MNC
                0806          IF (useMNC .AND. snapshot_mnc) THEN
b22b541fe9 Ed H*0807            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
cb356b4c5f Ed H*0808      &          offsets, myThid)
b22b541fe9 Ed H*0809            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
cb356b4c5f Ed H*0810      &          offsets, myThid)
ef1e764710 Ed H*0811          ENDIF
                0812 #endif /*  ALLOW_MNC  */
f7d48db11c Jean*0813        ENDIF
ef1e764710 Ed H*0814 
711329b07b Jean*0815 #ifdef ALLOW_DIAGNOSTICS
                0816        IF ( useDiagnostics ) THEN
                0817          CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
                0818          CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
                0819        ENDIF
                0820 #endif /* ALLOW_DIAGNOSTICS */
cab1a69b8a Jean*0821 
481f592dad Jean*0822 C--   Vertical shear terms (-w*du/dr & -w*dv/dr)
6d4b75df68 Jean*0823        IF ( .NOT. momImplVertAdv ) THEN
31fb0e0e6d Jean*0824         CALL MOM_VI_U_VERTSHEAR( bi, bj, k, deepFacAdv,
                0825      &                           uVel, wVel, uCf, myThid )
6d4b75df68 Jean*0826         DO j=jMin,jMax
                0827          DO i=iMin,iMax
                0828           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
                0829          ENDDO
481f592dad Jean*0830         ENDDO
31fb0e0e6d Jean*0831         CALL MOM_VI_V_VERTSHEAR( bi, bj, k, deepFacAdv,
                0832      &                           vVel, wVel, vCf, myThid )
6d4b75df68 Jean*0833         DO j=jMin,jMax
                0834          DO i=iMin,iMax
                0835           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
                0836          ENDDO
481f592dad Jean*0837         ENDDO
711329b07b Jean*0838 #ifdef ALLOW_DIAGNOSTICS
                0839         IF ( useDiagnostics ) THEN
                0840          CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
                0841          CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
                0842         ENDIF
                0843 #endif /* ALLOW_DIAGNOSTICS */
6d4b75df68 Jean*0844        ENDIF
aea29c8517 Alis*0845 
                0846 C--   Bernoulli term
f1a4cec01a Jean*0847        CALL MOM_VI_U_GRAD_KE(bi,bj,k,KE,uCf,myThid)
481f592dad Jean*0848        DO j=jMin,jMax
                0849         DO i=iMin,iMax
                0850          gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
                0851         ENDDO
aea29c8517 Alis*0852        ENDDO
f1a4cec01a Jean*0853        CALL MOM_VI_V_GRAD_KE(bi,bj,k,KE,vCf,myThid)
481f592dad Jean*0854        DO j=jMin,jMax
                0855         DO i=iMin,iMax
                0856          gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
                0857         ENDDO
                0858        ENDDO
f7d48db11c Jean*0859        IF ( writeDiag ) THEN
ef1e764710 Ed H*0860          IF (snapshot_mdsio) THEN
                0861            CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
                0862            CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
                0863          ENDIF
                0864 #ifdef ALLOW_MNC
                0865          IF (useMNC .AND. snapshot_mnc) THEN
b22b541fe9 Ed H*0866            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
cb356b4c5f Ed H*0867      &          offsets, myThid)
b22b541fe9 Ed H*0868            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
cb356b4c5f Ed H*0869      &          offsets, myThid)
df1ec59c8a Jean*0870          ENDIF
ef1e764710 Ed H*0871 #endif /*  ALLOW_MNC  */
f7d48db11c Jean*0872        ENDIF
                0873 
481f592dad Jean*0874 C--   end if momAdvection
                0875       ENDIF
                0876 
3daafce20b Jean*0877 C--   3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
427e24e121 Jean*0878       IF ( select3dCoriScheme.GE.1 ) THEN
039fe61d35 Jean*0879         CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
                0880         DO j=jMin,jMax
                0881          DO i=iMin,iMax
                0882           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
                0883          ENDDO
                0884         ENDDO
427e24e121 Jean*0885        IF ( usingCurvilinearGrid .OR. rotateGrid ) THEN
039fe61d35 Jean*0886 C-     presently, non zero angleSinC array only supported with Curvilinear-Grid
                0887         CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
                0888         DO j=jMin,jMax
                0889          DO i=iMin,iMax
                0890           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
                0891          ENDDO
                0892         ENDDO
                0893        ENDIF
                0894       ENDIF
                0895 
                0896 C--   Non-Hydrostatic (spherical) metric terms
31fb0e0e6d Jean*0897 #ifdef MOM_USE_OLD_DEEP_VERT_ADV
039fe61d35 Jean*0898       IF ( useNHMTerms ) THEN
31fb0e0e6d Jean*0899 #else
                0900       IF ( useNHMTerms .AND. .NOT.deepAtmosphere ) THEN
                0901 #endif
039fe61d35 Jean*0902        CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
                0903        CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
                0904        DO j=jMin,jMax
                0905         DO i=iMin,iMax
31fb0e0e6d Jean*0906          gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) + uCf(i,j)
                0907          gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) + vCf(i,j)
039fe61d35 Jean*0908         ENDDO
                0909        ENDDO
31fb0e0e6d Jean*0910 #ifdef ALLOW_DIAGNOSTICS
                0911        IF ( useDiagnostics ) THEN
                0912         CALL DIAGNOSTICS_FILL( uCf, 'Um_Metr ', k,1,2, bi,bj, myThid )
                0913         CALL DIAGNOSTICS_FILL( vCf, 'Vm_Metr ', k,1,2, bi,bj, myThid )
                0914        ENDIF
                0915 #endif /* ALLOW_DIAGNOSTICS */
039fe61d35 Jean*0916       ENDIF
df1ec59c8a Jean*0917 
481f592dad Jean*0918 C--   Set du/dt & dv/dt on boundaries to zero
aea29c8517 Alis*0919       DO j=jMin,jMax
                0920        DO i=iMin,iMax
481f592dad Jean*0921         gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
                0922         gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
aea29c8517 Alis*0923        ENDDO
                0924       ENDDO
481f592dad Jean*0925 
5751fc37e0 Jean*0926 #ifdef ALLOW_DEBUG
8830b8f970 Jean*0927       IF ( debugLevel .GE. debLevC
5751fc37e0 Jean*0928      &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0
                0929      &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
                0930      &   .AND. useCubedSphereExchange ) THEN
e85db1faec Jean*0931         CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
4667e4066d Jean*0932      &             guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
5751fc37e0 Jean*0933       ENDIF
                0934 #endif /* ALLOW_DEBUG */
aea29c8517 Alis*0935 
f7d48db11c Jean*0936       IF ( writeDiag ) THEN
ed2dbbe83d Jean*0937         IF (useBiharmonicVisc) THEN
                0938          CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
                0939      &                        bi,bj,k, myIter, myThid )
                0940          CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
                0941      &                        bi,bj,k, myIter, myThid )
                0942          CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
                0943      &                        bi,bj,k, myIter, myThid )
                0944          CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
                0945      &                        bi,bj,k, myIter, myThid )
                0946         ENDIF
ef1e764710 Ed H*0947         IF (snapshot_mdsio) THEN
df1ec59c8a Jean*0948          CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
ed2dbbe83d Jean*0949          CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
df1ec59c8a Jean*0950          CALL WRITE_LOCAL_RL('KE','I10',1,KE,     bi,bj,k,myIter,myThid)
                0951          CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv,   bi,bj,k,myIter,myThid)
                0952          CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
ed2dbbe83d Jean*0953          CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
                0954      &                        bi,bj,k, myIter, myThid )
df1ec59c8a Jean*0955          CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
                0956          CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
ef1e764710 Ed H*0957         ENDIF
                0958 #ifdef ALLOW_MNC
                0959         IF (useMNC .AND. snapshot_mnc) THEN
df1ec59c8a Jean*0960           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
ed2dbbe83d Jean*0961      &          offsets, myThid)
                0962           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
df1ec59c8a Jean*0963      &          offsets, myThid)
                0964           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
                0965      &          offsets, myThid)
                0966           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
cb356b4c5f Ed H*0967      &          offsets, myThid)
b22b541fe9 Ed H*0968           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
ed2dbbe83d Jean*0969      &          offsets, myThid)
                0970           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
cb356b4c5f Ed H*0971      &          offsets, myThid)
b22b541fe9 Ed H*0972           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
cb356b4c5f Ed H*0973      &          offsets, myThid)
b22b541fe9 Ed H*0974           CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
cb356b4c5f Ed H*0975      &          offsets, myThid)
ef1e764710 Ed H*0976         ENDIF
                0977 #endif /*  ALLOW_MNC  */
3a279374db Alis*0978       ENDIF
3add9696e1 Jean*0979 
711329b07b Jean*0980 #ifdef ALLOW_DIAGNOSTICS
                0981       IF ( useDiagnostics ) THEN
ed2dbbe83d Jean*0982         CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
df1ec59c8a Jean*0983         CALL DIAGNOSTICS_FILL(KE,     'momKE   ',k,1,2,bi,bj,myThid)
711329b07b Jean*0984        IF (momViscosity) THEN
df1ec59c8a Jean*0985         CALL DIAGNOSTICS_FILL(hDiv,   'momHDiv ',k,1,2,bi,bj,myThid)
ed2dbbe83d Jean*0986        ENDIF
                0987        IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
                0988         CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
                0989         CALL DIAGNOSTICS_FILL(strainBC,'Strain  ',k,1,2,bi,bj,myThid)
f59d76b0dd Ed D*0990 C         stretching will be zero, unless using QG Leith
a125ab7be7 Jean*0991         IF ( viscC2LeithQG.NE.zeroRL ) THEN
f59d76b0dd Ed D*0992           CALL DIAGNOSTICS_FILL(stretching,
                0993      I                          'Stretch ',k,1,2,bi,bj,myThid)
                0994         ENDIF
711329b07b Jean*0995        ENDIF
07e17fa184 Jean*0996         CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
df1ec59c8a Jean*0997      &                                'Um_Advec',k,1,2,bi,bj,myThid)
07e17fa184 Jean*0998         CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
df1ec59c8a Jean*0999      &                                'Vm_Advec',k,1,2,bi,bj,myThid)
711329b07b Jean*1000       ENDIF
                1001 #endif /* ALLOW_DIAGNOSTICS */
                1002 
99bc607d7a Ed H*1003 #endif /* ALLOW_MOM_VECINV */
cab1a69b8a Jean*1004 
aea29c8517 Alis*1005       RETURN
                1006       END