/[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.51 by baylor, Tue Sep 27 13:38:42 2005 UTC revision 1.78 by jmc, Wed Nov 30 00:11:22 2016 UTC
# Line 2  C $Header$ Line 2  C $Header$
2  C $Name$  C $Name$
3    
4  #include "MOM_VECINV_OPTIONS.h"  #include "MOM_VECINV_OPTIONS.h"
5    #ifdef ALLOW_AUTODIFF
6    # include "AUTODIFF_OPTIONS.h"
7    #endif
8    #ifdef ALLOW_MOM_COMMON
9    # include "MOM_COMMON_OPTIONS.h"
10    #endif
11    
12        SUBROUTINE MOM_VECINV(        SUBROUTINE MOM_VECINV(
13       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,       I        bi,bj,k,iMin,iMax,jMin,jMax,
14       I        KappaRU, KappaRV,       I        kappaRU, kappaRV,
15       U        fVerU, fVerV,       I        fVerUkm, fVerVkm,
16         O        fVerUkp, fVerVkp,
17       O        guDiss, gvDiss,       O        guDiss, gvDiss,
18       I        myTime, myIter, myThid)       I        myTime, myIter, myThid )
19  C     /==========================================================\  C     *==========================================================*
20  C     | S/R MOM_VECINV                                           |  C     | S/R MOM_VECINV                                           |
21  C     | o Form the right hand-side of the momentum equation.     |  C     | o Form the right hand-side of the momentum equation.     |
22  C     |==========================================================|  C     *==========================================================*
23  C     | Terms are evaluated one layer at a time working from     |  C     | Terms are evaluated one layer at a time working from     |
24  C     | the bottom to the top. The vertically integrated         |  C     | the bottom to the top. The vertically integrated         |
25  C     | barotropic flow tendency term is evluated by summing the |  C     | barotropic flow tendency term is evluated by summing the |
# Line 23  C     | for the diffusion equation bc wi Line 30  C     | for the diffusion equation bc wi
30  C     | form produces a diffusive flux that does not scale with  |  C     | form produces a diffusive flux that does not scale with  |
31  C     | open-area. Need to do something to solidfy this and to   |  C     | open-area. Need to do something to solidfy this and to   |
32  C     | deal "properly" with thin walls.                         |  C     | deal "properly" with thin walls.                         |
33  C     \==========================================================/  C     *==========================================================*
34        IMPLICIT NONE        IMPLICIT NONE
35    
36  C     == Global variables ==  C     == Global variables ==
37  #include "SIZE.h"  #include "SIZE.h"
 #include "DYNVARS.h"  
38  #include "EEPARAMS.h"  #include "EEPARAMS.h"
39  #include "PARAMS.h"  #include "PARAMS.h"
 #ifdef ALLOW_MNC  
 #include "MNC_PARAMS.h"  
 #endif  
40  #include "GRID.h"  #include "GRID.h"
41    #include "SURFACE.h"
42    #include "DYNVARS.h"
43    #ifdef ALLOW_MOM_COMMON
44    # include "MOM_VISC.h"
45    #endif
46  #ifdef ALLOW_TIMEAVE  #ifdef ALLOW_TIMEAVE
47  #include "TIMEAVE_STATV.h"  # include "TIMEAVE_STATV.h"
48    #endif
49    #ifdef ALLOW_MNC
50    # include "MNC_PARAMS.h"
51    #endif
52    #ifdef ALLOW_AUTODIFF_TAMC
53    # include "tamc.h"
54    # include "tamc_keys.h"
55  #endif  #endif
56    
57  C     == Routine arguments ==  C     == Routine arguments ==
58  C     fVerU  :: Flux of momentum in the vertical direction, out of the upper  C     bi,bj   :: current tile indices
59  C     fVerV  :: face of a cell K ( flux into the cell above ).  C     k       :: current vertical level
60  C     guDiss :: dissipation tendency (all explicit terms), u component  C     iMin,iMax,jMin,jMax :: loop ranges
61  C     gvDiss :: dissipation tendency (all explicit terms), v component  C     fVerU   :: Flux of momentum in the vertical direction, out of the upper
62  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 ).
63  C                                      results will be set.  C     fVerUkm :: vertical viscous flux of U, interface above (k-1/2)
64  C     kUp, kDown                     - Index for upper and lower layers.  C     fVerVkm :: vertical viscous flux of V, interface above (k-1/2)
65  C     myThid - Instance number for this innvocation of CALC_MOM_RHS  C     fVerUkp :: vertical viscous flux of U, interface below (k+1/2)
66        _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  C     fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
67        _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
68        _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  C     guDiss  :: dissipation tendency (all explicit terms), u component
69        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  C     gvDiss  :: dissipation tendency (all explicit terms), v component
70    C     myTime  :: current time
71    C     myIter  :: current time-step number
72    C     myThid  :: my Thread Id number
73          INTEGER bi,bj,k
74          INTEGER iMin,iMax,jMin,jMax
75          _RL kappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
76          _RL kappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
77          _RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78          _RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79          _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80          _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       INTEGER kUp,kDown  
83        _RL     myTime        _RL     myTime
84        INTEGER myIter        INTEGER myIter
85        INTEGER myThid        INTEGER myThid
       INTEGER bi,bj,iMin,iMax,jMin,jMax  
86    
87  #ifdef ALLOW_MOM_VECINV  #ifdef ALLOW_MOM_VECINV
88    
# Line 67  C     == Functions == Line 91  C     == Functions ==
91        EXTERNAL DIFFERENT_MULTIPLE        EXTERNAL DIFFERENT_MULTIPLE
92    
93  C     == Local variables ==  C     == Local variables ==
94    C     strainBC :: same as strain but account for no-slip BC
95    C     vort3BC  :: same as vort3  but account for no-slip BC
96        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL      vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL      uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99        _RL      vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100  c     _RL      mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS hFacZ   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101        _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS h0FacZ  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102        _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2u   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2v   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL zStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109        _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110        _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL strain  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111  C     I,J,K - Loop counters        _RL strainBC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112        INTEGER i,j,k        _RL KE      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113  C     xxxFac - On-off tracer parameters used for switching terms off.        _RL omega3  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114        _RL  ArDudrFac        _RL vort3   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115  c     _RL  mtFacU        _RL vort3BC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116        _RL  ArDvdrFac        _RL hDiv    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
 c     _RL  mtFacV  
       LOGICAL bottomDragTerms  
       LOGICAL writeDiag  
       _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
117        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121        LOGICAL harmonic,biharmonic,useVariableViscosity  C     i,j    :: Loop counters
122          INTEGER i,j
123    C     xxxFac :: On-off tracer parameters used for switching terms off.
124          _RL  ArDudrFac
125          _RL  ArDvdrFac
126          _RL  sideMaskFac
127          LOGICAL bottomDragTerms
128          LOGICAL writeDiag
129    #ifdef ALLOW_AUTODIFF_TAMC
130          INTEGER imomkey
131    #endif
132    
133  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
134        INTEGER offsets(9)        INTEGER offsets(9)
135          CHARACTER*(1) pf
136  #endif  #endif
137    
138  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF
139  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
140  C--   the kUp is still required  C--   the kUp is still required
141  C--   In the case of mom_fluxform Kup is set as well  C--   In the case of mom_fluxform kUp is set as well
142  C--   (at least in part)  C--   (at least in part)
143        fVerU(1,1,kUp) = fVerU(1,1,kUp)        fVerUkm(1,1) = fVerUkm(1,1)
144        fVerV(1,1,kUp) = fVerV(1,1,kUp)        fVerVkm(1,1) = fVerVkm(1,1)
145  #endif  #endif
146    
147    #ifdef ALLOW_AUTODIFF_TAMC
148              act0 = k - 1
149              max0 = Nr
150              act1 = bi - myBxLo(myThid)
151              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
152              act2 = bj - myByLo(myThid)
153              max2 = myByHi(myThid) - myByLo(myThid) + 1
154              act3 = myThid - 1
155              max3 = nTx*nTy
156              act4 = ikey_dynamics - 1
157              imomkey = (act0 + 1)
158         &                    + act1*max0
159         &                    + act2*max0*max1
160         &                    + act3*max0*max1*max2
161         &                    + act4*max0*max1*max2*max3
162    #endif /* ALLOW_AUTODIFF_TAMC */
163    
164        writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)        writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
165    
166  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
167        IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN        IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
168            IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
169              pf(1:1) = 'D'
170            ELSE
171              pf(1:1) = 'R'
172            ENDIF
173          IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN          IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
174            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
175            CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)            CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
# Line 128  C--   (at least in part) Line 180  C--   (at least in part)
180            offsets(i) = 0            offsets(i) = 0
181          ENDDO          ENDDO
182          offsets(3) = k          offsets(3) = k
183  C       write(*,*) 'offsets = ',(offsets(i),i=1,9)  c       write(*,*) 'offsets = ',(offsets(i),i=1,9)
184        ENDIF        ENDIF
185  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
186    
187  C     Initialise intermediate terms  C--   Initialise intermediate terms
188        DO J=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
189         DO I=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
190          vF(i,j)    = 0.          vF(i,j)    = 0.
191          vrF(i,j)   = 0.          vrF(i,j)   = 0.
192          uCf(i,j)   = 0.          uCf(i,j)   = 0.
193          vCf(i,j)   = 0.          vCf(i,j)   = 0.
 c       mT(i,j)    = 0.  
194          del2u(i,j) = 0.          del2u(i,j) = 0.
195          del2v(i,j) = 0.          del2v(i,j) = 0.
196          dStar(i,j) = 0.          dStar(i,j) = 0.
# Line 148  c       mT(i,j)    = 0. Line 199  c       mT(i,j)    = 0.
199          gvDiss(i,j)= 0.          gvDiss(i,j)= 0.
200          vort3(i,j) = 0.          vort3(i,j) = 0.
201          omega3(i,j)= 0.          omega3(i,j)= 0.
202          ke(i,j)    = 0.          KE(i,j)    = 0.
203          viscAh_Z(i,j) = 0.  C-    need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
204          viscAh_D(i,j) = 0.          hDiv(i,j)  = 0.
205          viscA4_Z(i,j) = 0.  c       viscAh_Z(i,j) = 0.
206          viscA4_D(i,j) = 0.  c       viscAh_D(i,j) = 0.
207    c       viscA4_Z(i,j) = 0.
208  #ifdef ALLOW_AUTODIFF_TAMC  c       viscA4_D(i,j) = 0.
209          strain(i,j)  = 0. _d 0          strain(i,j)  = 0. _d 0
210            strainBC(i,j)= 0. _d 0
211          tension(i,j) = 0. _d 0          tension(i,j) = 0. _d 0
212    #ifdef ALLOW_AUTODIFF
213            hFacZ(i,j)   = 0. _d 0
214  #endif  #endif
215         ENDDO         ENDDO
216        ENDDO        ENDDO
# Line 164  c       mT(i,j)    = 0. Line 218  c       mT(i,j)    = 0.
218  C--   Term by term tracer parmeters  C--   Term by term tracer parmeters
219  C     o U momentum equation  C     o U momentum equation
220        ArDudrFac    = vfFacMom*1.        ArDudrFac    = vfFacMom*1.
 c     mTFacU       = mtFacMom*1.  
221  C     o V momentum equation  C     o V momentum equation
222        ArDvdrFac    = vfFacMom*1.        ArDvdrFac    = vfFacMom*1.
 c     mTFacV       = mtFacMom*1.  
223    
224        IF (     no_slip_bottom  C note: using standard stencil (no mask) results in under-estimating
225       &    .OR. bottomDragQuadratic.NE.0.  C       vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
226       &    .OR. bottomDragLinear.NE.0.) THEN        IF ( no_slip_sides ) THEN
227            sideMaskFac = sideDragFactor
228          ELSE
229            sideMaskFac = 0. _d 0
230          ENDIF
231    
232          IF ( selectImplicitDrag.EQ.0 .AND.
233         &      (  no_slip_bottom
234         &    .OR. selectBotDragQuadr.GE.0
235         &    .OR. bottomDragLinear.NE.0. ) ) THEN
236         bottomDragTerms=.TRUE.         bottomDragTerms=.TRUE.
237        ELSE        ELSE
238         bottomDragTerms=.FALSE.         bottomDragTerms=.FALSE.
# Line 192  C note (jmc) : Dissipation and Vort3 adv Line 253  C note (jmc) : Dissipation and Vort3 adv
253  C              use the same maskZ (and hFacZ)  => needs 2 call(s)  C              use the same maskZ (and hFacZ)  => needs 2 call(s)
254  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)
255    
256        CALL MOM_CALC_KE(bi,bj,k,0,uFld,vFld,KE,myThid)        CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
   
       CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)  
257    
258        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
259    
260        CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)  C-    mask vort3 and account for no-slip / free-slip BC in vort3BC:
261          DO j=1-OLy,sNy+OLy
262           DO i=1-OLx,sNx+OLx
263             vort3BC(i,j) = vort3(i,j)
264             IF ( hFacZ(i,j).EQ.zeroRS ) THEN
265               vort3BC(i,j) = sideMaskFac*vort3BC(i,j)
266               vort3(i,j)   = 0.
267             ENDIF
268           ENDDO
269          ENDDO
270    
271        CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)        IF (momViscosity) THEN
272    C--    For viscous term, compute horizontal divergence, tension & strain
273    C      and mask relative vorticity (free-slip case):
274    
275        IF (useAbsVorticity)         DO j=1-OLy,sNy+OLy
276       & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)          DO i=1-OLx,sNx+OLx
277              h0FacZ(i,j) = hFacZ(i,j)
278            ENDDO
279           ENDDO
280    #ifdef NONLIN_FRSURF
281           IF ( no_slip_sides .AND. nonlinFreeSurf.GT.0 ) THEN
282            DO j=2-OLy,sNy+OLy
283             DO i=2-OLx,sNx+OLx
284              h0FacZ(i,j) = MIN(
285         &       MIN( h0FacW(i,j,k,bi,bj), h0FacW(i,j-1,k,bi,bj) ),
286         &       MIN( h0FacS(i,j,k,bi,bj), h0FacS(i-1,j,k,bi,bj) ) )
287             ENDDO
288            ENDDO
289           ENDIF
290    #endif /* NONLIN_FRSURF */
291    
292        IF (momViscosity) THEN         CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
293  C      Calculate Viscosities  
294        CALL MOM_CALC_VISC(         IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
295       I        bi,bj,k,          CALL MOM_CALC_TENSION( bi,bj,k,uFld,vFld,tension,myThid )
296       O        viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,          CALL MOM_CALC_STRAIN( bi,bj,k,uFld,vFld,hFacZ,strain,myThid )
297       O        harmonic,biharmonic,useVariableViscosity,  C-    mask strain and account for no-slip / free-slip BC in strainBC:
298       I        hDiv,vort3,tension,strain,KE,hfacZ,          DO j=1-OLy,sNy+OLy
299       I        myThid)           DO i=1-OLx,sNx+OLx
300               strainBC(i,j) = strain(i,j)
301               IF ( hFacZ(i,j).EQ.zeroRS ) THEN
302                 strainBC(i,j) = sideMaskFac*strainBC(i,j)
303                 strain(i,j)   = 0.
304               ENDIF
305             ENDDO
306            ENDDO
307           ENDIF
308    
309    C--    Calculate Lateral Viscosities
310           DO j=1-OLy,sNy+OLy
311            DO i=1-OLx,sNx+OLx
312             viscAh_D(i,j) = viscAhD
313             viscAh_Z(i,j) = viscAhZ
314             viscA4_D(i,j) = viscA4D
315             viscA4_Z(i,j) = viscA4Z
316            ENDDO
317           ENDDO
318           IF ( useVariableVisc ) THEN
319    C-     uses vort3BC & strainBC which account for no-slip / free-slip BC
320             CALL MOM_CALC_VISC( bi, bj, k,
321         O            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
322         I            hDiv, vort3BC, tension, strainBC, KE, hfacZ,
323         I            myThid )
324           ENDIF
325    
326  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
327         IF (biharmonic) THEN         IF (useBiharmonicVisc) THEN
328           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
329       O                      del2u,del2v,       O                      del2u,del2v,
330       &                      myThid)       I                      myThid)
331           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
332           CALL MOM_CALC_RELVORT3(bi,bj,k,           CALL MOM_CALC_RELVORT3(bi,bj,k,
333       &                          del2u,del2v,hFacZ,zStar,myThid)       &                          del2u,del2v,hFacZ,zStar,myThid)
334         ENDIF         ENDIF
335    
336  C      Calculate dissipation terms for U and V equations  C---   Calculate dissipation terms for U and V equations
337    
338  C      in terms of tension and strain  C-     in terms of tension and strain
339         IF (useStrainTensionVisc) THEN         IF (useStrainTensionVisc) THEN
340           CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,  C      use masked strain as if free-slip since side-drag is computed separately
341       I                    hFacZ,           CALL MOM_HDISSIP( bi, bj, k,
342       I                    viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       I            tension, strain, hFacZ,
343       I                    harmonic,biharmonic,useVariableViscosity,       I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
344       O                    guDiss,gvDiss,       I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
345       I                    myThid)       O            guDiss, gvDiss,
346         I            myThid )
347         ELSE         ELSE
348  C      in terms of vorticity and divergence  C-     in terms of vorticity and divergence
349           CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,           CALL MOM_VI_HDISSIP( bi, bj, k,
350       I                       hFacZ,dStar,zStar,       I            hDiv, vort3, dStar, zStar, hFacZ,
351       I                       viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
352       I                       harmonic,biharmonic,useVariableViscosity,       I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
353       O                       guDiss,gvDiss,       O            guDiss, gvDiss,
354       &                       myThid)               I            myThid )
355         ENDIF         ENDIF
       ENDIF  
   
 C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:  
 c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)  
356    
357  C---- Zonal momentum equation starts here  C---  Other dissipation terms in Zonal momentum equation
358    
359  C--   Vertical flux (fVer is at upper face of "u" cell)  C--   Vertical flux (fVer is at upper face of "u" cell)
   
360  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
361        IF (momViscosity.AND..NOT.implicitViscosity) THEN        IF ( .NOT.implicitViscosity ) THEN
362         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,kappaRU,vrF,myThid)
   
363  C     Combine fluxes  C     Combine fluxes
364         DO j=jMin,jMax         DO j=jMin,jMax
365          DO i=iMin,iMax          DO i=iMin,iMax
366           fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)           fVerUkp(i,j) = ArDudrFac*vrF(i,j)
367          ENDDO          ENDDO
368         ENDDO         ENDDO
   
369  C--   Tendency is minus divergence of the fluxes  C--   Tendency is minus divergence of the fluxes
370         DO j=2-Oly,sNy+Oly-1  C     vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
371          DO i=2-Olx,sNx+Olx-1         DO j=jMin,jMax
372            DO i=iMin,iMax
373           guDiss(i,j) = guDiss(i,j)           guDiss(i,j) = guDiss(i,j)
374       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
375       &   *recip_rAw(i,j,bi,bj)       &   *recip_rAw(i,j,bi,bj)
376       &  *(       &   *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
377       &    fVerU(i,j,kDown) - fVerU(i,j,kUp)       &   *recip_deepFac2C(k)*recip_rhoFacC(k)
      &   )*rkSign  
378          ENDDO          ENDDO
379         ENDDO         ENDDO
380        ENDIF        ENDIF
381    
382  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
383        IF (momViscosity.AND.no_slip_sides) THEN        IF ( no_slip_sides ) THEN
384  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
385         CALL MOM_U_SIDEDRAG(         CALL MOM_U_SIDEDRAG( bi, bj, k,
386       I        bi,bj,k,       I          uFld, del2u, h0FacZ,
387       I        uFld, del2u, hFacZ,       I          viscAh_Z, viscA4_Z,
388       I        viscAh_Z,viscA4_Z,       I          useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
389       I        harmonic,biharmonic,useVariableViscosity,       O          vF,
390       O        vF,       I          myThid )
      I        myThid)  
391         DO j=jMin,jMax         DO j=jMin,jMax
392          DO i=iMin,iMax          DO i=iMin,iMax
393           guDiss(i,j) = guDiss(i,j)+vF(i,j)           guDiss(i,j) = guDiss(i,j)+vF(i,j)
# Line 294  C-     No-slip BCs impose a drag at wall Line 396  C-     No-slip BCs impose a drag at wall
396        ENDIF        ENDIF
397    
398  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
399        IF (momViscosity.AND.bottomDragTerms) THEN        IF ( bottomDragTerms ) THEN
400         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)         CALL MOM_U_BOTTOMDRAG( bi, bj, k,
401         I            uFld, vFld, KE, kappaRU,
402         O            vF,
403         I            myThid )
404         DO j=jMin,jMax         DO j=jMin,jMax
405          DO i=iMin,iMax          DO i=iMin,iMax
406           guDiss(i,j) = guDiss(i,j)+vF(i,j)           guDiss(i,j) = guDiss(i,j)+vF(i,j)
407          ENDDO          ENDDO
408         ENDDO         ENDDO
409        ENDIF        ENDIF
410    #ifdef ALLOW_SHELFICE
411          IF ( useShelfIce ) THEN
412           CALL SHELFICE_U_DRAG( bi, bj, k,
413         I            uFld, vFld, KE, kappaRU,
414         O            vF,
415         I            myThid )
416           DO j=jMin,jMax
417            DO i=iMin,iMax
418             guDiss(i,j) = guDiss(i,j) + vF(i,j)
419            ENDDO
420           ENDDO
421          ENDIF
422    #endif /* ALLOW_SHELFICE */
423    
424  C--   Metric terms for curvilinear grid systems  C---  Other dissipation terms in Meridional momentum equation
 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  
   
 C---- Meridional momentum equation starts here  
425    
426  C--   Vertical flux (fVer is at upper face of "v" cell)  C--   Vertical flux (fVer is at upper face of "v" cell)
   
427  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
428        IF (momViscosity.AND..NOT.implicitViscosity) THEN        IF ( .NOT.implicitViscosity ) THEN
429         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,kappaRV,vrF,myThid)
   
430  C     Combine fluxes -> fVerV  C     Combine fluxes -> fVerV
431         DO j=jMin,jMax         DO j=jMin,jMax
432          DO i=iMin,iMax          DO i=iMin,iMax
433           fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)           fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
434          ENDDO          ENDDO
435         ENDDO         ENDDO
   
436  C--   Tendency is minus divergence of the fluxes  C--   Tendency is minus divergence of the fluxes
437    C     vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
438         DO j=jMin,jMax         DO j=jMin,jMax
439          DO i=iMin,iMax          DO i=iMin,iMax
440           gvDiss(i,j) = gvDiss(i,j)           gvDiss(i,j) = gvDiss(i,j)
441       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
442       &    *recip_rAs(i,j,bi,bj)       &   *recip_rAs(i,j,bi,bj)
443       &  *(       &   *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
444       &    fVerV(i,j,kDown) - fVerV(i,j,kUp)       &   *recip_deepFac2C(k)*recip_rhoFacC(k)
      &   )*rkSign  
445          ENDDO          ENDDO
446         ENDDO         ENDDO
447        ENDIF        ENDIF
448    
449  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
450        IF (momViscosity.AND.no_slip_sides) THEN        IF ( no_slip_sides ) THEN
451  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
452         CALL MOM_V_SIDEDRAG(         CALL MOM_V_SIDEDRAG( bi, bj, k,
453       I        bi,bj,k,       I          vFld, del2v, h0FacZ,
454       I        vFld, del2v, hFacZ,       I          viscAh_Z, viscA4_Z,
455       I        viscAh_Z,viscA4_Z,       I          useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
456       I        harmonic,biharmonic,useVariableViscosity,       O          vF,
457       O        vF,       I          myThid )
      I        myThid)  
458         DO j=jMin,jMax         DO j=jMin,jMax
459          DO i=iMin,iMax          DO i=iMin,iMax
460           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
461          ENDDO          ENDDO
462         ENDDO         ENDDO
463        ENDIF        ENDIF
464    
465  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
466        IF (momViscosity.AND.bottomDragTerms) THEN        IF ( bottomDragTerms ) THEN
467         CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)         CALL MOM_V_BOTTOMDRAG( bi, bj, k,
468         I            uFld, vFld, KE, kappaRV,
469         O            vF,
470         I            myThid )
471         DO j=jMin,jMax         DO j=jMin,jMax
472          DO i=iMin,iMax          DO i=iMin,iMax
473           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
474          ENDDO          ENDDO
475         ENDDO         ENDDO
476        ENDIF        ENDIF
477    #ifdef ALLOW_SHELFICE
478          IF ( useShelfIce ) THEN
479           CALL SHELFICE_V_DRAG( bi, bj, k,
480         I            uFld, vFld, KE, kappaRV,
481         O            vF,
482         I            myThid )
483           DO j=jMin,jMax
484            DO i=iMin,iMax
485             gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
486            ENDDO
487           ENDDO
488          ENDIF
489    #endif /* ALLOW_SHELFICE */
490    
491    C--   if (momViscosity) end of block.
492          ENDIF
493    
494  C--   Metric terms for curvilinear grid systems  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:
495  c     IF (usingSphericalPolarMTerms) THEN  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
496  C      o Spherical polar grid metric terms  
497  c      CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
498  c      DO j=jMin,jMax  
499  c       DO i=iMin,iMax  C---  Prepare for Advection & Coriolis terms:
500  c        gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)  C-    calculate absolute vorticity
501  c       ENDDO        IF (useAbsVorticity)
502  c      ENDDO       &  CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
 c     ENDIF  
503    
504  C--   Horizontal Coriolis terms  C--   Horizontal Coriolis terms
505  c     IF (useCoriolis .AND. .NOT.useCDscheme  c     IF (useCoriolis .AND. .NOT.useCDscheme
# Line 387  C- jmc: change it to keep the Coriolis t Line 509  C- jmc: change it to keep the Coriolis t
509       &     .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )       &     .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
510       &   ) THEN       &   ) THEN
511         IF (useAbsVorticity) THEN         IF (useAbsVorticity) THEN
512          CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,          CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,omega3,hFacZ,r_hFacZ,
513       &                         uCf,myThid)       &                         uCf,myThid)
514          CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,          CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,omega3,hFacZ,r_hFacZ,
515       &                         vCf,myThid)       &                         vCf,myThid)
516         ELSE         ELSE
517          CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,          CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
# Line 401  C- jmc: change it to keep the Coriolis t Line 523  C- jmc: change it to keep the Coriolis t
523           gV(i,j,k,bi,bj) = vCf(i,j)           gV(i,j,k,bi,bj) = vCf(i,j)
524          ENDDO          ENDDO
525         ENDDO         ENDDO
   
526         IF ( writeDiag ) THEN         IF ( writeDiag ) THEN
527           IF (snapshot_mdsio) THEN           IF (snapshot_mdsio) THEN
528             CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)             CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
# Line 409  C- jmc: change it to keep the Coriolis t Line 530  C- jmc: change it to keep the Coriolis t
530           ENDIF           ENDIF
531  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
532           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
533             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fV', uCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
534       &          offsets, myThid)       &          offsets, myThid)
535             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fU', vCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
536       &          offsets, myThid)       &          offsets, myThid)
537           ENDIF           ENDIF
538  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 422  C- jmc: change it to keep the Coriolis t Line 543  C- jmc: change it to keep the Coriolis t
543           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
544         ENDIF         ENDIF
545  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
   
546        ELSE        ELSE
547         DO j=jMin,jMax         DO j=jMin,jMax
548          DO i=iMin,iMax          DO i=iMin,iMax
# Line 434  C- jmc: change it to keep the Coriolis t Line 554  C- jmc: change it to keep the Coriolis t
554    
555        IF (momAdvection) THEN        IF (momAdvection) THEN
556  C--   Horizontal advection of relative (or absolute) vorticity  C--   Horizontal advection of relative (or absolute) vorticity
557         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
558         &     .AND.useAbsVorticity ) THEN
559          CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,          CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
560       &                         uCf,myThid)       &                         uCf,myThid)
561         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
562          CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,          CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
563       &                         uCf,myThid)       &                         uCf,myThid)
564         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
565          CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,          CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,omega3,hFacZ,r_hFacZ,
566       &                         uCf,myThid)       &                         uCf,myThid)
567         ELSE         ELSE
568          CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,          CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
# Line 452  C--   Horizontal advection of relative ( Line 573  C--   Horizontal advection of relative (
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         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
577          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,       &     .AND.useAbsVorticity ) THEN
578            CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,uFld,omega3,r_hFacZ,
579       &                         vCf,myThid)       &                         vCf,myThid)
580         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
581          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,uFld,vort3, r_hFacZ,
582       &                         vCf,myThid)       &                         vCf,myThid)
583         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
584          CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,          CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,omega3,hFacZ,r_hFacZ,
585       &                         vCf,myThid)       &                         vCf,myThid)
586         ELSE         ELSE
587          CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,          CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
# Line 478  C--   Horizontal advection of relative ( Line 600  C--   Horizontal advection of relative (
600           ENDIF           ENDIF
601  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
602           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
603             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zV', uCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
604       &          offsets, myThid)       &          offsets, myThid)
605             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zU', vCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
606       &          offsets, myThid)       &          offsets, myThid)
607           ENDIF           ENDIF
608  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 503  C--   Horizontal advection of relative ( Line 625  C--   Horizontal advection of relative (
625    
626  C--   Vertical shear terms (-w*du/dr & -w*dv/dr)  C--   Vertical shear terms (-w*du/dr & -w*dv/dr)
627         IF ( .NOT. momImplVertAdv ) THEN         IF ( .NOT. momImplVertAdv ) THEN
628          CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)          CALL MOM_VI_U_VERTSHEAR(bi,bj,k,uVel,wVel,uCf,myThid)
629          DO j=jMin,jMax          DO j=jMin,jMax
630           DO i=iMin,iMax           DO i=iMin,iMax
631            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)
632           ENDDO           ENDDO
633          ENDDO          ENDDO
634          CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)          CALL MOM_VI_V_VERTSHEAR(bi,bj,k,vVel,wVel,vCf,myThid)
635          DO j=jMin,jMax          DO j=jMin,jMax
636           DO i=iMin,iMax           DO i=iMin,iMax
637            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)
# Line 524  C--   Vertical shear terms (-w*du/dr & - Line 646  C--   Vertical shear terms (-w*du/dr & -
646         ENDIF         ENDIF
647    
648  C--   Bernoulli term  C--   Bernoulli term
649         CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)         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_GRAD_KE(bi,bj,K,KE,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)
# Line 543  C--   Bernoulli term Line 665  C--   Bernoulli term
665           ENDIF           ENDIF
666  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
667           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
668             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEx', uCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
669       &          offsets, myThid)       &          offsets, myThid)
670             CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEy', vCf,             CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
671       &          offsets, myThid)       &          offsets, myThid)
672          ENDIF           ENDIF
673  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
674         ENDIF         ENDIF
675    
676  C--   end if momAdvection  C--   end if momAdvection
677        ENDIF        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
702            DO i=iMin,iMax
703             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
704            ENDDO
705           ENDDO
706           CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
707           DO j=jMin,jMax
708            DO i=iMin,iMax
709             gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
710            ENDDO
711           ENDDO
712          ENDIF
713    
714  C--   Set du/dt & dv/dt on boundaries to zero  C--   Set du/dt & dv/dt on boundaries to zero
715        DO j=jMin,jMax        DO j=jMin,jMax
716         DO i=iMin,iMax         DO i=iMin,iMax
# Line 563  C--   Set du/dt & dv/dt on boundaries to Line 720  C--   Set du/dt & dv/dt on boundaries to
720        ENDDO        ENDDO
721    
722  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
723        IF ( debugLevel .GE. debLevB        IF ( debugLevel .GE. debLevC
724       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0
725       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
726       &   .AND. useCubedSphereExchange ) THEN       &   .AND. useCubedSphereExchange ) THEN
# Line 573  C--   Set du/dt & dv/dt on boundaries to Line 730  C--   Set du/dt & dv/dt on boundaries to
730  #endif /* ALLOW_DEBUG */  #endif /* ALLOW_DEBUG */
731    
732        IF ( writeDiag ) THEN        IF ( writeDiag ) THEN
733            IF (useBiharmonicVisc) THEN
734             CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
735         &                        bi,bj,k, myIter, myThid )
736             CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
737         &                        bi,bj,k, myIter, myThid )
738             CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
739         &                        bi,bj,k, myIter, myThid )
740             CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
741         &                        bi,bj,k, myIter, myThid )
742            ENDIF
743          IF (snapshot_mdsio) THEN          IF (snapshot_mdsio) THEN
744            CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
745            CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,           CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
746       &         myThid)           CALL WRITE_LOCAL_RL('KE','I10',1,KE,     bi,bj,k,myIter,myThid)
747            CALL WRITE_LOCAL_RL('Du','I10',1,guDiss,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv,   bi,bj,k,myIter,myThid)
748            CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
749            CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
750            CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)       &                        bi,bj,k, myIter, myThid )
751            CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
752            CALL WRITE_LOCAL_RL('D','I10',1,hDiv,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
753          ENDIF          ENDIF
754  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
755          IF (useMNC .AND. snapshot_mnc) THEN          IF (useMNC .AND. snapshot_mnc) THEN
756            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Ds',strain,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
757       &          offsets, myThid)       &          offsets, myThid)
758            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dt',tension,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
759       &          offsets, myThid)       &          offsets, myThid)
760            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Du',guDiss,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
761       &          offsets, myThid)       &          offsets, myThid)
762            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dv',gvDiss,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
763       &          offsets, myThid)       &          offsets, myThid)
764            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Z3',vort3,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
765       &          offsets, myThid)       &          offsets, myThid)
766            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'W3',omega3,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
767       &          offsets, myThid)       &          offsets, myThid)
768            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'KE',KE,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
769       &          offsets, myThid)       &          offsets, myThid)
770            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'D', hDiv,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
771       &          offsets, myThid)       &          offsets, myThid)
772          ENDIF          ENDIF
773  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 608  C--   Set du/dt & dv/dt on boundaries to Line 775  C--   Set du/dt & dv/dt on boundaries to
775    
776  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
777        IF ( useDiagnostics ) THEN        IF ( useDiagnostics ) THEN
778          CALL DIAGNOSTICS_FILL(KE,    'momKE   ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
779          CALL DIAGNOSTICS_FILL(hDiv,  'momHDiv ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(KE,     'momKE   ',k,1,2,bi,bj,myThid)
         CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid)  
         CALL DIAGNOSTICS_FILL(strain,'Strain  ',k,1,2,bi,bj,myThid)  
         CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)  
         CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),  
      &                               'Um_Advec',k,1,2,bi,bj,myThid)  
         CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),  
      &                               'Vm_Advec',k,1,2,bi,bj,myThid)  
780         IF (momViscosity) THEN         IF (momViscosity) THEN
781          CALL DIAGNOSTICS_FILL(guDiss,'Um_Diss ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(hDiv,   'momHDiv ',k,1,2,bi,bj,myThid)
         CALL DIAGNOSTICS_FILL(gvDiss,'Vm_Diss ',k,1,2,bi,bj,myThid)  
782         ENDIF         ENDIF
783           IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
784            CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
785            CALL DIAGNOSTICS_FILL(strainBC,'Strain  ',k,1,2,bi,bj,myThid)
786           ENDIF
787            CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
788         &                                'Um_Advec',k,1,2,bi,bj,myThid)
789            CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
790         &                                'Vm_Advec',k,1,2,bi,bj,myThid)
791        ENDIF        ENDIF
792  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
793    

Legend:
Removed from v.1.51  
changed lines
  Added in v.1.78

  ViewVC Help
Powered by ViewVC 1.1.22