/[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.26 - (hide annotations) (download)
Sun Oct 10 06:08:49 2004 UTC (19 years, 7 months ago) by edhill
Branch: MAIN
Changes since 1.25: +1 -4 lines
 o move useMNC and related runtime switches to PARAMS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22