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

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

  ViewVC Help
Powered by ViewVC 1.1.22