/[MITgcm]/MITgcm/pkg/mom_vecinv/mom_vecinv.F
ViewVC logotype

Annotation of /MITgcm/pkg/mom_vecinv/mom_vecinv.F

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


Revision 1.19 - (hide annotations) (download)
Wed May 26 14:50:10 2004 UTC (20 years ago) by adcroft
Branch: MAIN
Changes since 1.18: +9 -4 lines
Added variable viscosity for the vector invariant equations
based on Leith, 1968, Phys. Fluids (10) 1409-1416
 - the use of the variable viscosty in the no-slip boundary conditions
   has not been implemented (but should be)
 - new parameters viscC2leith and viscC4leith are non-dimensional
 - I decided to modulate the variable viscosuty with the same viscAhMax
   and viscA4max; ideally we should have another maximum based on dx^2/dt etc.

1 adcroft 1.19 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.18 2004/05/24 20:03:49 adcroft Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4 edhill 1.10 #include "PACKAGES_CONFIG.h"
5 adcroft 1.1 #include "CPP_OPTIONS.h"
6    
7     SUBROUTINE MOM_VECINV(
8     I bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
9 jmc 1.4 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
10 adcroft 1.1 U fVerU, fVerV,
11 jmc 1.15 I myTime, myIter, myThid)
12 adcroft 1.1 C /==========================================================\
13     C | S/R MOM_VECINV |
14     C | o Form the right hand-side of the momentum equation. |
15     C |==========================================================|
16     C | Terms are evaluated one layer at a time working from |
17     C | the bottom to the top. The vertically integrated |
18     C | barotropic flow tendency term is evluated by summing the |
19     C | tendencies. |
20     C | Notes: |
21     C | We have not sorted out an entirely satisfactory formula |
22     C | for the diffusion equation bc with lopping. The present |
23     C | form produces a diffusive flux that does not scale with |
24     C | open-area. Need to do something to solidfy this and to |
25     C | deal "properly" with thin walls. |
26     C \==========================================================/
27     IMPLICIT NONE
28    
29     C == Global variables ==
30     #include "SIZE.h"
31     #include "DYNVARS.h"
32     #include "EEPARAMS.h"
33     #include "PARAMS.h"
34     #include "GRID.h"
35 jmc 1.7 #ifdef ALLOW_TIMEAVE
36     #include "TIMEAVE_STATV.h"
37     #endif
38 adcroft 1.1
39     C == Routine arguments ==
40     C fVerU - Flux of momentum in the vertical
41     C fVerV direction out of the upper face of a cell K
42     C ( flux into the cell above ).
43 jmc 1.4 C dPhiHydX,Y :: Gradient (X & Y dir.) of Hydrostatic Potential
44 adcroft 1.1 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
45     C results will be set.
46     C kUp, kDown - Index for upper and lower layers.
47     C myThid - Instance number for this innvocation of CALC_MOM_RHS
48 jmc 1.4 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
49     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
50 adcroft 1.1 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
51     _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
52     _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
53     _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
54     INTEGER kUp,kDown
55 jmc 1.15 _RL myTime
56 adcroft 1.2 INTEGER myIter
57 adcroft 1.1 INTEGER myThid
58     INTEGER bi,bj,iMin,iMax,jMin,jMax
59    
60 edhill 1.11 #ifdef ALLOW_MOM_VECINV
61 jmc 1.7
62 adcroft 1.2 C == Functions ==
63     LOGICAL DIFFERENT_MULTIPLE
64     EXTERNAL DIFFERENT_MULTIPLE
65    
66 adcroft 1.1 C == Local variables ==
67     _RL aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68     _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
69     _RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70     _RL uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RL vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72     _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73     _RL pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74     _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75     _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 adcroft 1.3 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77     _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 adcroft 1.1 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79     _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     C I,J,K - Loop counters
89     INTEGER i,j,k
90     C rVelMaskOverride - Factor for imposing special surface boundary conditions
91     C ( set according to free-surface condition ).
92     C hFacROpen - Lopped cell factos used tohold fraction of open
93     C hFacRClosed and closed cell wall.
94     _RL rVelMaskOverride
95     C xxxFac - On-off tracer parameters used for switching terms off.
96     _RL uDudxFac
97     _RL AhDudxFac
98     _RL A4DuxxdxFac
99     _RL vDudyFac
100     _RL AhDudyFac
101     _RL A4DuyydyFac
102     _RL rVelDudrFac
103     _RL ArDudrFac
104     _RL fuFac
105     _RL phxFac
106     _RL mtFacU
107     _RL uDvdxFac
108     _RL AhDvdxFac
109     _RL A4DvxxdxFac
110     _RL vDvdyFac
111     _RL AhDvdyFac
112     _RL A4DvyydyFac
113     _RL rVelDvdrFac
114     _RL ArDvdrFac
115     _RL fvFac
116     _RL phyFac
117     _RL vForcFac
118     _RL mtFacV
119     _RL wVelBottomOverride
120     LOGICAL bottomDragTerms
121 jmc 1.15 LOGICAL writeDiag
122 adcroft 1.1 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123     _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124     _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125     _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126    
127 heimbach 1.9 #ifdef ALLOW_AUTODIFF_TAMC
128     C-- only the kDown part of fverU/V is set in this subroutine
129     C-- the kUp is still required
130     C-- In the case of mom_fluxform Kup is set as well
131     C-- (at least in part)
132     fVerU(1,1,kUp) = fVerU(1,1,kUp)
133     fVerV(1,1,kUp) = fVerV(1,1,kUp)
134     #endif
135    
136 adcroft 1.1 rVelMaskOverride=1.
137     IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
138     wVelBottomOverride=1.
139     IF (k.EQ.Nr) wVelBottomOverride=0.
140 jmc 1.15 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime,
141     & myTime-deltaTClock)
142 adcroft 1.1
143     C Initialise intermediate terms
144     DO J=1-OLy,sNy+OLy
145     DO I=1-OLx,sNx+OLx
146     aF(i,j) = 0.
147     vF(i,j) = 0.
148     vrF(i,j) = 0.
149     uCf(i,j) = 0.
150     vCf(i,j) = 0.
151     mT(i,j) = 0.
152     pF(i,j) = 0.
153     del2u(i,j) = 0.
154     del2v(i,j) = 0.
155     dStar(i,j) = 0.
156     zStar(i,j) = 0.
157     uDiss(i,j) = 0.
158     vDiss(i,j) = 0.
159     vort3(i,j) = 0.
160     omega3(i,j) = 0.
161     ke(i,j) = 0.
162 heimbach 1.8 #ifdef ALLOW_AUTODIFF_TAMC
163     strain(i,j) = 0. _d 0
164     tension(i,j) = 0. _d 0
165     #endif
166 adcroft 1.1 ENDDO
167     ENDDO
168    
169     C-- Term by term tracer parmeters
170     C o U momentum equation
171     uDudxFac = afFacMom*1.
172     AhDudxFac = vfFacMom*1.
173     A4DuxxdxFac = vfFacMom*1.
174     vDudyFac = afFacMom*1.
175     AhDudyFac = vfFacMom*1.
176     A4DuyydyFac = vfFacMom*1.
177     rVelDudrFac = afFacMom*1.
178     ArDudrFac = vfFacMom*1.
179     mTFacU = mtFacMom*1.
180     fuFac = cfFacMom*1.
181     phxFac = pfFacMom*1.
182     C o V momentum equation
183     uDvdxFac = afFacMom*1.
184     AhDvdxFac = vfFacMom*1.
185     A4DvxxdxFac = vfFacMom*1.
186     vDvdyFac = afFacMom*1.
187     AhDvdyFac = vfFacMom*1.
188     A4DvyydyFac = vfFacMom*1.
189     rVelDvdrFac = afFacMom*1.
190     ArDvdrFac = vfFacMom*1.
191     mTFacV = mtFacMom*1.
192     fvFac = cfFacMom*1.
193     phyFac = pfFacMom*1.
194     vForcFac = foFacMom*1.
195    
196     IF ( no_slip_bottom
197     & .OR. bottomDragQuadratic.NE.0.
198     & .OR. bottomDragLinear.NE.0.) THEN
199     bottomDragTerms=.TRUE.
200     ELSE
201     bottomDragTerms=.FALSE.
202     ENDIF
203    
204     C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
205     IF (staggerTimeStep) THEN
206     phxFac = 0.
207     phyFac = 0.
208     ENDIF
209    
210     C-- Calculate open water fraction at vorticity points
211     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
212    
213     C---- Calculate common quantities used in both U and V equations
214     C Calculate tracer cell face open areas
215     DO j=1-OLy,sNy+OLy
216     DO i=1-OLx,sNx+OLx
217     xA(i,j) = _dyG(i,j,bi,bj)
218     & *drF(k)*_hFacW(i,j,k,bi,bj)
219     yA(i,j) = _dxG(i,j,bi,bj)
220     & *drF(k)*_hFacS(i,j,k,bi,bj)
221     ENDDO
222     ENDDO
223    
224     C Make local copies of horizontal flow field
225     DO j=1-OLy,sNy+OLy
226     DO i=1-OLx,sNx+OLx
227     uFld(i,j) = uVel(i,j,k,bi,bj)
228     vFld(i,j) = vVel(i,j,k,bi,bj)
229     ENDDO
230     ENDDO
231    
232 jmc 1.7 C note (jmc) : Dissipation and Vort3 advection do not necesary
233     C use the same maskZ (and hFacZ) => needs 2 call(s)
234     c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
235    
236 adcroft 1.16 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
237 adcroft 1.1
238 adcroft 1.17 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
239 adcroft 1.1
240 adcroft 1.18 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
241 adcroft 1.1
242 adcroft 1.18 c CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
243 adcroft 1.1
244     IF (momViscosity) THEN
245     C Calculate del^2 u and del^2 v for bi-harmonic term
246 adcroft 1.19 IF (viscA4.NE.0.
247     & .OR. viscA4Grid.NE.0.
248     & .OR. viscC4leith.NE.0.
249     & ) THEN
250 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
251     O del2u,del2v,
252     & myThid)
253 adcroft 1.17 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
254 adcroft 1.18 CALL MOM_CALC_RELVORT3(
255 adcroft 1.2 & bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)
256     ENDIF
257 adcroft 1.1 C Calculate dissipation terms for U and V equations
258 adcroft 1.2 C in terms of vorticity and divergence
259 adcroft 1.19 IF (viscAh.NE.0. .OR. viscA4.NE.0.
260     & .OR. viscAhGrid.NE.0. .OR. viscA4Grid.NE.0.
261     & .OR. viscC2leith.NE.0. .OR. viscC4leith.NE.0.
262     & ) THEN
263 adcroft 1.2 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,
264     O uDiss,vDiss,
265     & myThid)
266     ENDIF
267 adcroft 1.3 C or in terms of tension and strain
268     IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
269     CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
270     O tension,
271     I myThid)
272     CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
273     O strain,
274     I myThid)
275     CALL MOM_HDISSIP(bi,bj,k,
276     I tension,strain,hFacZ,viscAtension,viscAstrain,
277     O uDiss,vDiss,
278     I myThid)
279     ENDIF
280 adcroft 1.1 ENDIF
281    
282 jmc 1.7 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
283     c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
284    
285 adcroft 1.1 C---- Zonal momentum equation starts here
286    
287     C-- Vertical flux (fVer is at upper face of "u" cell)
288    
289     C Eddy component of vertical flux (interior component only) -> vrF
290     IF (momViscosity.AND..NOT.implicitViscosity)
291     & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
292    
293     C Combine fluxes
294     DO j=jMin,jMax
295     DO i=iMin,iMax
296     fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
297     ENDDO
298     ENDDO
299    
300     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
301     DO j=2-Oly,sNy+Oly-1
302     DO i=2-Olx,sNx+Olx-1
303     gU(i,j,k,bi,bj) = uDiss(i,j)
304     & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
305     & *recip_rAw(i,j,bi,bj)
306     & *(
307     & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
308     & )
309 jmc 1.4 & - phxFac*dPhiHydX(i,j)
310 adcroft 1.1 ENDDO
311     ENDDO
312    
313     C-- No-slip and drag BCs appear as body forces in cell abutting topography
314     IF (momViscosity.AND.no_slip_sides) THEN
315     C- No-slip BCs impose a drag at walls...
316     CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
317     DO j=jMin,jMax
318     DO i=iMin,iMax
319     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
320     ENDDO
321     ENDDO
322     ENDIF
323 heimbach 1.8
324 adcroft 1.1 C- No-slip BCs impose a drag at bottom
325     IF (momViscosity.AND.bottomDragTerms) THEN
326     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
327     DO j=jMin,jMax
328     DO i=iMin,iMax
329     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
330     ENDDO
331     ENDDO
332     ENDIF
333    
334     C-- Metric terms for curvilinear grid systems
335     c IF (usingSphericalPolarMTerms) THEN
336     C o Spherical polar grid metric terms
337     c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
338     c DO j=jMin,jMax
339     c DO i=iMin,iMax
340     c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
341     c ENDDO
342     c ENDDO
343     c ENDIF
344    
345     C---- Meridional momentum equation starts here
346    
347     C-- Vertical flux (fVer is at upper face of "v" cell)
348    
349     C Eddy component of vertical flux (interior component only) -> vrF
350     IF (momViscosity.AND..NOT.implicitViscosity)
351     & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
352    
353     C Combine fluxes -> fVerV
354     DO j=jMin,jMax
355     DO i=iMin,iMax
356     fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
357     ENDDO
358     ENDDO
359    
360     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
361     DO j=jMin,jMax
362     DO i=iMin,iMax
363     gV(i,j,k,bi,bj) = vDiss(i,j)
364     & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
365     & *recip_rAs(i,j,bi,bj)
366     & *(
367     & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
368     & )
369 jmc 1.4 & - phyFac*dPhiHydY(i,j)
370 adcroft 1.1 ENDDO
371     ENDDO
372    
373     C-- No-slip and drag BCs appear as body forces in cell abutting topography
374     IF (momViscosity.AND.no_slip_sides) THEN
375     C- No-slip BCs impose a drag at walls...
376     CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
377     DO j=jMin,jMax
378     DO i=iMin,iMax
379     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
380     ENDDO
381     ENDDO
382     ENDIF
383     C- No-slip BCs impose a drag at bottom
384     IF (momViscosity.AND.bottomDragTerms) THEN
385     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
386     DO j=jMin,jMax
387     DO i=iMin,iMax
388     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
389     ENDDO
390     ENDDO
391     ENDIF
392    
393     C-- Metric terms for curvilinear grid systems
394     c IF (usingSphericalPolarMTerms) THEN
395     C o Spherical polar grid metric terms
396     c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
397     c DO j=jMin,jMax
398     c DO i=iMin,iMax
399     c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
400     c ENDDO
401     c ENDDO
402     c ENDIF
403    
404 jmc 1.5 C-- Horizontal Coriolis terms
405 jmc 1.6 IF (useCoriolis .AND. .NOT.useCDscheme) THEN
406 jmc 1.7 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,omega3,hFacZ,r_hFacZ,
407 jmc 1.5 & uCf,vCf,myThid)
408     DO j=jMin,jMax
409     DO i=iMin,iMax
410     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
411     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
412     ENDDO
413 adcroft 1.1 ENDDO
414 jmc 1.15 IF ( writeDiag ) THEN
415     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
416     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
417     ENDIF
418 jmc 1.5 ENDIF
419 adcroft 1.1
420 jmc 1.5 IF (momAdvection) THEN
421     C-- Horizontal advection of relative vorticity
422     c CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,r_hFacZ,uCf,myThid)
423 jmc 1.7 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3,hFacZ,r_hFacZ,
424     & uCf,myThid)
425 jmc 1.5 c CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid)
426     DO j=jMin,jMax
427     DO i=iMin,iMax
428     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
429     ENDDO
430 adcroft 1.1 ENDDO
431 jmc 1.5 c CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,r_hFacZ,vCf,myThid)
432 jmc 1.7 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3,hFacZ,r_hFacZ,
433     & vCf,myThid)
434 jmc 1.5 c CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid)
435     DO j=jMin,jMax
436     DO i=iMin,iMax
437     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
438     ENDDO
439 adcroft 1.1 ENDDO
440    
441 jmc 1.15 IF ( writeDiag ) THEN
442     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
443     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
444     ENDIF
445 jmc 1.7 #ifdef ALLOW_TIMEAVE
446 dimitri 1.13 #ifndef HRCUBE
447 jmc 1.7 IF (taveFreq.GT.0.) THEN
448     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
449     & Nr, k, bi, bj, myThid)
450     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
451     & Nr, k, bi, bj, myThid)
452     ENDIF
453 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
454     #endif /* ndef HRCUBE */
455 jmc 1.7
456 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
457 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
458     CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
459     DO j=jMin,jMax
460     DO i=iMin,iMax
461     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
462     ENDDO
463 jmc 1.5 ENDDO
464 jmc 1.12 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
465     DO j=jMin,jMax
466     DO i=iMin,iMax
467     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
468     ENDDO
469 jmc 1.5 ENDDO
470 jmc 1.12 ENDIF
471 adcroft 1.1
472     C-- Bernoulli term
473 jmc 1.5 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
474     DO j=jMin,jMax
475     DO i=iMin,iMax
476     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
477     ENDDO
478     ENDDO
479     CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
480     DO j=jMin,jMax
481     DO i=iMin,iMax
482     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
483     ENDDO
484 adcroft 1.1 ENDDO
485 jmc 1.15 IF ( writeDiag ) THEN
486     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
487     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
488     ENDIF
489    
490 jmc 1.5 C-- end if momAdvection
491     ENDIF
492    
493     C-- Set du/dt & dv/dt on boundaries to zero
494 adcroft 1.1 DO j=jMin,jMax
495     DO i=iMin,iMax
496 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
497     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
498 adcroft 1.1 ENDDO
499     ENDDO
500 jmc 1.5
501 adcroft 1.2
502 jmc 1.15 IF ( writeDiag ) THEN
503 adcroft 1.3 CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
504     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
505 adcroft 1.2 CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid)
506     CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid)
507 adcroft 1.3 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
508 jmc 1.5 c CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
509 adcroft 1.3 CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
510     CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid)
511 adcroft 1.1 ENDIF
512 jmc 1.7
513 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
514 adcroft 1.1
515     RETURN
516     END

  ViewVC Help
Powered by ViewVC 1.1.22