/[MITgcm]/MITgcm/model/src/calc_gw.F
ViewVC logotype

Diff of /MITgcm/model/src/calc_gw.F

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

revision 1.24 by jmc, Thu Dec 15 02:06:06 2005 UTC revision 1.25 by jmc, Thu Dec 15 21:08:59 2005 UTC
# Line 2  C $Header$ Line 2  C $Header$
2  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
3  C $Name$  C $Name$
4    
5  #include "PACKAGES_CONFIG.h"  c #include "PACKAGES_CONFIG.h"
6  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
7    
8  CBOP  CBOP
9  C     !ROUTINE: CALC_GW  C     !ROUTINE: CALC_GW
10  C     !INTERFACE:  C     !INTERFACE:
11        SUBROUTINE CALC_GW(            SUBROUTINE CALC_GW(
12       I           myTime, myIter, myThid )       I                    myTime, myIter, myThid )
13  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
14  C     *==========================================================*  C     *==========================================================*
15  C     | S/R CALC_GW                                                C     | S/R CALC_GW
16  C     | o Calculate vert. velocity tendency terms ( NH, QH only )  C     | o Calculate vert. velocity tendency terms ( NH, QH only )
17  C     *==========================================================*  C     *==========================================================*
18  C     | In NH and QH, the vertical momentum tendency must be  C     | In NH and QH, the vertical momentum tendency must be
19  C     | calculated explicitly and included as a source term  C     | calculated explicitly and included as a source term
20  C     | for a 3d pressure eqn. Calculate that term here.  C     | for a 3d pressure eqn. Calculate that term here.
21  C     | This routine is not used in HYD calculations.  C     | This routine is not used in HYD calculations.
22  C     *==========================================================*  C     *==========================================================*
23  C     \ev  C     \ev
24    
25  C     !USES:  C     !USES:
26        IMPLICIT NONE        IMPLICIT NONE
# Line 67  C     i,j,k - Loop counters Line 67  C     i,j,k - Loop counters
67        _RS  hFacStmp        _RS  hFacStmp
68        _RS  hFacCtmp        _RS  hFacCtmp
69        _RS  recip_hFacCtmp        _RS  recip_hFacCtmp
       _RL ab15,ab05  
70        _RL slipSideFac        _RL slipSideFac
71        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
72    
# Line 75  C     i,j,k - Loop counters Line 74  C     i,j,k - Loop counters
74        PARAMETER(Half=0.5D0)        PARAMETER(Half=0.5D0)
75  CEOP  CEOP
76    
77    C Catch barotropic mode
78          IF ( Nr .LT. 2 ) RETURN
79    
80        iMin = 1        iMin = 1
81        iMax = sNx        iMax = sNx
82        jMin = 1        jMin = 1
83        jMax = sNy        jMax = sNy
84    
 C     Adams-Bashforth timestepping weights  
       IF (myIter .EQ. 0) THEN  
        ab15 =  1.0 _d 0  
        ab05 =  0.0 _d 0  
       ELSE  
        ab15 =  1.5 _d 0 + abeps  
        ab05 = -0.5 _d 0 - abeps  
       ENDIF  
   
85  C     Lateral friction (no-slip, free slip, or half slip):  C     Lateral friction (no-slip, free slip, or half slip):
86        IF ( no_slip_sides ) THEN        IF ( no_slip_sides ) THEN
87          slipSideFac = -1. _d 0          slipSideFac = -1. _d 0
88        ELSE        ELSE
89          slipSideFac =  1. _d 0          slipSideFac =  1. _d 0
90        ENDIF        ENDIF
91  CML   half slip was used before ; keep the line for now, but half slip is  CML   half slip was used before ; keep the line for now, but half slip is
92  CML   not used anywhere in the code as far as I can see.  CML   not used anywhere in the code as far as I can see.
93  C        slipSideFac = 0. _d 0  C        slipSideFac = 0. _d 0
94    
95    C For each tile
96        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
97         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
98    
99    C Initialise gW to zero
100          DO K=1,Nr          DO K=1,Nr
101           DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
102            DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
# Line 108  C        slipSideFac = 0. _d 0 Line 104  C        slipSideFac = 0. _d 0
104            ENDDO            ENDDO
105           ENDDO           ENDDO
106          ENDDO          ENDDO
        ENDDO  
       ENDDO  
   
 C Catch barotropic mode  
       IF ( Nr .LT. 2 ) RETURN  
   
 C For each tile  
       DO bj=myByLo(myThid),myByHi(myThid)  
        DO bi=myBxLo(myThid),myBxHi(myThid)  
107    
108  C Boundaries condition at top  C Boundaries condition at top
109          DO J=jMin,jMax          DO J=jMin,jMax
# Line 166  C     Meridional flux d/dy W Line 153  C     Meridional flux d/dy W
153  #endif  #endif
154             ENDDO             ENDDO
155            ENDDO            ENDDO
156              
157  C     del^2 W  C     del^2 W
158  C     Difference of zonal fluxes ...  C     Difference of zonal fluxes ...
159            DO j=1-Oly,sNy+Oly            DO j=1-Oly,sNy+Oly
# Line 222  C     at the southern face Line 209  C     at the southern face
209       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
210             Flx_NS(I,J,bi,bj)=             Flx_NS(I,J,bi,bj)=
211       &     tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))       &     tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
212       &    -viscAhW*_recip_dyC(I,J,bi,bj)         &    -viscAhW*_recip_dyC(I,J,bi,bj)
213       &       *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))       &       *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
214       &        +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)       &        +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
215       &         *wVel(I,J,K,bi,bj))       &         *wVel(I,J,K,bi,bj))
# Line 261  C     at the western face Line 248  C     at the western face
248       &    +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))       &    +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
249  #ifdef COSINEMETH_III  #ifdef COSINEMETH_III
250       &                *sqCosFacU(j,bi,bj)       &                *sqCosFacU(j,bi,bj)
251  #else  #else
252       &                *CosFacU(j,bi,bj)       &                *CosFacU(j,bi,bj)
253  #endif  #endif
254  C     The last term is the weighted average of the viscous stress at the open  C     The last term is the weighted average of the viscous stress at the open
# Line 276  C Flux on Lower face Line 263  C Flux on Lower face
263           DO J=jMin,jMax           DO J=jMin,jMax
264            DO I=iMin,iMax            DO I=iMin,iMax
265             Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)             Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
266             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
267       &         +wOverRide*wVel(I,J,Kp1,bi,bj))       &         +wOverRide*wVel(I,J,Kp1,bi,bj))
268             Flx_Dn(I,J,bi,bj)=             Flx_Dn(I,J,bi,bj)=
269       &     tmp_WbarZ*tmp_WbarZ       &     tmp_WbarZ*tmp_WbarZ
270       &    -viscAr*recip_drF(K)       &    -viscAr*recip_drF(K)
271       &       *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )       &       *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
272            ENDDO            ENDDO
273           ENDDO           ENDDO
# Line 297  C        Divergence of fluxes Line 284  C        Divergence of fluxes
284       &              Flx_Up(I,J,bi,bj)  -Flx_Dn(I,J,bi,bj) )       &              Flx_Up(I,J,bi,bj)  -Flx_Dn(I,J,bi,bj) )
285       &       )       &       )
286  caja    * recip_hFacU(I,J,K,bi,bj)  caja    * recip_hFacU(I,J,K,bi,bj)
287  caja   NOTE:  This should be included    caja   NOTE:  This should be included
288  caja   but we need an hFacUW (above U points)  caja   but we need an hFacUW (above U points)
289  caja           and an hFacUS (above V points) too...  caja           and an hFacUS (above V points) too...
290            ENDDO            ENDDO
291           ENDDO           ENDDO
         ENDDO  
        ENDDO  
       ENDDO  
292    
293    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
294    
295        DO bj=myByLo(myThid),myByHi(myThid)  C-    Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights)
296         DO bi=myBxLo(myThid),myBxHi(myThid)  C     and save gW_[n] into gwNm1 for the next time step.
297          DO K=2,Nr  c#ifdef ALLOW_ADAMSBASHFORTH_3
298           DO j=jMin,jMax  c        CALL ADAMS_BASHFORTH3(
299            DO i=iMin,iMax  c    I                          bi, bj, k,
300             wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)  c    U                          gW, gwNm,
301       &     +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)  c    I                          momStartAB, myIter, myThid )
302       &                 +ab05*gwNm1(i,j,k,bi,bj) )  c#else /* ALLOW_ADAMSBASHFORTH_3 */
303             IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.           CALL ADAMS_BASHFORTH2(
304             gwNm1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)       I                          bi, bj, k,
305            ENDDO       U                          gW, gwNm1,
306           ENDDO       I                          myIter, myThid )
307          ENDDO  c#endif /* ALLOW_ADAMSBASHFORTH_3 */
        ENDDO  
       ENDDO  
308    
309  #ifdef ALLOW_OBCS  C-    end of the k loop
       IF (useOBCS) THEN  
 C--   This call is aesthetic: it makes the W field  
 C     consistent with the OBs but this has no algorithmic  
 C     impact. This is purely for diagnostic purposes.  
        DO bj=myByLo(myThid),myByHi(myThid)  
         DO bi=myBxLo(myThid),myBxHi(myThid)  
          DO K=1,Nr  
           CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )  
          ENDDO  
310          ENDDO          ENDDO
311    
312    C-    end of bi,bj loops
313         ENDDO         ENDDO
314        ENDIF        ENDDO
 #endif /* ALLOW_OBCS */  
315    
316  #endif /* ALLOW_NONHYDROSTATIC */  #endif /* ALLOW_NONHYDROSTATIC */
317    

Legend:
Removed from v.1.24  
changed lines
  Added in v.1.25

  ViewVC Help
Powered by ViewVC 1.1.22