/[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.39 - (show annotations) (download)
Tue Apr 22 22:20:31 2008 UTC (16 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59r, checkpoint61f, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61h, checkpoint61i
Changes since 1.38: +24 -2 lines
add diagnostics for vertical momentum tendencies

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

  ViewVC Help
Powered by ViewVC 1.1.22