/[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.54 - (hide annotations) (download)
Wed Oct 12 01:52:09 2005 UTC (18 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57v_post, checkpoint57y_pre, checkpoint57w_post, checkpoint57x_post
Changes since 1.53: +144 -85 lines
apply free-slip / no-slip BC on vorticity & strain.

1 jmc 1.54 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.53 2005/09/29 12:19:52 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.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 jmc 1.54 C myThid :: my Thread Id number
51 adcroft 1.1 _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 jmc 1.54 _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 jmc 1.54 _RS hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76     _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79     _RL del2u (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RL del2v (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL dStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL zStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL strain (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL KE (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RL omega3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL vort3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     _RL hDiv (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89     _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90     _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91     _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92     _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     C i,j,k :: Loop counters
94 adcroft 1.1 INTEGER i,j,k
95     C xxxFac - On-off tracer parameters used for switching terms off.
96     _RL ArDudrFac
97 jmc 1.29 c _RL mtFacU
98 adcroft 1.1 _RL ArDvdrFac
99 jmc 1.29 c _RL mtFacV
100 jmc 1.54 _RL sideMaskFac
101 adcroft 1.1 LOGICAL bottomDragTerms
102 jmc 1.15 LOGICAL writeDiag
103 baylor 1.50 LOGICAL harmonic,biharmonic,useVariableViscosity
104 adcroft 1.1
105 edhill 1.25 #ifdef ALLOW_MNC
106     INTEGER offsets(9)
107 edhill 1.53 CHARACTER*(1) pf
108 edhill 1.25 #endif
109    
110 heimbach 1.9 #ifdef ALLOW_AUTODIFF_TAMC
111     C-- only the kDown part of fverU/V is set in this subroutine
112     C-- the kUp is still required
113     C-- In the case of mom_fluxform Kup is set as well
114     C-- (at least in part)
115     fVerU(1,1,kUp) = fVerU(1,1,kUp)
116     fVerV(1,1,kUp) = fVerV(1,1,kUp)
117     #endif
118    
119 jmc 1.38 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
120 adcroft 1.1
121 edhill 1.24 #ifdef ALLOW_MNC
122     IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
123 edhill 1.53 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
124     pf(1:1) = 'D'
125     ELSE
126     pf(1:1) = 'R'
127     ENDIF
128 edhill 1.25 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
129     CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
130 edhill 1.39 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
131 edhill 1.25 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
132 edhill 1.39 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
133 edhill 1.25 ENDIF
134     DO i = 1,9
135     offsets(i) = 0
136     ENDDO
137     offsets(3) = k
138     C write(*,*) 'offsets = ',(offsets(i),i=1,9)
139 edhill 1.24 ENDIF
140     #endif /* ALLOW_MNC */
141    
142 adcroft 1.1 C Initialise intermediate terms
143     DO J=1-OLy,sNy+OLy
144     DO I=1-OLx,sNx+OLx
145 jmc 1.31 vF(i,j) = 0.
146     vrF(i,j) = 0.
147 adcroft 1.1 uCf(i,j) = 0.
148     vCf(i,j) = 0.
149 jmc 1.31 c mT(i,j) = 0.
150 adcroft 1.1 del2u(i,j) = 0.
151     del2v(i,j) = 0.
152     dStar(i,j) = 0.
153     zStar(i,j) = 0.
154 jmc 1.31 guDiss(i,j)= 0.
155     gvDiss(i,j)= 0.
156 adcroft 1.1 vort3(i,j) = 0.
157 jmc 1.31 omega3(i,j)= 0.
158 jmc 1.54 KE(i,j) = 0.
159 baylor 1.50 viscAh_Z(i,j) = 0.
160     viscAh_D(i,j) = 0.
161     viscA4_Z(i,j) = 0.
162     viscA4_D(i,j) = 0.
163    
164 heimbach 1.8 #ifdef ALLOW_AUTODIFF_TAMC
165     strain(i,j) = 0. _d 0
166     tension(i,j) = 0. _d 0
167     #endif
168 adcroft 1.1 ENDDO
169     ENDDO
170    
171     C-- Term by term tracer parmeters
172     C o U momentum equation
173     ArDudrFac = vfFacMom*1.
174 jmc 1.29 c mTFacU = mtFacMom*1.
175 adcroft 1.1 C o V momentum equation
176     ArDvdrFac = vfFacMom*1.
177 jmc 1.29 c mTFacV = mtFacMom*1.
178 adcroft 1.1
179 jmc 1.54 C note: using standard stencil (no mask) results in under-estimating
180     C vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
181     IF ( no_slip_sides ) THEN
182     sideMaskFac = sideDragFactor
183     ELSE
184     sideMaskFac = 0. _d 0
185     ENDIF
186    
187 adcroft 1.1 IF ( no_slip_bottom
188     & .OR. bottomDragQuadratic.NE.0.
189     & .OR. bottomDragLinear.NE.0.) THEN
190     bottomDragTerms=.TRUE.
191     ELSE
192     bottomDragTerms=.FALSE.
193     ENDIF
194    
195     C-- Calculate open water fraction at vorticity points
196     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
197    
198     C Make local copies of horizontal flow field
199     DO j=1-OLy,sNy+OLy
200     DO i=1-OLx,sNx+OLx
201     uFld(i,j) = uVel(i,j,k,bi,bj)
202     vFld(i,j) = vVel(i,j,k,bi,bj)
203     ENDDO
204     ENDDO
205    
206 jmc 1.7 C note (jmc) : Dissipation and Vort3 advection do not necesary
207     C use the same maskZ (and hFacZ) => needs 2 call(s)
208     c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
209    
210 jmc 1.52 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
211 adcroft 1.1
212 jmc 1.54 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
213 adcroft 1.1
214 jmc 1.54 IF (momViscosity) THEN
215     C-- For viscous term, compute horizontal divergence, tension & strain
216     C and mask relative vorticity (free-slip case):
217 adcroft 1.1
218 jmc 1.54 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
219 adcroft 1.1
220 jmc 1.52 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
221    
222     CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
223    
224 jmc 1.54 C- account for no-slip / free-slip BC:
225     DO j=1-Oly,sNy+Oly
226     DO i=1-Olx,sNx+Olx
227     IF ( hFacZ(i,j).EQ.0. ) THEN
228     vort3(i,j) = sideMaskFac*vort3(i,j)
229     strain(i,j) = sideMaskFac*strain(i,j)
230     ENDIF
231     ENDDO
232     ENDDO
233    
234     C-- Calculate Viscosities
235 jmc 1.52 CALL MOM_CALC_VISC(
236 baylor 1.50 I bi,bj,k,
237     O viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
238     O harmonic,biharmonic,useVariableViscosity,
239     I hDiv,vort3,tension,strain,KE,hfacZ,
240     I myThid)
241    
242 adcroft 1.1 C Calculate del^2 u and del^2 v for bi-harmonic term
243 baylor 1.50 IF (biharmonic) THEN
244 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
245     O del2u,del2v,
246     & myThid)
247 jmc 1.48 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
248     CALL MOM_CALC_RELVORT3(bi,bj,k,
249     & del2u,del2v,hFacZ,zStar,myThid)
250 adcroft 1.2 ENDIF
251 baylor 1.47
252 jmc 1.54 C- Strain diagnostics:
253     IF ( writeDiag ) THEN
254     IF (snapshot_mdsio) THEN
255     CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
256     ENDIF
257     #ifdef ALLOW_MNC
258     IF (useMNC .AND. snapshot_mnc) THEN
259     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
260     & offsets, myThid)
261     ENDIF
262     #endif /* ALLOW_MNC */
263     ENDIF
264     #ifdef ALLOW_DIAGNOSTICS
265     IF ( useDiagnostics ) THEN
266     CALL DIAGNOSTICS_FILL(strain, 'Strain ',k,1,2,bi,bj,myThid)
267     ENDIF
268     #endif /* ALLOW_DIAGNOSTICS */
269    
270     C--- Calculate dissipation terms for U and V equations
271 baylor 1.47
272     C in terms of tension and strain
273     IF (useStrainTensionVisc) THEN
274 jmc 1.54 C mask strain as if free-slip since side-drag is computed separately
275     DO j=1-Oly,sNy+Oly
276     DO i=1-Olx,sNx+Olx
277     IF ( hFacZ(i,j).EQ.0. ) strain(i,j) = 0. _d 0
278     ENDDO
279     ENDDO
280 baylor 1.47 CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
281     I hFacZ,
282 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
283     I harmonic,biharmonic,useVariableViscosity,
284 baylor 1.47 O guDiss,gvDiss,
285     I myThid)
286     ELSE
287 adcroft 1.2 C in terms of vorticity and divergence
288 baylor 1.47 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
289     I hFacZ,dStar,zStar,
290 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
291     I harmonic,biharmonic,useVariableViscosity,
292 jmc 1.31 O guDiss,gvDiss,
293 baylor 1.47 & myThid)
294 adcroft 1.3 ENDIF
295 jmc 1.54 C-- if (momViscosity) end of block.
296 adcroft 1.1 ENDIF
297    
298 jmc 1.7 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
299     c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
300    
301 jmc 1.54 C--- Other dissipation terms in Zonal momentum equation
302 adcroft 1.1
303     C-- Vertical flux (fVer is at upper face of "u" cell)
304    
305     C Eddy component of vertical flux (interior component only) -> vrF
306 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
307 jmc 1.44 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
308 adcroft 1.1
309     C Combine fluxes
310 jmc 1.31 DO j=jMin,jMax
311     DO i=iMin,iMax
312     fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
313     ENDDO
314 adcroft 1.1 ENDDO
315    
316 jmc 1.31 C-- Tendency is minus divergence of the fluxes
317     DO j=2-Oly,sNy+Oly-1
318     DO i=2-Olx,sNx+Olx-1
319     guDiss(i,j) = guDiss(i,j)
320 adcroft 1.1 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
321     & *recip_rAw(i,j,bi,bj)
322     & *(
323 jmc 1.42 & fVerU(i,j,kDown) - fVerU(i,j,kUp)
324     & )*rkSign
325 jmc 1.31 ENDDO
326 adcroft 1.1 ENDDO
327 jmc 1.31 ENDIF
328 adcroft 1.1
329     C-- No-slip and drag BCs appear as body forces in cell abutting topography
330     IF (momViscosity.AND.no_slip_sides) THEN
331     C- No-slip BCs impose a drag at walls...
332 baylor 1.50 CALL MOM_U_SIDEDRAG(
333     I bi,bj,k,
334     I uFld, del2u, hFacZ,
335     I viscAh_Z,viscA4_Z,
336     I harmonic,biharmonic,useVariableViscosity,
337     O vF,
338     I myThid)
339 adcroft 1.1 DO j=jMin,jMax
340     DO i=iMin,iMax
341 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
342 adcroft 1.1 ENDDO
343     ENDDO
344     ENDIF
345     C- No-slip BCs impose a drag at bottom
346     IF (momViscosity.AND.bottomDragTerms) THEN
347     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
348     DO j=jMin,jMax
349     DO i=iMin,iMax
350 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
351 adcroft 1.1 ENDDO
352     ENDDO
353     ENDIF
354    
355 jmc 1.54 C--- Other dissipation terms in Meridional momentum equation
356 adcroft 1.1
357     C-- Vertical flux (fVer is at upper face of "v" cell)
358    
359     C Eddy component of vertical flux (interior component only) -> vrF
360 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
361 jmc 1.44 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
362 adcroft 1.1
363     C Combine fluxes -> fVerV
364 jmc 1.31 DO j=jMin,jMax
365     DO i=iMin,iMax
366     fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
367     ENDDO
368 adcroft 1.1 ENDDO
369    
370 jmc 1.31 C-- Tendency is minus divergence of the fluxes
371     DO j=jMin,jMax
372     DO i=iMin,iMax
373     gvDiss(i,j) = gvDiss(i,j)
374 adcroft 1.1 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
375     & *recip_rAs(i,j,bi,bj)
376     & *(
377 jmc 1.42 & fVerV(i,j,kDown) - fVerV(i,j,kUp)
378     & )*rkSign
379 jmc 1.31 ENDDO
380 adcroft 1.1 ENDDO
381 jmc 1.31 ENDIF
382 adcroft 1.1
383     C-- No-slip and drag BCs appear as body forces in cell abutting topography
384     IF (momViscosity.AND.no_slip_sides) THEN
385     C- No-slip BCs impose a drag at walls...
386 baylor 1.50 CALL MOM_V_SIDEDRAG(
387     I bi,bj,k,
388     I vFld, del2v, hFacZ,
389     I viscAh_Z,viscA4_Z,
390     I harmonic,biharmonic,useVariableViscosity,
391     O vF,
392     I myThid)
393 adcroft 1.1 DO j=jMin,jMax
394     DO i=iMin,iMax
395 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
396 adcroft 1.1 ENDDO
397     ENDDO
398     ENDIF
399     C- No-slip BCs impose a drag at bottom
400     IF (momViscosity.AND.bottomDragTerms) THEN
401     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
402     DO j=jMin,jMax
403     DO i=iMin,iMax
404 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
405 adcroft 1.1 ENDDO
406     ENDDO
407     ENDIF
408    
409 jmc 1.54 C- Vorticity diagnostics:
410     IF ( writeDiag ) THEN
411     IF (snapshot_mdsio) THEN
412     CALL WRITE_LOCAL_RL('Z3','I10',1,vort3, bi,bj,k,myIter,myThid)
413     ENDIF
414     #ifdef ALLOW_MNC
415     IF (useMNC .AND. snapshot_mnc) THEN
416     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
417     & offsets, myThid)
418     ENDIF
419     #endif /* ALLOW_MNC */
420     ENDIF
421     #ifdef ALLOW_DIAGNOSTICS
422     IF ( useDiagnostics ) THEN
423     CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid)
424     ENDIF
425     #endif /* ALLOW_DIAGNOSTICS */
426    
427     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
428    
429     C--- Prepare for Advection & Coriolis terms:
430     C- Mask relative vorticity and calculate absolute vorticity
431     DO j=1-Oly,sNy+Oly
432     DO i=1-Olx,sNx+Olx
433     IF ( hFacZ(i,j).EQ.0. ) vort3(i,j) = 0.
434     ENDDO
435     ENDDO
436     IF (useAbsVorticity)
437     & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
438 adcroft 1.1
439 jmc 1.5 C-- Horizontal Coriolis terms
440 jmc 1.37 c IF (useCoriolis .AND. .NOT.useCDscheme
441     c & .AND. .NOT. useAbsVorticity) THEN
442     C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
443 jmc 1.46 IF ( useCoriolis .AND.
444 jmc 1.37 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
445     & ) THEN
446     IF (useAbsVorticity) THEN
447     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
448     & uCf,myThid)
449     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
450     & vCf,myThid)
451     ELSE
452     CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
453     & uCf,vCf,myThid)
454     ENDIF
455 jmc 1.5 DO j=jMin,jMax
456     DO i=iMin,iMax
457 jmc 1.43 gU(i,j,k,bi,bj) = uCf(i,j)
458     gV(i,j,k,bi,bj) = vCf(i,j)
459 jmc 1.5 ENDDO
460 adcroft 1.1 ENDDO
461 jmc 1.46
462 jmc 1.15 IF ( writeDiag ) THEN
463 edhill 1.24 IF (snapshot_mdsio) THEN
464     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
465     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
466     ENDIF
467     #ifdef ALLOW_MNC
468     IF (useMNC .AND. snapshot_mnc) THEN
469 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
470 edhill 1.25 & offsets, myThid)
471 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
472 edhill 1.25 & offsets, myThid)
473 edhill 1.24 ENDIF
474     #endif /* ALLOW_MNC */
475 jmc 1.15 ENDIF
476 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
477     IF ( useDiagnostics ) THEN
478     CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
479     CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
480     ENDIF
481     #endif /* ALLOW_DIAGNOSTICS */
482    
483 jmc 1.31 ELSE
484     DO j=jMin,jMax
485     DO i=iMin,iMax
486 jmc 1.43 gU(i,j,k,bi,bj) = 0. _d 0
487     gV(i,j,k,bi,bj) = 0. _d 0
488 jmc 1.31 ENDDO
489     ENDDO
490 jmc 1.5 ENDIF
491 adcroft 1.1
492 jmc 1.5 IF (momAdvection) THEN
493 jmc 1.41 C-- Horizontal advection of relative (or absolute) vorticity
494     IF (highOrderVorticity.AND.useAbsVorticity) THEN
495     CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
496 adcroft 1.20 & uCf,myThid)
497 jmc 1.40 ELSEIF (highOrderVorticity) THEN
498 jmc 1.41 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
499     & uCf,myThid)
500     ELSEIF (useAbsVorticity) THEN
501     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
502 jmc 1.40 & uCf,myThid)
503 adcroft 1.20 ELSE
504 jmc 1.41 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
505 adcroft 1.20 & uCf,myThid)
506     ENDIF
507 jmc 1.5 DO j=jMin,jMax
508     DO i=iMin,iMax
509     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
510     ENDDO
511 adcroft 1.1 ENDDO
512 jmc 1.41 IF (highOrderVorticity.AND.useAbsVorticity) THEN
513     CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
514 adcroft 1.20 & vCf,myThid)
515 jmc 1.40 ELSEIF (highOrderVorticity) THEN
516 jmc 1.41 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
517     & vCf,myThid)
518     ELSEIF (useAbsVorticity) THEN
519     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
520 jmc 1.40 & vCf,myThid)
521 adcroft 1.20 ELSE
522 jmc 1.41 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
523 adcroft 1.20 & vCf,myThid)
524     ENDIF
525 jmc 1.5 DO j=jMin,jMax
526     DO i=iMin,iMax
527     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
528     ENDDO
529 adcroft 1.1 ENDDO
530    
531 jmc 1.15 IF ( writeDiag ) THEN
532 edhill 1.24 IF (snapshot_mdsio) THEN
533     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
534     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
535     ENDIF
536     #ifdef ALLOW_MNC
537     IF (useMNC .AND. snapshot_mnc) THEN
538 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
539 edhill 1.25 & offsets, myThid)
540 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
541 edhill 1.25 & offsets, myThid)
542 edhill 1.24 ENDIF
543     #endif /* ALLOW_MNC */
544 jmc 1.15 ENDIF
545 edhill 1.24
546 jmc 1.7 #ifdef ALLOW_TIMEAVE
547     IF (taveFreq.GT.0.) THEN
548     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
549     & Nr, k, bi, bj, myThid)
550     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
551     & Nr, k, bi, bj, myThid)
552     ENDIF
553 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
554 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
555     IF ( useDiagnostics ) THEN
556     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
557     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
558     ENDIF
559     #endif /* ALLOW_DIAGNOSTICS */
560 jmc 1.7
561 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
562 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
563     CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
564     DO j=jMin,jMax
565     DO i=iMin,iMax
566     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
567     ENDDO
568 jmc 1.5 ENDDO
569 jmc 1.12 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
570     DO j=jMin,jMax
571     DO i=iMin,iMax
572     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
573     ENDDO
574 jmc 1.5 ENDDO
575 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
576     IF ( useDiagnostics ) THEN
577     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
578     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
579     ENDIF
580     #endif /* ALLOW_DIAGNOSTICS */
581 jmc 1.12 ENDIF
582 adcroft 1.1
583     C-- Bernoulli term
584 jmc 1.5 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
585     DO j=jMin,jMax
586     DO i=iMin,iMax
587     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
588     ENDDO
589     ENDDO
590     CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
591     DO j=jMin,jMax
592     DO i=iMin,iMax
593     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
594     ENDDO
595 adcroft 1.1 ENDDO
596 jmc 1.15 IF ( writeDiag ) THEN
597 edhill 1.24 IF (snapshot_mdsio) THEN
598     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
599     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
600     ENDIF
601     #ifdef ALLOW_MNC
602     IF (useMNC .AND. snapshot_mnc) THEN
603 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
604 edhill 1.25 & offsets, myThid)
605 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
606 edhill 1.25 & offsets, myThid)
607 jmc 1.54 ENDIF
608 edhill 1.24 #endif /* ALLOW_MNC */
609 jmc 1.15 ENDIF
610    
611 jmc 1.5 C-- end if momAdvection
612     ENDIF
613    
614 jmc 1.54 C-- Metric terms for curvilinear grid systems
615     c IF (usingSphericalPolarMTerms) THEN
616     C o Spherical polar grid metric terms
617     c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
618     c DO j=jMin,jMax
619     c DO i=iMin,iMax
620     c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
621     c ENDDO
622     c ENDDO
623     c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
624     c DO j=jMin,jMax
625     c DO i=iMin,iMax
626     c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
627     c ENDDO
628     c ENDDO
629     c ENDIF
630    
631 jmc 1.5 C-- Set du/dt & dv/dt on boundaries to zero
632 adcroft 1.1 DO j=jMin,jMax
633     DO i=iMin,iMax
634 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
635     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
636 adcroft 1.1 ENDDO
637     ENDDO
638 jmc 1.5
639 jmc 1.22 #ifdef ALLOW_DEBUG
640     IF ( debugLevel .GE. debLevB
641     & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
642     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
643     & .AND. useCubedSphereExchange ) THEN
644 jmc 1.23 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
645 jmc 1.31 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
646 jmc 1.22 ENDIF
647     #endif /* ALLOW_DEBUG */
648 adcroft 1.2
649 jmc 1.15 IF ( writeDiag ) THEN
650 edhill 1.24 IF (snapshot_mdsio) THEN
651 jmc 1.54 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
652     CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
653     CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
654     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
655     CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
656     CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
657 edhill 1.24 ENDIF
658     #ifdef ALLOW_MNC
659     IF (useMNC .AND. snapshot_mnc) THEN
660 jmc 1.54 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
661     & offsets, myThid)
662     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
663     & offsets, myThid)
664     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
665 edhill 1.25 & offsets, myThid)
666 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
667 edhill 1.25 & offsets, myThid)
668 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
669 edhill 1.25 & offsets, myThid)
670 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
671 edhill 1.25 & offsets, myThid)
672 edhill 1.24 ENDIF
673     #endif /* ALLOW_MNC */
674 adcroft 1.1 ENDIF
675 jmc 1.41
676 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
677     IF ( useDiagnostics ) THEN
678 jmc 1.54 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
679 jmc 1.46 IF (momViscosity) THEN
680 jmc 1.54 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
681 jmc 1.52 CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
682 jmc 1.54 CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
683     CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
684 jmc 1.46 ENDIF
685 jmc 1.54 CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),
686     & 'Um_Advec',k,1,2,bi,bj,myThid)
687     CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),
688     & 'Vm_Advec',k,1,2,bi,bj,myThid)
689 jmc 1.46 ENDIF
690     #endif /* ALLOW_DIAGNOSTICS */
691    
692 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
693 adcroft 1.1
694     RETURN
695     END

  ViewVC Help
Powered by ViewVC 1.1.22