/[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.40 - (hide annotations) (download)
Fri Feb 27 01:47:41 2009 UTC (15 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61n, checkpoint61q, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61x, checkpoint61y
Changes since 1.39: +54 -56 lines
- fix filling of diagnostics: Wm_Diss & Wm_Advec
- fix viscA4 on CS-grid : use rLow & rSurf @ W & S location (new)
   + apply FILL_CS_CORNER_TR_RL to wVel before del2w calculation

1 jmc 1.40 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.39 2008/04/22 22:20:31 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 adcroft 1.1 _RL wOverride
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.28 C-- Boundaries condition at top
147 jmc 1.30 DO j=1-OLy,sNy+OLy
148     DO i=1-OLx,sNx+OLx
149     flxAdvUp(i,j) = 0.
150     flxDisUp(i,j) = 0.
151 jmc 1.28 ENDDO
152     ENDDO
153 adcroft 1.1
154 jmc 1.28 C--- Sweep down column
155     DO k=2,Nr
156     kp1=k+1
157     wOverRide=1.
158     IF (k.EQ.Nr) THEN
159     kp1=Nr
160 adcroft 1.1 wOverRide=0.
161 jmc 1.28 ENDIF
162 jmc 1.30 C-- Compute grid factor arround a W-point:
163 jmc 1.36 #ifdef CALC_GW_NEW_THICK
164     DO j=1-Oly,sNy+Oly
165     DO i=1-Olx,sNx+Olx
166     IF ( maskC(i,j,k-1,bi,bj).EQ.0. .OR.
167     & maskC(i,j, k ,bi,bj).EQ.0. ) THEN
168     recip_rThickC(i,j) = 0.
169     ELSE
170     C- valid in z & p coord.; also accurate if Interface @ middle between 2 centers
171     recip_rThickC(i,j) = 1. _d 0 /
172     & ( MIN( Ro_surf(i,j,bi,bj),rC(k-1) )
173     & - MAX( R_low(i,j,bi,bj), rC(k) )
174     & )
175     ENDIF
176     ENDDO
177     ENDDO
178     IF (momViscosity) THEN
179     DO j=1-Oly,sNy+Oly
180     DO i=1-Olx,sNx+Olx
181     rThickC_C(i,j) = MAX( zeroRS,
182     & MIN( Ro_surf(i,j,bi,bj), rC(k-1) )
183     & -MAX( R_low(i,j,bi,bj), rC(k) )
184     & )
185     ENDDO
186     ENDDO
187     DO j=1-Oly,sNy+Oly
188     DO i=1-Olx+1,sNx+Olx
189     rThickC_W(i,j) = MAX( zeroRS,
190 jmc 1.40 & MIN( rSurfW(i,j,bi,bj), rC(k-1) )
191     & -MAX( rLowW(i,j,bi,bj), rC(k) )
192 jmc 1.36 & )
193     C W-Cell Western face area:
194     xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
195     c & *deepFacF(k)
196     ENDDO
197     ENDDO
198     DO j=1-Oly+1,sNy+Oly
199     DO i=1-Olx,sNx+Olx
200     rThickC_S(i,j) = MAX( zeroRS,
201 jmc 1.40 & MIN( rSurfS(i,j,bi,bj), rC(k-1) )
202     & -MAX( rLowS(i,j,bi,bj), rC(k) )
203 jmc 1.36 & )
204     C W-Cell Southern face area:
205     yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
206     c & *deepFacF(k)
207     C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC.
208     C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted)
209     ENDDO
210     ENDDO
211     ENDIF
212     #else /* CALC_GW_NEW_THICK */
213 jmc 1.30 DO j=1-Oly,sNy+Oly
214     DO i=1-Olx,sNx+Olx
215     C- note: assume fluid @ smaller k than bottom: does not work in p-coordinate !
216     IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN
217     recip_rThickC(i,j) = 0.
218     ELSE
219     recip_rThickC(i,j) = 1. _d 0 /
220 jmc 1.33 & ( drF(k-1)*halfRS
221 jmc 1.30 & + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )
222     & )
223     ENDIF
224 jmc 1.34 c IF (momViscosity) THEN
225     #ifdef NONLIN_FRSURF
226 jmc 1.35 rThickC_C(i,j) =
227 jmc 1.34 & drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
228     & + drF( k )*MIN( h0FacC(i,j,k ,bi,bj), halfRS )
229     #else
230 jmc 1.35 rThickC_C(i,j) =
231 jmc 1.34 & drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
232     & + drF( k )*MIN( _hFacC(i,j,k ,bi,bj), halfRS )
233     #endif
234     rThickC_W(i,j) =
235     & drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )
236     & + drF( k )*MIN( _hFacW(i,j,k ,bi,bj), halfRS )
237     rThickC_S(i,j) =
238     & drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )
239     & + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )
240 jmc 1.30 C W-Cell Western face area:
241     xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
242 jmc 1.35 c & *deepFacF(k)
243 jmc 1.30 C W-Cell Southern face area:
244     yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
245 jmc 1.35 c & *deepFacF(k)
246     C deep-model: xA,yA is only used for viscous flux, in terms like: xA/dxC,yA/dyC.
247     C this gives deepFacF*recip_deepFacF => cancel each other (and therefore omitted)
248 jmc 1.34 c ENDIF
249 jmc 1.30 ENDDO
250     ENDDO
251 jmc 1.36 #endif /* CALC_GW_NEW_THICK */
252 jmc 1.30
253 jmc 1.28 C-- horizontal bi-harmonic dissipation
254     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
255 jmc 1.40
256     C- local copy of wVel:
257     DO j=1-Oly,sNy+Oly
258     DO i=1-Olx,sNx+Olx
259     wFld(i,j) = wVel(i,j,k,bi,bj)
260     ENDDO
261     ENDDO
262 jmc 1.28 C- calculate the horizontal Laplacian of vertical flow
263 mlosch 1.18 C Zonal flux d/dx W
264 jmc 1.40 IF ( useCubedSphereExchange ) THEN
265     C to compute d/dx(W), fill corners with appropriate values:
266     CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
267     & wFld, bi,bj, myThid )
268     ENDIF
269 mlosch 1.18 DO j=1-Oly,sNy+Oly
270 jmc 1.30 flx_EW(1-Olx,j)=0.
271 mlosch 1.18 DO i=1-Olx+1,sNx+Olx
272 jmc 1.30 flx_EW(i,j) =
273 jmc 1.40 & ( wFld(i,j) - wFld(i-1,j) )
274 jmc 1.30 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
275 mlosch 1.18 #ifdef COSINEMETH_III
276 jmc 1.35 & *sqCosFacU(j,bi,bj)
277 mlosch 1.18 #endif
278     ENDDO
279 jmc 1.28 ENDDO
280 jmc 1.40
281 mlosch 1.18 C Meridional flux d/dy W
282 jmc 1.40 IF ( useCubedSphereExchange ) THEN
283     C to compute d/dy(W), fill corners with appropriate values:
284     CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
285     & wFld, bi,bj, myThid )
286     ENDIF
287 mlosch 1.18 DO i=1-Olx,sNx+Olx
288 jmc 1.30 flx_NS(i,1-Oly)=0.
289 mlosch 1.18 ENDDO
290     DO j=1-Oly+1,sNy+Oly
291     DO i=1-Olx,sNx+Olx
292 jmc 1.30 flx_NS(i,j) =
293 jmc 1.40 & ( wFld(i,j) - wFld(i,j-1) )
294 jmc 1.30 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
295 mlosch 1.18 #ifdef ISOTROPIC_COS_SCALING
296     #ifdef COSINEMETH_III
297 jmc 1.30 & *sqCosFacV(j,bi,bj)
298 mlosch 1.18 #endif
299     #endif
300     ENDDO
301     ENDDO
302 jmc 1.25
303 mlosch 1.18 C del^2 W
304 jmc 1.40 C Divergence of horizontal fluxes
305     DO j=1-Oly,sNy+Oly-1
306 mlosch 1.18 DO i=1-Olx,sNx+Olx-1
307 jmc 1.40 del2w(i,j) = ( ( flx_EW(i+1,j)-flx_EW(i,j) )
308     & +( flx_NS(i,j+1)-flx_NS(i,j) )
309 jmc 1.30 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
310 jmc 1.35 & *recip_deepFac2F(k)
311 mlosch 1.18 ENDDO
312     ENDDO
313 jmc 1.40 C end if biharmonic viscosity
314 jmc 1.28 ENDIF
315 mlosch 1.18
316 jmc 1.30 IF (momViscosity) THEN
317     C Viscous Flux on Western face
318     DO j=jMin,jMax
319     DO i=iMin,iMax+1
320     flx_EW(i,j)=
321     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i-1,j,k,bi,bj))*halfRL
322     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
323 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
324 jmc 1.34 & *cosFacU(j,bi,bj)
325 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL
326     & *(del2w(i,j)-del2w(i-1,j))
327 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
328 mlosch 1.18 #ifdef COSINEMETH_III
329 jmc 1.30 & *sqCosFacU(j,bi,bj)
330 mlosch 1.18 #else
331 jmc 1.34 & *cosFacU(j,bi,bj)
332 mlosch 1.18 #endif
333 jmc 1.30 ENDDO
334     ENDDO
335     C Viscous Flux on Southern face
336     DO j=jMin,jMax+1
337     DO i=iMin,iMax
338     flx_NS(i,j)=
339     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i,j-1,k,bi,bj))*halfRL
340     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
341 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
342 jmc 1.34 #ifdef ISOTROPIC_COS_SCALING
343     & *cosFacV(j,bi,bj)
344     #endif
345 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL
346     & *(del2w(i,j)-del2w(i,j-1))
347 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
348 jmc 1.30 #ifdef ISOTROPIC_COS_SCALING
349 mlosch 1.18 #ifdef COSINEMETH_III
350 jmc 1.30 & *sqCosFacV(j,bi,bj)
351 jmc 1.25 #else
352 jmc 1.34 & *cosFacV(j,bi,bj)
353 jmc 1.30 #endif
354 mlosch 1.18 #endif
355 jmc 1.30 ENDDO
356     ENDDO
357     C Viscous Flux on Lower face of W-Cell (= at tracer-cell center, level k)
358     DO j=jMin,jMax
359     DO i=iMin,iMax
360     C Interpolate vert viscosity to center of tracer-cell (level k):
361     viscLoc = ( KappaRU(i,j,k) +KappaRU(i+1,j,k)
362     & +KappaRU(i,j,kp1)+KappaRU(i+1,j,kp1)
363     & +KappaRV(i,j,k) +KappaRV(i,j+1,k)
364     & +KappaRV(i,j,kp1)+KappaRV(i,j+1,kp1)
365     & )*0.125 _d 0
366     flx_Dn(i,j) =
367     & - viscLoc*( wVel(i,j,kp1,bi,bj)*wOverRide
368     & -wVel(i,j, k ,bi,bj) )*rkSign
369 jmc 1.31 & *recip_drF(k)*rA(i,j,bi,bj)
370 jmc 1.35 & *deepFac2C(k)*rhoFacC(k)
371 jmc 1.30 ENDDO
372     ENDDO
373     C Tendency is minus divergence of viscous fluxes:
374 jmc 1.35 C anelastic: vert.visc.flx is scaled by rhoFac but hor.visc.fluxes are not
375 jmc 1.30 DO j=jMin,jMax
376     DO i=iMin,iMax
377     gwDiss(i,j) =
378 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
379     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
380     & + ( flx_Dn(i,j)-flxDisUp(i,j) )*rkSign
381 jmc 1.35 & *recip_rhoFacF(k)
382 jmc 1.31 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
383 jmc 1.35 & *recip_deepFac2F(k)
384 jmc 1.28 C-- prepare for next level (k+1)
385 jmc 1.30 flxDisUp(i,j)=flx_Dn(i,j)
386     ENDDO
387     ENDDO
388     ENDIF
389    
390 jmc 1.36 IF ( momViscosity .AND. no_slip_sides ) THEN
391 jmc 1.30 C- No-slip BCs impose a drag at walls...
392 jmc 1.34 CALL MOM_W_SIDEDRAG(
393     I bi,bj,k,
394     I wVel, del2w,
395     I rThickC_C, recip_rThickC,
396     I viscAh_W, viscA4_W,
397     O gwAdd,
398     I myThid )
399     DO j=jMin,jMax
400     DO i=iMin,iMax
401     gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)
402     ENDDO
403     ENDDO
404 jmc 1.30 ENDIF
405    
406     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
407    
408     IF ( momAdvection ) THEN
409     C Advective Flux on Western face
410     DO j=jMin,jMax
411     DO i=iMin,iMax+1
412     C transport through Western face area:
413     uTrans = (
414     & drF(k-1)*_hFacW(i,j,k-1,bi,bj)*uVel(i,j,k-1,bi,bj)
415 jmc 1.35 & *rhoFacC(k-1)
416 jmc 1.30 & + drF( k )*_hFacW(i,j, k ,bi,bj)*uVel(i,j, k ,bi,bj)
417 jmc 1.35 & *rhoFacC(k)
418     & )*halfRL*_dyG(i,j,bi,bj)*deepFacF(k)
419 jmc 1.30 flx_EW(i,j)=
420     & uTrans*(wVel(i,j,k,bi,bj)+wVel(i-1,j,k,bi,bj))*halfRL
421     ENDDO
422     ENDDO
423     C Advective Flux on Southern face
424     DO j=jMin,jMax+1
425     DO i=iMin,iMax
426     C transport through Southern face area:
427     vTrans = (
428     & drF(k-1)*_hFacS(i,j,k-1,bi,bj)*vVel(i,j,k-1,bi,bj)
429 jmc 1.35 & *rhoFacC(k-1)
430 jmc 1.30 & +drF( k )*_hFacS(i,j, k ,bi,bj)*vVel(i,j, k ,bi,bj)
431 jmc 1.35 & *rhoFacC(k)
432     & )*halfRL*_dxG(i,j,bi,bj)*deepFacF(k)
433 jmc 1.30 flx_NS(i,j)=
434     & vTrans*(wVel(i,j,k,bi,bj)+wVel(i,j-1,k,bi,bj))*halfRL
435     ENDDO
436     ENDDO
437     C Advective Flux on Lower face of W-Cell (= at tracer-cell center, level k)
438     DO j=jMin,jMax
439     DO i=iMin,iMax
440 jmc 1.36 C NH in p-coord.: advect wSpeed [m/s] with rTrans
441     tmp_WbarZ = halfRL*
442     & ( wVel(i,j, k ,bi,bj)*rVel2wUnit(k)
443     & +wVel(i,j,kp1,bi,bj)*rVel2wUnit(kp1)*wOverRide )
444 jmc 1.30 C transport through Lower face area:
445 jmc 1.35 rTrans = halfRL*
446     & ( wVel(i,j, k ,bi,bj)*deepFac2F( k )*rhoFacF( k )
447     & +wVel(i,j,kp1,bi,bj)*deepFac2F(kp1)*rhoFacF(kp1)
448     & *wOverRide
449     & )*rA(i,j,bi,bj)
450 jmc 1.31 flx_Dn(i,j) = rTrans*tmp_WbarZ
451 jmc 1.30 ENDDO
452     ENDDO
453     C Tendency is minus divergence of advective fluxes:
454 jmc 1.35 C anelastic: all transports & advect. fluxes are scaled by rhoFac
455 jmc 1.30 DO j=jMin,jMax
456     DO i=iMin,iMax
457     gW(i,j,k,bi,bj) =
458 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
459     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
460 jmc 1.36 & + ( flx_Dn(i,j)-flxAdvUp(i,j) )*rkSign*wUnit2rVel(k)
461 jmc 1.31 & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
462 jmc 1.35 & *recip_deepFac2F(k)*recip_rhoFacF(k)
463 jmc 1.30 C-- prepare for next level (k+1)
464     flxAdvUp(i,j)=flx_Dn(i,j)
465     ENDDO
466     ENDDO
467     ENDIF
468    
469     IF ( useNHMTerms ) THEN
470     CALL MOM_W_METRIC_NH(
471     I bi,bj,k,
472     I uVel, vVel,
473     O gwAdd,
474     I myThid )
475     DO j=jMin,jMax
476     DO i=iMin,iMax
477     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
478     ENDDO
479     ENDDO
480     ENDIF
481 jmc 1.32 IF ( use3dCoriolis ) THEN
482 jmc 1.30 CALL MOM_W_CORIOLIS_NH(
483     I bi,bj,k,
484     I uVel, vVel,
485     O gwAdd,
486     I myThid )
487     DO j=jMin,jMax
488     DO i=iMin,iMax
489     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
490     ENDDO
491     ENDDO
492     ENDIF
493 adcroft 1.1
494 jmc 1.25 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
495    
496 jmc 1.39 #ifdef ALLOW_DIAGNOSTICS
497 jmc 1.40 IF ( diagDiss ) THEN
498     CALL DIAGNOSTICS_FILL( gwDiss, 'Wm_Diss ',
499     & k, 1, 2, bi,bj, myThid )
500     C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
501     C does it only if k=1 (never the case here)
502     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Diss ',bi,bj,myThid)
503     ENDIF
504     IF ( diagAdvec ) THEN
505     CALL DIAGNOSTICS_FILL( gW, 'Wm_Advec',
506     & k,Nr, 1, bi,bj, myThid )
507     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('Wm_Advec',bi,bj,myThid)
508     ENDIF
509 jmc 1.39 #endif /* ALLOW_DIAGNOSTICS */
510    
511 jmc 1.30 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 jmc 1.25 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 jmc 1.28 c CALL ADAMS_BASHFORTH3(
524     c I bi, bj, k,
525     c U gW, gwNm,
526 jmc 1.37 c I nHydStartAB, myIter, myThid )
527 jmc 1.25 c#else /* ALLOW_ADAMSBASHFORTH_3 */
528 jmc 1.28 CALL ADAMS_BASHFORTH2(
529     I bi, bj, k,
530     U gW, gwNm1,
531 jmc 1.37 I nHydStartAB, myIter, myThid )
532 jmc 1.25 c#endif /* ALLOW_ADAMSBASHFORTH_3 */
533 jmc 1.21
534 jmc 1.30 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 jmc 1.25 C- end of the k loop
544 adcroft 1.4 ENDDO
545    
546 jmc 1.38 #ifdef ALLOW_DIAGNOSTICS
547     IF (useDiagnostics) THEN
548     CALL DIAGNOSTICS_FILL(viscAh_W,'VISCAHW ',0,Nr,1,bi,bj,myThid)
549     CALL DIAGNOSTICS_FILL(viscA4_W,'VISCA4W ',0,Nr,1,bi,bj,myThid)
550     ENDIF
551     #endif /* ALLOW_DIAGNOSTICS */
552    
553 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
554    
555     RETURN
556     END

  ViewVC Help
Powered by ViewVC 1.1.22