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

Legend:
Removed from v.1.14  
changed lines
  Added in v.1.73

  ViewVC Help
Powered by ViewVC 1.1.22