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

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

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


Revision 1.14 - (show annotations) (download)
Wed May 19 08:42:00 2004 UTC (20 years ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint53c_post, checkpoint53d_pre
Changes since 1.13: +31 -24 lines
o changes to calc_gw.F
 - bug fix for w^2 term near the bottom boundary
 - (hopefully) improve the lateral slip boundary condtions for use with
   partial/looped cells
 - because the bug fix changes two verification experiments anyway (exp5
   and plume on slope), change the lateral boundary condition from half slip
   to the value of no_slip_sides
o update the verification experiments exp5 and plume_on_slope

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.13 2004/04/06 00:31:54 jmc Exp $
2 C !DESCRIPTION: \bv
3 C $Name: $
4
5 #include "PACKAGES_CONFIG.h"
6 #include "CPP_OPTIONS.h"
7
8 CBOP
9 C !ROUTINE: CALC_GW
10 C !INTERFACE:
11 SUBROUTINE CALC_GW(
12 I myThid)
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | S/R CALC_GW
16 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 ==
28 #include "SIZE.h"
29 #include "DYNVARS.h"
30 #include "EEPARAMS.h"
31 #include "PARAMS.h"
32 #include "GRID.h"
33 #include "GW.h"
34 #include "CG3D.h"
35
36 C !INPUT/OUTPUT PARAMETERS:
37 C == Routine arguments ==
38 C myThid - Instance number for this innvocation of CALC_GW
39 INTEGER myThid
40
41 #ifdef ALLOW_NONHYDROSTATIC
42
43 C !LOCAL VARIABLES:
44 C == Local variables ==
45 C bi, bj, :: Loop counters
46 C iMin, iMax,
47 C jMin, jMax
48 C flx_NS :: Temp. used for fVol meridional terms.
49 C flx_EW :: Temp. used for fVol zonal terms.
50 C flx_Up :: Temp. used for fVol vertical terms.
51 C flx_Dn :: Temp. used for fVol vertical terms.
52 INTEGER bi,bj,iMin,iMax,jMin,jMax
53 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
54 _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55 _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
56 _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57 C I,J,K - Loop counters
58 INTEGER i,j,k, kP1, kUp
59 _RL wOverride
60 _RS hFacROpen
61 _RS hFacRClosed
62 _RL ab15,ab05
63 _RL slipSideFac
64 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
65
66 _RL Half
67 _RL One
68 PARAMETER(Half=0.5D0)
69 PARAMETER(One=0.5D0)
70 CEOP
71
72 ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
73
74 iMin = 1
75 iMax = sNx
76 jMin = 1
77 jMax = sNy
78
79 C Adams-Bashforth timestepping weights
80 ab15 = 1.5 _d 0 + abeps
81 ab05 = -0.5 _d 0 - abeps
82
83 C Lateral friction (no-slip, free slip, or half slip):
84 IF ( no_slip_sides ) THEN
85 slipSideFac = -One
86 ELSE
87 slipSideFac = One
88 ENDIF
89 C- half slip was used before ; keep it for now.
90 C slipSideFac = 0. _d 0
91
92 DO bj=myByLo(myThid),myByHi(myThid)
93 DO bi=myBxLo(myThid),myBxHi(myThid)
94 DO K=1,Nr
95 DO j=1-OLy,sNy+OLy
96 DO i=1-OLx,sNx+OLx
97 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
98 gW(i,j,k,bi,bj) = 0.
99 ENDDO
100 ENDDO
101 ENDDO
102 ENDDO
103 ENDDO
104
105 C Catch barotropic mode
106 IF ( Nr .LT. 2 ) RETURN
107
108 C For each tile
109 DO bj=myByLo(myThid),myByHi(myThid)
110 DO bi=myBxLo(myThid),myBxHi(myThid)
111
112 C Boundaries condition at top
113 DO J=jMin,jMax
114 DO I=iMin,iMax
115 Flx_Dn(I,J,bi,bj)=0.
116 ENDDO
117 ENDDO
118
119 C Sweep down column
120 DO K=2,Nr
121 Kp1=K+1
122 wOverRide=1.
123 if (K.EQ.Nr) then
124 Kp1=Nr
125 wOverRide=0.
126 endif
127 C Flux on Southern face
128 DO J=jMin,jMax+1
129 DO I=iMin,iMax
130 tmp_VbarZ=Half*(
131 & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
132 & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
133 Flx_NS(I,J,bi,bj)=
134 & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
135 & -viscAh*_recip_dyC(I,J,bi,bj)
136 & *(1. _d 0 + slipSideFac*
137 & (maskS(I,J,K-1,bi,bj)*maskS(I,J,K,bi,bj)-2. _d 0))
138 & *(max(hFacS(I,J,K-1,bi,bj)-Half,0)
139 & +min(hFacS(I,J,K,bi,bj),Half))
140 & *(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
141 ENDDO
142 ENDDO
143 C Flux on Western face
144 DO J=jMin,jMax
145 DO I=iMin,iMax+1
146 tmp_UbarZ=Half*(
147 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
148 & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
149 Flx_EW(I,J,bi,bj)=
150 & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
151 & -viscAh*_recip_dxC(I,J,bi,bj)
152 & *(1. _d 0 + slipSideFac*
153 & (maskW(I,J,K-1,bi,bj)*maskW(I,J,K,bi,bj)-1. _d 0))
154 & *(max(hFacW(I,J,K-1,bi,bj)-Half,0)
155 & +min(hFacW(I,J,K,bi,bj),Half))
156 & *(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
157 ENDDO
158 ENDDO
159 C Flux on Lower face
160 DO J=jMin,jMax
161 DO I=iMin,iMax
162 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
163 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
164 & +wOverRide*wVel(I,J,Kp1,bi,bj))
165 Flx_Dn(I,J,bi,bj)=
166 & tmp_WbarZ*tmp_WbarZ
167 & -viscAr*recip_drF(K)
168 & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
169 ENDDO
170 ENDDO
171 C Divergence of fluxes
172 DO J=jMin,jMax
173 DO I=iMin,iMax
174 gW(I,J,K,bi,bj) = 0.
175 & -(
176 & +_recip_dxF(I,J,bi,bj)*(
177 & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
178 & +_recip_dyF(I,J,bi,bj)*(
179 & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
180 & +recip_drC(K) *(
181 & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
182 & )
183 caja * recip_hFacU(I,J,K,bi,bj)
184 caja NOTE: This should be included
185 caja but we need an hFacUW (above U points)
186 caja and an hFacUS (above V points) too...
187 ENDDO
188 ENDDO
189 ENDDO
190 ENDDO
191 ENDDO
192
193
194 DO bj=myByLo(myThid),myByHi(myThid)
195 DO bi=myBxLo(myThid),myBxHi(myThid)
196 DO K=2,Nr
197 DO j=jMin,jMax
198 DO i=iMin,iMax
199 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
200 & +deltatMom*( ab15*gW(i,j,k,bi,bj)
201 & +ab05*gWNM1(i,j,k,bi,bj) )
202 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
203 ENDDO
204 ENDDO
205 ENDDO
206 ENDDO
207 ENDDO
208
209 #ifdef ALLOW_OBCS
210 IF (useOBCS) THEN
211 C-- This call is aesthetic: it makes the W field
212 C consistent with the OBs but this has no algorithmic
213 C impact. This is purely for diagnostic purposes.
214 DO bj=myByLo(myThid),myByHi(myThid)
215 DO bi=myBxLo(myThid),myBxHi(myThid)
216 DO K=1,Nr
217 CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
218 ENDDO
219 ENDDO
220 ENDDO
221 ENDIF
222 #endif /* ALLOW_OBCS */
223
224 #endif /* ALLOW_NONHYDROSTATIC */
225
226 RETURN
227 END

  ViewVC Help
Powered by ViewVC 1.1.22