/[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.38 - (hide annotations) (download)
Sun May 15 03:04:56 2005 UTC (19 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57h_done
Changes since 1.37: +4 -5 lines
remove "baseTime" (no used) from arg. list of DIFF_BASE_MULTIPLE
and rename it: DIFFERENT_MULTIPLE

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

  ViewVC Help
Powered by ViewVC 1.1.22