/[MITgcm]/MITgcm/pkg/mom_vecinv/mom_vecinv.F
ViewVC logotype

Diff of /MITgcm/pkg/mom_vecinv/mom_vecinv.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.11 by edhill, Tue Nov 4 19:51:54 2003 UTC revision 1.63 by jmc, Tue Mar 16 00:16:50 2010 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2  C $Name$  C $Name$
3    
4  #include "PACKAGES_CONFIG.h"  #include "MOM_VECINV_OPTIONS.h"
 #include "CPP_OPTIONS.h"  
5    
6        SUBROUTINE MOM_VECINV(        SUBROUTINE MOM_VECINV(
7       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
8       I        dPhiHydX,dPhiHydY,KappaRU,KappaRV,       I        KappaRU, KappaRV,
9       U        fVerU, fVerV,       U        fVerU, fVerV,
10       I        myCurrentTime, myIter, myThid)       O        guDiss, gvDiss,
11         I        myTime, myIter, myThid)
12  C     /==========================================================\  C     /==========================================================\
13  C     | S/R MOM_VECINV                                           |  C     | S/R MOM_VECINV                                           |
14  C     | o Form the right hand-side of the momentum equation.     |  C     | o Form the right hand-side of the momentum equation.     |
# Line 31  C     == Global variables == Line 31  C     == Global variables ==
31  #include "DYNVARS.h"  #include "DYNVARS.h"
32  #include "EEPARAMS.h"  #include "EEPARAMS.h"
33  #include "PARAMS.h"  #include "PARAMS.h"
34    #ifdef ALLOW_MNC
35    #include "MNC_PARAMS.h"
36    #endif
37  #include "GRID.h"  #include "GRID.h"
38  #ifdef ALLOW_TIMEAVE  #ifdef ALLOW_TIMEAVE
39  #include "TIMEAVE_STATV.h"  #include "TIMEAVE_STATV.h"
40  #endif  #endif
41    #ifdef ALLOW_AUTODIFF_TAMC
42    # include "tamc.h"
43    # include "tamc_keys.h"
44    #endif
45    
46  C     == Routine arguments ==  C     == Routine arguments ==
47  C     fVerU   - Flux of momentum in the vertical  C     fVerU  :: Flux of momentum in the vertical direction, out of the upper
48  C     fVerV     direction out of the upper face of a cell K  C     fVerV  :: face of a cell K ( flux into the cell above ).
49  C               ( flux into the cell above ).  C     guDiss :: dissipation tendency (all explicit terms), u component
50  C     dPhiHydX,Y :: Gradient (X & Y dir.) of Hydrostatic Potential  C     gvDiss :: dissipation tendency (all explicit terms), v component
51  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
52  C                                      results will be set.  C                                      results will be set.
53  C     kUp, kDown                     - Index for upper and lower layers.  C     kUp, kDown                     - Index for upper and lower layers.
54  C     myThid - Instance number for this innvocation of CALC_MOM_RHS  C     myThid :: my Thread Id number
       _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)  
       _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)  
55        _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
56        _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
57        _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
58        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
59          _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
60          _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61        INTEGER kUp,kDown        INTEGER kUp,kDown
62        _RL     myCurrentTime        _RL     myTime
63        INTEGER myIter        INTEGER myIter
64        INTEGER myThid        INTEGER myThid
65        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
# Line 64  C     == Functions == Line 71  C     == Functions ==
71        EXTERNAL DIFFERENT_MULTIPLE        EXTERNAL DIFFERENT_MULTIPLE
72    
73  C     == Local variables ==  C     == Local variables ==
       _RL      aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
74        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75        _RL      vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76        _RL      uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77        _RL      vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78        _RL      mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS hFacZ   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79        _RL      pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80        _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81        _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82        _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2u   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83        _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2v   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84        _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85        _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL zStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86        _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87        _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL strain  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88        _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KE      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL omega3  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90        _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vort3   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL hDiv    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94  C     I,J,K - Loop counters        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95          _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96    C     i,j,k  :: Loop counters
97        INTEGER i,j,k        INTEGER i,j,k
 C     rVelMaskOverride - Factor for imposing special surface boundary conditions  
 C                        ( set according to free-surface condition ).  
 C     hFacROpen        - Lopped cell factos used tohold fraction of open  
 C     hFacRClosed        and closed cell wall.  
       _RL  rVelMaskOverride  
98  C     xxxFac - On-off tracer parameters used for switching terms off.  C     xxxFac - On-off tracer parameters used for switching terms off.
       _RL  uDudxFac  
       _RL  AhDudxFac  
       _RL  A4DuxxdxFac  
       _RL  vDudyFac  
       _RL  AhDudyFac  
       _RL  A4DuyydyFac  
       _RL  rVelDudrFac  
99        _RL  ArDudrFac        _RL  ArDudrFac
       _RL  fuFac  
       _RL  phxFac  
       _RL  mtFacU  
       _RL  uDvdxFac  
       _RL  AhDvdxFac  
       _RL  A4DvxxdxFac  
       _RL  vDvdyFac  
       _RL  AhDvdyFac  
       _RL  A4DvyydyFac  
       _RL  rVelDvdrFac  
100        _RL  ArDvdrFac        _RL  ArDvdrFac
101        _RL  fvFac        _RL  sideMaskFac
       _RL  phyFac  
       _RL  vForcFac  
       _RL  mtFacV  
       _RL wVelBottomOverride  
102        LOGICAL bottomDragTerms        LOGICAL bottomDragTerms
103        _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        LOGICAL writeDiag
104        _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        LOGICAL harmonic,biharmonic,useVariableViscosity
105        _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  #ifdef ALLOW_AUTODIFF_TAMC
106        _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        INTEGER imomkey
107    #endif
108    
109    #ifdef ALLOW_MNC
110          INTEGER offsets(9)
111          CHARACTER*(1) pf
112    #endif
113    
114  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
115  C--   only the kDown part of fverU/V is set in this subroutine  C--   only the kDown part of fverU/V is set in this subroutine
# Line 132  C--   (at least in part) Line 120  C--   (at least in part)
120        fVerV(1,1,kUp) = fVerV(1,1,kUp)        fVerV(1,1,kUp) = fVerV(1,1,kUp)
121  #endif  #endif
122    
123        rVelMaskOverride=1.  #ifdef ALLOW_AUTODIFF_TAMC
124        IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac            act0 = k - 1
125        wVelBottomOverride=1.            max0 = Nr
126        IF (k.EQ.Nr) wVelBottomOverride=0.            act1 = bi - myBxLo(myThid)
127              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
128  C     Initialise intermediate terms            act2 = bj - myByLo(myThid)
129        DO J=1-OLy,sNy+OLy            max2 = myByHi(myThid) - myByLo(myThid) + 1
130         DO I=1-OLx,sNx+OLx            act3 = myThid - 1
131          aF(i,j)   = 0.            max3 = nTx*nTy
132          vF(i,j)   = 0.            act4 = ikey_dynamics - 1
133          vrF(i,j)  = 0.            imomkey = (act0 + 1)
134         &                    + act1*max0
135         &                    + act2*max0*max1
136         &                    + act3*max0*max1*max2
137         &                    + act4*max0*max1*max2*max3
138    #endif /* ALLOW_AUTODIFF_TAMC */
139    
140          writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
141    
142    #ifdef ALLOW_MNC
143          IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
144            IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
145              pf(1:1) = 'D'
146            ELSE
147              pf(1:1) = 'R'
148            ENDIF
149            IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
150              CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
151              CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
152              CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
153              CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
154            ENDIF
155            DO i = 1,9
156              offsets(i) = 0
157            ENDDO
158            offsets(3) = k
159    c       write(*,*) 'offsets = ',(offsets(i),i=1,9)
160          ENDIF
161    #endif /*  ALLOW_MNC  */
162    
163    C--   Initialise intermediate terms
164          DO j=1-OLy,sNy+OLy
165           DO i=1-OLx,sNx+OLx
166            vF(i,j)    = 0.
167            vrF(i,j)   = 0.
168          uCf(i,j)   = 0.          uCf(i,j)   = 0.
169          vCf(i,j)   = 0.          vCf(i,j)   = 0.
         mT(i,j)   = 0.  
         pF(i,j)   = 0.  
170          del2u(i,j) = 0.          del2u(i,j) = 0.
171          del2v(i,j) = 0.          del2v(i,j) = 0.
172          dStar(i,j) = 0.          dStar(i,j) = 0.
173          zStar(i,j) = 0.          zStar(i,j) = 0.
174          uDiss(i,j) = 0.          guDiss(i,j)= 0.
175          vDiss(i,j) = 0.          gvDiss(i,j)= 0.
176          vort3(i,j) = 0.          vort3(i,j) = 0.
177          omega3(i,j) = 0.          omega3(i,j)= 0.
178          ke(i,j) = 0.          KE(i,j)    = 0.
179  #ifdef ALLOW_AUTODIFF_TAMC  C-    need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
180            hDiv(i,j)  = 0.
181            viscAh_Z(i,j) = 0.
182            viscAh_D(i,j) = 0.
183            viscA4_Z(i,j) = 0.
184            viscA4_D(i,j) = 0.
185    
186          strain(i,j)  = 0. _d 0          strain(i,j)  = 0. _d 0
187          tension(i,j) = 0. _d 0          tension(i,j) = 0. _d 0
188    #ifdef ALLOW_AUTODIFF_TAMC
189            hFacZ(i,j)   = 0. _d 0
190  #endif  #endif
191         ENDDO         ENDDO
192        ENDDO        ENDDO
193    
194  C--   Term by term tracer parmeters  C--   Term by term tracer parmeters
195  C     o U momentum equation  C     o U momentum equation
       uDudxFac     = afFacMom*1.  
       AhDudxFac    = vfFacMom*1.  
       A4DuxxdxFac  = vfFacMom*1.  
       vDudyFac     = afFacMom*1.  
       AhDudyFac    = vfFacMom*1.  
       A4DuyydyFac  = vfFacMom*1.  
       rVelDudrFac  = afFacMom*1.  
196        ArDudrFac    = vfFacMom*1.        ArDudrFac    = vfFacMom*1.
       mTFacU       = mtFacMom*1.  
       fuFac        = cfFacMom*1.  
       phxFac       = pfFacMom*1.  
197  C     o V momentum equation  C     o V momentum equation
       uDvdxFac     = afFacMom*1.  
       AhDvdxFac    = vfFacMom*1.  
       A4DvxxdxFac  = vfFacMom*1.  
       vDvdyFac     = afFacMom*1.  
       AhDvdyFac    = vfFacMom*1.  
       A4DvyydyFac  = vfFacMom*1.  
       rVelDvdrFac  = afFacMom*1.  
198        ArDvdrFac    = vfFacMom*1.        ArDvdrFac    = vfFacMom*1.
199        mTFacV       = mtFacMom*1.  
200        fvFac        = cfFacMom*1.  C note: using standard stencil (no mask) results in under-estimating
201        phyFac       = pfFacMom*1.  C       vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
202        vForcFac     = foFacMom*1.        IF ( no_slip_sides ) THEN
203            sideMaskFac = sideDragFactor
204          ELSE
205            sideMaskFac = 0. _d 0
206          ENDIF
207    
208        IF (     no_slip_bottom        IF (     no_slip_bottom
209       &    .OR. bottomDragQuadratic.NE.0.       &    .OR. bottomDragQuadratic.NE.0.
# Line 198  C     o V momentum equation Line 213  C     o V momentum equation
213         bottomDragTerms=.FALSE.         bottomDragTerms=.FALSE.
214        ENDIF        ENDIF
215    
 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP  
       IF (staggerTimeStep) THEN  
         phxFac = 0.  
         phyFac = 0.  
       ENDIF  
   
216  C--   Calculate open water fraction at vorticity points  C--   Calculate open water fraction at vorticity points
217        CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)        CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
218    
 C---- Calculate common quantities used in both U and V equations  
 C     Calculate tracer cell face open areas  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         xA(i,j) = _dyG(i,j,bi,bj)  
      &   *drF(k)*_hFacW(i,j,k,bi,bj)  
         yA(i,j) = _dxG(i,j,bi,bj)  
      &   *drF(k)*_hFacS(i,j,k,bi,bj)  
        ENDDO  
       ENDDO  
   
219  C     Make local copies of horizontal flow field  C     Make local copies of horizontal flow field
220        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
221         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
# Line 230  C note (jmc) : Dissipation and Vort3 adv Line 228  C note (jmc) : Dissipation and Vort3 adv
228  C              use the same maskZ (and hFacZ)  => needs 2 call(s)  C              use the same maskZ (and hFacZ)  => needs 2 call(s)
229  c     CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)  c     CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
230    
231        CALL MOM_VI_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)        CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
232    
233          CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
234    
235          IF (momViscosity) THEN
236    C--    For viscous term, compute horizontal divergence, tension & strain
237    C      and mask relative vorticity (free-slip case):
238    
239    #ifdef ALLOW_AUTODIFF_TAMC
240    CADJ STORE vort3(:,:) =
241    CADJ &     comlev1_bibj_k, key = imomkey, byte = isbyte
242    #endif
243    
244        CALL MOM_VI_CALC_HDIV(bi,bj,k,uFld,vFld,hDiv,myThid)         CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
245    
246        CALL MOM_VI_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)         CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
247    
248  c     CALL MOM_VI_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)         CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
249    
250    C-     account for no-slip / free-slip BC:
251           DO j=1-Oly,sNy+Oly
252            DO i=1-Olx,sNx+Olx
253              IF ( hFacZ(i,j).EQ.0. ) THEN
254                vort3(i,j)  = sideMaskFac*vort3(i,j)
255                strain(i,j) = sideMaskFac*strain(i,j)
256              ENDIF
257            ENDDO
258           ENDDO
259    
260    C--    Calculate Viscosities
261           CALL MOM_CALC_VISC(
262         I        bi,bj,k,
263         O        viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
264         O        harmonic,biharmonic,useVariableViscosity,
265         I        hDiv,vort3,tension,strain,KE,hfacZ,
266         I        myThid)
267    
       IF (momViscosity) THEN  
268  C      Calculate del^2 u and del^2 v for bi-harmonic term  C      Calculate del^2 u and del^2 v for bi-harmonic term
269         IF (viscA4.NE.0.) THEN         IF (biharmonic) THEN
270           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
271       O                      del2u,del2v,       O                      del2u,del2v,
272       &                      myThid)       &                      myThid)
273           CALL MOM_VI_CALC_HDIV(bi,bj,k,del2u,del2v,dStar,myThid)           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
274           CALL MOM_VI_CALC_RELVORT3(           CALL MOM_CALC_RELVORT3(bi,bj,k,
275       &                         bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)       &                          del2u,del2v,hFacZ,zStar,myThid)
276         ENDIF         ENDIF
277  C      Calculate dissipation terms for U and V equations  
278  C      in terms of vorticity and divergence  C-    Strain diagnostics:
279         IF (viscAh.NE.0. .OR. viscA4.NE.0.) THEN         IF ( writeDiag ) THEN
280           CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,          IF (snapshot_mdsio) THEN
281       O                       uDiss,vDiss,            CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
282       &                       myThid)          ENDIF
283    #ifdef ALLOW_MNC
284            IF (useMNC .AND. snapshot_mnc) THEN
285              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
286         &          offsets, myThid)
287            ENDIF
288    #endif /*  ALLOW_MNC  */
289           ENDIF
290    #ifdef ALLOW_DIAGNOSTICS
291           IF ( useDiagnostics ) THEN
292            CALL DIAGNOSTICS_FILL(strain, 'Strain  ',k,1,2,bi,bj,myThid)
293         ENDIF         ENDIF
294  C      or in terms of tension and strain  #endif /* ALLOW_DIAGNOSTICS */
295         IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN  
296           CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,  C---   Calculate dissipation terms for U and V equations
297       O                         tension,  
298       I                         myThid)  C      in terms of tension and strain
299           CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,         IF (useStrainTensionVisc) THEN
300       O                        strain,  C        mask strain as if free-slip since side-drag is computed separately
301       I                        myThid)           DO j=1-Oly,sNy+Oly
302           CALL MOM_HDISSIP(bi,bj,k,            DO i=1-Olx,sNx+Olx
303       I                    tension,strain,hFacZ,viscAtension,viscAstrain,              IF ( hFacZ(i,j).EQ.0. ) strain(i,j) = 0. _d 0
304       O                    uDiss,vDiss,            ENDDO
305             ENDDO
306             CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
307         I                    hFacZ,
308         I                    viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
309         I                    harmonic,biharmonic,useVariableViscosity,
310         O                    guDiss,gvDiss,
311       I                    myThid)       I                    myThid)
312           ELSE
313    C      in terms of vorticity and divergence
314             CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
315         I                       hFacZ,dStar,zStar,
316         I                       viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
317         I                       harmonic,biharmonic,useVariableViscosity,
318         O                       guDiss,gvDiss,
319         &                       myThid)
320         ENDIF         ENDIF
321    C--   if (momViscosity) end of block.
322        ENDIF        ENDIF
323    
324  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:
325  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
326    
327  C---- Zonal momentum equation starts here  C---  Other dissipation terms in Zonal momentum equation
328    
329  C--   Vertical flux (fVer is at upper face of "u" cell)  C--   Vertical flux (fVer is at upper face of "u" cell)
330    
331  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
332        IF (momViscosity.AND..NOT.implicitViscosity)        IF (momViscosity.AND..NOT.implicitViscosity) THEN
333       & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
334    
335  C     Combine fluxes  C     Combine fluxes
336        DO j=jMin,jMax         DO j=jMin,jMax
337         DO i=iMin,iMax          DO i=iMin,iMax
338          fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)           fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
339            ENDDO
340         ENDDO         ENDDO
       ENDDO  
341    
342  C--   Tendency is minus divergence of the fluxes + coriolis + pressure term  C--   Tendency is minus divergence of the fluxes
343        DO j=2-Oly,sNy+Oly-1         DO j=2-Oly,sNy+Oly-1
344         DO i=2-Olx,sNx+Olx-1          DO i=2-Olx,sNx+Olx-1
345          gU(i,j,k,bi,bj) = uDiss(i,j)           guDiss(i,j) = guDiss(i,j)
346       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
347       &   *recip_rAw(i,j,bi,bj)       &   *recip_rAw(i,j,bi,bj)
348       &  *(       &  *(
349       &   +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac       &    fVerU(i,j,kDown) - fVerU(i,j,kUp)
350       &   )       &   )*rkSign
351       &  - phxFac*dPhiHydX(i,j)          ENDDO
352         ENDDO         ENDDO
353        ENDDO        ENDIF
354    
355  C-- No-slip and drag BCs appear as body forces in cell abutting topography  C-- No-slip and drag BCs appear as body forces in cell abutting topography
356        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
357  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
358         CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)         CALL MOM_U_SIDEDRAG(
359         I        bi,bj,k,
360         I        uFld, del2u, hFacZ,
361         I        viscAh_Z,viscA4_Z,
362         I        harmonic,biharmonic,useVariableViscosity,
363         O        vF,
364         I        myThid)
365         DO j=jMin,jMax         DO j=jMin,jMax
366          DO i=iMin,iMax          DO i=iMin,iMax
367           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)           guDiss(i,j) = guDiss(i,j)+vF(i,j)
368          ENDDO          ENDDO
369         ENDDO         ENDDO
370        ENDIF        ENDIF
   
371  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
372        IF (momViscosity.AND.bottomDragTerms) THEN        IF (momViscosity.AND.bottomDragTerms) THEN
373         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
374         DO j=jMin,jMax         DO j=jMin,jMax
375          DO i=iMin,iMax          DO i=iMin,iMax
376           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)           guDiss(i,j) = guDiss(i,j)+vF(i,j)
377          ENDDO          ENDDO
378         ENDDO         ENDDO
379        ENDIF        ENDIF
380    #ifdef ALLOW_SHELFICE
381          IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
382           CALL SHELFICE_U_DRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
383           DO j=jMin,jMax
384            DO i=iMin,iMax
385             guDiss(i,j) = guDiss(i,j) + vF(i,j)
386            ENDDO
387           ENDDO
388          ENDIF
389    #endif /* ALLOW_SHELFICE */
390    
 C--   Metric terms for curvilinear grid systems  
 c     IF (usingSphericalPolarMTerms) THEN  
 C      o Spherical polar grid metric terms  
 c      CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)  
 c      DO j=jMin,jMax  
 c       DO i=iMin,iMax  
 c        gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)  
 c       ENDDO  
 c      ENDDO  
 c     ENDIF  
391    
392  C---- Meridional momentum equation starts here  C---  Other dissipation terms in Meridional momentum equation
393    
394  C--   Vertical flux (fVer is at upper face of "v" cell)  C--   Vertical flux (fVer is at upper face of "v" cell)
395    
396  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
397        IF (momViscosity.AND..NOT.implicitViscosity)        IF (momViscosity.AND..NOT.implicitViscosity) THEN
398       & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
399    
400  C     Combine fluxes -> fVerV  C     Combine fluxes -> fVerV
401        DO j=jMin,jMax         DO j=jMin,jMax
402         DO i=iMin,iMax          DO i=iMin,iMax
403          fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)           fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
404            ENDDO
405         ENDDO         ENDDO
       ENDDO  
406    
407  C--   Tendency is minus divergence of the fluxes + coriolis + pressure term  C--   Tendency is minus divergence of the fluxes
408        DO j=jMin,jMax         DO j=jMin,jMax
409         DO i=iMin,iMax          DO i=iMin,iMax
410          gV(i,j,k,bi,bj) = vDiss(i,j)           gvDiss(i,j) = gvDiss(i,j)
411       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
412       &    *recip_rAs(i,j,bi,bj)       &    *recip_rAs(i,j,bi,bj)
413       &  *(       &  *(
414       &   +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac       &    fVerV(i,j,kDown) - fVerV(i,j,kUp)
415       &   )       &   )*rkSign
416       &  - phyFac*dPhiHydY(i,j)          ENDDO
417         ENDDO         ENDDO
418        ENDDO        ENDIF
419    
420  C-- No-slip and drag BCs appear as body forces in cell abutting topography  C-- No-slip and drag BCs appear as body forces in cell abutting topography
421        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
422  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
423         CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)         CALL MOM_V_SIDEDRAG(
424         I        bi,bj,k,
425         I        vFld, del2v, hFacZ,
426         I        viscAh_Z,viscA4_Z,
427         I        harmonic,biharmonic,useVariableViscosity,
428         O        vF,
429         I        myThid)
430         DO j=jMin,jMax         DO j=jMin,jMax
431          DO i=iMin,iMax          DO i=iMin,iMax
432           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
433          ENDDO          ENDDO
434         ENDDO         ENDDO
435        ENDIF        ENDIF
# Line 376  C-    No-slip BCs impose a drag at botto Line 438  C-    No-slip BCs impose a drag at botto
438         CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)         CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
439         DO j=jMin,jMax         DO j=jMin,jMax
440          DO i=iMin,iMax          DO i=iMin,iMax
441           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
442          ENDDO          ENDDO
443         ENDDO         ENDDO
444        ENDIF        ENDIF
445    #ifdef ALLOW_SHELFICE
446          IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
447             CALL SHELFICE_V_DRAG(bi,bj,k,vFld,KE,KappaRU,vF,myThid)
448             DO j=jMin,jMax
449              DO i=iMin,iMax
450               gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
451              ENDDO
452             ENDDO
453            ENDIF
454    #endif /* ALLOW_SHELFICE */
455    
456    
457    C-    Vorticity diagnostics:
458          IF ( writeDiag ) THEN
459            IF (snapshot_mdsio) THEN
460              CALL WRITE_LOCAL_RL('Z3','I10',1,vort3, bi,bj,k,myIter,myThid)
461            ENDIF
462    #ifdef ALLOW_MNC
463            IF (useMNC .AND. snapshot_mnc) THEN
464              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
465         &          offsets, myThid)
466            ENDIF
467    #endif /*  ALLOW_MNC  */
468          ENDIF
469    #ifdef ALLOW_DIAGNOSTICS
470          IF ( useDiagnostics ) THEN
471            CALL DIAGNOSTICS_FILL(vort3,  'momVort3',k,1,2,bi,bj,myThid)
472          ENDIF
473    #endif /* ALLOW_DIAGNOSTICS */
474    
475  C--   Metric terms for curvilinear grid systems  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
476  c     IF (usingSphericalPolarMTerms) THEN  
477  C      o Spherical polar grid metric terms  C---  Prepare for Advection & Coriolis terms:
478  c      CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)  C-    Mask relative vorticity and calculate absolute vorticity
479  c      DO j=jMin,jMax        DO j=1-Oly,sNy+Oly
480  c       DO i=iMin,iMax         DO i=1-Olx,sNx+Olx
481  c        gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)           IF ( hFacZ(i,j).EQ.0. ) vort3(i,j) = 0.
482  c       ENDDO         ENDDO
483  c      ENDDO        ENDDO
484  c     ENDIF        IF (useAbsVorticity)
485         &  CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
486    
487  C--   Horizontal Coriolis terms  C--   Horizontal Coriolis terms
488        IF (useCoriolis .AND. .NOT.useCDscheme) THEN  c     IF (useCoriolis .AND. .NOT.useCDscheme
489         CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,omega3,hFacZ,r_hFacZ,  c    &    .AND. .NOT. useAbsVorticity) THEN
490       &                      uCf,vCf,myThid)  C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
491          IF ( useCoriolis .AND.
492         &     .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
493         &   ) THEN
494           IF (useAbsVorticity) THEN
495            CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
496         &                         uCf,myThid)
497            CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
498         &                         vCf,myThid)
499           ELSE
500            CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
501         &                       uCf,vCf,myThid)
502           ENDIF
503         DO j=jMin,jMax         DO j=jMin,jMax
504          DO i=iMin,iMax          DO i=iMin,iMax
505           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)           gU(i,j,k,bi,bj) = uCf(i,j)
506           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)           gV(i,j,k,bi,bj) = vCf(i,j)
507            ENDDO
508           ENDDO
509           IF ( writeDiag ) THEN
510             IF (snapshot_mdsio) THEN
511               CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
512               CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
513             ENDIF
514    #ifdef ALLOW_MNC
515             IF (useMNC .AND. snapshot_mnc) THEN
516               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
517         &          offsets, myThid)
518               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
519         &          offsets, myThid)
520             ENDIF
521    #endif /*  ALLOW_MNC  */
522           ENDIF
523    #ifdef ALLOW_DIAGNOSTICS
524           IF ( useDiagnostics ) THEN
525             CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
526             CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
527           ENDIF
528    #endif /* ALLOW_DIAGNOSTICS */
529          ELSE
530           DO j=jMin,jMax
531            DO i=iMin,iMax
532             gU(i,j,k,bi,bj) = 0. _d 0
533             gV(i,j,k,bi,bj) = 0. _d 0
534          ENDDO          ENDDO
535         ENDDO         ENDDO
536        ENDIF        ENDIF
537    
538        IF (momAdvection) THEN        IF (momAdvection) THEN
539  C--   Horizontal advection of relative vorticity  C--   Horizontal advection of relative (or absolute) vorticity
540  c      CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,r_hFacZ,uCf,myThid)         IF ( (highOrderVorticity.OR.upwindVorticity)
541         CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3,hFacZ,r_hFacZ,       &     .AND.useAbsVorticity ) THEN
542       &                        uCf,myThid)          CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
543  c      CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid)       &                         uCf,myThid)
544           ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
545            CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
546         &                         uCf,myThid)
547           ELSEIF ( useAbsVorticity ) THEN
548            CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
549         &                         uCf,myThid)
550           ELSE
551            CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
552         &                         uCf,myThid)
553           ENDIF
554         DO j=jMin,jMax         DO j=jMin,jMax
555          DO i=iMin,iMax          DO i=iMin,iMax
556           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
557          ENDDO          ENDDO
558         ENDDO         ENDDO
559  c      CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,r_hFacZ,vCf,myThid)         IF ( (highOrderVorticity.OR.upwindVorticity)
560         CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3,hFacZ,r_hFacZ,       &     .AND.useAbsVorticity ) THEN
561       &                        vCf,myThid)          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
562  c      CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid)       &                         vCf,myThid)
563           ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
564            CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
565         &                         vCf,myThid)
566           ELSEIF ( useAbsVorticity ) THEN
567            CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
568         &                         vCf,myThid)
569           ELSE
570            CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
571         &                         vCf,myThid)
572           ENDIF
573         DO j=jMin,jMax         DO j=jMin,jMax
574          DO i=iMin,iMax          DO i=iMin,iMax
575           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
576          ENDDO          ENDDO
577         ENDDO         ENDDO
578    
579           IF ( writeDiag ) THEN
580             IF (snapshot_mdsio) THEN
581               CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
582               CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
583             ENDIF
584    #ifdef ALLOW_MNC
585             IF (useMNC .AND. snapshot_mnc) THEN
586               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
587         &          offsets, myThid)
588               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
589         &          offsets, myThid)
590             ENDIF
591    #endif /*  ALLOW_MNC  */
592           ENDIF
593    
594  #ifdef ALLOW_TIMEAVE  #ifdef ALLOW_TIMEAVE
595         IF (taveFreq.GT.0.) THEN         IF (taveFreq.GT.0.) THEN
596           CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,           CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
# Line 432  c      CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K Line 598  c      CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K
598           CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,           CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
599       &                           Nr, k, bi, bj, myThid)       &                           Nr, k, bi, bj, myThid)
600         ENDIF         ENDIF
601  #endif  #endif /* ALLOW_TIMEAVE */
602    #ifdef ALLOW_DIAGNOSTICS
603           IF ( useDiagnostics ) THEN
604             CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
605             CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
606           ENDIF
607    #endif /* ALLOW_DIAGNOSTICS */
608    
609  C--   Vertical shear terms (-w*du/dr & -w*dv/dr)  C--   Vertical shear terms (-w*du/dr & -w*dv/dr)
610         CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)         IF ( .NOT. momImplVertAdv ) THEN
611            CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
612            DO j=jMin,jMax
613             DO i=iMin,iMax
614              gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
615             ENDDO
616            ENDDO
617            CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
618            DO j=jMin,jMax
619             DO i=iMin,iMax
620              gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
621             ENDDO
622            ENDDO
623    #ifdef ALLOW_DIAGNOSTICS
624            IF ( useDiagnostics ) THEN
625             CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
626             CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
627            ENDIF
628    #endif /* ALLOW_DIAGNOSTICS */
629           ENDIF
630    
631    C--   Bernoulli term
632           CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
633         DO j=jMin,jMax         DO j=jMin,jMax
634          DO i=iMin,iMax          DO i=iMin,iMax
635           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
636          ENDDO          ENDDO
637         ENDDO         ENDDO
638         CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)         CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
639         DO j=jMin,jMax         DO j=jMin,jMax
640          DO i=iMin,iMax          DO i=iMin,iMax
641           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
642          ENDDO          ENDDO
643         ENDDO         ENDDO
644           IF ( writeDiag ) THEN
645             IF (snapshot_mdsio) THEN
646               CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
647               CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
648             ENDIF
649    #ifdef ALLOW_MNC
650             IF (useMNC .AND. snapshot_mnc) THEN
651               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
652         &          offsets, myThid)
653               CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
654         &          offsets, myThid)
655             ENDIF
656    #endif /*  ALLOW_MNC  */
657           ENDIF
658    
659  C--   Bernoulli term  C--   end if momAdvection
660         CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)        ENDIF
661    
662    C--   3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
663          IF ( use3dCoriolis ) THEN
664            CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
665            DO j=jMin,jMax
666             DO i=iMin,iMax
667              gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
668             ENDDO
669            ENDDO
670           IF ( usingCurvilinearGrid ) THEN
671    C-     presently, non zero angleSinC array only supported with Curvilinear-Grid
672            CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
673            DO j=jMin,jMax
674             DO i=iMin,iMax
675              gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
676             ENDDO
677            ENDDO
678           ENDIF
679          ENDIF
680    
681    C--   Non-Hydrostatic (spherical) metric terms
682          IF ( useNHMTerms ) THEN
683           CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
684         DO j=jMin,jMax         DO j=jMin,jMax
685          DO i=iMin,iMax          DO i=iMin,iMax
686           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
687          ENDDO          ENDDO
688         ENDDO         ENDDO
689         CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)         CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
690         DO j=jMin,jMax         DO j=jMin,jMax
691          DO i=iMin,iMax          DO i=iMin,iMax
692           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
693          ENDDO          ENDDO
694         ENDDO         ENDDO
 C--   end if momAdvection  
695        ENDIF        ENDIF
696    
697  C--   Set du/dt & dv/dt on boundaries to zero  C--   Set du/dt & dv/dt on boundaries to zero
# Line 472  C--   Set du/dt & dv/dt on boundaries to Line 702  C--   Set du/dt & dv/dt on boundaries to
702         ENDDO         ENDDO
703        ENDDO        ENDDO
704    
705    #ifdef ALLOW_DEBUG
706          IF ( debugLevel .GE. debLevB
707         &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0
708         &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
709         &   .AND. useCubedSphereExchange ) THEN
710            CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
711         &             guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
712          ENDIF
713    #endif /* ALLOW_DEBUG */
714    
715        IF (        IF ( writeDiag ) THEN
716       &  DIFFERENT_MULTIPLE(diagFreq,myCurrentTime,          IF (snapshot_mdsio) THEN
717       &                     myCurrentTime-deltaTClock)           CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
718       & ) THEN           CALL WRITE_LOCAL_RL('KE','I10',1,KE,     bi,bj,k,myIter,myThid)
719         CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv,   bi,bj,k,myIter,myThid)
720         CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
721         CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
722         CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
723         CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid)          ENDIF
724         CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid)  #ifdef ALLOW_MNC
725         CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)          IF (useMNC .AND. snapshot_mnc) THEN
726  c      CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
727         CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)       &          offsets, myThid)
728         CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid)            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
729         &          offsets, myThid)
730              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
731         &          offsets, myThid)
732              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
733         &          offsets, myThid)
734              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
735         &          offsets, myThid)
736              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
737         &          offsets, myThid)
738            ENDIF
739    #endif /*  ALLOW_MNC  */
740          ENDIF
741    
742    #ifdef ALLOW_DIAGNOSTICS
743          IF ( useDiagnostics ) THEN
744            CALL DIAGNOSTICS_FILL(KE,     'momKE   ',k,1,2,bi,bj,myThid)
745           IF (momViscosity) THEN
746            CALL DIAGNOSTICS_FILL(hDiv,   'momHDiv ',k,1,2,bi,bj,myThid)
747            CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
748            CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
749            CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
750           ENDIF
751            CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),
752         &                                'Um_Advec',k,1,2,bi,bj,myThid)
753            CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),
754         &                                'Vm_Advec',k,1,2,bi,bj,myThid)
755        ENDIF        ENDIF
756    #endif /* ALLOW_DIAGNOSTICS */
757    
758  #endif /* ALLOW_MOM_VECINV */  #endif /* ALLOW_MOM_VECINV */
759    

Legend:
Removed from v.1.11  
changed lines
  Added in v.1.63

  ViewVC Help
Powered by ViewVC 1.1.22