/[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.41 - (hide annotations) (download)
Mon Nov 30 19:22:45 2009 UTC (14 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61z
Changes since 1.40: +45 -9 lines
Fix missing vertical flux of vert. momentum near surface (k=1). This fix
  a spurious source of energy for simple baroclinic adjusment test case.

1 jmc 1.41 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.40 2009/02/27 01:47:41 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     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 jmc 1.34 C rThickC_C :: thickness (in r-units) of W-Cell (centered on W pt)
65 jmc 1.30 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 jmc 1.40 C gwAdd :: other tendencies (Coriolis, Metric-terms)
73     C del2w :: laplacian of wVel
74     C wFld :: local copy of wVel
75 jmc 1.30 C i,j,k :: Loop counters
76 jmc 1.28 INTEGER iMin,iMax,jMin,jMax
77 jmc 1.30 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79     _RL rThickC_W (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RL rThickC_S (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 jmc 1.34 _RL rThickC_C (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 jmc 1.30 _RL recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 jmc 1.28 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 jmc 1.30 _RL flxAdvUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL flxDisUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     _RL gwDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89     _RL gwAdd (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 jmc 1.28 _RL del2w (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 jmc 1.40 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92 jmc 1.28 INTEGER i,j,k, kp1
93 jmc 1.41 _RL wOverride, surfFac
94 jmc 1.30 _RL tmp_WbarZ
95     _RL uTrans, vTrans, rTrans
96 jmc 1.28 _RL viscLoc
97 jmc 1.30 _RL halfRL
98     _RS halfRS, zeroRS
99     PARAMETER( halfRL = 0.5D0 )
100     PARAMETER( halfRS = 0.5 , zeroRS = 0. )
101 jmc 1.36 PARAMETER( iMin = 1 , iMax = sNx )
102     PARAMETER( jMin = 1 , jMax = sNy )
103 cnh 1.9 CEOP
104 jmc 1.39 #ifdef ALLOW_DIAGNOSTICS
105     LOGICAL diagDiss, diagAdvec
106     LOGICAL DIAGNOSTICS_IS_ON
107     EXTERNAL DIAGNOSTICS_IS_ON
108     #endif /* ALLOW_DIAGNOSTICS */
109 edhill 1.10
110 jmc 1.39 C-- Catch barotropic mode
111 jmc 1.25 IF ( Nr .LT. 2 ) RETURN
112    
113 jmc 1.39 #ifdef ALLOW_DIAGNOSTICS
114     IF ( useDiagnostics ) THEN
115     diagDiss = DIAGNOSTICS_IS_ON( 'Wm_Diss ', myThid )
116     diagAdvec = DIAGNOSTICS_IS_ON( 'Wm_Advec', myThid )
117     ELSE
118     diagDiss = .FALSE.
119     diagAdvec = .FALSE.
120     ENDIF
121     #endif /* ALLOW_DIAGNOSTICS */
122    
123 jmc 1.28 C-- Initialise gW to zero
124     DO k=1,Nr
125     DO j=1-OLy,sNy+OLy
126     DO i=1-OLx,sNx+OLx
127 adcroft 1.1 gW(i,j,k,bi,bj) = 0.
128     ENDDO
129     ENDDO
130 jmc 1.28 ENDDO
131 jmc 1.30 C- Initialise gwDiss to zero
132     DO j=1-OLy,sNy+OLy
133     DO i=1-OLx,sNx+OLx
134     gwDiss(i,j) = 0.
135     ENDDO
136     ENDDO
137 jmc 1.40 IF (momViscosity) THEN
138     C- Initialize del2w to zero:
139     DO j=1-Oly,sNy+Oly
140     DO i=1-Olx,sNx+Olx
141     del2w(i,j) = 0. _d 0
142     ENDDO
143     ENDDO
144     ENDIF
145 adcroft 1.1
146 jmc 1.41 C-- Boundaries condition at top (vertical advection of vertical momentum):
147     C Setting surfFac=0 ignores surface wVel (as if wVel_k=1 = 0)
148     surfFac = 1. _d 0
149     c surfFac = 0. _d 0
150 adcroft 1.1
151 jmc 1.28 C--- Sweep down column
152     DO k=2,Nr
153     kp1=k+1
154     wOverRide=1.
155     IF (k.EQ.Nr) THEN
156     kp1=Nr
157 adcroft 1.1 wOverRide=0.
158 jmc 1.28 ENDIF
159 jmc 1.30 C-- Compute grid factor arround a W-point:
160 jmc 1.36 #ifdef CALC_GW_NEW_THICK
161     DO j=1-Oly,sNy+Oly
162     DO i=1-Olx,sNx+Olx
163     IF ( maskC(i,j,k-1,bi,bj).EQ.0. .OR.
164     & maskC(i,j, k ,bi,bj).EQ.0. ) THEN
165     recip_rThickC(i,j) = 0.
166     ELSE
167     C- valid in z & p coord.; also accurate if Interface @ middle between 2 centers
168     recip_rThickC(i,j) = 1. _d 0 /
169     & ( MIN( Ro_surf(i,j,bi,bj),rC(k-1) )
170     & - MAX( R_low(i,j,bi,bj), rC(k) )
171     & )
172     ENDIF
173     ENDDO
174     ENDDO
175     IF (momViscosity) THEN
176     DO j=1-Oly,sNy+Oly
177     DO i=1-Olx,sNx+Olx
178     rThickC_C(i,j) = MAX( zeroRS,
179     & MIN( Ro_surf(i,j,bi,bj), rC(k-1) )
180     & -MAX( R_low(i,j,bi,bj), rC(k) )
181     & )
182     ENDDO
183     ENDDO
184     DO j=1-Oly,sNy+Oly
185     DO i=1-Olx+1,sNx+Olx
186     rThickC_W(i,j) = MAX( zeroRS,
187 jmc 1.40 & MIN( rSurfW(i,j,bi,bj), rC(k-1) )
188     & -MAX( rLowW(i,j,bi,bj), rC(k) )
189 jmc 1.36 & )
190     C W-Cell Western face area:
191     xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
192     c & *deepFacF(k)
193     ENDDO
194     ENDDO
195     DO j=1-Oly+1,sNy+Oly
196     DO i=1-Olx,sNx+Olx
197     rThickC_S(i,j) = MAX( zeroRS,
198 jmc 1.40 & MIN( rSurfS(i,j,bi,bj), rC(k-1) )
199     & -MAX( rLowS(i,j,bi,bj), rC(k) )
200 jmc 1.36 & )
201     C W-Cell Southern face area:
202     yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
203     c & *deepFacF(k)
204     C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC.
205     C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted)
206     ENDDO
207     ENDDO
208     ENDIF
209     #else /* CALC_GW_NEW_THICK */
210 jmc 1.30 DO j=1-Oly,sNy+Oly
211     DO i=1-Olx,sNx+Olx
212     C- note: assume fluid @ smaller k than bottom: does not work in p-coordinate !
213     IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN
214     recip_rThickC(i,j) = 0.
215     ELSE
216     recip_rThickC(i,j) = 1. _d 0 /
217 jmc 1.33 & ( drF(k-1)*halfRS
218 jmc 1.30 & + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )
219     & )
220     ENDIF
221 jmc 1.34 c IF (momViscosity) THEN
222     #ifdef NONLIN_FRSURF
223 jmc 1.35 rThickC_C(i,j) =
224 jmc 1.34 & drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
225     & + drF( k )*MIN( h0FacC(i,j,k ,bi,bj), halfRS )
226     #else
227 jmc 1.35 rThickC_C(i,j) =
228 jmc 1.34 & drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
229     & + drF( k )*MIN( _hFacC(i,j,k ,bi,bj), halfRS )
230     #endif
231     rThickC_W(i,j) =
232     & drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )
233     & + drF( k )*MIN( _hFacW(i,j,k ,bi,bj), halfRS )
234     rThickC_S(i,j) =
235     & drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )
236     & + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )
237 jmc 1.30 C W-Cell Western face area:
238     xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
239 jmc 1.35 c & *deepFacF(k)
240 jmc 1.30 C W-Cell Southern face area:
241     yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
242 jmc 1.35 c & *deepFacF(k)
243     C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC.
244     C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted)
245 jmc 1.34 c ENDIF
246 jmc 1.30 ENDDO
247     ENDDO
248 jmc 1.36 #endif /* CALC_GW_NEW_THICK */
249 jmc 1.30
250 jmc 1.28 C-- horizontal bi-harmonic dissipation
251     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
252 jmc 1.40
253     C- local copy of wVel:
254     DO j=1-Oly,sNy+Oly
255     DO i=1-Olx,sNx+Olx
256     wFld(i,j) = wVel(i,j,k,bi,bj)
257     ENDDO
258     ENDDO
259 jmc 1.28 C- calculate the horizontal Laplacian of vertical flow
260 mlosch 1.18 C Zonal flux d/dx W
261 jmc 1.40 IF ( useCubedSphereExchange ) THEN
262     C to compute d/dx(W), fill corners with appropriate values:
263     CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
264     & wFld, bi,bj, myThid )
265     ENDIF
266 mlosch 1.18 DO j=1-Oly,sNy+Oly
267 jmc 1.30 flx_EW(1-Olx,j)=0.
268 mlosch 1.18 DO i=1-Olx+1,sNx+Olx
269 jmc 1.30 flx_EW(i,j) =
270 jmc 1.40 & ( wFld(i,j) - wFld(i-1,j) )
271 jmc 1.30 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
272 mlosch 1.18 #ifdef COSINEMETH_III
273 jmc 1.35 & *sqCosFacU(j,bi,bj)
274 mlosch 1.18 #endif
275     ENDDO
276 jmc 1.28 ENDDO
277 jmc 1.40
278 mlosch 1.18 C Meridional flux d/dy W
279 jmc 1.40 IF ( useCubedSphereExchange ) THEN
280     C to compute d/dy(W), fill corners with appropriate values:
281     CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
282     & wFld, bi,bj, myThid )
283     ENDIF
284 mlosch 1.18 DO i=1-Olx,sNx+Olx
285 jmc 1.30 flx_NS(i,1-Oly)=0.
286 mlosch 1.18 ENDDO
287     DO j=1-Oly+1,sNy+Oly
288     DO i=1-Olx,sNx+Olx
289 jmc 1.30 flx_NS(i,j) =
290 jmc 1.40 & ( wFld(i,j) - wFld(i,j-1) )
291 jmc 1.30 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
292 mlosch 1.18 #ifdef ISOTROPIC_COS_SCALING
293     #ifdef COSINEMETH_III
294 jmc 1.30 & *sqCosFacV(j,bi,bj)
295 mlosch 1.18 #endif
296     #endif
297     ENDDO
298     ENDDO
299 jmc 1.25
300 mlosch 1.18 C del^2 W
301 jmc 1.40 C Divergence of horizontal fluxes
302     DO j=1-Oly,sNy+Oly-1
303 mlosch 1.18 DO i=1-Olx,sNx+Olx-1
304 jmc 1.40 del2w(i,j) = ( ( flx_EW(i+1,j)-flx_EW(i,j) )
305     & +( flx_NS(i,j+1)-flx_NS(i,j) )
306 jmc 1.30 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
307 jmc 1.35 & *recip_deepFac2F(k)
308 mlosch 1.18 ENDDO
309     ENDDO
310 jmc 1.40 C end if biharmonic viscosity
311 jmc 1.28 ENDIF
312 mlosch 1.18
313 jmc 1.30 IF (momViscosity) THEN
314     C Viscous Flux on Western face
315     DO j=jMin,jMax
316     DO i=iMin,iMax+1
317     flx_EW(i,j)=
318     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i-1,j,k,bi,bj))*halfRL
319     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
320 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
321 jmc 1.34 & *cosFacU(j,bi,bj)
322 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL
323     & *(del2w(i,j)-del2w(i-1,j))
324 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
325 mlosch 1.18 #ifdef COSINEMETH_III
326 jmc 1.30 & *sqCosFacU(j,bi,bj)
327 mlosch 1.18 #else
328 jmc 1.34 & *cosFacU(j,bi,bj)
329 mlosch 1.18 #endif
330 jmc 1.30 ENDDO
331     ENDDO
332     C Viscous Flux on Southern face
333     DO j=jMin,jMax+1
334     DO i=iMin,iMax
335     flx_NS(i,j)=
336     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i,j-1,k,bi,bj))*halfRL
337     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
338 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
339 jmc 1.34 #ifdef ISOTROPIC_COS_SCALING
340     & *cosFacV(j,bi,bj)
341     #endif
342 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL
343     & *(del2w(i,j)-del2w(i,j-1))
344 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
345 jmc 1.30 #ifdef ISOTROPIC_COS_SCALING
346 mlosch 1.18 #ifdef COSINEMETH_III
347 jmc 1.30 & *sqCosFacV(j,bi,bj)
348 jmc 1.25 #else
349 jmc 1.34 & *cosFacV(j,bi,bj)
350 jmc 1.30 #endif
351 mlosch 1.18 #endif
352 jmc 1.30 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 jmc 1.31 & *recip_drF(k)*rA(i,j,bi,bj)
367 jmc 1.35 & *deepFac2C(k)*rhoFacC(k)
368 jmc 1.30 ENDDO
369     ENDDO
370 jmc 1.41 IF ( k.EQ.2 ) THEN
371     C Viscous Flux on Upper face of W-Cell (= at tracer-cell center, level k-1)
372     DO j=jMin,jMax
373     DO i=iMin,iMax
374     C Interpolate horizontally (but not vertically) vert viscosity to center:
375     C Although background visc. might be defined at k=1, this is not
376     C generally true when using variable visc. (from vertical mixing scheme).
377     C Therefore, no vert. interp. and only horizontal interpolation.
378     viscLoc = ( KappaRU(i,j,k) + KappaRU(i+1,j,k)
379     & +KappaRV(i,j,k) + KappaRV(i,j+1,k)
380     & )*0.25 _d 0
381     flxDisUp(i,j) =
382     & - viscLoc*( wVel(i,j, k ,bi,bj)
383     & -wVel(i,j,k-1,bi,bj) )*rkSign
384     & *recip_drF(k-1)*rA(i,j,bi,bj)
385     & *deepFac2C(k-1)*rhoFacC(k-1)
386     C to recover old (before 2009/11/30) results (since flxDisUp(k=2) was zero)
387     c flxDisUp(i,j) = 0.
388     ENDDO
389     ENDDO
390     ENDIF
391 jmc 1.30 C Tendency is minus divergence of viscous fluxes:
392 jmc 1.35 C anelastic: vert.visc.flx is scaled by rhoFac but hor.visc.fluxes are not
393 jmc 1.30 DO j=jMin,jMax
394     DO i=iMin,iMax
395     gwDiss(i,j) =
396 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
397     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
398     & + ( flx_Dn(i,j)-flxDisUp(i,j) )*rkSign
399 jmc 1.35 & *recip_rhoFacF(k)
400 jmc 1.31 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
401 jmc 1.35 & *recip_deepFac2F(k)
402 jmc 1.28 C-- prepare for next level (k+1)
403 jmc 1.30 flxDisUp(i,j)=flx_Dn(i,j)
404     ENDDO
405     ENDDO
406     ENDIF
407    
408 jmc 1.36 IF ( momViscosity .AND. no_slip_sides ) THEN
409 jmc 1.30 C- No-slip BCs impose a drag at walls...
410 jmc 1.34 CALL MOM_W_SIDEDRAG(
411     I bi,bj,k,
412     I wVel, del2w,
413     I rThickC_C, recip_rThickC,
414     I viscAh_W, viscA4_W,
415     O gwAdd,
416     I myThid )
417     DO j=jMin,jMax
418     DO i=iMin,iMax
419     gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)
420     ENDDO
421     ENDDO
422 jmc 1.30 ENDIF
423    
424     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
425    
426     IF ( momAdvection ) THEN
427     C Advective Flux on Western face
428     DO j=jMin,jMax
429     DO i=iMin,iMax+1
430     C transport through Western face area:
431     uTrans = (
432     & drF(k-1)*_hFacW(i,j,k-1,bi,bj)*uVel(i,j,k-1,bi,bj)
433 jmc 1.35 & *rhoFacC(k-1)
434 jmc 1.30 & + drF( k )*_hFacW(i,j, k ,bi,bj)*uVel(i,j, k ,bi,bj)
435 jmc 1.35 & *rhoFacC(k)
436     & )*halfRL*_dyG(i,j,bi,bj)*deepFacF(k)
437 jmc 1.30 flx_EW(i,j)=
438     & uTrans*(wVel(i,j,k,bi,bj)+wVel(i-1,j,k,bi,bj))*halfRL
439     ENDDO
440     ENDDO
441     C Advective Flux on Southern face
442     DO j=jMin,jMax+1
443     DO i=iMin,iMax
444     C transport through Southern face area:
445     vTrans = (
446     & drF(k-1)*_hFacS(i,j,k-1,bi,bj)*vVel(i,j,k-1,bi,bj)
447 jmc 1.35 & *rhoFacC(k-1)
448 jmc 1.30 & +drF( k )*_hFacS(i,j, k ,bi,bj)*vVel(i,j, k ,bi,bj)
449 jmc 1.35 & *rhoFacC(k)
450     & )*halfRL*_dxG(i,j,bi,bj)*deepFacF(k)
451 jmc 1.30 flx_NS(i,j)=
452     & vTrans*(wVel(i,j,k,bi,bj)+wVel(i,j-1,k,bi,bj))*halfRL
453     ENDDO
454     ENDDO
455     C Advective Flux on Lower face of W-Cell (= at tracer-cell center, level k)
456     DO j=jMin,jMax
457     DO i=iMin,iMax
458 jmc 1.36 C NH in p-coord.: advect wSpeed [m/s] with rTrans
459     tmp_WbarZ = halfRL*
460     & ( wVel(i,j, k ,bi,bj)*rVel2wUnit(k)
461     & +wVel(i,j,kp1,bi,bj)*rVel2wUnit(kp1)*wOverRide )
462 jmc 1.30 C transport through Lower face area:
463 jmc 1.35 rTrans = halfRL*
464     & ( wVel(i,j, k ,bi,bj)*deepFac2F( k )*rhoFacF( k )
465     & +wVel(i,j,kp1,bi,bj)*deepFac2F(kp1)*rhoFacF(kp1)
466     & *wOverRide
467     & )*rA(i,j,bi,bj)
468 jmc 1.31 flx_Dn(i,j) = rTrans*tmp_WbarZ
469 jmc 1.30 ENDDO
470     ENDDO
471 jmc 1.41 IF ( k.EQ.2 ) THEN
472     C Advective Flux on Upper face of W-Cell (= at tracer-cell center, level k-1)
473     DO j=jMin,jMax
474     DO i=iMin,iMax
475     tmp_WbarZ = halfRL*
476     & ( wVel(i,j,k-1,bi,bj)*rVel2wUnit(k-1)*surfFac
477     & +wVel(i,j, k, bi,bj)*rVel2wUnit( k ) )
478     rTrans = halfRL*
479     & ( wVel(i,j,k-1,bi,bj)*deepFac2F(k-1)*rhoFacF(k-1)
480     & *surfFac
481     & +wVel(i,j, k ,bi,bj)*deepFac2F( k )*rhoFacF( k )
482     & )*rA(i,j,bi,bj)
483     flxAdvUp(i,j) = rTrans*tmp_WbarZ
484     C to recover old (before 2009/11/30) results (since flxAdvUp(k=2) was zero)
485     c flxAdvUp(i,j) = 0.
486     ENDDO
487     ENDDO
488     ENDIF
489 jmc 1.30 C Tendency is minus divergence of advective fluxes:
490 jmc 1.35 C anelastic: all transports & advect. fluxes are scaled by rhoFac
491 jmc 1.30 DO j=jMin,jMax
492     DO i=iMin,iMax
493     gW(i,j,k,bi,bj) =
494 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
495     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
496 jmc 1.36 & + ( flx_Dn(i,j)-flxAdvUp(i,j) )*rkSign*wUnit2rVel(k)
497 jmc 1.31 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
498 jmc 1.35 & *recip_deepFac2F(k)*recip_rhoFacF(k)
499 jmc 1.30 C-- prepare for next level (k+1)
500     flxAdvUp(i,j)=flx_Dn(i,j)
501     ENDDO
502     ENDDO
503     ENDIF
504    
505     IF ( useNHMTerms ) THEN
506     CALL MOM_W_METRIC_NH(
507     I bi,bj,k,
508     I uVel, vVel,
509     O gwAdd,
510     I myThid )
511     DO j=jMin,jMax
512     DO i=iMin,iMax
513     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
514     ENDDO
515     ENDDO
516     ENDIF
517 jmc 1.32 IF ( use3dCoriolis ) THEN
518 jmc 1.30 CALL MOM_W_CORIOLIS_NH(
519     I bi,bj,k,
520     I uVel, vVel,
521     O gwAdd,
522     I myThid )
523     DO j=jMin,jMax
524     DO i=iMin,iMax
525     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
526     ENDDO
527     ENDDO
528     ENDIF
529 adcroft 1.1
530 jmc 1.25 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
531    
532 jmc 1.39 #ifdef ALLOW_DIAGNOSTICS
533 jmc 1.40 IF ( diagDiss ) THEN
534     CALL DIAGNOSTICS_FILL( gwDiss, 'Wm_Diss ',
535     & k, 1, 2, bi,bj, myThid )
536     C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
537     C does it only if k=1 (never the case here)
538     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Diss ',bi,bj,myThid)
539     ENDIF
540     IF ( diagAdvec ) THEN
541     CALL DIAGNOSTICS_FILL( gW, 'Wm_Advec',
542     & k,Nr, 1, bi,bj, myThid )
543     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Advec',bi,bj,myThid)
544     ENDIF
545 jmc 1.39 #endif /* ALLOW_DIAGNOSTICS */
546    
547 jmc 1.30 C-- Dissipation term inside the Adams-Bashforth:
548     IF ( momViscosity .AND. momDissip_In_AB) THEN
549     DO j=jMin,jMax
550     DO i=iMin,iMax
551     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
552     ENDDO
553     ENDDO
554     ENDIF
555    
556 jmc 1.25 C- Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights)
557     C and save gW_[n] into gwNm1 for the next time step.
558     c#ifdef ALLOW_ADAMSBASHFORTH_3
559 jmc 1.28 c CALL ADAMS_BASHFORTH3(
560     c I bi, bj, k,
561     c U gW, gwNm,
562 jmc 1.37 c I nHydStartAB, myIter, myThid )
563 jmc 1.25 c#else /* ALLOW_ADAMSBASHFORTH_3 */
564 jmc 1.28 CALL ADAMS_BASHFORTH2(
565     I bi, bj, k,
566     U gW, gwNm1,
567 jmc 1.37 I nHydStartAB, myIter, myThid )
568 jmc 1.25 c#endif /* ALLOW_ADAMSBASHFORTH_3 */
569 jmc 1.21
570 jmc 1.30 C-- Dissipation term outside the Adams-Bashforth:
571     IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN
572     DO j=jMin,jMax
573     DO i=iMin,iMax
574     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
575     ENDDO
576     ENDDO
577     ENDIF
578    
579 jmc 1.25 C- end of the k loop
580 adcroft 1.4 ENDDO
581    
582 jmc 1.38 #ifdef ALLOW_DIAGNOSTICS
583     IF (useDiagnostics) THEN
584     CALL DIAGNOSTICS_FILL(viscAh_W,'VISCAHW ',0,Nr,1,bi,bj,myThid)
585     CALL DIAGNOSTICS_FILL(viscA4_W,'VISCA4W ',0,Nr,1,bi,bj,myThid)
586     ENDIF
587     #endif /* ALLOW_DIAGNOSTICS */
588    
589 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
590    
591     RETURN
592     END

  ViewVC Help
Powered by ViewVC 1.1.22