/[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.19 - (show annotations) (download)
Tue May 31 14:49:38 2005 UTC (18 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57i_post, checkpoint57n_post, checkpoint57l_post, checkpoint57j_post, checkpoint57o_post, checkpoint57k_post
Changes since 1.18: +2 -2 lines
Added run-time parameters nh_Am2 which scales the non-hydrostatic terms and
changes internal scales (i.e. allows convection at different Rayleigh numbers).

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.18 2004/12/09 09:01:07 mlosch 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 _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
58 _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
59 _RL del2w(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
60 C I,J,K - Loop counters
61 INTEGER i,j,k, kP1, kUp
62 _RL wOverride
63 _RS hFacWtmp
64 _RS hFacStmp
65 _RS hFacCtmp
66 _RS recip_hFacCtmp
67 _RL ab15,ab05
68 _RL slipSideFac
69 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
70
71 _RL Half
72 PARAMETER(Half=0.5D0)
73 CEOP
74
75 ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
76
77 iMin = 1
78 iMax = sNx
79 jMin = 1
80 jMax = sNy
81
82 C Adams-Bashforth timestepping weights
83 ab15 = 1.5 _d 0 + abeps
84 ab05 = -0.5 _d 0 - abeps
85
86 C Lateral friction (no-slip, free slip, or half slip):
87 IF ( no_slip_sides ) THEN
88 slipSideFac = -1. _d 0
89 ELSE
90 slipSideFac = 1. _d 0
91 ENDIF
92 CML half slip was used before ; keep the line for now, but half slip is
93 CML not used anywhere in the code as far as I can see.
94 C slipSideFac = 0. _d 0
95
96 DO bj=myByLo(myThid),myByHi(myThid)
97 DO bi=myBxLo(myThid),myBxHi(myThid)
98 DO K=1,Nr
99 DO j=1-OLy,sNy+OLy
100 DO i=1-OLx,sNx+OLx
101 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
102 gW(i,j,k,bi,bj) = 0.
103 ENDDO
104 ENDDO
105 ENDDO
106 ENDDO
107 ENDDO
108
109 C Catch barotropic mode
110 IF ( Nr .LT. 2 ) RETURN
111
112 C For each tile
113 DO bj=myByLo(myThid),myByHi(myThid)
114 DO bi=myBxLo(myThid),myBxHi(myThid)
115
116 C Boundaries condition at top
117 DO J=jMin,jMax
118 DO I=iMin,iMax
119 Flx_Dn(I,J,bi,bj)=0.
120 ENDDO
121 ENDDO
122
123 C Sweep down column
124 DO K=2,Nr
125 Kp1=K+1
126 wOverRide=1.
127 if (K.EQ.Nr) then
128 Kp1=Nr
129 wOverRide=0.
130 endif
131 C horizontal bi-harmonic dissipation
132 IF (momViscosity .AND. viscA4W.NE.0. ) THEN
133 C calculate the horizontal Laplacian of vertical flow
134 C Zonal flux d/dx W
135 DO j=1-Oly,sNy+Oly
136 fZon(1-Olx,j)=0.
137 DO i=1-Olx+1,sNx+Olx
138 fZon(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
139 & *_dyG(i,j,bi,bj)
140 & *_recip_dxC(i,j,bi,bj)
141 & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
142 #ifdef COSINEMETH_III
143 & *sqcosFacU(J,bi,bj)
144 #endif
145 ENDDO
146 ENDDO
147 C Meridional flux d/dy W
148 DO i=1-Olx,sNx+Olx
149 fMer(I,1-Oly)=0.
150 ENDDO
151 DO j=1-Oly+1,sNy+Oly
152 DO i=1-Olx,sNx+Olx
153 fMer(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
154 & *_dxG(i,j,bi,bj)
155 & *_recip_dyC(i,j,bi,bj)
156 & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
157 #ifdef ISOTROPIC_COS_SCALING
158 #ifdef COSINEMETH_III
159 & *sqCosFacV(j,bi,bj)
160 #endif
161 #endif
162 ENDDO
163 ENDDO
164
165 C del^2 W
166 C Difference of zonal fluxes ...
167 DO j=1-Oly,sNy+Oly
168 DO i=1-Olx,sNx+Olx-1
169 del2w(i,j)=fZon(i+1,j)-fZon(i,j)
170 ENDDO
171 del2w(sNx+Olx,j)=0.
172 ENDDO
173
174 C ... add difference of meridional fluxes and divide by volume
175 DO j=1-Oly,sNy+Oly-1
176 DO i=1-Olx,sNx+Olx
177 C First compute the fraction of open water for the w-control volume
178 C at the southern face
179 hFacCtmp=max(hFacC(I,J,K-1,bi,bj)-Half,0. _d 0)
180 & + min(hFacC(I,J,K ,bi,bj),Half)
181 recip_hFacCtmp = 0. _d 0
182 IF (hFacCtmp .GT. 0.) THEN
183 recip_hFacCtmp = 1./hFacCtmp
184 ELSE
185 recip_hFacCtmp = 0. _d 0
186 ENDIF
187 del2w(i,j)=recip_rA(i,j,bi,bj)
188 & *recip_drC(k)*recip_hFacCtmp
189 & *(
190 & del2w(i,j)
191 & +( fMer(i,j+1)-fMer(i,j) )
192 & )
193 ENDDO
194 ENDDO
195 C-- No-slip BCs impose a drag at walls...
196 CML ************************************************************
197 CML No-slip Boundary conditions for bi-harmonic dissipation
198 CML need to be implemented here!
199 CML ************************************************************
200 ENDIF
201
202 C Flux on Southern face
203 DO J=jMin,jMax+1
204 DO I=iMin,iMax
205 C First compute the fraction of open water for the w-control volume
206 C at the southern face
207 hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
208 & + min(hFacS(I,J,K ,bi,bj),Half)
209 tmp_VbarZ=Half*(
210 & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
211 & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
212 Flx_NS(I,J,bi,bj)=
213 & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
214 & -viscAhW*_recip_dyC(I,J,bi,bj)
215 & *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
216 & +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
217 & *wVel(I,J,K,bi,bj))
218 & +viscA4W*_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
219 #ifdef ISOTROPIC_COS_SCALING
220 #ifdef COSINEMETH_III
221 & *sqCosFacV(j,bi,bj)
222 #else
223 & *CosFacV(j,bi,bj)
224 #endif
225 #endif
226 C The last term is the weighted average of the viscous stress at the open
227 C fraction of the w control volume and at the closed fraction of the
228 C the control volume. A more compact but less intelligible version
229 C of the last three lines is:
230 CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
231 CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
232 ENDDO
233 ENDDO
234 C Flux on Western face
235 DO J=jMin,jMax
236 DO I=iMin,iMax+1
237 C First compute the fraction of open water for the w-control volume
238 C at the western face
239 hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
240 & + min(hFacW(I,J,K ,bi,bj),Half)
241 tmp_UbarZ=Half*(
242 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
243 & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
244 Flx_EW(I,J,bi,bj)=
245 & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
246 & -viscAhW*_recip_dxC(I,J,bi,bj)
247 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
248 & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
249 & *wVel(I,J,K,bi,bj) )
250 & +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
251 #ifdef COSINEMETH_III
252 & *sqCosFacU(j,bi,bj)
253 #else
254 & *CosFacU(j,bi,bj)
255 #endif
256 C The last term is the weighted average of the viscous stress at the open
257 C fraction of the w control volume and at the closed fraction of the
258 C the control volume. A more compact but less intelligible version
259 C of the last three lines is:
260 CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
261 CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
262 ENDDO
263 ENDDO
264 C Flux on Lower face
265 DO J=jMin,jMax
266 DO I=iMin,iMax
267 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
268 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
269 & +wOverRide*wVel(I,J,Kp1,bi,bj))
270 Flx_Dn(I,J,bi,bj)=
271 & tmp_WbarZ*tmp_WbarZ
272 & -viscAr*recip_drF(K)
273 & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
274 ENDDO
275 ENDDO
276 C Divergence of fluxes
277 DO J=jMin,jMax
278 DO I=iMin,iMax
279 gW(I,J,K,bi,bj) = 0.
280 & -(
281 & +_recip_dxF(I,J,bi,bj)*(
282 & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
283 & +_recip_dyF(I,J,bi,bj)*(
284 & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
285 & +recip_drC(K) *(
286 & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
287 & )
288 caja * recip_hFacU(I,J,K,bi,bj)
289 caja NOTE: This should be included
290 caja but we need an hFacUW (above U points)
291 caja and an hFacUS (above V points) too...
292 ENDDO
293 ENDDO
294 ENDDO
295 ENDDO
296 ENDDO
297
298
299 DO bj=myByLo(myThid),myByHi(myThid)
300 DO bi=myBxLo(myThid),myBxHi(myThid)
301 DO K=2,Nr
302 DO j=jMin,jMax
303 DO i=iMin,iMax
304 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
305 & +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)
306 & +ab05*gWNM1(i,j,k,bi,bj) )
307 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
308 ENDDO
309 ENDDO
310 ENDDO
311 ENDDO
312 ENDDO
313
314 #ifdef ALLOW_OBCS
315 IF (useOBCS) THEN
316 C-- This call is aesthetic: it makes the W field
317 C consistent with the OBs but this has no algorithmic
318 C impact. This is purely for diagnostic purposes.
319 DO bj=myByLo(myThid),myByHi(myThid)
320 DO bi=myBxLo(myThid),myBxHi(myThid)
321 DO K=1,Nr
322 CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
323 ENDDO
324 ENDDO
325 ENDDO
326 ENDIF
327 #endif /* ALLOW_OBCS */
328
329 #endif /* ALLOW_NONHYDROSTATIC */
330
331 RETURN
332 END

  ViewVC Help
Powered by ViewVC 1.1.22