/[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.28 - (hide annotations) (download)
Tue Nov 2 01:04:08 2004 UTC (19 years, 7 months ago) by jmc
Branch: MAIN
Changes since 1.27: +5 -4 lines
updated (horiz. viscosity for Divergence and Vorticity)

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

  ViewVC Help
Powered by ViewVC 1.1.22