/[MITgcm]/MITgcm/pkg/mom_fluxform/mom_fluxform.F
ViewVC logotype

Diff of /MITgcm/pkg/mom_fluxform/mom_fluxform.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.2 by adcroft, Fri Aug 17 18:40:30 2001 UTC revision 1.8 by jmc, Sun Jan 26 21:18:50 2003 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2  C $Name$  C $Name$
3    
4    CBOI
5    C !TITLE: pkg/mom\_advdiff
6    C !AUTHORS: adcroft@mit.edu
7    C !INTRODUCTION: Flux-form Momentum Equations Package
8    C
9    C Package "mom\_fluxform" provides methods for calculating explicit terms
10    C in the momentum equation cast in flux-form:
11    C \begin{eqnarray*}
12    C G^u & = & -\frac{1}{\rho} \partial_x \phi_h
13    C           -\nabla \cdot {\bf v} u
14    C           -fv
15    C           +\frac{1}{\rho} \nabla \cdot {\bf \tau}^x
16    C           + \mbox{metrics}
17    C \\
18    C G^v & = & -\frac{1}{\rho} \partial_y \phi_h
19    C           -\nabla \cdot {\bf v} v
20    C           +fu
21    C           +\frac{1}{\rho} \nabla \cdot {\bf \tau}^y
22    C           + \mbox{metrics}
23    C \end{eqnarray*}
24    C where ${\bf v}=(u,v,w)$ and $\tau$, the stress tensor, includes surface
25    C stresses as well as internal viscous stresses.
26    CEOI
27    
28  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
29    
30    CBOP
31    C !ROUTINE: MOM_FLUXFORM
32    
33    C !INTERFACE: ==========================================================
34        SUBROUTINE MOM_FLUXFORM(        SUBROUTINE MOM_FLUXFORM(
35       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
36       I        phi_hyd,KappaRU,KappaRV,       I        phi_hyd,KappaRU,KappaRV,
37       U        fVerU, fVerV,       U        fVerU, fVerV,
38       I        myCurrentTime, myIter, myThid)       I        myTime,myIter,myThid)
 C     /==========================================================\  
 C     | S/R MOM_FLUXFORM                                         |  
 C     | o Form the right hand-side of the momentum equation.     |  
 C     |==========================================================|  
 C     | Terms are evaluated one layer at a time working from     |  
 C     | the bottom to the top. The vertically integrated         |  
 C     | barotropic flow tendency term is evluated by summing the |  
 C     | tendencies.                                              |  
 C     | Notes:                                                   |  
 C     | We have not sorted out an entirely satisfactory formula  |  
 C     | for the diffusion equation bc with lopping. The present  |  
 C     | form produces a diffusive flux that does not scale with  |  
 C     | open-area. Need to do something to solidfy this and to   |  
 C     | deal "properly" with thin walls.                         |  
 C     \==========================================================/  
       IMPLICIT NONE  
39    
40    C !DESCRIPTION:
41    C Calculates all the horizontal accelerations except for the implicit surface
42    C pressure gradient and implciit vertical viscosity.
43    
44    C !USES: ===============================================================
45  C     == Global variables ==  C     == Global variables ==
46          IMPLICIT NONE
47  #include "SIZE.h"  #include "SIZE.h"
48  #include "DYNVARS.h"  #include "DYNVARS.h"
49  #include "FFIELDS.h"  #include "FFIELDS.h"
# Line 34  C     == Global variables == Line 52  C     == Global variables ==
52  #include "GRID.h"  #include "GRID.h"
53  #include "SURFACE.h"  #include "SURFACE.h"
54    
55  C     == Routine arguments ==  C !INPUT PARAMETERS: ===================================================
56  C     fZon    - Work array for flux of momentum in the east-west  C  bi,bj                :: tile indices
57  C               direction at the west face of a cell.  C  iMin,iMax,jMin,jMAx  :: loop ranges
58  C     fMer    - Work array for flux of momentum in the north-south  C  k                    :: vertical level
59  C               direction at the south face of a cell.  C  kUp                  :: =1 or 2 for consecutive k
60  C     fVerU   - Flux of momentum in the vertical  C  kDown                :: =2 or 1 for consecutive k
61  C     fVerV     direction out of the upper face of a cell K  C  phi_hyd              :: hydrostatic pressure (perturbation)
62  C               ( flux into the cell above ).  C  KappaRU              :: vertical viscosity
63  C     phi_hyd - Hydrostatic pressure  C  KappaRV              :: vertical viscosity
64  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation  C  fVerU                :: vertical flux of U, 2 1/2 dim for pipe-lining
65  C                                      results will be set.  C  fVerV                :: vertical flux of V, 2 1/2 dim for pipe-lining
66  C     kUp, kDown                     - Index for upper and lower layers.  C  myTime               :: current time
67  C     myThid - Instance number for this innvocation of CALC_MOM_RHS  C  myIter               :: current time-step number
68    C  myThid               :: thread number
69          INTEGER bi,bj,iMin,iMax,jMin,jMax
70          INTEGER k,kUp,kDown
71        _RL phi_hyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL phi_hyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
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 fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
75        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
76        INTEGER kUp,kDown        _RL     myTime
       _RL     myCurrentTime  
77        INTEGER myIter        INTEGER myIter
78        INTEGER myThid        INTEGER myThid
       INTEGER bi,bj,iMin,iMax,jMin,jMax  
79    
80  C     == Local variables ==  C !OUTPUT PARAMETERS: ==================================================
81  C     ab15, ab05    - Weights for Adams-Bashforth time stepping scheme.  C None - updates gU() and gV() in common blocks
82  C     i,j,k         - Loop counters  
83    C !LOCAL VARIABLES: ====================================================
84    C  i,j                  :: loop indices
85    C  aF                   :: advective flux
86    C  vF                   :: viscous flux
87    C  v4F                  :: bi-harmonic viscous flux
88    C  vrF                  :: vertical viscous flux
89    C  cF                   :: Coriolis acceleration
90    C  mT                   :: Metric terms
91    C  pF                   :: Pressure gradient
92    C  fZon                 :: zonal fluxes
93    C  fMer                 :: meridional fluxes
94          INTEGER i,j
95          _RL aF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96          _RL vF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97          _RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98          _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99          _RL cF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100          _RL mT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101          _RL pF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102          _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103          _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104  C     wMaskOverride - Land sea flag override for top layer.  C     wMaskOverride - Land sea flag override for top layer.
105  C     afFacMom      - Tracer parameters for turning terms  C     afFacMom      - Tracer parameters for turning terms
106  C     vfFacMom        on and off.  C     vfFacMom        on and off.
# Line 70  C     mTFacMom        pfFacMom - Pressur Line 110  C     mTFacMom        pfFacMom - Pressur
110  C                     cfFacMom - Coriolis terms  C                     cfFacMom - Coriolis terms
111  C                     foFacMom - Forcing  C                     foFacMom - Forcing
112  C                     mTFacMom - Metric term  C                     mTFacMom - Metric term
 C     vF            - Temporary holding viscous term (Laplacian)  
 C     v4F           - Temporary holding viscous term (Biharmonic)  
 C     cF            - Temporary holding coriolis term.  
 C     mT            - Temporary holding metric terms(s).  
 C     pF            - Temporary holding pressure|potential gradient terms.  
113  C     uDudxFac, AhDudxFac, etc ... individual term tracer parameters  C     uDudxFac, AhDudxFac, etc ... individual term tracer parameters
       _RL      aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      cF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL      pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL    fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL    fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
114        _RS    hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS    hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115        _RS  r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS  r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116        _RS      xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS      xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
# Line 93  C     uDudxFac, AhDudxFac, etc ... indiv Line 119  C     uDudxFac, AhDudxFac, etc ... indiv
119        _RL  vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL  vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120        _RL  uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL  uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121        _RL  vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL  vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122          _RL  rTransU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123          _RL  rTransV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124  C     I,J,K - Loop counters  C     I,J,K - Loop counters
       INTEGER i,j,k  
125  C     rVelMaskOverride - Factor for imposing special surface boundary conditions  C     rVelMaskOverride - Factor for imposing special surface boundary conditions
126  C                        ( set according to free-surface condition ).  C                        ( set according to free-surface condition ).
127  C     hFacROpen        - Lopped cell factos used tohold fraction of open  C     hFacROpen        - Lopped cell factos used tohold fraction of open
# Line 124  C     xxxFac - On-off tracer parameters Line 151  C     xxxFac - On-off tracer parameters
151        _RL  phyFac        _RL  phyFac
152        _RL  vForcFac        _RL  vForcFac
153        _RL  mtFacV        _RL  mtFacV
 C     ab05, ab15 - Adams-Bashforth time-stepping weights.  
       _RL  ab05, ab15  
154        INTEGER km1,kp1        INTEGER km1,kp1
155        _RL wVelBottomOverride        _RL wVelBottomOverride
156        LOGICAL bottomDragTerms        LOGICAL bottomDragTerms
157        _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
158    CEOP
159    
160        km1=MAX(1,k-1)        km1=MAX(1,k-1)
161        kp1=MIN(Nr,k+1)        kp1=MIN(Nr,k+1)
# Line 150  C     Initialise intermediate terms Line 176  C     Initialise intermediate terms
176          pF(i,j)   = 0.          pF(i,j)   = 0.
177          fZon(i,j) = 0.          fZon(i,j) = 0.
178          fMer(i,j) = 0.          fMer(i,j) = 0.
179            rTransU(i,j) = 0.
180            rTransV(i,j) = 0.
181         ENDDO         ENDDO
182        ENDDO        ENDDO
183    
# Line 194  C-- with stagger time stepping, grad Phi Line 222  C-- with stagger time stepping, grad Phi
222          phyFac = 0.          phyFac = 0.
223        ENDIF        ENDIF
224    
 C--   Adams-Bashforth weighting factors  
       ab15   =  1.5 _d 0 + abEps  
       ab05   = -0.5 _d 0 - abEps  
     
225  C--   Calculate open water fraction at vorticity points  C--   Calculate open water fraction at vorticity points
226        CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)        CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
227    
# Line 230  C     Calculate velocity field "volume t Line 254  C     Calculate velocity field "volume t
254    
255        CALL MOM_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)        CALL MOM_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)
256    
257    C---  First call (k=1): compute vertical adv. flux fVerU(kUp) & fVerV(kUp)
258          IF (momAdvection.AND.k.EQ.1) THEN
259    
260    C-    Calculate vertical transports above U & V points (West & South face):
261           CALL MOM_CALC_RTRANS( k, bi, bj,
262         O                       rTransU, rTransV,
263         I                       myTime, myIter, myThid)
264    
265    C-    Free surface correction term (flux at k=1)
266           CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,rTransU,af,myThid)
267           DO j=jMin,jMax
268            DO i=iMin,iMax
269             fVerU(i,j,kUp) = af(i,j)
270            ENDDO
271           ENDDO
272    
273           CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,rTransV,af,myThid)
274           DO j=jMin,jMax
275            DO i=iMin,iMax
276             fVerV(i,j,kUp) = af(i,j)
277            ENDDO
278           ENDDO
279    
280    C---  endif momAdvection & k=1
281          ENDIF
282    
283    
284    C---  Calculate vertical transports (at k+1) below U & V points :
285          IF (momAdvection) THEN
286           CALL MOM_CALC_RTRANS( k+1, bi, bj,
287         O                       rTransU, rTransV,
288         I                       myTime, myIter, myThid)
289          ENDIF
290    
291    
292  C---- Zonal momentum equation starts here  C---- Zonal momentum equation starts here
293    
294  C     Bi-harmonic term del^2 U -> v4F  C     Bi-harmonic term del^2 U -> v4F
# Line 274  C     Combine fluxes -> fMer Line 333  C     Combine fluxes -> fMer
333    
334  C--   Vertical flux (fVer is at upper face of "u" cell)  C--   Vertical flux (fVer is at upper face of "u" cell)
335    
 C--   Free surface correction term (flux at k=1)  
       IF (momAdvection.AND.k.EQ.1) THEN  
        CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,af,myThid)  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          fVerU(i,j,kUp) = af(i,j)  
         ENDDO  
        ENDDO  
       ENDIF  
336  C     Mean flow component of vertical flux (at k+1) -> aF  C     Mean flow component of vertical flux (at k+1) -> aF
337        IF (momAdvection)        IF (momAdvection)
338       & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,af,myThid)       & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,rTransU,af,myThid)
339    
340  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
341        IF (momViscosity.AND..NOT.implicitViscosity)        IF (momViscosity.AND..NOT.implicitViscosity)
# Line 327  C--   Tendency is minus divergence of th Line 377  C--   Tendency is minus divergence of th
377         ENDDO         ENDDO
378        ENDDO        ENDDO
379    
380    #ifdef NONLIN_FRSURF
381    C-- account for 3.D divergence of the flow in rStar coordinate:
382          IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
383           DO j=jMin,jMax
384            DO i=iMin,iMax
385             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
386         &     - (rStarExpW(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
387         &       *uVel(i,j,k,bi,bj)
388            ENDDO
389           ENDDO
390          ENDIF
391          IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
392           DO j=jMin,jMax
393            DO i=iMin,iMax
394             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
395         &     - rStarDhWDt(i,j,bi,bj)*uVel(i,j,k,bi,bj)
396            ENDDO
397           ENDDO
398          ENDIF
399    #endif /* NONLIN_FRSURF */
400    
401  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
402        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
403  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
# Line 351  C--   Forcing term Line 422  C--   Forcing term
422        IF (momForcing)        IF (momForcing)
423       &  CALL EXTERNAL_FORCING_U(       &  CALL EXTERNAL_FORCING_U(
424       I     iMin,iMax,jMin,jMax,bi,bj,k,       I     iMin,iMax,jMin,jMax,bi,bj,k,
425       I     myCurrentTime,myThid)       I     myTime,myThid)
426    
427  C--   Metric terms for curvilinear grid systems  C--   Metric terms for curvilinear grid systems
428        IF (usingSphericalPolarMTerms) THEN        IF (useNHMTerms) THEN
429  C      o Spherical polar grid metric terms  C      o Non-hydrosatic metric terms
430         CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)         CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
431         DO j=jMin,jMax         DO j=jMin,jMax
432          DO i=iMin,iMax          DO i=iMin,iMax
433           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)           gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
434          ENDDO          ENDDO
435         ENDDO         ENDDO
436          ENDIF
437          IF (usingSphericalPolarMTerms) THEN
438         CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)         CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)
439         DO j=jMin,jMax         DO j=jMin,jMax
440          DO i=iMin,iMax          DO i=iMin,iMax
# Line 422  C     Combine fluxes -> fMer Line 495  C     Combine fluxes -> fMer
495    
496  C--   Vertical flux (fVer is at upper face of "v" cell)  C--   Vertical flux (fVer is at upper face of "v" cell)
497    
 C--   Free surface correction term (flux at k=1)  
       IF (momAdvection.AND.k.EQ.1) THEN  
        CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,af,myThid)  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          fVerV(i,j,kUp) = af(i,j)  
         ENDDO  
        ENDDO  
       ENDIF  
498  C     o Mean flow component of vertical flux  C     o Mean flow component of vertical flux
499        IF (momAdvection)        IF (momAdvection)
500       & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,af,myThid)       & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,rTransV,af,myThid)
501    
502  C     Eddy component of vertical flux (interior component only) -> vrF  C     Eddy component of vertical flux (interior component only) -> vrF
503        IF (momViscosity.AND..NOT.implicitViscosity)        IF (momViscosity.AND..NOT.implicitViscosity)
# Line 475  C--   Tendency is minus divergence of th Line 539  C--   Tendency is minus divergence of th
539         ENDDO         ENDDO
540        ENDDO        ENDDO
541    
542    #ifdef NONLIN_FRSURF
543    C-- account for 3.D divergence of the flow in rStar coordinate:
544          IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
545           DO j=jMin,jMax
546            DO i=iMin,iMax
547             gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
548         &     - (rStarExpS(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
549         &       *vVel(i,j,k,bi,bj)
550            ENDDO
551           ENDDO
552          ENDIF
553          IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
554           DO j=jMin,jMax
555            DO i=iMin,iMax
556             gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
557         &     - rStarDhSDt(i,j,bi,bj)*vVel(i,j,k,bi,bj)
558            ENDDO
559           ENDDO
560          ENDIF
561    #endif /* NONLIN_FRSURF */
562    
563  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
564        IF (momViscosity.AND.no_slip_sides) THEN        IF (momViscosity.AND.no_slip_sides) THEN
565  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
# Line 499  C--   Forcing term Line 584  C--   Forcing term
584        IF (momForcing)        IF (momForcing)
585       & CALL EXTERNAL_FORCING_V(       & CALL EXTERNAL_FORCING_V(
586       I     iMin,iMax,jMin,jMax,bi,bj,k,       I     iMin,iMax,jMin,jMax,bi,bj,k,
587       I     myCurrentTime,myThid)       I     myTime,myThid)
588    
589  C--   Metric terms for curvilinear grid systems  C--   Metric terms for curvilinear grid systems
590        IF (usingSphericalPolarMTerms) THEN        IF (useNHMTerms) THEN
591  C      o Spherical polar grid metric terms  C      o Spherical polar grid metric terms
592         CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)         CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
593         DO j=jMin,jMax         DO j=jMin,jMax
# Line 510  C      o Spherical polar grid metric ter Line 595  C      o Spherical polar grid metric ter
595           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)           gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
596          ENDDO          ENDDO
597         ENDDO         ENDDO
598          ENDIF
599          IF (usingSphericalPolarMTerms) THEN
600         CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)         CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)
601         DO j=jMin,jMax         DO j=jMin,jMax
602          DO i=iMin,iMax          DO i=iMin,iMax
# Line 543  C     Note. As coded here, coriolis will Line 630  C     Note. As coded here, coriolis will
630         ENDDO         ENDDO
631        ENDDO        ENDDO
632  #endif /* INCLUDE_CD_CODE */  #endif /* INCLUDE_CD_CODE */
633          IF (nonHydrostatic.OR.quasiHydrostatic) THEN
634           CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,cf,myThid)
635           DO j=jMin,jMax
636            DO i=iMin,iMax
637             gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
638            ENDDO
639           ENDDO
640          ENDIF
641    
642        RETURN        RETURN
643        END        END

Legend:
Removed from v.1.2  
changed lines
  Added in v.1.8

  ViewVC Help
Powered by ViewVC 1.1.22