/[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.3 by adcroft, Fri Mar 24 18:31:31 2000 UTC revision 1.12 by jmc, Mon Apr 5 21:46:18 2004 UTC
# Line 1  Line 1 
1    C $Header$
2    C     !DESCRIPTION: \bv
3    C $Name$
4    
5    #include "PACKAGES_CONFIG.h"
6  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
7    
8    CBOP
9    C     !ROUTINE: CALC_GW
10    C     !INTERFACE:
11        SUBROUTINE CALC_GW(            SUBROUTINE CALC_GW(    
12       I        myThid)       I        myThid)
13  C     /==========================================================\  C     !DESCRIPTION: \bv
14  C     | S/R CALC_GW                                              |  C     *==========================================================*
15  C     \==========================================================/  C     | S/R CALC_GW                                              
16        IMPLICIT NONE  C     | o Calculate vert. velocity tendency terms ( NH, QH only )
17    C     *==========================================================*
18    C     | In NH and QH, the vertical momentum tendency must be
19    C     | calculated explicitly and included as a source term
20    C     | for a 3d pressure eqn. Calculate that term here.
21    C     | This routine is not used in HYD calculations.
22    C     *==========================================================*
23    C     \ev
24    
25    C     !USES:
26          IMPLICIT NONE
27  C     == Global variables ==  C     == Global variables ==
28  #include "SIZE.h"  #include "SIZE.h"
29  #include "DYNVARS.h"  #include "DYNVARS.h"
# Line 14  C     == Global variables == Line 31  C     == Global variables ==
31  #include "EEPARAMS.h"  #include "EEPARAMS.h"
32  #include "PARAMS.h"  #include "PARAMS.h"
33  #include "GRID.h"  #include "GRID.h"
 #include "CG2D.h"  
34  #include "GW.h"  #include "GW.h"
35  #include "CG3D.h"  #include "CG3D.h"
36    
37    C     !INPUT/OUTPUT PARAMETERS:
38  C     == Routine arguments ==  C     == Routine arguments ==
39  C     myThid - Instance number for this innvocation of CALC_GW  C     myThid - Instance number for this innvocation of CALC_GW
40        INTEGER myThid        INTEGER myThid
41    
42  #ifdef ALLOW_NONHYDROSTATIC  #ifdef ALLOW_NONHYDROSTATIC
43    
44    C     !LOCAL VARIABLES:
45  C     == Local variables ==  C     == Local variables ==
46    C     bi, bj,      :: Loop counters
47    C     iMin, iMax,
48    C     jMin, jMax
49    C     flx_NS       :: Temp. used for fVol meridional terms.
50    C     flx_EW       :: Temp. used for fVol zonal terms.
51    C     flx_Up       :: Temp. used for fVol vertical terms.
52    C     flx_Dn       :: Temp. used for fVol vertical terms.
53        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
       _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      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)  
   
54        _RL    flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL    flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55        _RL    flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL    flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
56        _RL    flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL    flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57        _RL    flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL    flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58  C     I,J,K - Loop counters  C     I,J,K - Loop counters
59        INTEGER i,j,k, kP1        INTEGER i,j,k, kP1, kUp
60        _RL  wOverride        _RL  wOverride
61        _RS  hFacROpen        _RS  hFacROpen
62        _RS  hFacRClosed        _RS  hFacRClosed
63        _RL ab15,ab05        _RL ab15,ab05
64          _RL slipSideFac
65        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
66    
67        _RL  Half        _RL  Half
# Line 54  C     I,J,K - Loop counters Line 71  C     I,J,K - Loop counters
71  #define In sNx  #define In sNx
72  #define J0 1  #define J0 1
73  #define Jn sNy  #define Jn sNy
74    CEOP
75    
76    ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
77    
78  C     Adams-Bashforth timestepping weights  C     Adams-Bashforth timestepping weights
79        ab15=1.5+abeps        ab15 =  1.5 _d 0 + abeps
80        ab05=-0.5-abeps        ab05 = -0.5 _d 0 - abeps
81    
82    C     Lateral friction (no-slip, free slip, or half slip):
83          IF ( no_slip_sides ) THEN
84            slipSideFac = -Half
85          ELSE
86            slipSideFac =  Half
87          ENDIF
88    C-    half slip was used before ; keep it for now.
89            slipSideFac = 0. _d 0
90    
91        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
92         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
93          DO K=1,Nr          DO K=1,Nr
94           DO j=1-OLy,sNy+OLy           DO j=1-OLy,sNy+OLy
95            DO i=1-OLx,sNx+OLx            DO i=1-OLx,sNx+OLx
96               gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
97             gW(i,j,k,bi,bj) = 0.             gW(i,j,k,bi,bj) = 0.
98            ENDDO            ENDDO
99           ENDDO           ENDDO
# Line 101  C Flux on Southern face Line 131  C Flux on Southern face
131       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
132             Flx_NS(I,J,bi,bj)=             Flx_NS(I,J,bi,bj)=
133       &     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))
134       &    -viscAh*_recip_dyC(I,J,bi,bj)*(       &    -viscAh*_recip_dyC(I,J,bi,bj)  
135       &                   wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj) )       &      *(1. _d 0 + slipSideFac*
136         &         (maskS(I,J,K-1,bi,bj)+maskS(I,J,K,bi,bj)-2. _d 0))
137         &                   *(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
138            ENDDO            ENDDO
139           ENDDO           ENDDO
140  C Flux on Western face  C Flux on Western face
# Line 113  C Flux on Western face Line 145  C Flux on Western face
145       &        +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))       &        +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
146             Flx_EW(I,J,bi,bj)=             Flx_EW(I,J,bi,bj)=
147       &     tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))       &     tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
148       &    -viscAh*_recip_dxC(I,J,bi,bj)*(       &    -viscAh*_recip_dxC(I,J,bi,bj)
149       &                   wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj) )       &      *(1. _d 0 + slipSideFac*
150         &         (maskW(I,J,K-1,bi,bj)+maskW(I,J,K,bi,bj)-2. _d 0))
151         &                   *(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
152            ENDDO            ENDDO
153           ENDDO           ENDDO
154  C Flux on Lower face  C Flux on Lower face
# Line 124  C Flux on Lower face Line 158  C Flux on Lower face
158             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj))             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj))
159             Flx_Dn(I,J,bi,bj)=             Flx_Dn(I,J,bi,bj)=
160       &     tmp_WbarZ*tmp_WbarZ       &     tmp_WbarZ*tmp_WbarZ
161       &    -viscAr*recip_drF(K)*( wVel(I,J,K,bi,bj)       &    -viscAr*recip_drF(K)
162       &                -wOverRide*wVel(I,J,Kp1,bi,bj) )       &       *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
163            ENDDO            ENDDO
164           ENDDO           ENDDO
165  C        Divergence of fluxes  C        Divergence of fluxes
# Line 165  caja           and an hFacUS (above V po Line 199  caja           and an hFacUS (above V po
199          ENDDO          ENDDO
200         ENDDO         ENDDO
201        ENDDO        ENDDO
202        DO bj=myByLo(myThid),myByHi(myThid)  
203         DO bi=myBxLo(myThid),myBxHi(myThid)  #ifdef ALLOW_OBCS
204          DO K=1,Nr        IF (useOBCS) THEN
205           DO j=J0,Jn  C--   This call is aesthetic: it makes the W field
206            DO i=I0,In  C     consistent with the OBs but this has no algorithmic
207             gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)  C     impact. This is purely for diagnostic purposes.
208            ENDDO         DO bj=myByLo(myThid),myByHi(myThid)
209            DO bi=myBxLo(myThid),myBxHi(myThid)
210             DO K=1,Nr
211              CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
212           ENDDO           ENDDO
213          ENDDO          ENDDO
214         ENDDO         ENDDO
215        ENDDO        ENDIF
216    #endif /* ALLOW_OBCS */
217    
218  #endif /* ALLOW_NONHYDROSTATIC */  #endif /* ALLOW_NONHYDROSTATIC */
219    
220        RETURN        RETURN
221        END        END
   

Legend:
Removed from v.1.3  
changed lines
  Added in v.1.12

  ViewVC Help
Powered by ViewVC 1.1.22