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

Legend:
Removed from v.1.6  
changed lines
  Added in v.1.69

  ViewVC Help
Powered by ViewVC 1.1.22