/[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.50 by baylor, Mon Sep 26 15:27:11 2005 UTC revision 1.70 by jmc, Thu Aug 1 20:12:42 2013 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_MOM_COMMON
6    # include "MOM_COMMON_OPTIONS.h"
7    #endif
8    
9        SUBROUTINE MOM_VECINV(        SUBROUTINE MOM_VECINV(
10       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,       I        bi,bj,k,iMin,iMax,jMin,jMax,
11       I        KappaRU, KappaRV,       I        KappaRU, KappaRV,
12       U        fVerU, fVerV,       I        fVerUkm, fVerVkm,
13         O        fVerUkp, fVerVkp,
14       O        guDiss, gvDiss,       O        guDiss, gvDiss,
15       I        myTime, myIter, myThid)       I        myTime, myIter, myThid )
16  C     /==========================================================\  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 23  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"
 #ifdef ALLOW_MNC  
 #include "MNC_PARAMS.h"  
 #endif  
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  #ifdef ALLOW_TIMEAVE
43  #include "TIMEAVE_STATV.h"  # 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  #endif
52    
53  C     == Routine arguments ==  C     == Routine arguments ==
54  C     fVerU  :: Flux of momentum in the vertical direction, out of the upper  C     bi,bj   :: current tile indices
55  C     fVerV  :: face of a cell K ( flux into the cell above ).  C     k       :: current vertical level
56  C     guDiss :: dissipation tendency (all explicit terms), u component  C     iMin,iMax,jMin,jMax :: loop ranges
57  C     gvDiss :: dissipation tendency (all explicit terms), v component  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    C     fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
63    
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          _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76          _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       INTEGER kUp,kDown  
79        _RL     myTime        _RL     myTime
80        INTEGER myIter        INTEGER myIter
81        INTEGER myThid        INTEGER myThid
       INTEGER bi,bj,iMin,iMax,jMin,jMax  
82    
83  #ifdef ALLOW_MOM_VECINV  #ifdef ALLOW_MOM_VECINV
84    
# Line 68  C     == Functions == Line 88  C     == Functions ==
88    
89  C     == Local variables ==  C     == Local variables ==
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  c     _RL      mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS hFacZ   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vFld    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2u   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99        _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL del2v   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100        _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101        _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL zStar   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102        _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103        _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL strain  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104        _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KE      (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105  C     I,J,K - Loop counters        _RL omega3  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106        INTEGER i,j,k        _RL vort3   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107  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)  
108        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112        LOGICAL harmonic,biharmonic,useVariableViscosity  C     i,j    :: Loop counters
113          INTEGER i,j
114    C     xxxFac :: On-off tracer parameters used for switching terms off.
115          _RL  ArDudrFac
116          _RL  ArDvdrFac
117          _RL  sideMaskFac
118          LOGICAL bottomDragTerms
119          LOGICAL writeDiag
120    #ifdef ALLOW_AUTODIFF_TAMC
121          INTEGER imomkey
122    #endif
123    
124  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
125        INTEGER offsets(9)        INTEGER offsets(9)
126          CHARACTER*(1) pf
127  #endif  #endif
128    
129  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 110  C--   only the kDown part of fverU/V is Line 131  C--   only the kDown part of fverU/V is
131  C--   the kUp is still required  C--   the kUp is still required
132  C--   In the case of mom_fluxform Kup is set as well  C--   In the case of mom_fluxform Kup is set as well
133  C--   (at least in part)  C--   (at least in part)
134        fVerU(1,1,kUp) = fVerU(1,1,kUp)        fVerUkm(1,1) = fVerUkm(1,1)
135        fVerV(1,1,kUp) = fVerV(1,1,kUp)        fVerVkm(1,1) = fVerVkm(1,1)
136  #endif  #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)        writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
156    
157  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
158        IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN        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          IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
165            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)            CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
166            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 171  C--   (at least in part)
171            offsets(i) = 0            offsets(i) = 0
172          ENDDO          ENDDO
173          offsets(3) = k          offsets(3) = k
174  C       write(*,*) 'offsets = ',(offsets(i),i=1,9)  c       write(*,*) 'offsets = ',(offsets(i),i=1,9)
175        ENDIF        ENDIF
176  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
177    
178  C     Initialise intermediate terms  C--   Initialise intermediate terms
179        DO J=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
180         DO I=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
181          vF(i,j)    = 0.          vF(i,j)    = 0.
182          vrF(i,j)   = 0.          vrF(i,j)   = 0.
183          uCf(i,j)   = 0.          uCf(i,j)   = 0.
184          vCf(i,j)   = 0.          vCf(i,j)   = 0.
 c       mT(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.
# Line 148  c       mT(i,j)    = 0. Line 190  c       mT(i,j)    = 0.
190          gvDiss(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          viscAh_Z(i,j) = 0.  C-    need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
195          viscAh_D(i,j) = 0.          hDiv(i,j)  = 0.
196          viscA4_Z(i,j) = 0.  c       viscAh_Z(i,j) = 0.
197          viscA4_D(i,j) = 0.  c       viscAh_D(i,j) = 0.
198    c       viscA4_Z(i,j) = 0.
199  #ifdef ALLOW_AUTODIFF_TAMC  c       viscA4_D(i,j) = 0.
200          strain(i,j)  = 0. _d 0          strain(i,j)  = 0. _d 0
201          tension(i,j) = 0. _d 0          tension(i,j) = 0. _d 0
202    #ifdef ALLOW_AUTODIFF_TAMC
203            hFacZ(i,j)   = 0. _d 0
204  #endif  #endif
205         ENDDO         ENDDO
206        ENDDO        ENDDO
# Line 164  c       mT(i,j)    = 0. Line 208  c       mT(i,j)    = 0.
208  C--   Term by term tracer parmeters  C--   Term by term tracer parmeters
209  C     o U momentum equation  C     o U momentum equation
210        ArDudrFac    = vfFacMom*1.        ArDudrFac    = vfFacMom*1.
 c     mTFacU       = mtFacMom*1.  
211  C     o V momentum equation  C     o V momentum equation
212        ArDvdrFac    = vfFacMom*1.        ArDvdrFac    = vfFacMom*1.
213  c     mTFacV       = mtFacMom*1.  
214    C note: using standard stencil (no mask) results in under-estimating
215    C       vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
216          IF ( no_slip_sides ) THEN
217            sideMaskFac = sideDragFactor
218          ELSE
219            sideMaskFac = 0. _d 0
220          ENDIF
221    
222        IF (     no_slip_bottom        IF (     no_slip_bottom
223       &    .OR. bottomDragQuadratic.NE.0.       &    .OR. bottomDragQuadratic.NE.0.
# Line 192  C note (jmc) : Dissipation and Vort3 adv Line 242  C note (jmc) : Dissipation and Vort3 adv
242  C              use the same maskZ (and hFacZ)  => needs 2 call(s)  C              use the same maskZ (and hFacZ)  => needs 2 call(s)
243  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)
244    
245        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)  
246    
247        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)        CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
248    
249        CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)        IF (momViscosity) THEN
250    C--    For viscous term, compute horizontal divergence, tension & strain
251    C      and mask relative vorticity (free-slip case):
252    
253        CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)  #ifdef ALLOW_AUTODIFF_TAMC
254    CADJ STORE vort3(:,:) =
255    CADJ &     comlev1_bibj_k, key = imomkey, byte = isbyte
256    #endif
257    
258        IF (useAbsVorticity)         CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
      & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)  
259    
260        IF (momViscosity) THEN         CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
261  C      Calculate Viscosities  
262        CALL MOM_CALC_VISC(         CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
263       I        bi,bj,k,  
264       O        viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,  C-     account for no-slip / free-slip BC:
265       O        harmonic,biharmonic,useVariableViscosity,         DO j=1-OLy,sNy+OLy
266       I        hDiv,vort3,tension,strain,KE,hfacZ,          DO i=1-OLx,sNx+OLx
267       I        myThid)            IF ( hFacZ(i,j).EQ.0. ) THEN
268                vort3(i,j)  = sideMaskFac*vort3(i,j)
269                strain(i,j) = sideMaskFac*strain(i,j)
270              ENDIF
271            ENDDO
272           ENDDO
273    
274    C--    Calculate Lateral Viscosities
275           DO j=1-OLy,sNy+OLy
276            DO i=1-OLx,sNx+OLx
277             viscAh_D(i,j) = viscAhD
278             viscAh_Z(i,j) = viscAhZ
279             viscA4_D(i,j) = viscA4D
280             viscA4_Z(i,j) = viscA4Z
281            ENDDO
282           ENDDO
283           IF ( useVariableVisc ) THEN
284             CALL MOM_CALC_VISC( bi, bj, k,
285         O            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
286         I            hDiv, vort3, tension, strain, KE, hfacZ,
287         I            myThid )
288           ENDIF
289    
290  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
291         IF (biharmonic) THEN         IF (useBiharmonicVisc) THEN
292           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,           CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
293       O                      del2u,del2v,       O                      del2u,del2v,
294       &                      myThid)       &                      myThid)
295           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)           CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
296           CALL MOM_CALC_RELVORT3(bi,bj,k,           CALL MOM_CALC_RELVORT3(bi,bj,k,
297       &                          del2u,del2v,hFacZ,zStar,myThid)       &                          del2u,del2v,hFacZ,zStar,myThid)
298             IF ( writeDiag ) THEN
299               CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
300         &                           bi,bj,k, myIter, myThid )
301               CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
302         &                           bi,bj,k, myIter, myThid )
303               CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
304         &                           bi,bj,k, myIter, myThid )
305               CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
306         &                           bi,bj,k, myIter, myThid )
307             ENDIF
308         ENDIF         ENDIF
309    
310  C      Calculate dissipation terms for U and V equations  C-    Strain diagnostics:
311           IF ( writeDiag ) THEN
312            IF (snapshot_mdsio) THEN
313              CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
314            ENDIF
315    #ifdef ALLOW_MNC
316            IF (useMNC .AND. snapshot_mnc) THEN
317              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
318         &          offsets, myThid)
319            ENDIF
320    #endif /*  ALLOW_MNC  */
321           ENDIF
322    #ifdef ALLOW_DIAGNOSTICS
323           IF ( useDiagnostics ) THEN
324            CALL DIAGNOSTICS_FILL(strain, 'Strain  ',k,1,2,bi,bj,myThid)
325           ENDIF
326    #endif /* ALLOW_DIAGNOSTICS */
327    
328    C---   Calculate dissipation terms for U and V equations
329    
330  C      in terms of tension and strain  C      in terms of tension and strain
331         IF (useStrainTensionVisc) THEN         IF (useStrainTensionVisc) THEN
332           CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,  C        mask strain as if free-slip since side-drag is computed separately
333       I                    hFacZ,           DO j=1-OLy,sNy+OLy
334       I                    viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,            DO i=1-OLx,sNx+OLx
335       I                    harmonic,biharmonic,useVariableViscosity,              IF ( hFacZ(i,j).EQ.0. ) strain(i,j) = 0. _d 0
336       O                    guDiss,gvDiss,            ENDDO
337       I                    myThid)           ENDDO
338             CALL MOM_HDISSIP( bi, bj, k,
339         I            hDiv, vort3, tension, strain, KE, hFacZ,
340         I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
341         I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
342         O            guDiss, gvDiss,
343         I            myThid )
344         ELSE         ELSE
345  C      in terms of vorticity and divergence  C      in terms of vorticity and divergence
346           CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,           CALL MOM_VI_HDISSIP( bi, bj, k,
347       I                       hFacZ,dStar,zStar,       I            hDiv, vort3, tension, strain, KE, hFacZ,dStar,zStar,
348       I                       viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,       I            viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
349       I                       harmonic,biharmonic,useVariableViscosity,       I            useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
350       O                       guDiss,gvDiss,       O            guDiss, gvDiss,
351       &                       myThid)               &            myThid )
352         ENDIF         ENDIF
       ENDIF  
353    
354  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:  C---  Other dissipation terms in Zonal momentum equation
 c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)  
   
 C---- Zonal momentum equation starts here  
355    
356  C--   Vertical flux (fVer is at upper face of "u" cell)  C--   Vertical flux (fVer is at upper face of "u" cell)
   
357  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
358        IF (momViscosity.AND..NOT.implicitViscosity) THEN        IF ( .NOT.implicitViscosity ) THEN
359         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)         CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
   
360  C     Combine fluxes  C     Combine fluxes
361         DO j=jMin,jMax         DO j=jMin,jMax
362          DO i=iMin,iMax          DO i=iMin,iMax
363           fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)           fVerUkp(i,j) = ArDudrFac*vrF(i,j)
364          ENDDO          ENDDO
365         ENDDO         ENDDO
   
366  C--   Tendency is minus divergence of the fluxes  C--   Tendency is minus divergence of the fluxes
367         DO j=2-Oly,sNy+Oly-1         DO j=jMin,jMax
368          DO i=2-Olx,sNx+Olx-1          DO i=iMin,iMax
369           guDiss(i,j) = guDiss(i,j)           guDiss(i,j) = guDiss(i,j)
370       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
371       &   *recip_rAw(i,j,bi,bj)       &   *recip_rAw(i,j,bi,bj)
372       &  *(       &   *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
      &    fVerU(i,j,kDown) - fVerU(i,j,kUp)  
      &   )*rkSign  
373          ENDDO          ENDDO
374         ENDDO         ENDDO
375        ENDIF        ENDIF
376    
377  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
378        IF (momViscosity.AND.no_slip_sides) THEN        IF ( no_slip_sides ) THEN
379  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
380         CALL MOM_U_SIDEDRAG(         CALL MOM_U_SIDEDRAG( bi, bj, k,
381       I        bi,bj,k,       I          uFld, del2u, hFacZ,
382       I        uFld, del2u, hFacZ,       I          viscAh_Z, viscA4_Z,
383       I        viscAh_Z,viscA4_Z,       I          useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
384       I        harmonic,biharmonic,useVariableViscosity,       O          vF,
385       O        vF,       I          myThid )
      I        myThid)  
386         DO j=jMin,jMax         DO j=jMin,jMax
387          DO i=iMin,iMax          DO i=iMin,iMax
388           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 391  C-     No-slip BCs impose a drag at wall
391        ENDIF        ENDIF
392    
393  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
394        IF (momViscosity.AND.bottomDragTerms) THEN        IF ( bottomDragTerms ) THEN
395         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)         CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
396         DO j=jMin,jMax         DO j=jMin,jMax
397          DO i=iMin,iMax          DO i=iMin,iMax
# Line 302  C-    No-slip BCs impose a drag at botto Line 399  C-    No-slip BCs impose a drag at botto
399          ENDDO          ENDDO
400         ENDDO         ENDDO
401        ENDIF        ENDIF
402    #ifdef ALLOW_SHELFICE
403          IF ( useShelfIce.AND.bottomDragTerms ) THEN
404           CALL SHELFICE_U_DRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
405           DO j=jMin,jMax
406            DO i=iMin,iMax
407             guDiss(i,j) = guDiss(i,j) + vF(i,j)
408            ENDDO
409           ENDDO
410          ENDIF
411    #endif /* ALLOW_SHELFICE */
412    
413  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  
414    
415  C--   Vertical flux (fVer is at upper face of "v" cell)  C--   Vertical flux (fVer is at upper face of "v" cell)
   
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) THEN        IF ( .NOT.implicitViscosity ) THEN
418         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)         CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
   
419  C     Combine fluxes -> fVerV  C     Combine fluxes -> fVerV
420         DO j=jMin,jMax         DO j=jMin,jMax
421          DO i=iMin,iMax          DO i=iMin,iMax
422           fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)           fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
423          ENDDO          ENDDO
424         ENDDO         ENDDO
   
425  C--   Tendency is minus divergence of the fluxes  C--   Tendency is minus divergence of the fluxes
426         DO j=jMin,jMax         DO j=jMin,jMax
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 ( 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( bi, bj, k,
440       I        bi,bj,k,       I          vFld, del2v, hFacZ,
441       I        vFld, del2v, hFacZ,       I          viscAh_Z, viscA4_Z,
442       I        viscAh_Z,viscA4_Z,       I          useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
443       I        harmonic,biharmonic,useVariableViscosity,       O          vF,
444       O        vF,       I          myThid )
      I        myThid)  
445         DO j=jMin,jMax         DO j=jMin,jMax
446          DO i=iMin,iMax          DO i=iMin,iMax
447           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)           gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
448          ENDDO          ENDDO
449         ENDDO         ENDDO
450        ENDIF        ENDIF
451    
452  C-    No-slip BCs impose a drag at bottom  C-    No-slip BCs impose a drag at bottom
453        IF (momViscosity.AND.bottomDragTerms) THEN        IF ( bottomDragTerms ) THEN
454         CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)         CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
455         DO j=jMin,jMax         DO j=jMin,jMax
456          DO i=iMin,iMax          DO i=iMin,iMax
# 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.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--   if (momViscosity) end of block.
473          ENDIF
474    
475  C--   Metric terms for curvilinear grid systems  C-    Return to standard hfacZ (min-4) and mask vort3 accordingly:
476  c     IF (usingSphericalPolarMTerms) THEN  c     CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
477  C      o Spherical polar grid metric terms  
478  c      CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)  C-    Vorticity diagnostics:
479  c      DO j=jMin,jMax        IF ( writeDiag ) THEN
480  c       DO i=iMin,iMax          IF (snapshot_mdsio) THEN
481  c        gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)            CALL WRITE_LOCAL_RL('Z3','I10',1,vort3, bi,bj,k,myIter,myThid)
482  c       ENDDO          ENDIF
483  c      ENDDO  #ifdef ALLOW_MNC
484  c     ENDIF          IF (useMNC .AND. snapshot_mnc) THEN
485              CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
486         &          offsets, myThid)
487            ENDIF
488    #endif /*  ALLOW_MNC  */
489          ENDIF
490    #ifdef ALLOW_DIAGNOSTICS
491          IF ( useDiagnostics ) THEN
492            CALL DIAGNOSTICS_FILL(vort3,  'momVort3',k,1,2,bi,bj,myThid)
493          ENDIF
494    #endif /* ALLOW_DIAGNOSTICS */
495    
496    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
497    
498    C---  Prepare for Advection & Coriolis terms:
499    C-    Mask relative vorticity and calculate absolute vorticity
500          DO j=1-OLy,sNy+OLy
501           DO i=1-OLx,sNx+OLx
502             IF ( hFacZ(i,j).EQ.0. ) vort3(i,j) = 0.
503           ENDDO
504          ENDDO
505          IF (useAbsVorticity)
506         &  CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
507    
508  C--   Horizontal Coriolis terms  C--   Horizontal Coriolis terms
509  c     IF (useCoriolis .AND. .NOT.useCDscheme  c     IF (useCoriolis .AND. .NOT.useCDscheme
# Line 401  C- jmc: change it to keep the Coriolis t Line 527  C- jmc: change it to keep the Coriolis t
527           gV(i,j,k,bi,bj) = vCf(i,j)           gV(i,j,k,bi,bj) = vCf(i,j)
528          ENDDO          ENDDO
529         ENDDO         ENDDO
   
530         IF ( writeDiag ) THEN         IF ( writeDiag ) THEN
531           IF (snapshot_mdsio) THEN           IF (snapshot_mdsio) THEN
532             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 534  C- jmc: change it to keep the Coriolis t
534           ENDIF           ENDIF
535  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
536           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
537             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,
538       &          offsets, myThid)       &          offsets, myThid)
539             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,
540       &          offsets, myThid)       &          offsets, myThid)
541           ENDIF           ENDIF
542  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 422  C- jmc: change it to keep the Coriolis t Line 547  C- jmc: change it to keep the Coriolis t
547           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)           CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
548         ENDIF         ENDIF
549  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
   
550        ELSE        ELSE
551         DO j=jMin,jMax         DO j=jMin,jMax
552          DO i=iMin,iMax          DO i=iMin,iMax
# Line 434  C- jmc: change it to keep the Coriolis t Line 558  C- jmc: change it to keep the Coriolis t
558    
559        IF (momAdvection) THEN        IF (momAdvection) THEN
560  C--   Horizontal advection of relative (or absolute) vorticity  C--   Horizontal advection of relative (or absolute) vorticity
561         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
562         &     .AND.useAbsVorticity ) THEN
563          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,
564       &                         uCf,myThid)       &                         uCf,myThid)
565         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
566          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,
567       &                         uCf,myThid)       &                         uCf,myThid)
568         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
569          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,
570       &                         uCf,myThid)       &                         uCf,myThid)
571         ELSE         ELSE
# Line 452  C--   Horizontal advection of relative ( Line 577  C--   Horizontal advection of relative (
577           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)
578          ENDDO          ENDDO
579         ENDDO         ENDDO
580         IF (highOrderVorticity.AND.useAbsVorticity) THEN         IF ( (highOrderVorticity.OR.upwindVorticity)
581         &     .AND.useAbsVorticity ) THEN
582          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,
583       &                         vCf,myThid)       &                         vCf,myThid)
584         ELSEIF (highOrderVorticity) THEN         ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
585          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,
586       &                         vCf,myThid)       &                         vCf,myThid)
587         ELSEIF (useAbsVorticity) THEN         ELSEIF ( useAbsVorticity ) THEN
588          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,
589       &                         vCf,myThid)       &                         vCf,myThid)
590         ELSE         ELSE
# Line 478  C--   Horizontal advection of relative ( Line 604  C--   Horizontal advection of relative (
604           ENDIF           ENDIF
605  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
606           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
607             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,
608       &          offsets, myThid)       &          offsets, myThid)
609             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,
610       &          offsets, myThid)       &          offsets, myThid)
611           ENDIF           ENDIF
612  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 543  C--   Bernoulli term Line 669  C--   Bernoulli term
669           ENDIF           ENDIF
670  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
671           IF (useMNC .AND. snapshot_mnc) THEN           IF (useMNC .AND. snapshot_mnc) THEN
672             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,
673       &          offsets, myThid)       &          offsets, myThid)
674             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,
675       &          offsets, myThid)       &          offsets, myThid)
676          ENDIF           ENDIF
677  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
678         ENDIF         ENDIF
679    
680  C--   end if momAdvection  C--   end if momAdvection
681        ENDIF        ENDIF
682    
683    C--   3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
684          IF ( use3dCoriolis ) THEN
685            CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
686            DO j=jMin,jMax
687             DO i=iMin,iMax
688              gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
689             ENDDO
690            ENDDO
691           IF ( usingCurvilinearGrid ) THEN
692    C-     presently, non zero angleSinC array only supported with Curvilinear-Grid
693            CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
694            DO j=jMin,jMax
695             DO i=iMin,iMax
696              gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
697             ENDDO
698            ENDDO
699           ENDIF
700          ENDIF
701    
702    C--   Non-Hydrostatic (spherical) metric terms
703          IF ( useNHMTerms ) THEN
704           CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
705           DO j=jMin,jMax
706            DO i=iMin,iMax
707             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
708            ENDDO
709           ENDDO
710           CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
711           DO j=jMin,jMax
712            DO i=iMin,iMax
713             gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
714            ENDDO
715           ENDDO
716          ENDIF
717    
718  C--   Set du/dt & dv/dt on boundaries to zero  C--   Set du/dt & dv/dt on boundaries to zero
719        DO j=jMin,jMax        DO j=jMin,jMax
720         DO i=iMin,iMax         DO i=iMin,iMax
# Line 563  C--   Set du/dt & dv/dt on boundaries to Line 724  C--   Set du/dt & dv/dt on boundaries to
724        ENDDO        ENDDO
725    
726  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
727        IF ( debugLevel .GE. debLevB        IF ( debugLevel .GE. debLevC
728       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0       &   .AND. k.EQ.4 .AND. myIter.EQ.nIter0
729       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
730       &   .AND. useCubedSphereExchange ) THEN       &   .AND. useCubedSphereExchange ) THEN
# Line 574  C--   Set du/dt & dv/dt on boundaries to Line 735  C--   Set du/dt & dv/dt on boundaries to
735    
736        IF ( writeDiag ) THEN        IF ( writeDiag ) THEN
737          IF (snapshot_mdsio) THEN          IF (snapshot_mdsio) THEN
738            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)
739            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)
740       &         myThid)           CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv,   bi,bj,k,myIter,myThid)
741            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)
742            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)
743            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)  
744          ENDIF          ENDIF
745  #ifdef ALLOW_MNC  #ifdef ALLOW_MNC
746          IF (useMNC .AND. snapshot_mnc) THEN          IF (useMNC .AND. snapshot_mnc) THEN
747            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,  
748       &          offsets, myThid)       &          offsets, myThid)
749            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,
750       &          offsets, myThid)       &          offsets, myThid)
751            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,
752       &          offsets, myThid)       &          offsets, myThid)
753            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,
754       &          offsets, myThid)       &          offsets, myThid)
755            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,
756       &          offsets, myThid)       &          offsets, myThid)
757            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,  
758       &          offsets, myThid)       &          offsets, myThid)
759          ENDIF          ENDIF
760  #endif /*  ALLOW_MNC  */  #endif /*  ALLOW_MNC  */
# Line 608  C--   Set du/dt & dv/dt on boundaries to Line 762  C--   Set du/dt & dv/dt on boundaries to
762    
763  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
764        IF ( useDiagnostics ) THEN        IF ( useDiagnostics ) THEN
765          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(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)  
766         IF (momViscosity) THEN         IF (momViscosity) THEN
767          CALL DIAGNOSTICS_FILL(guDiss,'Um_Diss ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(hDiv,   'momHDiv ',k,1,2,bi,bj,myThid)
768          CALL DIAGNOSTICS_FILL(gvDiss,'Vm_Diss ',k,1,2,bi,bj,myThid)          CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
769         ENDIF          CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
770            CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
771           ENDIF
772            CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
773         &                                'Um_Advec',k,1,2,bi,bj,myThid)
774            CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
775         &                                'Vm_Advec',k,1,2,bi,bj,myThid)
776        ENDIF        ENDIF
777  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
778    

Legend:
Removed from v.1.50  
changed lines
  Added in v.1.70

  ViewVC Help
Powered by ViewVC 1.1.22