/[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.32 - (hide annotations) (download)
Sat Dec 4 05:59:50 2004 UTC (19 years, 5 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint57, checkpoint57a_post, checkpoint57a_pre
Changes since 1.31: +3 -3 lines
Added CPP option MINIMAL_TAVE_OUTPUT for minimal time-averaged output:
S, T, U, V, W, ETA, and phiHydLow.

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

  ViewVC Help
Powered by ViewVC 1.1.22