/[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.32 by jmc, Thu Jul 13 03:01:21 2006 UTC revision 1.34 by jmc, Thu Jul 13 21:32:38 2006 UTC
# Line 30  C     == Global variables == Line 30  C     == Global variables ==
30  #include "EEPARAMS.h"  #include "EEPARAMS.h"
31  #include "PARAMS.h"  #include "PARAMS.h"
32  #include "GRID.h"  #include "GRID.h"
33    #include "SURFACE.h"
34  #include "DYNVARS.h"  #include "DYNVARS.h"
35  #include "NH_VARS.h"  #include "NH_VARS.h"
36    
# Line 58  C     xA            :: W-Cell face area Line 59  C     xA            :: W-Cell face area
59  C     yA            :: W-Cell face area normal to Y  C     yA            :: W-Cell face area normal to Y
60  C     rThickC_W     :: thickness (in r-units) of W-Cell at Western Edge  C     rThickC_W     :: thickness (in r-units) of W-Cell at Western Edge
61  C     rThickC_S     :: thickness (in r-units) of W-Cell at Southern Edge  C     rThickC_S     :: thickness (in r-units) of W-Cell at Southern Edge
62    C     rThickC_C     :: thickness (in r-units) of W-Cell (centered on W pt)
63  C     recip_rThickC :: reciprol thickness of W-Cell (centered on W-point)  C     recip_rThickC :: reciprol thickness of W-Cell (centered on W-point)
64  C     flx_NS        :: vertical momentum flux, meridional direction  C     flx_NS        :: vertical momentum flux, meridional direction
65  C     flx_EW        :: vertical momentum flux, zonal direction  C     flx_EW        :: vertical momentum flux, zonal direction
# Line 71  C     i,j,k         :: Loop counters Line 73  C     i,j,k         :: Loop counters
73        _RS    yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS    yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74        _RL    rThickC_W    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL    rThickC_W    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75        _RL    rThickC_S    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL    rThickC_S    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76          _RL    rThickC_C    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77        _RL    recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL    recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78        _RL    flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL    flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79        _RL    flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL    flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
# Line 134  C--   Compute grid factor arround a W-po Line 137  C--   Compute grid factor arround a W-po
137          DO j=1-Oly,sNy+Oly          DO j=1-Oly,sNy+Oly
138           DO i=1-Olx,sNx+Olx           DO i=1-Olx,sNx+Olx
139  C-    note: assume fluid @ smaller k than bottom: does not work in p-coordinate !  C-    note: assume fluid @ smaller k than bottom: does not work in p-coordinate !
            rThickC_W(i,j) =  
      &          drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )  
      &        + drF( k )*MIN( _hFacW(i,j,k  ,bi,bj), halfRS )  
            rThickC_S(i,j) =  
      &          drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )  
      &        + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )  
140             IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN             IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN
141               recip_rThickC(i,j) = 0.               recip_rThickC(i,j) = 0.
142             ELSE             ELSE
143               recip_rThickC(i,j) = 1. _d 0 /               recip_rThickC(i,j) = 1. _d 0 /
144       &        ( drF(k-1)*halfRS +       &        ( drF(k-1)*halfRS
145       &        + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )       &        + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )
146       &        )       &        )
147             ENDIF             ENDIF
148    c          IF (momViscosity) THEN
149    #ifdef NONLIN_FRSURF
150               rThickC_C(i,j) =
151         &          drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
152         &        + drF( k )*MIN( h0FacC(i,j,k  ,bi,bj), halfRS )
153    #else
154               rThickC_C(i,j) =
155         &          drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
156         &        + drF( k )*MIN( _hFacC(i,j,k  ,bi,bj), halfRS )
157    #endif
158               rThickC_W(i,j) =
159         &          drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )
160         &        + drF( k )*MIN( _hFacW(i,j,k  ,bi,bj), halfRS )
161               rThickC_S(i,j) =
162         &          drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )
163         &        + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )
164  C     W-Cell Western face area:  C     W-Cell Western face area:
165             xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)             xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
166  C     W-Cell Southern face area:  C     W-Cell Southern face area:
167             yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)             yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
168    c          ENDIF
169           ENDDO           ENDDO
170          ENDDO          ENDDO
171    
# Line 227  C Viscous Flux on Western face Line 241  C Viscous Flux on Western face
241       &              *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))       &              *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
242       &              *_recip_dxC(i,j,bi,bj)*xA(i,j)       &              *_recip_dxC(i,j,bi,bj)*xA(i,j)
243  cOld &              *_recip_dxC(i,j,bi,bj)*rThickC_W(i,j)  cOld &              *_recip_dxC(i,j,bi,bj)*rThickC_W(i,j)
244         &              *cosFacU(j,bi,bj)
245       &       + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL       &       + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL
246       &              *(del2w(i,j)-del2w(i-1,j))       &              *(del2w(i,j)-del2w(i-1,j))
247       &              *_recip_dxC(i,j,bi,bj)*xA(i,j)       &              *_recip_dxC(i,j,bi,bj)*xA(i,j)
# Line 234  cOld &              *_recip_dxC(i,j,bi,b Line 249  cOld &              *_recip_dxC(i,j,bi,b
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             ENDDO             ENDDO
255            ENDDO            ENDDO
# Line 246  C Viscous Flux on Southern face Line 261  C Viscous Flux on Southern face
261       &              *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))       &              *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
262       &              *_recip_dyC(i,j,bi,bj)*yA(i,j)       &              *_recip_dyC(i,j,bi,bj)*yA(i,j)
263  cOld &              *_recip_dyC(i,j,bi,bj)*rThickC_S(i,j)  cOld &              *_recip_dyC(i,j,bi,bj)*rThickC_S(i,j)
264    #ifdef ISOTROPIC_COS_SCALING
265         &              *cosFacV(j,bi,bj)
266    #endif
267       &       + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL       &       + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL
268       &              *(del2w(i,j)-del2w(i,j-1))       &              *(del2w(i,j)-del2w(i,j-1))
269       &              *_recip_dyC(i,j,bi,bj)*yA(i,j)       &              *_recip_dyC(i,j,bi,bj)*yA(i,j)
# Line 254  cOld &              *_recip_dyC(i,j,bi,b Line 272  cOld &              *_recip_dyC(i,j,bi,b
272  #ifdef COSINEMETH_III  #ifdef COSINEMETH_III
273       &              *sqCosFacV(j,bi,bj)       &              *sqCosFacV(j,bi,bj)
274  #else  #else
275       &              *CosFacV(j,bi,bj)       &              *cosFacV(j,bi,bj)
276  #endif  #endif
277  #endif  #endif
278             ENDDO             ENDDO
# Line 298  C--        prepare for next level (k+1) Line 316  C--        prepare for next level (k+1)
316    
317          IF (no_slip_sides) THEN          IF (no_slip_sides) THEN
318  C-     No-slip BCs impose a drag at walls...  C-     No-slip BCs impose a drag at walls...
319  c         CALL MOM_W_SIDEDRAG(            CALL MOM_W_SIDEDRAG(
320  c    I        bi,bj,k,       I               bi,bj,k,
321  c    O        gwAdd,       I               wVel, del2w,
322  c    I        myThid)       I               rThickC_C, recip_rThickC,
323  c         DO j=jMin,jMax       I               viscAh_W, viscA4_W,
324  c          DO i=iMin,iMax       O               gwAdd,
325  c           gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)       I               myThid )
326  c          ENDDO            DO j=jMin,jMax
327  c         ENDDO             DO i=iMin,iMax
328                gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)
329               ENDDO
330              ENDDO
331          ENDIF          ENDIF
332    
333  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|

Legend:
Removed from v.1.32  
changed lines
  Added in v.1.34

  ViewVC Help
Powered by ViewVC 1.1.22