/[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.30 - (hide annotations) (download)
Fri Nov 5 19:23:06 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
Changes since 1.29: +3 -12 lines
forgot 1 viscA4 last time (horiz. viscosity for Divergence and Vorticity)

1 jmc 1.30 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.29 2004/11/05 18:39:15 jmc 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 vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70     _RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RL uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72     _RL vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 jmc 1.29 c _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 adcroft 1.1 _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     _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     C I,J,K - Loop counters
87     INTEGER i,j,k
88     C xxxFac - On-off tracer parameters used for switching terms off.
89     _RL ArDudrFac
90     _RL phxFac
91 jmc 1.29 c _RL mtFacU
92 adcroft 1.1 _RL ArDvdrFac
93     _RL phyFac
94 jmc 1.29 c _RL mtFacV
95 adcroft 1.1 LOGICAL bottomDragTerms
96 jmc 1.15 LOGICAL writeDiag
97 adcroft 1.1 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98     _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99     _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100     _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101    
102 edhill 1.25 #ifdef ALLOW_MNC
103     INTEGER offsets(9)
104     #endif
105    
106 heimbach 1.9 #ifdef ALLOW_AUTODIFF_TAMC
107     C-- only the kDown part of fverU/V is set in this subroutine
108     C-- the kUp is still required
109     C-- In the case of mom_fluxform Kup is set as well
110     C-- (at least in part)
111     fVerU(1,1,kUp) = fVerU(1,1,kUp)
112     fVerV(1,1,kUp) = fVerV(1,1,kUp)
113     #endif
114    
115 jmc 1.15 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime,
116     & myTime-deltaTClock)
117 adcroft 1.1
118 edhill 1.24 #ifdef ALLOW_MNC
119     IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
120 edhill 1.25 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
121     CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
122     CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
123     CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
124     ENDIF
125     DO i = 1,9
126     offsets(i) = 0
127     ENDDO
128     offsets(3) = k
129     C write(*,*) 'offsets = ',(offsets(i),i=1,9)
130 edhill 1.24 ENDIF
131     #endif /* ALLOW_MNC */
132    
133 adcroft 1.1 C Initialise intermediate terms
134     DO J=1-OLy,sNy+OLy
135     DO I=1-OLx,sNx+OLx
136     vF(i,j) = 0.
137     vrF(i,j) = 0.
138     uCf(i,j) = 0.
139     vCf(i,j) = 0.
140 jmc 1.29 c mT(i,j) = 0.
141 adcroft 1.1 del2u(i,j) = 0.
142     del2v(i,j) = 0.
143     dStar(i,j) = 0.
144     zStar(i,j) = 0.
145     uDiss(i,j) = 0.
146     vDiss(i,j) = 0.
147     vort3(i,j) = 0.
148     omega3(i,j) = 0.
149     ke(i,j) = 0.
150 heimbach 1.8 #ifdef ALLOW_AUTODIFF_TAMC
151     strain(i,j) = 0. _d 0
152     tension(i,j) = 0. _d 0
153     #endif
154 adcroft 1.1 ENDDO
155     ENDDO
156    
157     C-- Term by term tracer parmeters
158     C o U momentum equation
159     ArDudrFac = vfFacMom*1.
160 jmc 1.29 c mTFacU = mtFacMom*1.
161 adcroft 1.1 phxFac = pfFacMom*1.
162     C o V momentum equation
163     ArDvdrFac = vfFacMom*1.
164 jmc 1.29 c mTFacV = mtFacMom*1.
165 adcroft 1.1 phyFac = pfFacMom*1.
166    
167     IF ( no_slip_bottom
168     & .OR. bottomDragQuadratic.NE.0.
169     & .OR. bottomDragLinear.NE.0.) THEN
170     bottomDragTerms=.TRUE.
171     ELSE
172     bottomDragTerms=.FALSE.
173     ENDIF
174    
175     C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
176     IF (staggerTimeStep) THEN
177     phxFac = 0.
178     phyFac = 0.
179     ENDIF
180    
181     C-- Calculate open water fraction at vorticity points
182     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
183    
184     C Make local copies of horizontal flow field
185     DO j=1-OLy,sNy+OLy
186     DO i=1-OLx,sNx+OLx
187     uFld(i,j) = uVel(i,j,k,bi,bj)
188     vFld(i,j) = vVel(i,j,k,bi,bj)
189     ENDDO
190     ENDDO
191    
192 jmc 1.7 C note (jmc) : Dissipation and Vort3 advection do not necesary
193     C use the same maskZ (and hFacZ) => needs 2 call(s)
194     c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
195    
196 adcroft 1.16 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
197 adcroft 1.1
198 adcroft 1.17 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
199 adcroft 1.1
200 adcroft 1.18 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
201 adcroft 1.1
202 adcroft 1.20 IF (useAbsVorticity)
203     & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
204 adcroft 1.1
205     IF (momViscosity) THEN
206     C Calculate del^2 u and del^2 v for bi-harmonic term
207 jmc 1.30 IF ( (viscA4.NE.0. .AND. no_slip_sides)
208     & .OR. viscA4D.NE.0. .OR. viscA4Z.NE.0.
209 adcroft 1.19 & .OR. viscA4Grid.NE.0.
210     & .OR. viscC4leith.NE.0.
211     & ) THEN
212 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
213     O del2u,del2v,
214     & myThid)
215 adcroft 1.17 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
216 adcroft 1.18 CALL MOM_CALC_RELVORT3(
217 adcroft 1.2 & bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)
218     ENDIF
219 adcroft 1.1 C Calculate dissipation terms for U and V equations
220 adcroft 1.2 C in terms of vorticity and divergence
221 jmc 1.28 IF ( viscAhD.NE.0. .OR. viscAhZ.NE.0.
222     & .OR. viscA4D.NE.0. .OR. viscA4Z.NE.0.
223     & .OR. viscAhGrid.NE.0. .OR. viscA4Grid.NE.0.
224     & .OR. viscC2leith.NE.0. .OR. viscC4leith.NE.0.
225 adcroft 1.19 & ) THEN
226 adcroft 1.2 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,
227     O uDiss,vDiss,
228     & myThid)
229     ENDIF
230 adcroft 1.3 C or in terms of tension and strain
231     IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
232     CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
233     O tension,
234     I myThid)
235     CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
236     O strain,
237     I myThid)
238     CALL MOM_HDISSIP(bi,bj,k,
239     I tension,strain,hFacZ,viscAtension,viscAstrain,
240     O uDiss,vDiss,
241     I myThid)
242     ENDIF
243 adcroft 1.1 ENDIF
244    
245 jmc 1.7 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
246     c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
247    
248 adcroft 1.1 C---- Zonal momentum equation starts here
249    
250     C-- Vertical flux (fVer is at upper face of "u" cell)
251    
252     C Eddy component of vertical flux (interior component only) -> vrF
253     IF (momViscosity.AND..NOT.implicitViscosity)
254     & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
255    
256     C Combine fluxes
257     DO j=jMin,jMax
258     DO i=iMin,iMax
259     fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
260     ENDDO
261     ENDDO
262    
263     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
264     DO j=2-Oly,sNy+Oly-1
265     DO i=2-Olx,sNx+Olx-1
266     gU(i,j,k,bi,bj) = uDiss(i,j)
267     & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
268     & *recip_rAw(i,j,bi,bj)
269     & *(
270     & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
271     & )
272 jmc 1.4 & - phxFac*dPhiHydX(i,j)
273 adcroft 1.1 ENDDO
274     ENDDO
275    
276     C-- No-slip and drag BCs appear as body forces in cell abutting topography
277     IF (momViscosity.AND.no_slip_sides) THEN
278     C- No-slip BCs impose a drag at walls...
279     CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
280     DO j=jMin,jMax
281     DO i=iMin,iMax
282     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
283     ENDDO
284     ENDDO
285     ENDIF
286 heimbach 1.8
287 adcroft 1.1 C- No-slip BCs impose a drag at bottom
288     IF (momViscosity.AND.bottomDragTerms) THEN
289     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
290     DO j=jMin,jMax
291     DO i=iMin,iMax
292     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
293     ENDDO
294     ENDDO
295     ENDIF
296    
297     C-- Metric terms for curvilinear grid systems
298     c IF (usingSphericalPolarMTerms) THEN
299     C o Spherical polar grid metric terms
300     c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
301     c DO j=jMin,jMax
302     c DO i=iMin,iMax
303     c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
304     c ENDDO
305     c ENDDO
306     c ENDIF
307    
308     C---- Meridional momentum equation starts here
309    
310     C-- Vertical flux (fVer is at upper face of "v" cell)
311    
312     C Eddy component of vertical flux (interior component only) -> vrF
313     IF (momViscosity.AND..NOT.implicitViscosity)
314     & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
315    
316     C Combine fluxes -> fVerV
317     DO j=jMin,jMax
318     DO i=iMin,iMax
319     fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
320     ENDDO
321     ENDDO
322    
323     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
324     DO j=jMin,jMax
325     DO i=iMin,iMax
326     gV(i,j,k,bi,bj) = vDiss(i,j)
327     & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
328     & *recip_rAs(i,j,bi,bj)
329     & *(
330     & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
331     & )
332 jmc 1.4 & - phyFac*dPhiHydY(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_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
340     DO j=jMin,jMax
341     DO i=iMin,iMax
342     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
343     ENDDO
344     ENDDO
345     ENDIF
346     C- No-slip BCs impose a drag at bottom
347     IF (momViscosity.AND.bottomDragTerms) THEN
348     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
349     DO j=jMin,jMax
350     DO i=iMin,iMax
351     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
352     ENDDO
353     ENDDO
354     ENDIF
355    
356     C-- Metric terms for curvilinear grid systems
357     c IF (usingSphericalPolarMTerms) THEN
358     C o Spherical polar grid metric terms
359     c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
360     c DO j=jMin,jMax
361     c DO i=iMin,iMax
362     c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
363     c ENDDO
364     c ENDDO
365     c ENDIF
366    
367 jmc 1.5 C-- Horizontal Coriolis terms
368 adcroft 1.20 IF (useCoriolis .AND. .NOT.useCDscheme
369     & .AND. .NOT. useAbsVorticity) THEN
370     CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
371 jmc 1.5 & uCf,vCf,myThid)
372     DO j=jMin,jMax
373     DO i=iMin,iMax
374     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
375     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
376     ENDDO
377 adcroft 1.1 ENDDO
378 jmc 1.15 IF ( writeDiag ) THEN
379 edhill 1.24 IF (snapshot_mdsio) THEN
380     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
381     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
382     ENDIF
383     #ifdef ALLOW_MNC
384     IF (useMNC .AND. snapshot_mnc) THEN
385 edhill 1.25 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fV', uCf,
386     & offsets, myThid)
387     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fU', vCf,
388     & offsets, myThid)
389 edhill 1.24 ENDIF
390     #endif /* ALLOW_MNC */
391 jmc 1.15 ENDIF
392 jmc 1.5 ENDIF
393 adcroft 1.1
394 jmc 1.5 IF (momAdvection) THEN
395     C-- Horizontal advection of relative vorticity
396 adcroft 1.20 IF (useAbsVorticity) THEN
397     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
398     & uCf,myThid)
399     ELSE
400     CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3,hFacZ,r_hFacZ,
401     & uCf,myThid)
402     ENDIF
403 jmc 1.5 c CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid)
404     DO j=jMin,jMax
405     DO i=iMin,iMax
406     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
407     ENDDO
408 adcroft 1.1 ENDDO
409 adcroft 1.20 IF (useAbsVorticity) THEN
410     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
411     & vCf,myThid)
412     ELSE
413     CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3,hFacZ,r_hFacZ,
414     & vCf,myThid)
415     ENDIF
416 jmc 1.5 c CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid)
417     DO j=jMin,jMax
418     DO i=iMin,iMax
419     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
420     ENDDO
421 adcroft 1.1 ENDDO
422    
423 jmc 1.15 IF ( writeDiag ) THEN
424 edhill 1.24 IF (snapshot_mdsio) THEN
425     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
426     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
427     ENDIF
428     #ifdef ALLOW_MNC
429     IF (useMNC .AND. snapshot_mnc) THEN
430 edhill 1.25 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zV', uCf,
431     & offsets, myThid)
432     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zU', vCf,
433     & offsets, myThid)
434 edhill 1.24 ENDIF
435     #endif /* ALLOW_MNC */
436 jmc 1.15 ENDIF
437 edhill 1.24
438 jmc 1.7 #ifdef ALLOW_TIMEAVE
439 dimitri 1.13 #ifndef HRCUBE
440 jmc 1.7 IF (taveFreq.GT.0.) THEN
441     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
442     & Nr, k, bi, bj, myThid)
443     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
444     & Nr, k, bi, bj, myThid)
445     ENDIF
446 jmc 1.22 #endif /* ndef HRCUBE */
447 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
448 jmc 1.7
449 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
450 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
451     CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
452     DO j=jMin,jMax
453     DO i=iMin,iMax
454     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
455     ENDDO
456 jmc 1.5 ENDDO
457 jmc 1.12 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
458     DO j=jMin,jMax
459     DO i=iMin,iMax
460     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
461     ENDDO
462 jmc 1.5 ENDDO
463 jmc 1.12 ENDIF
464 adcroft 1.1
465     C-- Bernoulli term
466 jmc 1.5 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
467     DO j=jMin,jMax
468     DO i=iMin,iMax
469     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
470     ENDDO
471     ENDDO
472     CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
473     DO j=jMin,jMax
474     DO i=iMin,iMax
475     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
476     ENDDO
477 adcroft 1.1 ENDDO
478 jmc 1.15 IF ( writeDiag ) THEN
479 edhill 1.24 IF (snapshot_mdsio) THEN
480     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
481     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
482     ENDIF
483     #ifdef ALLOW_MNC
484     IF (useMNC .AND. snapshot_mnc) THEN
485 edhill 1.25 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEx', uCf,
486     & offsets, myThid)
487     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEy', vCf,
488     & offsets, myThid)
489     ENDIF
490 edhill 1.24 #endif /* ALLOW_MNC */
491 jmc 1.15 ENDIF
492    
493 jmc 1.5 C-- end if momAdvection
494     ENDIF
495    
496     C-- Set du/dt & dv/dt on boundaries to zero
497 adcroft 1.1 DO j=jMin,jMax
498     DO i=iMin,iMax
499 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
500     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
501 adcroft 1.1 ENDDO
502     ENDDO
503 jmc 1.5
504 jmc 1.22 #ifdef ALLOW_DEBUG
505     IF ( debugLevel .GE. debLevB
506     & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
507     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
508     & .AND. useCubedSphereExchange ) THEN
509 jmc 1.23 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
510     & uDiss,vDiss, k, standardMessageUnit,bi,bj,myThid )
511 jmc 1.22 ENDIF
512     #endif /* ALLOW_DEBUG */
513 adcroft 1.2
514 jmc 1.15 IF ( writeDiag ) THEN
515 edhill 1.24 IF (snapshot_mdsio) THEN
516     CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
517     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,
518     & myThid)
519     CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid)
520     CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid)
521     CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
522     CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
523     CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
524     CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid)
525     ENDIF
526     #ifdef ALLOW_MNC
527     IF (useMNC .AND. snapshot_mnc) THEN
528 edhill 1.25 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Ds',strain,
529     & offsets, myThid)
530     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dt',tension,
531     & offsets, myThid)
532     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Du',uDiss,
533     & offsets, myThid)
534     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dv',vDiss,
535     & offsets, myThid)
536     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Z3',vort3,
537     & offsets, myThid)
538     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'W3',omega3,
539     & offsets, myThid)
540     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'KE',KE,
541     & offsets, myThid)
542     CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'D', hdiv,
543     & offsets, myThid)
544 edhill 1.24 ENDIF
545     #endif /* ALLOW_MNC */
546 adcroft 1.1 ENDIF
547 edhill 1.24
548 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
549 adcroft 1.1
550     RETURN
551     END

  ViewVC Help
Powered by ViewVC 1.1.22