/[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.35 - (show annotations) (download)
Tue Dec 5 05:25:08 2006 UTC (17 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58t_post, mitgcm_mapl_00, checkpoint58v_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.34: +26 -8 lines
start to implement deep-atmosphere and/or anelastic formulation

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.34 2006/07/13 21:32:38 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: CALC_GW
8 C !INTERFACE:
9 SUBROUTINE CALC_GW(
10 I bi, bj, KappaRU, KappaRV,
11 I myTime, myIter, myThid )
12 C !DESCRIPTION: \bv
13 C *==========================================================*
14 C | S/R CALC_GW
15 C | o Calculate vertical velocity tendency terms
16 C | ( Non-Hydrostatic only )
17 C *==========================================================*
18 C | In NH, 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 "EEPARAMS.h"
30 #include "PARAMS.h"
31 #include "GRID.h"
32 #include "SURFACE.h"
33 #include "DYNVARS.h"
34 #include "NH_VARS.h"
35
36 C !INPUT/OUTPUT PARAMETERS:
37 C == Routine arguments ==
38 C bi,bj :: current tile indices
39 C KappaRU :: vertical viscosity at U points
40 C KappaRV :: vertical viscosity at V points
41 C myTime :: Current time in simulation
42 C myIter :: Current iteration number in simulation
43 C myThid :: Thread number for this instance of the routine.
44 INTEGER bi,bj
45 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
46 _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
47 _RL myTime
48 INTEGER myIter
49 INTEGER myThid
50
51 #ifdef ALLOW_NONHYDROSTATIC
52
53 C !LOCAL VARIABLES:
54 C == Local variables ==
55 C iMin,iMax
56 C jMin,jMax
57 C xA :: W-Cell face area normal to X
58 C yA :: W-Cell face area normal to Y
59 C rThickC_W :: thickness (in r-units) of W-Cell at Western Edge
60 C rThickC_S :: thickness (in r-units) of W-Cell at Southern Edge
61 C rThickC_C :: thickness (in r-units) of W-Cell (centered on W pt)
62 C recip_rThickC :: reciprol thickness of W-Cell (centered on W-point)
63 C flx_NS :: vertical momentum flux, meridional direction
64 C flx_EW :: vertical momentum flux, zonal direction
65 C flxAdvUp :: vertical mom. advective flux, vertical direction (@ level k-1)
66 C flxDisUp :: vertical mom. dissipation flux, vertical direction (@ level k-1)
67 C flx_Dn :: vertical momentum flux, vertical direction (@ level k)
68 C gwDiss :: vertical momentum dissipation tendency
69 C i,j,k :: Loop counters
70 INTEGER iMin,iMax,jMin,jMax
71 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 _RL rThickC_W (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RL rThickC_S (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL rThickC_C (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL flxAdvUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL flxDisUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL gwDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL gwAdd (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL del2w (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 INTEGER i,j,k, kp1
86 _RL wOverride
87 _RL tmp_WbarZ
88 _RL uTrans, vTrans, rTrans
89 _RL viscLoc
90 _RL halfRL
91 _RS halfRS, zeroRS
92 PARAMETER( halfRL = 0.5D0 )
93 PARAMETER( halfRS = 0.5 , zeroRS = 0. )
94 CEOP
95
96 C Catch barotropic mode
97 IF ( Nr .LT. 2 ) RETURN
98
99 iMin = 1
100 iMax = sNx
101 jMin = 1
102 jMax = sNy
103
104 C-- Initialise gW to zero
105 DO k=1,Nr
106 DO j=1-OLy,sNy+OLy
107 DO i=1-OLx,sNx+OLx
108 gW(i,j,k,bi,bj) = 0.
109 ENDDO
110 ENDDO
111 ENDDO
112 C- Initialise gwDiss to zero
113 DO j=1-OLy,sNy+OLy
114 DO i=1-OLx,sNx+OLx
115 gwDiss(i,j) = 0.
116 ENDDO
117 ENDDO
118
119 C-- Boundaries condition at top
120 DO j=1-OLy,sNy+OLy
121 DO i=1-OLx,sNx+OLx
122 flxAdvUp(i,j) = 0.
123 flxDisUp(i,j) = 0.
124 ENDDO
125 ENDDO
126
127 C--- Sweep down column
128 DO k=2,Nr
129 kp1=k+1
130 wOverRide=1.
131 IF (k.EQ.Nr) THEN
132 kp1=Nr
133 wOverRide=0.
134 ENDIF
135 C-- Compute grid factor arround a W-point:
136 DO j=1-Oly,sNy+Oly
137 DO i=1-Olx,sNx+Olx
138 C- note: assume fluid @ smaller k than bottom: does not work in p-coordinate !
139 IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN
140 recip_rThickC(i,j) = 0.
141 ELSE
142 recip_rThickC(i,j) = 1. _d 0 /
143 & ( drF(k-1)*halfRS
144 & + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )
145 & )
146 ENDIF
147 c IF (momViscosity) THEN
148 #ifdef NONLIN_FRSURF
149 rThickC_C(i,j) =
150 & drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
151 & + drF( k )*MIN( h0FacC(i,j,k ,bi,bj), halfRS )
152 #else
153 rThickC_C(i,j) =
154 & drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
155 & + drF( k )*MIN( _hFacC(i,j,k ,bi,bj), halfRS )
156 #endif
157 rThickC_W(i,j) =
158 & drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )
159 & + drF( k )*MIN( _hFacW(i,j,k ,bi,bj), halfRS )
160 rThickC_S(i,j) =
161 & drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )
162 & + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )
163 C W-Cell Western face area:
164 xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
165 c & *deepFacF(k)
166 C W-Cell Southern face area:
167 yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
168 c & *deepFacF(k)
169 C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC.
170 C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted)
171 c ENDIF
172 ENDDO
173 ENDDO
174
175 C-- horizontal bi-harmonic dissipation
176 IF (momViscosity .AND. viscA4W.NE.0. ) THEN
177 C- calculate the horizontal Laplacian of vertical flow
178 C Zonal flux d/dx W
179 DO j=1-Oly,sNy+Oly
180 flx_EW(1-Olx,j)=0.
181 DO i=1-Olx+1,sNx+Olx
182 flx_EW(i,j) =
183 & (wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
184 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
185 #ifdef COSINEMETH_III
186 & *sqCosFacU(j,bi,bj)
187 #endif
188 ENDDO
189 ENDDO
190 C Meridional flux d/dy W
191 DO i=1-Olx,sNx+Olx
192 flx_NS(i,1-Oly)=0.
193 ENDDO
194 DO j=1-Oly+1,sNy+Oly
195 DO i=1-Olx,sNx+Olx
196 flx_NS(i,j) =
197 & (wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
198 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
199 #ifdef ISOTROPIC_COS_SCALING
200 #ifdef COSINEMETH_III
201 & *sqCosFacV(j,bi,bj)
202 #endif
203 #endif
204 ENDDO
205 ENDDO
206
207 C del^2 W
208 C Difference of zonal fluxes ...
209 DO j=1-Oly,sNy+Oly
210 DO i=1-Olx,sNx+Olx-1
211 del2w(i,j)=flx_EW(i+1,j)-flx_EW(i,j)
212 ENDDO
213 del2w(sNx+Olx,j)=0.
214 ENDDO
215
216 C ... add difference of meridional fluxes and divide by volume
217 DO j=1-Oly,sNy+Oly-1
218 DO i=1-Olx,sNx+Olx
219 del2w(i,j) = ( del2w(i,j)
220 & +(flx_NS(i,j+1)-flx_NS(i,j))
221 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
222 & *recip_deepFac2F(k)
223 ENDDO
224 ENDDO
225 C-- No-slip BCs impose a drag at walls...
226 CML ************************************************************
227 CML No-slip Boundary conditions for bi-harmonic dissipation
228 CML need to be implemented here!
229 CML ************************************************************
230 ELSE
231 C- Initialize del2w to zero:
232 DO j=1-Oly,sNy+Oly
233 DO i=1-Olx,sNx+Olx
234 del2w(i,j) = 0. _d 0
235 ENDDO
236 ENDDO
237 ENDIF
238
239 IF (momViscosity) THEN
240 C Viscous Flux on Western face
241 DO j=jMin,jMax
242 DO i=iMin,iMax+1
243 flx_EW(i,j)=
244 & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i-1,j,k,bi,bj))*halfRL
245 & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
246 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
247 cOld & *_recip_dxC(i,j,bi,bj)*rThickC_W(i,j)
248 & *cosFacU(j,bi,bj)
249 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL
250 & *(del2w(i,j)-del2w(i-1,j))
251 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
252 cOld & *_recip_dxC(i,j,bi,bj)*drC(k)
253 #ifdef COSINEMETH_III
254 & *sqCosFacU(j,bi,bj)
255 #else
256 & *cosFacU(j,bi,bj)
257 #endif
258 ENDDO
259 ENDDO
260 C Viscous Flux on Southern face
261 DO j=jMin,jMax+1
262 DO i=iMin,iMax
263 flx_NS(i,j)=
264 & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i,j-1,k,bi,bj))*halfRL
265 & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
266 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
267 cOld & *_recip_dyC(i,j,bi,bj)*rThickC_S(i,j)
268 #ifdef ISOTROPIC_COS_SCALING
269 & *cosFacV(j,bi,bj)
270 #endif
271 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL
272 & *(del2w(i,j)-del2w(i,j-1))
273 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
274 cOld & *_recip_dyC(i,j,bi,bj)*drC(k)
275 #ifdef ISOTROPIC_COS_SCALING
276 #ifdef COSINEMETH_III
277 & *sqCosFacV(j,bi,bj)
278 #else
279 & *cosFacV(j,bi,bj)
280 #endif
281 #endif
282 ENDDO
283 ENDDO
284 C Viscous Flux on Lower face of W-Cell (= at tracer-cell center, level k)
285 DO j=jMin,jMax
286 DO i=iMin,iMax
287 C Interpolate vert viscosity to center of tracer-cell (level k):
288 viscLoc = ( KappaRU(i,j,k) +KappaRU(i+1,j,k)
289 & +KappaRU(i,j,kp1)+KappaRU(i+1,j,kp1)
290 & +KappaRV(i,j,k) +KappaRV(i,j+1,k)
291 & +KappaRV(i,j,kp1)+KappaRV(i,j+1,kp1)
292 & )*0.125 _d 0
293 flx_Dn(i,j) =
294 & - viscLoc*( wVel(i,j,kp1,bi,bj)*wOverRide
295 & -wVel(i,j, k ,bi,bj) )*rkSign
296 & *recip_drF(k)*rA(i,j,bi,bj)
297 & *deepFac2C(k)*rhoFacC(k)
298 cOld & *recip_drF(k)
299 ENDDO
300 ENDDO
301 C Tendency is minus divergence of viscous fluxes:
302 C anelastic: vert.visc.flx is scaled by rhoFac but hor.visc.fluxes are not
303 DO j=jMin,jMax
304 DO i=iMin,iMax
305 gwDiss(i,j) =
306 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
307 & + ( flx_NS(i,j+1)-flx_NS(i,j) )
308 & + ( flx_Dn(i,j)-flxDisUp(i,j) )*rkSign
309 & *recip_rhoFacF(k)
310 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
311 & *recip_deepFac2F(k)
312 cOld gwDiss(i,j) =
313 cOld & -(
314 cOld & +_recip_dxF(i,j,bi,bj)*( flx_EW(i+1,j)-flx_EW(i,j) )
315 cOld & +_recip_dyF(i,j,bi,bj)*( flx_NS(i,j+1)-flx_NS(i,j) )
316 cOld & + ( flxDisUp(i,j)-flx_Dn(i,j) )
317 c & )*recip_rThickC(i,j)
318 cOld & )*recip_drC(k)
319 C-- prepare for next level (k+1)
320 flxDisUp(i,j)=flx_Dn(i,j)
321 ENDDO
322 ENDDO
323 ENDIF
324
325 IF (no_slip_sides) THEN
326 C- No-slip BCs impose a drag at walls...
327 CALL MOM_W_SIDEDRAG(
328 I bi,bj,k,
329 I wVel, del2w,
330 I rThickC_C, recip_rThickC,
331 I viscAh_W, viscA4_W,
332 O gwAdd,
333 I myThid )
334 DO j=jMin,jMax
335 DO i=iMin,iMax
336 gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)
337 ENDDO
338 ENDDO
339 ENDIF
340
341 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
342
343 IF ( momAdvection ) THEN
344 C Advective Flux on Western face
345 DO j=jMin,jMax
346 DO i=iMin,iMax+1
347 C transport through Western face area:
348 uTrans = (
349 & drF(k-1)*_hFacW(i,j,k-1,bi,bj)*uVel(i,j,k-1,bi,bj)
350 & *rhoFacC(k-1)
351 & + drF( k )*_hFacW(i,j, k ,bi,bj)*uVel(i,j, k ,bi,bj)
352 & *rhoFacC(k)
353 & )*halfRL*_dyG(i,j,bi,bj)*deepFacF(k)
354 cOld & )*halfRL
355 flx_EW(i,j)=
356 & uTrans*(wVel(i,j,k,bi,bj)+wVel(i-1,j,k,bi,bj))*halfRL
357 ENDDO
358 ENDDO
359 C Advective Flux on Southern face
360 DO j=jMin,jMax+1
361 DO i=iMin,iMax
362 C transport through Southern face area:
363 vTrans = (
364 & drF(k-1)*_hFacS(i,j,k-1,bi,bj)*vVel(i,j,k-1,bi,bj)
365 & *rhoFacC(k-1)
366 & +drF( k )*_hFacS(i,j, k ,bi,bj)*vVel(i,j, k ,bi,bj)
367 & *rhoFacC(k)
368 & )*halfRL*_dxG(i,j,bi,bj)*deepFacF(k)
369 cOld & )*halfRL
370 flx_NS(i,j)=
371 & vTrans*(wVel(i,j,k,bi,bj)+wVel(i,j-1,k,bi,bj))*halfRL
372 ENDDO
373 ENDDO
374 C Advective Flux on Lower face of W-Cell (= at tracer-cell center, level k)
375 DO j=jMin,jMax
376 DO i=iMin,iMax
377 tmp_WbarZ = halfRL*( wVel(i,j, k ,bi,bj)
378 & +wVel(i,j,kp1,bi,bj)*wOverRide )
379 C transport through Lower face area:
380 rTrans = halfRL*
381 & ( wVel(i,j, k ,bi,bj)*deepFac2F( k )*rhoFacF( k )
382 & +wVel(i,j,kp1,bi,bj)*deepFac2F(kp1)*rhoFacF(kp1)
383 & *wOverRide
384 & )*rA(i,j,bi,bj)
385 flx_Dn(i,j) = rTrans*tmp_WbarZ
386 cOld flx_Dn(i,j) = tmp_WbarZ*tmp_WbarZ
387 ENDDO
388 ENDDO
389 C Tendency is minus divergence of advective fluxes:
390 C anelastic: all transports & advect. fluxes are scaled by rhoFac
391 DO j=jMin,jMax
392 DO i=iMin,iMax
393 gW(i,j,k,bi,bj) =
394 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
395 & + ( flx_NS(i,j+1)-flx_NS(i,j) )
396 & + ( flx_Dn(i,j)-flxAdvUp(i,j) )*rkSign
397 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
398 & *recip_deepFac2F(k)*recip_rhoFacF(k)
399 cOld gW(i,j,k,bi,bj) =
400 cOld & -(
401 cOld & +_recip_dxF(i,j,bi,bj)*( flx_EW(i+1,j)-flx_EW(i,j) )
402 cOld & +_recip_dyF(i,j,bi,bj)*( flx_NS(i,j+1)-flx_NS(i,j) )
403 cOld & + ( flxAdvUp(i,j)-flx_Dn(i,j) )
404 c & )*recip_rThickC(i,j)
405 cOld & )*recip_drC(k)
406 C-- prepare for next level (k+1)
407 flxAdvUp(i,j)=flx_Dn(i,j)
408 ENDDO
409 ENDDO
410 ENDIF
411
412 IF ( useNHMTerms ) THEN
413 CALL MOM_W_METRIC_NH(
414 I bi,bj,k,
415 I uVel, vVel,
416 O gwAdd,
417 I myThid )
418 DO j=jMin,jMax
419 DO i=iMin,iMax
420 gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
421 ENDDO
422 ENDDO
423 ENDIF
424 IF ( use3dCoriolis ) THEN
425 CALL MOM_W_CORIOLIS_NH(
426 I bi,bj,k,
427 I uVel, vVel,
428 O gwAdd,
429 I myThid )
430 DO j=jMin,jMax
431 DO i=iMin,iMax
432 gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
433 ENDDO
434 ENDDO
435 ENDIF
436
437 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
438
439 C-- Dissipation term inside the Adams-Bashforth:
440 IF ( momViscosity .AND. momDissip_In_AB) THEN
441 DO j=jMin,jMax
442 DO i=iMin,iMax
443 gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
444 ENDDO
445 ENDDO
446 ENDIF
447
448 C- Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights)
449 C and save gW_[n] into gwNm1 for the next time step.
450 c#ifdef ALLOW_ADAMSBASHFORTH_3
451 c CALL ADAMS_BASHFORTH3(
452 c I bi, bj, k,
453 c U gW, gwNm,
454 c I momStartAB, myIter, myThid )
455 c#else /* ALLOW_ADAMSBASHFORTH_3 */
456 CALL ADAMS_BASHFORTH2(
457 I bi, bj, k,
458 U gW, gwNm1,
459 I myIter, myThid )
460 c#endif /* ALLOW_ADAMSBASHFORTH_3 */
461
462 C-- Dissipation term outside the Adams-Bashforth:
463 IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN
464 DO j=jMin,jMax
465 DO i=iMin,iMax
466 gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
467 ENDDO
468 ENDDO
469 ENDIF
470
471 C- end of the k loop
472 ENDDO
473
474 #endif /* ALLOW_NONHYDROSTATIC */
475
476 RETURN
477 END

  ViewVC Help
Powered by ViewVC 1.1.22