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

Legend:
Removed from v.1.17  
changed lines
  Added in v.1.55

  ViewVC Help
Powered by ViewVC 1.1.22