/[MITgcm]/MITgcm/model/src/calc_gw.F
ViewVC logotype

Annotation of /MITgcm/model/src/calc_gw.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.38 - (hide annotations) (download)
Mon Mar 24 00:32:53 2008 UTC (16 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59p
Changes since 1.37: +9 -23 lines
fix diagnostics VISCAHW & VISCA4W.

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

  ViewVC Help
Powered by ViewVC 1.1.22