1 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.49 2012/08/06 16:54:40 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 |
#ifdef ALLOW_ADDFLUID |
39 |
# include "FFIELDS.h" |
40 |
#endif |
41 |
|
42 |
C !INPUT/OUTPUT PARAMETERS: |
43 |
C == Routine arguments == |
44 |
C bi,bj :: current tile indices |
45 |
C KappaRU :: vertical viscosity at U points |
46 |
C KappaRV :: vertical viscosity at V points |
47 |
C myTime :: Current time in simulation |
48 |
C myIter :: Current iteration number in simulation |
49 |
C myThid :: Thread number for this instance of the routine. |
50 |
INTEGER bi,bj |
51 |
_RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
52 |
_RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
53 |
_RL myTime |
54 |
INTEGER myIter |
55 |
INTEGER myThid |
56 |
|
57 |
#ifdef ALLOW_NONHYDROSTATIC |
58 |
|
59 |
C !LOCAL VARIABLES: |
60 |
C == Local variables == |
61 |
C iMin,iMax |
62 |
C jMin,jMax |
63 |
C xA :: W-Cell face area normal to X |
64 |
C yA :: W-Cell face area normal to Y |
65 |
C rThickC_W :: thickness (in r-units) of W-Cell at Western Edge |
66 |
C rThickC_S :: thickness (in r-units) of W-Cell at Southern Edge |
67 |
C rThickC_C :: thickness (in r-units) of W-Cell (centered on W pt) |
68 |
C recip_rThickC :: reciprol thickness of W-Cell (centered on W-point) |
69 |
C flx_NS :: vertical momentum flux, meridional direction |
70 |
C flx_EW :: vertical momentum flux, zonal direction |
71 |
C flxAdvUp :: vertical mom. advective flux, vertical direction (@ level k-1) |
72 |
C flxDisUp :: vertical mom. dissipation flux, vertical direction (@ level k-1) |
73 |
C flx_Dn :: vertical momentum flux, vertical direction (@ level k) |
74 |
C gwDiss :: vertical momentum dissipation tendency |
75 |
C gwAdd :: other tendencies (Coriolis, Metric-terms) |
76 |
C gw_AB :: tendency increment from Adams-Bashforth |
77 |
C del2w :: laplacian of wVel |
78 |
C wFld :: local copy of wVel |
79 |
C i,j,k :: Loop counters |
80 |
INTEGER iMin,iMax,jMin,jMax |
81 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
82 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
83 |
_RL rThickC_W (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
84 |
_RL rThickC_S (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
85 |
_RL rThickC_C (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
86 |
_RL recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
87 |
_RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
88 |
_RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
89 |
_RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
90 |
_RL flxAdvUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
_RL flxDisUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
_RL gwDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
93 |
_RL gwAdd (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
94 |
_RL gw_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
95 |
_RL del2w (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
96 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
97 |
INTEGER i,j,k, km1, kp1 |
98 |
_RL mskM1, mskP1 |
99 |
_RL tmp_WbarZ |
100 |
_RL uTrans, vTrans, rTrans |
101 |
_RL viscLoc |
102 |
PARAMETER( iMin = 1 , iMax = sNx ) |
103 |
PARAMETER( jMin = 1 , jMax = sNy ) |
104 |
CEOP |
105 |
#ifdef ALLOW_DIAGNOSTICS |
106 |
LOGICAL diagDiss, diagAdvec, diag_AB |
107 |
LOGICAL DIAGNOSTICS_IS_ON |
108 |
EXTERNAL DIAGNOSTICS_IS_ON |
109 |
#endif /* ALLOW_DIAGNOSTICS */ |
110 |
|
111 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
112 |
|
113 |
#ifdef ALLOW_DIAGNOSTICS |
114 |
IF ( useDiagnostics ) THEN |
115 |
diagDiss = DIAGNOSTICS_IS_ON( 'Wm_Diss ', myThid ) |
116 |
diagAdvec = DIAGNOSTICS_IS_ON( 'Wm_Advec', myThid ) |
117 |
diag_AB = DIAGNOSTICS_IS_ON( 'AB_gW ', myThid ) |
118 |
ELSE |
119 |
diagDiss = .FALSE. |
120 |
diagAdvec = .FALSE. |
121 |
diag_AB = .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 |
IF (momViscosity) THEN |
140 |
C- Initialize del2w to zero: |
141 |
DO j=1-OLy,sNy+OLy |
142 |
DO i=1-OLx,sNx+OLx |
143 |
del2w(i,j) = 0. _d 0 |
144 |
ENDDO |
145 |
ENDDO |
146 |
ENDIF |
147 |
|
148 |
C-- Boundaries condition at top (vertical advection of vertical momentum): |
149 |
DO j=1-OLy,sNy+OLy |
150 |
DO i=1-OLx,sNx+OLx |
151 |
flxAdvUp(i,j) = 0. |
152 |
c flxDisUp(i,j) = 0. |
153 |
ENDDO |
154 |
ENDDO |
155 |
|
156 |
|
157 |
C--- Sweep down column |
158 |
DO k=1,Nr |
159 |
km1 = MAX( k-1, 1 ) |
160 |
kp1 = MIN( k+1,Nr ) |
161 |
mskM1 = 1. |
162 |
mskP1 = 1. |
163 |
IF ( k.EQ. 1 ) mskM1 = 0. |
164 |
IF ( k.EQ.Nr ) mskP1 = 0. |
165 |
IF ( k.GT.1 ) THEN |
166 |
C-- Compute grid factor arround a W-point: |
167 |
#ifdef CALC_GW_NEW_THICK |
168 |
DO j=1-OLy,sNy+OLy |
169 |
DO i=1-OLx,sNx+OLx |
170 |
IF ( maskC(i,j,k-1,bi,bj).EQ.0. .OR. |
171 |
& maskC(i,j, k ,bi,bj).EQ.0. ) THEN |
172 |
recip_rThickC(i,j) = 0. |
173 |
ELSE |
174 |
C- valid in z & p coord.; also accurate if Interface @ middle between 2 centers |
175 |
recip_rThickC(i,j) = 1. _d 0 / |
176 |
& ( MIN( Ro_surf(i,j,bi,bj),rC(k-1) ) |
177 |
& - MAX( R_low(i,j,bi,bj), rC(k) ) |
178 |
& ) |
179 |
ENDIF |
180 |
ENDDO |
181 |
ENDDO |
182 |
IF (momViscosity) THEN |
183 |
DO j=1-OLy,sNy+OLy |
184 |
DO i=1-OLx,sNx+OLx |
185 |
rThickC_C(i,j) = MAX( zeroRS, |
186 |
& MIN( Ro_surf(i,j,bi,bj), rC(k-1) ) |
187 |
& -MAX( R_low(i,j,bi,bj), rC(k) ) |
188 |
& ) |
189 |
ENDDO |
190 |
ENDDO |
191 |
DO j=1-OLy,sNy+OLy |
192 |
DO i=1-OLx+1,sNx+OLx |
193 |
rThickC_W(i,j) = MAX( zeroRS, |
194 |
& MIN( rSurfW(i,j,bi,bj), rC(k-1) ) |
195 |
& -MAX( rLowW(i,j,bi,bj), rC(k) ) |
196 |
& ) |
197 |
C W-Cell Western face area: |
198 |
xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j) |
199 |
c & *deepFacF(k) |
200 |
ENDDO |
201 |
ENDDO |
202 |
DO j=1-OLy+1,sNy+OLy |
203 |
DO i=1-OLx,sNx+OLx |
204 |
rThickC_S(i,j) = MAX( zeroRS, |
205 |
& MIN( rSurfS(i,j,bi,bj), rC(k-1) ) |
206 |
& -MAX( rLowS(i,j,bi,bj), rC(k) ) |
207 |
& ) |
208 |
C W-Cell Southern face area: |
209 |
yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j) |
210 |
c & *deepFacF(k) |
211 |
C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC. |
212 |
C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted) |
213 |
ENDDO |
214 |
ENDDO |
215 |
ENDIF |
216 |
#else /* CALC_GW_NEW_THICK */ |
217 |
DO j=1-OLy,sNy+OLy |
218 |
DO i=1-OLx,sNx+OLx |
219 |
C- note: assume fluid @ smaller k than bottom: does not work in p-coordinate ! |
220 |
IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN |
221 |
recip_rThickC(i,j) = 0. |
222 |
ELSE |
223 |
recip_rThickC(i,j) = 1. _d 0 / |
224 |
& ( drF(k-1)*halfRS |
225 |
& + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS ) |
226 |
& ) |
227 |
ENDIF |
228 |
c IF (momViscosity) THEN |
229 |
#ifdef NONLIN_FRSURF |
230 |
rThickC_C(i,j) = |
231 |
& drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS ) |
232 |
& + drF( k )*MIN( h0FacC(i,j,k ,bi,bj), halfRS ) |
233 |
#else |
234 |
rThickC_C(i,j) = |
235 |
& drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS ) |
236 |
& + drF( k )*MIN( _hFacC(i,j,k ,bi,bj), halfRS ) |
237 |
#endif |
238 |
rThickC_W(i,j) = |
239 |
& drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS ) |
240 |
& + drF( k )*MIN( _hFacW(i,j,k ,bi,bj), halfRS ) |
241 |
rThickC_S(i,j) = |
242 |
& drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS ) |
243 |
& + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS ) |
244 |
C W-Cell Western face area: |
245 |
xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j) |
246 |
c & *deepFacF(k) |
247 |
C W-Cell Southern face area: |
248 |
yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j) |
249 |
c & *deepFacF(k) |
250 |
C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC. |
251 |
C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted) |
252 |
c ENDIF |
253 |
ENDDO |
254 |
ENDDO |
255 |
#endif /* CALC_GW_NEW_THICK */ |
256 |
ELSEIF ( selectNHfreeSurf.GE.1 ) THEN |
257 |
DO j=1-OLy,sNy+OLy |
258 |
DO i=1-OLx,sNx+OLx |
259 |
recip_rThickC(i,j) = recip_drC(k) |
260 |
c rThickC_C(i,j) = drC(k) |
261 |
c rThickC_W(i,j) = drC(k) |
262 |
c rThickC_S(i,j) = drC(k) |
263 |
c xA(i,j) = _dyG(i,j,bi,bj)*drC(k) |
264 |
c yA(i,j) = _dxG(i,j,bi,bj)*drC(k) |
265 |
ENDDO |
266 |
ENDDO |
267 |
ENDIF |
268 |
|
269 |
C-- local copy of wVel: |
270 |
DO j=1-OLy,sNy+OLy |
271 |
DO i=1-OLx,sNx+OLx |
272 |
wFld(i,j) = wVel(i,j,k,bi,bj) |
273 |
ENDDO |
274 |
ENDDO |
275 |
|
276 |
C-- horizontal bi-harmonic dissipation |
277 |
IF ( momViscosity .AND. k.GT.1 .AND. viscA4W.NE.0. ) THEN |
278 |
|
279 |
C- calculate the horizontal Laplacian of vertical flow |
280 |
C Zonal flux d/dx W |
281 |
IF ( useCubedSphereExchange ) THEN |
282 |
C to compute d/dx(W), fill corners with appropriate values: |
283 |
CALL FILL_CS_CORNER_TR_RL( 1, .FALSE., |
284 |
& wFld, bi,bj, myThid ) |
285 |
ENDIF |
286 |
DO j=1-OLy,sNy+OLy |
287 |
flx_EW(1-OLx,j)=0. |
288 |
DO i=1-OLx+1,sNx+OLx |
289 |
flx_EW(i,j) = |
290 |
& ( wFld(i,j) - wFld(i-1,j) ) |
291 |
& *_recip_dxC(i,j,bi,bj)*xA(i,j) |
292 |
#ifdef COSINEMETH_III |
293 |
& *sqCosFacU(j,bi,bj) |
294 |
#endif |
295 |
#ifdef ALLOW_OBCS |
296 |
& *maskInW(i,j,bi,bj) |
297 |
#endif |
298 |
ENDDO |
299 |
ENDDO |
300 |
|
301 |
C Meridional flux d/dy W |
302 |
IF ( useCubedSphereExchange ) THEN |
303 |
C to compute d/dy(W), fill corners with appropriate values: |
304 |
CALL FILL_CS_CORNER_TR_RL( 2, .FALSE., |
305 |
& wFld, bi,bj, myThid ) |
306 |
ENDIF |
307 |
DO i=1-OLx,sNx+OLx |
308 |
flx_NS(i,1-OLy)=0. |
309 |
ENDDO |
310 |
DO j=1-OLy+1,sNy+OLy |
311 |
DO i=1-OLx,sNx+OLx |
312 |
flx_NS(i,j) = |
313 |
& ( wFld(i,j) - wFld(i,j-1) ) |
314 |
& *_recip_dyC(i,j,bi,bj)*yA(i,j) |
315 |
#ifdef ISOTROPIC_COS_SCALING |
316 |
#ifdef COSINEMETH_III |
317 |
& *sqCosFacV(j,bi,bj) |
318 |
#endif |
319 |
#endif |
320 |
#ifdef ALLOW_OBCS |
321 |
& *maskInS(i,j,bi,bj) |
322 |
#endif |
323 |
ENDDO |
324 |
ENDDO |
325 |
|
326 |
C del^2 W |
327 |
C Divergence of horizontal fluxes |
328 |
DO j=1-OLy,sNy+OLy-1 |
329 |
DO i=1-OLx,sNx+OLx-1 |
330 |
del2w(i,j) = ( ( flx_EW(i+1,j)-flx_EW(i,j) ) |
331 |
& +( flx_NS(i,j+1)-flx_NS(i,j) ) |
332 |
& )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j) |
333 |
& *recip_deepFac2F(k) |
334 |
ENDDO |
335 |
ENDDO |
336 |
C end if biharmonic viscosity |
337 |
ENDIF |
338 |
|
339 |
IF ( momViscosity .AND. k.GT.1 ) THEN |
340 |
C Viscous Flux on Western face |
341 |
DO j=jMin,jMax |
342 |
DO i=iMin,iMax+1 |
343 |
flx_EW(i,j)= |
344 |
& - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i-1,j,k,bi,bj))*halfRL |
345 |
& *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj)) |
346 |
& *_recip_dxC(i,j,bi,bj)*xA(i,j) |
347 |
& *cosFacU(j,bi,bj) |
348 |
& + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL |
349 |
& *(del2w(i,j)-del2w(i-1,j)) |
350 |
& *_recip_dxC(i,j,bi,bj)*xA(i,j) |
351 |
#ifdef COSINEMETH_III |
352 |
& *sqCosFacU(j,bi,bj) |
353 |
#else |
354 |
& *cosFacU(j,bi,bj) |
355 |
#endif |
356 |
ENDDO |
357 |
ENDDO |
358 |
C Viscous Flux on Southern face |
359 |
DO j=jMin,jMax+1 |
360 |
DO i=iMin,iMax |
361 |
flx_NS(i,j)= |
362 |
& - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i,j-1,k,bi,bj))*halfRL |
363 |
& *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj)) |
364 |
& *_recip_dyC(i,j,bi,bj)*yA(i,j) |
365 |
#ifdef ISOTROPIC_COS_SCALING |
366 |
& *cosFacV(j,bi,bj) |
367 |
#endif |
368 |
& + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL |
369 |
& *(del2w(i,j)-del2w(i,j-1)) |
370 |
& *_recip_dyC(i,j,bi,bj)*yA(i,j) |
371 |
#ifdef ISOTROPIC_COS_SCALING |
372 |
#ifdef COSINEMETH_III |
373 |
& *sqCosFacV(j,bi,bj) |
374 |
#else |
375 |
& *cosFacV(j,bi,bj) |
376 |
#endif |
377 |
#endif |
378 |
ENDDO |
379 |
ENDDO |
380 |
C Viscous Flux on Lower face of W-Cell (= at tracer-cell center, level k) |
381 |
DO j=jMin,jMax |
382 |
DO i=iMin,iMax |
383 |
C Interpolate vert viscosity to center of tracer-cell (level k): |
384 |
viscLoc = ( KappaRU(i,j,k) +KappaRU(i+1,j,k) |
385 |
& +KappaRU(i,j,kp1)+KappaRU(i+1,j,kp1) |
386 |
& +KappaRV(i,j,k) +KappaRV(i,j+1,k) |
387 |
& +KappaRV(i,j,kp1)+KappaRV(i,j+1,kp1) |
388 |
& )*0.125 _d 0 |
389 |
flx_Dn(i,j) = |
390 |
& - viscLoc*( wVel(i,j,kp1,bi,bj)*mskP1 |
391 |
& -wVel(i,j, k ,bi,bj) )*rkSign |
392 |
& *recip_drF(k)*rA(i,j,bi,bj) |
393 |
& *deepFac2C(k)*rhoFacC(k) |
394 |
ENDDO |
395 |
ENDDO |
396 |
IF ( k.EQ.2 ) THEN |
397 |
C Viscous Flux on Upper face of W-Cell (= at tracer-cell center, level k-1) |
398 |
DO j=jMin,jMax |
399 |
DO i=iMin,iMax |
400 |
C Interpolate horizontally (but not vertically) vert viscosity to center: |
401 |
C Although background visc. might be defined at k=1, this is not |
402 |
C generally true when using variable visc. (from vertical mixing scheme). |
403 |
C Therefore, no vert. interp. and only horizontal interpolation. |
404 |
viscLoc = ( KappaRU(i,j,k) + KappaRU(i+1,j,k) |
405 |
& +KappaRV(i,j,k) + KappaRV(i,j+1,k) |
406 |
& )*0.25 _d 0 |
407 |
flxDisUp(i,j) = |
408 |
& - viscLoc*( wVel(i,j, k ,bi,bj) |
409 |
& -wVel(i,j,k-1,bi,bj) )*rkSign |
410 |
& *recip_drF(k-1)*rA(i,j,bi,bj) |
411 |
& *deepFac2C(k-1)*rhoFacC(k-1) |
412 |
C to recover old (before 2009/11/30) results (since flxDisUp(k=2) was zero) |
413 |
c flxDisUp(i,j) = 0. |
414 |
ENDDO |
415 |
ENDDO |
416 |
ENDIF |
417 |
C Tendency is minus divergence of viscous fluxes: |
418 |
C anelastic: vert.visc.flx is scaled by rhoFac but hor.visc.fluxes are not |
419 |
DO j=jMin,jMax |
420 |
DO i=iMin,iMax |
421 |
gwDiss(i,j) = |
422 |
& -( ( flx_EW(i+1,j)-flx_EW(i,j) ) |
423 |
& + ( flx_NS(i,j+1)-flx_NS(i,j) ) |
424 |
& + ( flx_Dn(i,j)-flxDisUp(i,j) )*rkSign |
425 |
& *recip_rhoFacF(k) |
426 |
& )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j) |
427 |
& *recip_deepFac2F(k) |
428 |
C-- prepare for next level (k+1) |
429 |
flxDisUp(i,j)=flx_Dn(i,j) |
430 |
ENDDO |
431 |
ENDDO |
432 |
ENDIF |
433 |
|
434 |
IF ( momViscosity .AND. k.GT.1 .AND. no_slip_sides ) THEN |
435 |
C- No-slip BCs impose a drag at walls... |
436 |
CALL MOM_W_SIDEDRAG( |
437 |
I bi,bj,k, |
438 |
I wVel, del2w, |
439 |
I rThickC_C, recip_rThickC, |
440 |
I viscAh_W, viscA4_W, |
441 |
O gwAdd, |
442 |
I myThid ) |
443 |
DO j=jMin,jMax |
444 |
DO i=iMin,iMax |
445 |
gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j) |
446 |
ENDDO |
447 |
ENDDO |
448 |
ENDIF |
449 |
|
450 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
451 |
|
452 |
IF ( momAdvection ) THEN |
453 |
|
454 |
IF ( k.GT.1 .OR. selectNHfreeSurf.GE.1 ) THEN |
455 |
C Advective Flux on Western face |
456 |
DO j=jMin,jMax |
457 |
DO i=iMin,iMax+1 |
458 |
C transport through Western face area: |
459 |
uTrans = ( |
460 |
& drF(km1)*_hFacW(i,j,km1,bi,bj)*uVel(i,j,km1,bi,bj) |
461 |
& *rhoFacC(km1)*mskM1 |
462 |
& + drF( k )*_hFacW(i,j, k ,bi,bj)*uVel(i,j, k ,bi,bj) |
463 |
& *rhoFacC(k) |
464 |
& )*halfRL*_dyG(i,j,bi,bj)*deepFacF(k) |
465 |
flx_EW(i,j) = uTrans*(wFld(i,j)+wFld(i-1,j))*halfRL |
466 |
c flx_EW(i,j)= |
467 |
c & uTrans*(wVel(i,j,k,bi,bj)+wVel(i-1,j,k,bi,bj))*halfRL |
468 |
ENDDO |
469 |
ENDDO |
470 |
C Advective Flux on Southern face |
471 |
DO j=jMin,jMax+1 |
472 |
DO i=iMin,iMax |
473 |
C transport through Southern face area: |
474 |
vTrans = ( |
475 |
& drF(km1)*_hFacS(i,j,km1,bi,bj)*vVel(i,j,km1,bi,bj) |
476 |
& *rhoFacC(km1)*mskM1 |
477 |
& +drF( k )*_hFacS(i,j, k ,bi,bj)*vVel(i,j, k ,bi,bj) |
478 |
& *rhoFacC(k) |
479 |
& )*halfRL*_dxG(i,j,bi,bj)*deepFacF(k) |
480 |
flx_NS(i,j) = vTrans*(wFld(i,j)+wFld(i,j-1))*halfRL |
481 |
c flx_NS(i,j)= |
482 |
c & vTrans*(wVel(i,j,k,bi,bj)+wVel(i,j-1,k,bi,bj))*halfRL |
483 |
ENDDO |
484 |
ENDDO |
485 |
ENDIF |
486 |
C Advective Flux on Lower face of W-Cell (= at tracer-cell center, level k) |
487 |
c IF (.TRUE.) THEN |
488 |
DO j=jMin,jMax |
489 |
DO i=iMin,iMax |
490 |
C NH in p-coord.: advect wSpeed [m/s] with rTrans |
491 |
tmp_WbarZ = halfRL* |
492 |
& ( wVel(i,j, k ,bi,bj)*rVel2wUnit( k ) |
493 |
& +wVel(i,j,kp1,bi,bj)*rVel2wUnit(kp1)*mskP1 ) |
494 |
C transport through Lower face area: |
495 |
rTrans = halfRL* |
496 |
& ( wVel(i,j, k ,bi,bj)*deepFac2F( k )*rhoFacF( k ) |
497 |
& +wVel(i,j,kp1,bi,bj)*deepFac2F(kp1)*rhoFacF(kp1) |
498 |
& *mskP1 |
499 |
& )*rA(i,j,bi,bj) |
500 |
flx_Dn(i,j) = rTrans*tmp_WbarZ |
501 |
ENDDO |
502 |
ENDDO |
503 |
c ENDIF |
504 |
IF ( k.EQ.1 .AND. selectNHfreeSurf.GE.1 ) THEN |
505 |
C Advective Flux on Upper face of W-Cell (= at surface) |
506 |
DO j=jMin,jMax |
507 |
DO i=iMin,iMax |
508 |
tmp_WbarZ = wVel(i,j,k,bi,bj)*rVel2wUnit(k) |
509 |
rTrans = wVel(i,j,k,bi,bj)*deepFac2F(k)*rhoFacF(k) |
510 |
& *rA(i,j,bi,bj) |
511 |
flxAdvUp(i,j) = rTrans*tmp_WbarZ |
512 |
c flxAdvUp(i,j) = 0. |
513 |
ENDDO |
514 |
ENDDO |
515 |
ENDIF |
516 |
|
517 |
IF ( k.GT.1 .OR. selectNHfreeSurf.GE.1 ) THEN |
518 |
C Tendency is minus divergence of advective fluxes: |
519 |
C anelastic: all transports & advect. fluxes are scaled by rhoFac |
520 |
DO j=jMin,jMax |
521 |
DO i=iMin,iMax |
522 |
C to recover old (before 2009/11/30) results (since flxAdvUp(k=2) was zero) |
523 |
c IF (k.EQ.2) flxAdvUp(i,j) = 0. |
524 |
gW(i,j,k,bi,bj) = |
525 |
& -( ( flx_EW(i+1,j)-flx_EW(i,j) ) |
526 |
& + ( flx_NS(i,j+1)-flx_NS(i,j) ) |
527 |
& + ( flx_Dn(i,j)-flxAdvUp(i,j) )*rkSign*wUnit2rVel(k) |
528 |
& )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j) |
529 |
& *recip_deepFac2F(k)*recip_rhoFacF(k) |
530 |
ENDDO |
531 |
ENDDO |
532 |
#ifdef ALLOW_ADDFLUID |
533 |
IF ( selectAddFluid.GE.1 ) THEN |
534 |
DO j=jMin,jMax |
535 |
DO i=iMin,iMax |
536 |
gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj) |
537 |
& + wVel(i,j,k,bi,bj)*mass2rUnit*0.5 _d 0 |
538 |
& *( addMass(i,j,k,bi,bj) |
539 |
& +addMass(i,j,km1,bi,bj)*mskM1 ) |
540 |
& *recip_rA(i,j,bi,bj)*recip_rThickC(i,j) |
541 |
& *recip_deepFac2F(k)*recip_rhoFacF(k) |
542 |
ENDDO |
543 |
ENDDO |
544 |
ENDIF |
545 |
#endif /* ALLOW_ADDFLUID */ |
546 |
ENDIF |
547 |
|
548 |
DO j=jMin,jMax |
549 |
DO i=iMin,iMax |
550 |
C-- prepare for next level (k+1) |
551 |
flxAdvUp(i,j)=flx_Dn(i,j) |
552 |
ENDDO |
553 |
ENDDO |
554 |
|
555 |
c ELSE |
556 |
C- if momAdvection / else |
557 |
c DO j=1-OLy,sNy+OLy |
558 |
c DO i=1-OLx,sNx+OLx |
559 |
c gW(i,j,k,bi,bj) = 0. _d 0 |
560 |
c ENDDO |
561 |
c ENDDO |
562 |
|
563 |
C- endif momAdvection. |
564 |
ENDIF |
565 |
|
566 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
567 |
|
568 |
IF ( useNHMTerms .AND. k.GT.1 ) THEN |
569 |
CALL MOM_W_METRIC_NH( |
570 |
I bi,bj,k, |
571 |
I uVel, vVel, |
572 |
O gwAdd, |
573 |
I myThid ) |
574 |
DO j=jMin,jMax |
575 |
DO i=iMin,iMax |
576 |
gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j) |
577 |
ENDDO |
578 |
ENDDO |
579 |
ENDIF |
580 |
IF ( use3dCoriolis .AND. k.GT.1 ) THEN |
581 |
CALL MOM_W_CORIOLIS_NH( |
582 |
I bi,bj,k, |
583 |
I uVel, vVel, |
584 |
O gwAdd, |
585 |
I myThid ) |
586 |
DO j=jMin,jMax |
587 |
DO i=iMin,iMax |
588 |
gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j) |
589 |
ENDDO |
590 |
ENDDO |
591 |
ENDIF |
592 |
|
593 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
594 |
|
595 |
#ifdef ALLOW_DIAGNOSTICS |
596 |
IF ( diagDiss ) THEN |
597 |
CALL DIAGNOSTICS_FILL( gwDiss, 'Wm_Diss ', |
598 |
& k, 1, 2, bi,bj, myThid ) |
599 |
C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL |
600 |
C does it only if k=1 (never the case here) |
601 |
c IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Diss ',bi,bj,myThid) |
602 |
ENDIF |
603 |
IF ( diagAdvec ) THEN |
604 |
CALL DIAGNOSTICS_FILL( gW, 'Wm_Advec', |
605 |
& k,Nr, 1, bi,bj, myThid ) |
606 |
c IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Advec',bi,bj,myThid) |
607 |
ENDIF |
608 |
#endif /* ALLOW_DIAGNOSTICS */ |
609 |
|
610 |
C-- Dissipation term inside the Adams-Bashforth: |
611 |
IF ( momViscosity .AND. momDissip_In_AB) THEN |
612 |
DO j=jMin,jMax |
613 |
DO i=iMin,iMax |
614 |
gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j) |
615 |
ENDDO |
616 |
ENDDO |
617 |
ENDIF |
618 |
|
619 |
C- Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights) |
620 |
C and save gW_[n] into gwNm1 for the next time step. |
621 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
622 |
CALL ADAMS_BASHFORTH3( |
623 |
I bi, bj, k, |
624 |
U gW, gwNm, gw_AB, |
625 |
I nHydStartAB, myIter, myThid ) |
626 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
627 |
CALL ADAMS_BASHFORTH2( |
628 |
I bi, bj, k, |
629 |
U gW, gwNm1, gw_AB, |
630 |
I nHydStartAB, myIter, myThid ) |
631 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
632 |
#ifdef ALLOW_DIAGNOSTICS |
633 |
IF ( diag_AB ) THEN |
634 |
CALL DIAGNOSTICS_FILL(gw_AB,'AB_gW ',k,1,2,bi,bj,myThid) |
635 |
ENDIF |
636 |
#endif /* ALLOW_DIAGNOSTICS */ |
637 |
|
638 |
C-- Dissipation term outside the Adams-Bashforth: |
639 |
IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN |
640 |
DO j=jMin,jMax |
641 |
DO i=iMin,iMax |
642 |
gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j) |
643 |
ENDDO |
644 |
ENDDO |
645 |
ENDIF |
646 |
|
647 |
C- end of the k loop |
648 |
ENDDO |
649 |
|
650 |
#ifdef ALLOW_DIAGNOSTICS |
651 |
IF (useDiagnostics) THEN |
652 |
CALL DIAGNOSTICS_FILL(viscAh_W,'VISCAHW ',0,Nr,1,bi,bj,myThid) |
653 |
CALL DIAGNOSTICS_FILL(viscA4_W,'VISCA4W ',0,Nr,1,bi,bj,myThid) |
654 |
ENDIF |
655 |
#endif /* ALLOW_DIAGNOSTICS */ |
656 |
|
657 |
#endif /* ALLOW_NONHYDROSTATIC */ |
658 |
|
659 |
RETURN |
660 |
END |