/[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.68 by heimbach, Mon Oct 1 15:46:33 2012 UTC
# Line 3  C $Name$ Line 3  C $Name$
3    
4  #include "MOM_VECINV_OPTIONS.h"  #include "MOM_VECINV_OPTIONS.h"
5    
6        SUBROUTINE MOM_VECINV(        SUBROUTINE MOM_VECINV(
7       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,       I        bi,bj,k,iMin,iMax,jMin,jMax,
8       I        KappaRU, KappaRV,       I        KappaRU, KappaRV,
9       U        fVerU, fVerV,       I        fVerUkm, fVerVkm,
10         O        fVerUkp, fVerVkp,
11       O        guDiss, gvDiss,       O        guDiss, gvDiss,
12       I        myTime, myIter, myThid)       I        myTime, myIter, myThid )
13  C     /==========================================================\  C     *==========================================================*
14  C     | S/R MOM_VECINV                                           |  C     | S/R MOM_VECINV                                           |
15  C     | o Form the right hand-side of the momentum equation.     |  C     | o Form the right hand-side of the momentum equation.     |
16  C     |==========================================================|  C     *==========================================================*
17  C     | Terms are evaluated one layer at a time working from     |  C     | Terms are evaluated one layer at a time working from     |
18  C     | the bottom to the top. The vertically integrated         |  C     | the bottom to the top. The vertically integrated         |
19  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 24  C     | for the diffusion equation bc wi
24  C     | form produces a diffusive flux that does not scale with  |  C     | form produces a diffusive flux that does not scale with  |
25  C     | open-area. Need to do something to solidfy this and to   |  C     | open-area. Need to do something to solidfy this and to   |
26  C     | deal "properly" with thin walls.                         |  C     | deal "properly" with thin walls.                         |
27  C     \==========================================================/  C     *==========================================================*
28        IMPLICIT NONE        IMPLICIT NONE
29    
30  C     == Global variables ==  C     == Global variables ==
# Line 38  C     == Global variables == Line 39  C     == Global variables ==
39  #ifdef ALLOW_TIMEAVE  #ifdef ALLOW_TIMEAVE
40  #include "TIMEAVE_STATV.h"  #include "TIMEAVE_STATV.h"
41  #endif  #endif
42    #ifdef ALLOW_AUTODIFF_TAMC
43    # include "tamc.h"
44    # include "tamc_keys.h"
45    #endif
46    
47  C     == Routine arguments ==  C     == Routine arguments ==
48  C     fVerU  :: Flux of momentum in the vertical direction, out of the upper  C     bi,bj   :: current tile indices
49  C     fVerV  :: face of a cell K ( flux into the cell above ).  C     k       :: current vertical level
50  C     guDiss :: dissipation tendency (all explicit terms), u component  C     iMin,iMax,jMin,jMax :: loop ranges
51  C     gvDiss :: dissipation tendency (all explicit terms), v component  C     fVerU   :: Flux of momentum in the vertical direction, out of the upper
52  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 ).
53  C                                      results will be set.  C     fVerUkm :: vertical viscous flux of U, interface above (k-1/2)
54  C     kUp, kDown                     - Index for upper and lower layers.  C     fVerVkm :: vertical viscous flux of V, interface above (k-1/2)
55  C     myThid - Instance number for this innvocation of CALC_MOM_RHS  C     fVerUkp :: vertical viscous flux of U, interface below (k+1/2)
56    C     fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
57    
58    C     guDiss  :: dissipation tendency (all explicit terms), u component
59    C     gvDiss  :: dissipation tendency (all explicit terms), v component
60    C     myTime  :: current time
61    C     myIter  :: current time-step number
62    C     myThid  :: my Thread Id number
63          INTEGER bi,bj,k
64          INTEGER iMin,iMax,jMin,jMax
65        _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
66        _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
67        _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69          _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70          _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       INTEGER kUp,kDown  
73        _RL     myTime        _RL     myTime
74        INTEGER myIter        INTEGER myIter
75        INTEGER myThid        INTEGER myThid
       INTEGER bi,bj,iMin,iMax,jMin,jMax  
76    
77  #ifdef ALLOW_MOM_VECINV  #ifdef ALLOW_MOM_VECINV
78    
# Line 68  C     == Functions == Line 82  C     == Functions ==
82    
83  C     == Local variables ==  C     == Local variables ==
84        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85        _RL      vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86        _RL      uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87        _RL      vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL      vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88  c     _RL      mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS hFacZ   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90        _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2u   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2v   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL zStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL strain  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KE      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99  C     I,J,K - Loop counters        _RL omega3  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100        INTEGER i,j,k        _RL vort3   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101  C     xxxFac - On-off tracer parameters used for switching terms off.        _RL hDiv    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL  ArDudrFac  
 c     _RL  mtFacU  
       _RL  ArDvdrFac  
 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)  
102        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106    C     i,j    :: Loop counters
107          INTEGER i,j
108    C     xxxFac :: On-off tracer parameters used for switching terms off.
109          _RL  ArDudrFac
110          _RL  ArDvdrFac
111          _RL  sideMaskFac
112          LOGICAL bottomDragTerms
113          LOGICAL writeDiag
114        LOGICAL harmonic,biharmonic,useVariableViscosity        LOGICAL harmonic,biharmonic,useVariableViscosity
115    #ifdef ALLOW_AUTODIFF_TAMC
116          INTEGER imomkey
117    #endif
118    
119  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
120        INTEGER offsets(9)        INTEGER offsets(9)
121          CHARACTER*(1) pf
122  #endif  #endif
123    
124  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 110  C--   only the kDown part of fverU/V is Line 126  C--   only the kDown part of fverU/V is
126  C--   the kUp is still required  C--   the kUp is still required
127  C--   In the case of mom_fluxform Kup is set as well  C--   In the case of mom_fluxform Kup is set as well
128  C--   (at least in part)  C--   (at least in part)
129        fVerU(1,1,kUp) = fVerU(1,1,kUp)        fVerUkm(1,1) = fVerUkm(1,1)
130        fVerV(1,1,kUp) = fVerV(1,1,kUp)        fVerVkm(1,1) = fVerVkm(1,1)
131  #endif  #endif
132    
133    #ifdef ALLOW_AUTODIFF_TAMC
134              act0 = k - 1
135              max0 = Nr
136              act1 = bi - myBxLo(myThid)
137              max1 = myBxHi(myThid) - myBxLo(myThid) + 1
138              act2 = bj - myByLo(myThid)
139              max2 = myByHi(myThid) - myByLo(myThid) + 1
140              act3 = myThid - 1
141              max3 = nTx*nTy
142              act4 = ikey_dynamics - 1
143              imomkey = (act0 + 1)
144         &                    + act1*max0
145         &                    + act2*max0*max1
146         &                    + act3*max0*max1*max2
147         &                    + act4*max0*max1*max2*max3
148    #endif /* ALLOW_AUTODIFF_TAMC */
149    
150        writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)        writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
151    
152  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
153        IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN        IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
154            IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
155              pf(1:1) = 'D'
156            ELSE
157              pf(1:1) = 'R'
158            ENDIF
159          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
160            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
161            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 166  C--   (at least in part)
166            offsets(i) = 0            offsets(i) = 0
167          ENDDO          ENDDO
168          offsets(3) = k          offsets(3) = k
169  C       write(*,*) 'offsets = ',(offsets(i),i=1,9)  c       write(*,*) 'offsets = ',(offsets(i),i=1,9)
170        ENDIF        ENDIF
171  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
172    
173  C     Initialise intermediate terms  C--   Initialise intermediate terms
174        DO J=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
175         DO I=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
176          vF(i,j)    = 0.          vF(i,j)    = 0.
177          vrF(i,j)   = 0.          vrF(i,j)   = 0.
178          uCf(i,j)   = 0.          uCf(i,j)   = 0.
179          vCf(i,j)   = 0.          vCf(i,j)   = 0.
 c       mT(i,j)    = 0.  
180          del2u(i,j) = 0.          del2u(i,j) = 0.
181          del2v(i,j) = 0.          del2v(i,j) = 0.
182          dStar(i,j) = 0.          dStar(i,j) = 0.
# Line 148  c       mT(i,j)    = 0. Line 185  c       mT(i,j)    = 0.
185          gvDiss(i,j)= 0.          gvDiss(i,j)= 0.
186          vort3(i,j) = 0.          vort3(i,j) = 0.
187          omega3(i,j)= 0.          omega3(i,j)= 0.
188          ke(i,j)    = 0.          KE(i,j)    = 0.
189    C-    need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
190            hDiv(i,j)  = 0.
191          viscAh_Z(i,j) = 0.          viscAh_Z(i,j) = 0.
192          viscAh_D(i,j) = 0.          viscAh_D(i,j) = 0.
193          viscA4_Z(i,j) = 0.          viscA4_Z(i,j) = 0.
194          viscA4_D(i,j) = 0.          viscA4_D(i,j) = 0.
195    
 #ifdef ALLOW_AUTODIFF_TAMC  
196          strain(i,j)  = 0. _d 0          strain(i,j)  = 0. _d 0
197          tension(i,j) = 0. _d 0          tension(i,j) = 0. _d 0
198    #ifdef ALLOW_AUTODIFF_TAMC
199            hFacZ(i,j)   = 0. _d 0
200  #endif  #endif
201         ENDDO         ENDDO
202        ENDDO        ENDDO
# Line 164  c       mT(i,j)    = 0. Line 204  c       mT(i,j)    = 0.
204  C--   Term by term tracer parmeters  C--   Term by term tracer parmeters
205  C     o U momentum equation  C     o U momentum equation
206        ArDudrFac    = vfFacMom*1.        ArDudrFac    = vfFacMom*1.
 c     mTFacU       = mtFacMom*1.  
207  C     o V momentum equation  C     o V momentum equation
208        ArDvdrFac    = vfFacMom*1.        ArDvdrFac    = vfFacMom*1.
209  c     mTFacV       = mtFacMom*1.  
210    C note: using standard stencil (no mask) results in under-estimating
211    C       vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
212          IF ( no_slip_sides ) THEN
213            sideMaskFac = sideDragFactor
214          ELSE
215            sideMaskFac = 0. _d 0
216          ENDIF
217    
218        IF (     no_slip_bottom        IF (     no_slip_bottom
219       &    .OR. bottomDragQuadratic.NE.0.       &    .OR. bottomDragQuadratic.NE.0.
# Line 192  C note (jmc) : Dissipation and Vort3 adv Line 238  C note (jmc) : Dissipation and Vort3 adv
238  C              use the same maskZ (and hFacZ)  => needs 2 call(s)  C              use the same maskZ (and hFacZ)  => needs 2 call(s)
239  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)
240    
241        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)  
242    
243        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
244    
245        CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)        IF (momViscosity) THEN
246    C--    For viscous term, compute horizontal divergence, tension & strain
247    C      and mask relative vorticity (free-slip case):
248    
249        CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)  #ifdef ALLOW_AUTODIFF_TAMC
250    CADJ STORE vort3(:,:) =
251    CADJ &     comlev1_bibj_k, key = imomkey, byte = isbyte
252    #endif
253    
254        IF (useAbsVorticity)         CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
      & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)  
255    
256        IF (momViscosity) THEN         CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
257  C      Calculate Viscosities  
258        CALL MOM_CALC_VISC(         CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
259    
260    C-     account for no-slip / free-slip BC:
261           DO j=1-OLy,sNy+OLy
262            DO i=1-OLx,sNx+OLx
263              IF ( hFacZ(i,j).EQ.0. ) THEN
264                vort3(i,j)  = sideMaskFac*vort3(i,j)
265                strain(i,j) = sideMaskFac*strain(i,j)
266              ENDIF
267            ENDDO
268           ENDDO
269    
270    C--    Calculate Viscosities
271           CALL MOM_CALC_VISC(
272       I        bi,bj,k,       I        bi,bj,k,
273       O        viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       O        viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
274       O        harmonic,biharmonic,useVariableViscosity,       O        harmonic,biharmonic,useVariableViscosity,
# Line 222  C      Calculate del^2 u and del^2 v for Line 283  C      Calculate del^2 u and del^2 v for
283           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
284           CALL MOM_CALC_RELVORT3(bi,bj,k,           CALL MOM_CALC_RELVORT3(bi,bj,k,
285       &                          del2u,del2v,hFacZ,zStar,myThid)       &                          del2u,del2v,hFacZ,zStar,myThid)
286             IF ( writeDiag ) THEN
287               CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
288         &                           bi,bj,k, myIter, myThid )
289               CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
290         &                           bi,bj,k, myIter, myThid )
291               CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
292         &                           bi,bj,k, myIter, myThid )
293               CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
294         &                           bi,bj,k, myIter, myThid )
295             ENDIF
296         ENDIF         ENDIF
297    
298  C      Calculate dissipation terms for U and V equations  C-    Strain diagnostics:
299           IF ( writeDiag ) THEN
300            IF (snapshot_mdsio) THEN
301              CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
302            ENDIF
303    #ifdef ALLOW_MNC
304            IF (useMNC .AND. snapshot_mnc) THEN
305              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
306         &          offsets, myThid)
307            ENDIF
308    #endif /*  ALLOW_MNC  */
309           ENDIF
310    #ifdef ALLOW_DIAGNOSTICS
311           IF ( useDiagnostics ) THEN
312            CALL DIAGNOSTICS_FILL(strain, 'Strain  ',k,1,2,bi,bj,myThid)
313           ENDIF
314    #endif /* ALLOW_DIAGNOSTICS */
315    
316    C---   Calculate dissipation terms for U and V equations
317    
318  C      in terms of tension and strain  C      in terms of tension and strain
319         IF (useStrainTensionVisc) THEN         IF (useStrainTensionVisc) THEN
320    C        mask strain as if free-slip since side-drag is computed separately
321             DO j=1-OLy,sNy+OLy
322              DO i=1-OLx,sNx+OLx
323                IF ( hFacZ(i,j).EQ.0. ) strain(i,j) = 0. _d 0
324              ENDDO
325             ENDDO
326           CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,           CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
327       I                    hFacZ,       I                    hFacZ,
328       I                    viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       I                    viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
# Line 241  C      in terms of vorticity and diverge Line 336  C      in terms of vorticity and diverge
336       I                       viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       I                       viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
337       I                       harmonic,biharmonic,useVariableViscosity,       I                       harmonic,biharmonic,useVariableViscosity,
338       O                       guDiss,gvDiss,       O                       guDiss,gvDiss,
339       &                       myThid)               &                       myThid)
340         ENDIF         ENDIF
341    C--   if (momViscosity) end of block.
342        ENDIF        ENDIF
343    
344  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:
345  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
346    
347  C---- Zonal momentum equation starts here  C---  Other dissipation terms in Zonal momentum equation
348    
349  C--   Vertical flux (fVer is at upper face of "u" cell)  C--   Vertical flux (fVer is at upper face of "u" cell)
350    
# Line 259  C     Eddy component of vertical flux (i Line 355  C     Eddy component of vertical flux (i
355  C     Combine fluxes  C     Combine fluxes
356         DO j=jMin,jMax         DO j=jMin,jMax
357          DO i=iMin,iMax          DO i=iMin,iMax
358           fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)           fVerUkp(i,j) = ArDudrFac*vrF(i,j)
359          ENDDO          ENDDO
360         ENDDO         ENDDO
361    
362  C--   Tendency is minus divergence of the fluxes  C--   Tendency is minus divergence of the fluxes
363         DO j=2-Oly,sNy+Oly-1         DO j=jMin,jMax
364          DO i=2-Olx,sNx+Olx-1          DO i=iMin,iMax
365           guDiss(i,j) = guDiss(i,j)           guDiss(i,j) = guDiss(i,j)
366       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
367       &   *recip_rAw(i,j,bi,bj)       &   *recip_rAw(i,j,bi,bj)
368       &  *(       &   *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
      &    fVerU(i,j,kDown) - fVerU(i,j,kUp)  
      &   )*rkSign  
369          ENDDO          ENDDO
370         ENDDO         ENDDO
371        ENDIF        ENDIF
372    
373  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
374        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
375  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
376         CALL MOM_U_SIDEDRAG(         CALL MOM_U_SIDEDRAG(
# Line 292  C-     No-slip BCs impose a drag at wall Line 386  C-     No-slip BCs impose a drag at wall
386          ENDDO          ENDDO
387         ENDDO         ENDDO
388        ENDIF        ENDIF
   
389  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
390        IF (momViscosity.AND.bottomDragTerms) THEN        IF (momViscosity.AND.bottomDragTerms) THEN
391         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
# Line 302  C-    No-slip BCs impose a drag at botto Line 395  C-    No-slip BCs impose a drag at botto
395          ENDDO          ENDDO
396         ENDDO         ENDDO
397        ENDIF        ENDIF
398    #ifdef ALLOW_SHELFICE
399          IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
400           CALL SHELFICE_U_DRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
401           DO j=jMin,jMax
402            DO i=iMin,iMax
403             guDiss(i,j) = guDiss(i,j) + vF(i,j)
404            ENDDO
405           ENDDO
406          ENDIF
407    #endif /* ALLOW_SHELFICE */
408    
 C--   Metric terms for curvilinear grid systems  
 c     IF (usingSphericalPolarMTerms) THEN  
 C      o Spherical polar grid metric terms  
 c      CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)  
 c      DO j=jMin,jMax  
 c       DO i=iMin,iMax  
 c        gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)  
 c       ENDDO  
 c      ENDDO  
 c     ENDIF  
409    
410  C---- Meridional momentum equation starts here  C---  Other dissipation terms in Meridional momentum equation
411    
412  C--   Vertical flux (fVer is at upper face of "v" cell)  C--   Vertical flux (fVer is at upper face of "v" cell)
413    
# Line 325  C     Eddy component of vertical flux (i Line 418  C     Eddy component of vertical flux (i
418  C     Combine fluxes -> fVerV  C     Combine fluxes -> fVerV
419         DO j=jMin,jMax         DO j=jMin,jMax
420          DO i=iMin,iMax          DO i=iMin,iMax
421           fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)           fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
422          ENDDO          ENDDO
423         ENDDO         ENDDO
424    
# Line 334  C--   Tendency is minus divergence of th Line 427  C--   Tendency is minus divergence of th
427          DO i=iMin,iMax          DO i=iMin,iMax
428           gvDiss(i,j) = gvDiss(i,j)           gvDiss(i,j) = gvDiss(i,j)
429       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
430       &    *recip_rAs(i,j,bi,bj)       &   *recip_rAs(i,j,bi,bj)
431       &  *(       &   *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
      &    fVerV(i,j,kDown) - fVerV(i,j,kUp)  
      &   )*rkSign  
432          ENDDO          ENDDO
433         ENDDO         ENDDO
434        ENDIF        ENDIF
435    
436  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
437        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
438  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
439         CALL MOM_V_SIDEDRAG(         CALL MOM_V_SIDEDRAG(
# Line 367  C-    No-slip BCs impose a drag at botto Line 458  C-    No-slip BCs impose a drag at botto
458          ENDDO          ENDDO
459         ENDDO         ENDDO
460        ENDIF        ENDIF
461    #ifdef ALLOW_SHELFICE
462          IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
463             CALL SHELFICE_V_DRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
464             DO j=jMin,jMax
465              DO i=iMin,iMax
466               gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
467              ENDDO
468             ENDDO
469            ENDIF
470    #endif /* ALLOW_SHELFICE */
471    
472  C--   Metric terms for curvilinear grid systems  
473  c     IF (usingSphericalPolarMTerms) THEN  C-    Vorticity diagnostics:
474  C      o Spherical polar grid metric terms        IF ( writeDiag ) THEN
475  c      CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)          IF (snapshot_mdsio) THEN
476  c      DO j=jMin,jMax            CALL WRITE_LOCAL_RL('Z3','I10',1,vort3, bi,bj,k,myIter,myThid)
477  c       DO i=iMin,iMax          ENDIF
478  c        gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)  #ifdef ALLOW_MNC
479  c       ENDDO          IF (useMNC .AND. snapshot_mnc) THEN
480  c      ENDDO            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
481  c     ENDIF       &          offsets, myThid)
482            ENDIF
483    #endif /*  ALLOW_MNC  */
484          ENDIF
485    #ifdef ALLOW_DIAGNOSTICS
486          IF ( useDiagnostics ) THEN
487            CALL DIAGNOSTICS_FILL(vort3,  'momVort3',k,1,2,bi,bj,myThid)
488          ENDIF
489    #endif /* ALLOW_DIAGNOSTICS */
490    
491    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
492    
493    C---  Prepare for Advection & Coriolis terms:
494    C-    Mask relative vorticity and calculate absolute vorticity
495          DO j=1-OLy,sNy+OLy
496           DO i=1-OLx,sNx+OLx
497             IF ( hFacZ(i,j).EQ.0. ) vort3(i,j) = 0.
498           ENDDO
499          ENDDO
500          IF (useAbsVorticity)
501         &  CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
502    
503  C--   Horizontal Coriolis terms  C--   Horizontal Coriolis terms
504  c     IF (useCoriolis .AND. .NOT.useCDscheme  c     IF (useCoriolis .AND. .NOT.useCDscheme
# Line 401  C- jmc: change it to keep the Coriolis t Line 522  C- jmc: change it to keep the Coriolis t
522           gV(i,j,k,bi,bj) = vCf(i,j)           gV(i,j,k,bi,bj) = vCf(i,j)
523          ENDDO          ENDDO
524         ENDDO         ENDDO
   
525         IF ( writeDiag ) THEN         IF ( writeDiag ) THEN
526           IF (snapshot_mdsio) THEN           IF (snapshot_mdsio) THEN
527             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 529  C- jmc: change it to keep the Coriolis t
529           ENDIF           ENDIF
530  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
531           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
532             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,
533       &          offsets, myThid)       &          offsets, myThid)
534             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,
535       &          offsets, myThid)       &          offsets, myThid)
536           ENDIF           ENDIF
537  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 422  C- jmc: change it to keep the Coriolis t Line 542  C- jmc: change it to keep the Coriolis t
542           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
543         ENDIF         ENDIF
544  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
   
545        ELSE        ELSE
546         DO j=jMin,jMax         DO j=jMin,jMax
547          DO i=iMin,iMax          DO i=iMin,iMax
# Line 434  C- jmc: change it to keep the Coriolis t Line 553  C- jmc: change it to keep the Coriolis t
553    
554        IF (momAdvection) THEN        IF (momAdvection) THEN
555  C--   Horizontal advection of relative (or absolute) vorticity  C--   Horizontal advection of relative (or absolute) vorticity
556         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
557         &     .AND.useAbsVorticity ) THEN
558          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,
559       &                         uCf,myThid)       &                         uCf,myThid)
560         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
561          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,
562       &                         uCf,myThid)       &                         uCf,myThid)
563         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
564          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,
565       &                         uCf,myThid)       &                         uCf,myThid)
566         ELSE         ELSE
# Line 452  C--   Horizontal advection of relative ( Line 572  C--   Horizontal advection of relative (
572           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)
573          ENDDO          ENDDO
574         ENDDO         ENDDO
575         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
576         &     .AND.useAbsVorticity ) THEN
577          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,          CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
578       &                         vCf,myThid)       &                         vCf,myThid)
579         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
580          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,
581       &                         vCf,myThid)       &                         vCf,myThid)
582         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
583          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,
584       &                         vCf,myThid)       &                         vCf,myThid)
585         ELSE         ELSE
# Line 478  C--   Horizontal advection of relative ( Line 599  C--   Horizontal advection of relative (
599           ENDIF           ENDIF
600  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
601           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
602             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,
603       &          offsets, myThid)       &          offsets, myThid)
604             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,
605       &          offsets, myThid)       &          offsets, myThid)
606           ENDIF           ENDIF
607  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 543  C--   Bernoulli term Line 664  C--   Bernoulli term
664           ENDIF           ENDIF
665  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
666           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
667             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,
668       &          offsets, myThid)       &          offsets, myThid)
669             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,
670       &          offsets, myThid)       &          offsets, myThid)
671          ENDIF           ENDIF
672  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
673         ENDIF         ENDIF
674    
675  C--   end if momAdvection  C--   end if momAdvection
676        ENDIF        ENDIF
677    
678    C--   3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
679          IF ( use3dCoriolis ) THEN
680            CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
681            DO j=jMin,jMax
682             DO i=iMin,iMax
683              gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
684             ENDDO
685            ENDDO
686           IF ( usingCurvilinearGrid ) THEN
687    C-     presently, non zero angleSinC array only supported with Curvilinear-Grid
688            CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
689            DO j=jMin,jMax
690             DO i=iMin,iMax
691              gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
692             ENDDO
693            ENDDO
694           ENDIF
695          ENDIF
696    
697    C--   Non-Hydrostatic (spherical) metric terms
698          IF ( useNHMTerms ) THEN
699           CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
700           DO j=jMin,jMax
701            DO i=iMin,iMax
702             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
703            ENDDO
704           ENDDO
705           CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
706           DO j=jMin,jMax
707            DO i=iMin,iMax
708             gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
709            ENDDO
710           ENDDO
711          ENDIF
712    
713  C--   Set du/dt & dv/dt on boundaries to zero  C--   Set du/dt & dv/dt on boundaries to zero
714        DO j=jMin,jMax        DO j=jMin,jMax
715         DO i=iMin,iMax         DO i=iMin,iMax
# Line 563  C--   Set du/dt & dv/dt on boundaries to Line 719  C--   Set du/dt & dv/dt on boundaries to
719        ENDDO        ENDDO
720    
721  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
722        IF ( debugLevel .GE. debLevB        IF ( debugLevel .GE. debLevC
723       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0
724       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
725       &   .AND. useCubedSphereExchange ) THEN       &   .AND. useCubedSphereExchange ) THEN
# Line 574  C--   Set du/dt & dv/dt on boundaries to Line 730  C--   Set du/dt & dv/dt on boundaries to
730    
731        IF ( writeDiag ) THEN        IF ( writeDiag ) THEN
732          IF (snapshot_mdsio) THEN          IF (snapshot_mdsio) THEN
733            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)
734            CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,           CALL WRITE_LOCAL_RL('KE','I10',1,KE,     bi,bj,k,myIter,myThid)
735       &         myThid)           CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv,   bi,bj,k,myIter,myThid)
736            CALL WRITE_LOCAL_RL('Du','I10',1,guDiss,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
737            CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
738            CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)           CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
           CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)  
           CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)  
           CALL WRITE_LOCAL_RL('D','I10',1,hDiv,bi,bj,k,myIter,myThid)  
739          ENDIF          ENDIF
740  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
741          IF (useMNC .AND. snapshot_mnc) THEN          IF (useMNC .AND. snapshot_mnc) THEN
742            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,
      &          offsets, myThid)  
           CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dt',tension,  
743       &          offsets, myThid)       &          offsets, myThid)
744            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,
745       &          offsets, myThid)       &          offsets, myThid)
746            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,
747       &          offsets, myThid)       &          offsets, myThid)
748            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,
749       &          offsets, myThid)       &          offsets, myThid)
750            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'W3',omega3,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
751       &          offsets, myThid)       &          offsets, myThid)
752            CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'KE',KE,            CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
      &          offsets, myThid)  
           CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'D', hDiv,  
753       &          offsets, myThid)       &          offsets, myThid)
754          ENDIF          ENDIF
755  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 608  C--   Set du/dt & dv/dt on boundaries to Line 757  C--   Set du/dt & dv/dt on boundaries to
757    
758  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
759        IF ( useDiagnostics ) THEN        IF ( useDiagnostics ) THEN
760          CALL DIAGNOSTICS_FILL(KE,    'momKE   ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(KE,     'momKE   ',k,1,2,bi,bj,myThid)
         CALL DIAGNOSTICS_FILL(hDiv,  'momHDiv ',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)  
761         IF (momViscosity) THEN         IF (momViscosity) THEN
762          CALL DIAGNOSTICS_FILL(guDiss,'Um_Diss ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(hDiv,   'momHDiv ',k,1,2,bi,bj,myThid)
763          CALL DIAGNOSTICS_FILL(gvDiss,'Vm_Diss ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
764            CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
765            CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
766         ENDIF         ENDIF
767            CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
768         &                                'Um_Advec',k,1,2,bi,bj,myThid)
769            CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
770         &                                'Vm_Advec',k,1,2,bi,bj,myThid)
771        ENDIF        ENDIF
772  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
773    

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

  ViewVC Help
Powered by ViewVC 1.1.22