/[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.44 - (hide annotations) (download)
Sat Jul 30 23:52:10 2005 UTC (18 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57p_post, checkpoint57q_post
Changes since 1.43: +3 -7 lines
does not rely on vrF initialization (safer).

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

  ViewVC Help
Powered by ViewVC 1.1.22