/[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.36 - (hide annotations) (download)
Mon Mar 12 23:48:24 2007 UTC (17 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint58w_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint58x_post
Changes since 1.35: +83 -11 lines
- use reference profile for w <-> omega conversion in NH-code
- change grid factor arround a W-point to work in p-coord. (keep
  the previous version within #ifndef CALC_GW_NEW_THICK for now)

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

  ViewVC Help
Powered by ViewVC 1.1.22