/[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.44 - (hide annotations) (download)
Tue May 11 20:54:11 2010 UTC (14 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w
Changes since 1.43: +3 -1 lines
add commented line to reproduce old results (with bug in vertical flux)

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

  ViewVC Help
Powered by ViewVC 1.1.22