/[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.36 - (show annotations) (download)
Mon Mar 12 23:48:24 2007 UTC (17 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint58w_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint58x_post
Changes since 1.35: +83 -11 lines
- use reference profile for w <-> omega conversion in NH-code
- change grid factor arround a W-point to work in p-coord. (keep
  the previous version within #ifndef CALC_GW_NEW_THICK for now)

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

  ViewVC Help
Powered by ViewVC 1.1.22