/[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.7 by jmc, Tue Mar 6 17:10:29 2001 UTC revision 1.24 by jmc, Thu Dec 15 02:06:06 2005 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C     !DESCRIPTION: \bv
3  C $Name$  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           myTime, myIter, 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"
 #include "DYNVARS.h"  
 #include "FFIELDS.h"  
29  #include "EEPARAMS.h"  #include "EEPARAMS.h"
30  #include "PARAMS.h"  #include "PARAMS.h"
31  #include "GRID.h"  #include "GRID.h"
32  #include "GW.h"  #include "DYNVARS.h"
33  #include "CG3D.h"  #include "NH_VARS.h"
34    
35    C     !INPUT/OUTPUT PARAMETERS:
36  C     == Routine arguments ==  C     == Routine arguments ==
37  C     myThid - Instance number for this innvocation of CALC_GW  C     myTime :: Current time in simulation
38    C     myIter :: Current iteration number in simulation
39    C     myThid :: Thread number for this instance of the routine.
40          _RL     myTime
41          INTEGER myIter
42        INTEGER myThid        INTEGER myThid
43    
44  #ifdef ALLOW_NONHYDROSTATIC  #ifdef ALLOW_NONHYDROSTATIC
45    
46    C     !LOCAL VARIABLES:
47  C     == Local variables ==  C     == Local variables ==
48    C     bi, bj,      :: Loop counters
49    C     iMin, iMax,
50    C     jMin, jMax
51    C     flx_NS       :: Temp. used for fVol meridional terms.
52    C     flx_EW       :: Temp. used for fVol zonal terms.
53    C     flx_Up       :: Temp. used for fVol vertical terms.
54    C     flx_Dn       :: Temp. used for fVol vertical terms.
55        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)  
   
56        _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)
57        _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)
58        _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)
59        _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)
60  C     I,J,K - Loop counters        _RL    fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61        INTEGER i,j,k, kP1, kUp        _RL    fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62          _RL    del2w(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63    C     i,j,k - Loop counters
64          INTEGER i,j,k, kP1
65        _RL  wOverride        _RL  wOverride
66        _RS  hFacROpen        _RS  hFacWtmp
67        _RS  hFacRClosed        _RS  hFacStmp
68          _RS  hFacCtmp
69          _RS  recip_hFacCtmp
70        _RL ab15,ab05        _RL ab15,ab05
71          _RL slipSideFac
72        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ        _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
73    
74        _RL  Half        _RL  Half
75        PARAMETER(Half=0.5D0)        PARAMETER(Half=0.5D0)
76    CEOP
77    
78  #define I0 1        iMin = 1
79  #define In sNx        iMax = sNx
80  #define J0 1        jMin = 1
81  #define Jn sNy        jMax = sNy
82    
83  C     Adams-Bashforth timestepping weights  C     Adams-Bashforth timestepping weights
84        ab15=1.5+abeps        IF (myIter .EQ. 0) THEN
85        ab05=-0.5-abeps         ab15 =  1.0 _d 0
86           ab05 =  0.0 _d 0
87          ELSE
88           ab15 =  1.5 _d 0 + abeps
89           ab05 = -0.5 _d 0 - abeps
90          ENDIF
91    
92    C     Lateral friction (no-slip, free slip, or half slip):
93          IF ( no_slip_sides ) THEN
94            slipSideFac = -1. _d 0
95          ELSE
96            slipSideFac =  1. _d 0
97          ENDIF
98    CML   half slip was used before ; keep the line for now, but half slip is
99    CML   not used anywhere in the code as far as I can see.
100    C        slipSideFac = 0. _d 0
101    
102        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
103         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
# Line 80  C For each tile Line 119  C For each tile
119         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
120    
121  C Boundaries condition at top  C Boundaries condition at top
122          DO J=J0,Jn          DO J=jMin,jMax
123           DO I=I0,In           DO I=iMin,iMax
124            Flx_Dn(I,J,bi,bj)=0.            Flx_Dn(I,J,bi,bj)=0.
125           ENDDO           ENDDO
126          ENDDO          ENDDO
# Line 94  C Sweep down column Line 133  C Sweep down column
133            Kp1=Nr            Kp1=Nr
134            wOverRide=0.            wOverRide=0.
135           endif           endif
136    C     horizontal bi-harmonic dissipation
137             IF (momViscosity .AND. viscA4W.NE.0. ) THEN
138    C     calculate the horizontal Laplacian of vertical flow
139    C     Zonal flux d/dx W
140              DO j=1-Oly,sNy+Oly
141               fZon(1-Olx,j)=0.
142               DO i=1-Olx+1,sNx+Olx
143                fZon(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
144         &           *_dyG(i,j,bi,bj)
145         &           *_recip_dxC(i,j,bi,bj)
146         &           *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
147    #ifdef COSINEMETH_III
148         &           *sqcosFacU(J,bi,bj)
149    #endif
150               ENDDO
151              ENDDO
152    C     Meridional flux d/dy W
153              DO i=1-Olx,sNx+Olx
154               fMer(I,1-Oly)=0.
155              ENDDO
156              DO j=1-Oly+1,sNy+Oly
157               DO i=1-Olx,sNx+Olx
158                fMer(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
159         &           *_dxG(i,j,bi,bj)
160         &           *_recip_dyC(i,j,bi,bj)
161         &           *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
162    #ifdef ISOTROPIC_COS_SCALING
163    #ifdef COSINEMETH_III
164         &           *sqCosFacV(j,bi,bj)
165    #endif
166    #endif
167               ENDDO
168              ENDDO
169              
170    C     del^2 W
171    C     Difference of zonal fluxes ...
172              DO j=1-Oly,sNy+Oly
173               DO i=1-Olx,sNx+Olx-1
174                del2w(i,j)=fZon(i+1,j)-fZon(i,j)
175               ENDDO
176               del2w(sNx+Olx,j)=0.
177              ENDDO
178    
179    C     ... add difference of meridional fluxes and divide by volume
180              DO j=1-Oly,sNy+Oly-1
181               DO i=1-Olx,sNx+Olx
182    C     First compute the fraction of open water for the w-control volume
183    C     at the southern face
184                hFacCtmp=max(hFacC(I,J,K-1,bi,bj)-Half,0. _d 0)
185         &           +   min(hFacC(I,J,K  ,bi,bj),Half)
186                IF (hFacCtmp .GT. 0.) THEN
187                 recip_hFacCtmp = 1./hFacCtmp
188                ELSE
189                 recip_hFacCtmp = 0. _d 0
190                ENDIF
191                del2w(i,j)=recip_rA(i,j,bi,bj)
192         &           *recip_drC(k)*recip_hFacCtmp
193         &           *(
194         &           del2w(i,j)
195         &           +( fMer(i,j+1)-fMer(i,j) )
196         &           )
197               ENDDO
198              ENDDO
199    C-- No-slip BCs impose a drag at walls...
200    CML ************************************************************
201    CML   No-slip Boundary conditions for bi-harmonic dissipation
202    CML   need to be implemented here!
203    CML ************************************************************
204             ELSE
205    C-    Initialize del2w to zero:
206              DO j=1-Oly,sNy+Oly
207               DO i=1-Olx,sNx+Olx
208                del2w(i,j) = 0. _d 0
209               ENDDO
210              ENDDO
211             ENDIF
212    
213  C Flux on Southern face  C Flux on Southern face
214           DO J=J0,Jn+1           DO J=jMin,jMax+1
215            DO I=I0,In            DO I=iMin,iMax
216    C     First compute the fraction of open water for the w-control volume
217    C     at the southern face
218               hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
219         &         +    min(hFacS(I,J,K  ,bi,bj),Half)
220             tmp_VbarZ=Half*(             tmp_VbarZ=Half*(
221       &          _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)       &          _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
222       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))       &         +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
223             Flx_NS(I,J,bi,bj)=             Flx_NS(I,J,bi,bj)=
224       &     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))
225       &    -viscAh*_recip_dyC(I,J,bi,bj)*(       &    -viscAhW*_recip_dyC(I,J,bi,bj)  
226       &                   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))
227         &        +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
228         &         *wVel(I,J,K,bi,bj))
229         &    +viscA4W*_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
230    #ifdef ISOTROPIC_COS_SCALING
231    #ifdef COSINEMETH_III
232         &    *sqCosFacV(j,bi,bj)
233    #else
234         &    *CosFacV(j,bi,bj)
235    #endif
236    #endif
237    C     The last term is the weighted average of the viscous stress at the open
238    C     fraction of the w control volume and at the closed fraction of the
239    C     the control volume. A more compact but less intelligible version
240    C     of the last three lines is:
241    CML     &       *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
242    CML     &       *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
243            ENDDO            ENDDO
244           ENDDO           ENDDO
245  C Flux on Western face  C Flux on Western face
246           DO J=J0,Jn           DO J=jMin,jMax
247            DO I=I0,In+1            DO I=iMin,iMax+1
248    C     First compute the fraction of open water for the w-control volume
249    C     at the western face
250               hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
251         &         +    min(hFacW(I,J,K  ,bi,bj),Half)
252             tmp_UbarZ=Half*(             tmp_UbarZ=Half*(
253       &         _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)       &         _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
254       &        +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))       &        +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
255             Flx_EW(I,J,bi,bj)=             Flx_EW(I,J,bi,bj)=
256       &     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))
257       &    -viscAh*_recip_dxC(I,J,bi,bj)*(       &    -viscAhW*_recip_dxC(I,J,bi,bj)
258       &                   wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj) )       &      *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
259         &       +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
260         &        *wVel(I,J,K,bi,bj) )
261         &    +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
262    #ifdef COSINEMETH_III
263         &                *sqCosFacU(j,bi,bj)
264    #else
265         &                *CosFacU(j,bi,bj)
266    #endif
267    C     The last term is the weighted average of the viscous stress at the open
268    C     fraction of the w control volume and at the closed fraction of the
269    C     the control volume. A more compact but less intelligible version
270    C     of the last three lines is:
271    CML     &       *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
272    CML     &       *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
273            ENDDO            ENDDO
274           ENDDO           ENDDO
275  C Flux on Lower face  C Flux on Lower face
276           DO J=J0,Jn           DO J=jMin,jMax
277            DO I=I0,In            DO I=iMin,iMax
278             Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)             Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
279             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj))             tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
280         &         +wOverRide*wVel(I,J,Kp1,bi,bj))
281             Flx_Dn(I,J,bi,bj)=             Flx_Dn(I,J,bi,bj)=
282       &     tmp_WbarZ*tmp_WbarZ       &     tmp_WbarZ*tmp_WbarZ
283       &    -viscAr*recip_drF(K)*( wVel(I,J,K,bi,bj)       &    -viscAr*recip_drF(K)
284       &                -wOverRide*wVel(I,J,Kp1,bi,bj) )       &       *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
285            ENDDO            ENDDO
286           ENDDO           ENDDO
287  C        Divergence of fluxes  C        Divergence of fluxes
288           DO J=J0,Jn           DO J=jMin,jMax
289            DO I=I0,In            DO I=iMin,iMax
290             gW(I,J,K,bi,bj) = 0.             gW(I,J,K,bi,bj) = 0.
291       &      -(       &      -(
292       &        +_recip_dxF(I,J,bi,bj)*(       &        +_recip_dxF(I,J,bi,bj)*(
# Line 151  caja           and an hFacUS (above V po Line 306  caja           and an hFacUS (above V po
306         ENDDO         ENDDO
307        ENDDO        ENDDO
308    
309        
310        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
311         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
312          DO K=2,Nr          DO K=2,Nr
313           DO j=J0,Jn           DO j=jMin,jMax
314            DO i=I0,In            DO i=iMin,iMax
315             wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)             wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
316       &     +deltatMom*( ab15*gW(i,j,k,bi,bj)       &     +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)
317       &                 +ab05*gWNM1(i,j,k,bi,bj) )       &                 +ab05*gwNm1(i,j,k,bi,bj) )
318             IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.             IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
319            ENDDO             gwNm1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
          ENDDO  
         ENDDO  
        ENDDO  
       ENDDO  
       DO bj=myByLo(myThid),myByHi(myThid)  
        DO bi=myBxLo(myThid),myBxHi(myThid)  
         DO K=1,Nr  
          DO j=J0,Jn  
           DO i=I0,In  
            gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)  
320            ENDDO            ENDDO
321           ENDDO           ENDDO
322          ENDDO          ENDDO
# Line 197  C     impact. This is purely for diagnos Line 342  C     impact. This is purely for diagnos
342    
343        RETURN        RETURN
344        END        END
   

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

  ViewVC Help
Powered by ViewVC 1.1.22