/[MITgcm]/MITgcm/pkg/mom_vecinv/mom_vecinv.F
ViewVC logotype

Contents of /MITgcm/pkg/mom_vecinv/mom_vecinv.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.72 - (show annotations) (download)
Tue Feb 11 20:24:06 2014 UTC (10 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64u
Changes since 1.71: +62 -84 lines
- remove unused arguments from S/R MOM_VI_HDISSIP & MOM_HDISSIP (as it used
  to be before mom_calc_visc.F)
- skip the call to MOM_CALC_TENSION & MOM_CALC_STRAIN if not needed.
- add 2nd copy of vort3 & strain that knows about lateral BC (free/no slip):
  this is simpler for the adjoint and for diagnostics.

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.71 2014/02/09 18:46:46 jmc Exp $
2 C $Name: $
3
4 #include "MOM_VECINV_OPTIONS.h"
5 #ifdef ALLOW_MOM_COMMON
6 # include "MOM_COMMON_OPTIONS.h"
7 #endif
8
9 SUBROUTINE MOM_VECINV(
10 I bi,bj,k,iMin,iMax,jMin,jMax,
11 I KappaRU, KappaRV,
12 I fVerUkm, fVerVkm,
13 O fVerUkp, fVerVkp,
14 O guDiss, gvDiss,
15 I myTime, myIter, myThid )
16 C *==========================================================*
17 C | S/R MOM_VECINV |
18 C | o Form the right hand-side of the momentum equation. |
19 C *==========================================================*
20 C | Terms are evaluated one layer at a time working from |
21 C | the bottom to the top. The vertically integrated |
22 C | barotropic flow tendency term is evluated by summing the |
23 C | tendencies. |
24 C | Notes: |
25 C | We have not sorted out an entirely satisfactory formula |
26 C | for the diffusion equation bc with lopping. The present |
27 C | form produces a diffusive flux that does not scale with |
28 C | open-area. Need to do something to solidfy this and to |
29 C | deal "properly" with thin walls. |
30 C *==========================================================*
31 IMPLICIT NONE
32
33 C == Global variables ==
34 #include "SIZE.h"
35 #include "EEPARAMS.h"
36 #include "PARAMS.h"
37 #include "GRID.h"
38 #include "SURFACE.h"
39 #include "DYNVARS.h"
40 #ifdef ALLOW_MOM_COMMON
41 # include "MOM_VISC.h"
42 #endif
43 #ifdef ALLOW_TIMEAVE
44 # include "TIMEAVE_STATV.h"
45 #endif
46 #ifdef ALLOW_MNC
47 # include "MNC_PARAMS.h"
48 #endif
49 #ifdef ALLOW_AUTODIFF_TAMC
50 # include "tamc.h"
51 # include "tamc_keys.h"
52 #endif
53
54 C == Routine arguments ==
55 C bi,bj :: current tile indices
56 C k :: current vertical level
57 C iMin,iMax,jMin,jMax :: loop ranges
58 C fVerU :: Flux of momentum in the vertical direction, out of the upper
59 C fVerV :: face of a cell K ( flux into the cell above ).
60 C fVerUkm :: vertical viscous flux of U, interface above (k-1/2)
61 C fVerVkm :: vertical viscous flux of V, interface above (k-1/2)
62 C fVerUkp :: vertical viscous flux of U, interface below (k+1/2)
63 C fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
64
65 C guDiss :: dissipation tendency (all explicit terms), u component
66 C gvDiss :: dissipation tendency (all explicit terms), v component
67 C myTime :: current time
68 C myIter :: current time-step number
69 C myThid :: my Thread Id number
70 INTEGER bi,bj,k
71 INTEGER iMin,iMax,jMin,jMax
72 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
73 _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
74 _RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL myTime
81 INTEGER myIter
82 INTEGER myThid
83
84 #ifdef ALLOW_MOM_VECINV
85
86 C == Functions ==
87 LOGICAL DIFFERENT_MULTIPLE
88 EXTERNAL DIFFERENT_MULTIPLE
89
90 C == Local variables ==
91 C strainBC :: same as strain but account for no-slip BC
92 C vort3BC :: same as vort3 but account for no-slip BC
93 _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RS hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RS h0FacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL del2u (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL del2v (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL dStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL zStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL strain (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 _RL strainBC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 _RL KE (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110 _RL omega3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111 _RL vort3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112 _RL vort3BC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113 _RL hDiv (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114 _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115 _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 C i,j :: Loop counters
119 INTEGER i,j
120 C xxxFac :: On-off tracer parameters used for switching terms off.
121 _RL ArDudrFac
122 _RL ArDvdrFac
123 _RL sideMaskFac
124 LOGICAL bottomDragTerms
125 LOGICAL writeDiag
126 #ifdef ALLOW_AUTODIFF_TAMC
127 INTEGER imomkey
128 #endif
129
130 #ifdef ALLOW_MNC
131 INTEGER offsets(9)
132 CHARACTER*(1) pf
133 #endif
134
135 #ifdef ALLOW_AUTODIFF_TAMC
136 C-- only the kDown part of fverU/V is set in this subroutine
137 C-- the kUp is still required
138 C-- In the case of mom_fluxform Kup is set as well
139 C-- (at least in part)
140 fVerUkm(1,1) = fVerUkm(1,1)
141 fVerVkm(1,1) = fVerVkm(1,1)
142 #endif
143
144 #ifdef ALLOW_AUTODIFF_TAMC
145 act0 = k - 1
146 max0 = Nr
147 act1 = bi - myBxLo(myThid)
148 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
149 act2 = bj - myByLo(myThid)
150 max2 = myByHi(myThid) - myByLo(myThid) + 1
151 act3 = myThid - 1
152 max3 = nTx*nTy
153 act4 = ikey_dynamics - 1
154 imomkey = (act0 + 1)
155 & + act1*max0
156 & + act2*max0*max1
157 & + act3*max0*max1*max2
158 & + act4*max0*max1*max2*max3
159 #endif /* ALLOW_AUTODIFF_TAMC */
160
161 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
162
163 #ifdef ALLOW_MNC
164 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
165 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
166 pf(1:1) = 'D'
167 ELSE
168 pf(1:1) = 'R'
169 ENDIF
170 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
171 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
172 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
173 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
174 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
175 ENDIF
176 DO i = 1,9
177 offsets(i) = 0
178 ENDDO
179 offsets(3) = k
180 c write(*,*) 'offsets = ',(offsets(i),i=1,9)
181 ENDIF
182 #endif /* ALLOW_MNC */
183
184 C-- Initialise intermediate terms
185 DO j=1-OLy,sNy+OLy
186 DO i=1-OLx,sNx+OLx
187 vF(i,j) = 0.
188 vrF(i,j) = 0.
189 uCf(i,j) = 0.
190 vCf(i,j) = 0.
191 del2u(i,j) = 0.
192 del2v(i,j) = 0.
193 dStar(i,j) = 0.
194 zStar(i,j) = 0.
195 guDiss(i,j)= 0.
196 gvDiss(i,j)= 0.
197 vort3(i,j) = 0.
198 omega3(i,j)= 0.
199 KE(i,j) = 0.
200 C- need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
201 hDiv(i,j) = 0.
202 c viscAh_Z(i,j) = 0.
203 c viscAh_D(i,j) = 0.
204 c viscA4_Z(i,j) = 0.
205 c viscA4_D(i,j) = 0.
206 strain(i,j) = 0. _d 0
207 strainBC(i,j)= 0. _d 0
208 tension(i,j) = 0. _d 0
209 #ifdef ALLOW_AUTODIFF_TAMC
210 hFacZ(i,j) = 0. _d 0
211 #endif
212 ENDDO
213 ENDDO
214
215 C-- Term by term tracer parmeters
216 C o U momentum equation
217 ArDudrFac = vfFacMom*1.
218 C o V momentum equation
219 ArDvdrFac = vfFacMom*1.
220
221 C note: using standard stencil (no mask) results in under-estimating
222 C vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
223 IF ( no_slip_sides ) THEN
224 sideMaskFac = sideDragFactor
225 ELSE
226 sideMaskFac = 0. _d 0
227 ENDIF
228
229 IF ( no_slip_bottom
230 & .OR. bottomDragQuadratic.NE.0.
231 & .OR. bottomDragLinear.NE.0.) THEN
232 bottomDragTerms=.TRUE.
233 ELSE
234 bottomDragTerms=.FALSE.
235 ENDIF
236
237 C-- Calculate open water fraction at vorticity points
238 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
239
240 C Make local copies of horizontal flow field
241 DO j=1-OLy,sNy+OLy
242 DO i=1-OLx,sNx+OLx
243 uFld(i,j) = uVel(i,j,k,bi,bj)
244 vFld(i,j) = vVel(i,j,k,bi,bj)
245 ENDDO
246 ENDDO
247
248 C note (jmc) : Dissipation and Vort3 advection do not necesary
249 C use the same maskZ (and hFacZ) => needs 2 call(s)
250 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
251
252 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
253
254 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
255
256 C- mask vort3 and account for no-slip / free-slip BC in vort3BC:
257 DO j=1-OLy,sNy+OLy
258 DO i=1-OLx,sNx+OLx
259 vort3BC(i,j) = vort3(i,j)
260 IF ( hFacZ(i,j).EQ.zeroRS ) THEN
261 vort3BC(i,j) = sideMaskFac*vort3BC(i,j)
262 vort3(i,j) = 0.
263 ENDIF
264 ENDDO
265 ENDDO
266
267 IF (momViscosity) THEN
268 C-- For viscous term, compute horizontal divergence, tension & strain
269 C and mask relative vorticity (free-slip case):
270
271 DO j=1-OLy,sNy+OLy
272 DO i=1-OLx,sNx+OLx
273 h0FacZ(i,j) = hFacZ(i,j)
274 ENDDO
275 ENDDO
276 #ifdef NONLIN_FRSURF
277 IF ( no_slip_sides .AND. nonlinFreeSurf.GT.0 ) THEN
278 DO j=2-OLy,sNy+OLy
279 DO i=2-OLx,sNx+OLx
280 h0FacZ(i,j) = MIN(
281 & MIN( h0FacW(i,j,k,bi,bj), h0FacW(i,j-1,k,bi,bj) ),
282 & MIN( h0FacS(i,j,k,bi,bj), h0FacS(i-1,j,k,bi,bj) ) )
283 ENDDO
284 ENDDO
285 ENDIF
286 #endif /* NONLIN_FRSURF */
287
288 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
289
290 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
291 CALL MOM_CALC_TENSION( bi,bj,k,uFld,vFld,tension,myThid )
292 CALL MOM_CALC_STRAIN( bi,bj,k,uFld,vFld,hFacZ,strain,myThid )
293 C- mask strain and account for no-slip / free-slip BC in strainBC:
294 DO j=1-OLy,sNy+OLy
295 DO i=1-OLx,sNx+OLx
296 strainBC(i,j) = strain(i,j)
297 IF ( hFacZ(i,j).EQ.zeroRS ) THEN
298 strainBC(i,j) = sideMaskFac*strainBC(i,j)
299 strain(i,j) = 0.
300 ENDIF
301 ENDDO
302 ENDDO
303 ENDIF
304
305 C-- Calculate Lateral Viscosities
306 DO j=1-OLy,sNy+OLy
307 DO i=1-OLx,sNx+OLx
308 viscAh_D(i,j) = viscAhD
309 viscAh_Z(i,j) = viscAhZ
310 viscA4_D(i,j) = viscA4D
311 viscA4_Z(i,j) = viscA4Z
312 ENDDO
313 ENDDO
314 IF ( useVariableVisc ) THEN
315 C- uses vort3BC & strainBC which account for no-slip / free-slip BC
316 CALL MOM_CALC_VISC( bi, bj, k,
317 O viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
318 I hDiv, vort3BC, tension, strainBC, KE, hfacZ,
319 I myThid )
320 ENDIF
321
322 C Calculate del^2 u and del^2 v for bi-harmonic term
323 IF (useBiharmonicVisc) THEN
324 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
325 O del2u,del2v,
326 I myThid)
327 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
328 CALL MOM_CALC_RELVORT3(bi,bj,k,
329 & del2u,del2v,hFacZ,zStar,myThid)
330 ENDIF
331
332 C--- Calculate dissipation terms for U and V equations
333
334 C- in terms of tension and strain
335 IF (useStrainTensionVisc) THEN
336 C use masked strain as if free-slip since side-drag is computed separately
337 CALL MOM_HDISSIP( bi, bj, k,
338 I tension, strain, hFacZ,
339 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
340 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
341 O guDiss, gvDiss,
342 I myThid )
343 ELSE
344 C- in terms of vorticity and divergence
345 CALL MOM_VI_HDISSIP( bi, bj, k,
346 I hDiv, vort3, dStar, zStar, hFacZ,
347 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
348 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
349 O guDiss, gvDiss,
350 I myThid )
351 ENDIF
352
353 C--- Other dissipation terms in Zonal momentum equation
354
355 C-- Vertical flux (fVer is at upper face of "u" cell)
356 C Eddy component of vertical flux (interior component only) -> vrF
357 IF ( .NOT.implicitViscosity ) THEN
358 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
359 C Combine fluxes
360 DO j=jMin,jMax
361 DO i=iMin,iMax
362 fVerUkp(i,j) = ArDudrFac*vrF(i,j)
363 ENDDO
364 ENDDO
365 C-- Tendency is minus divergence of the fluxes
366 DO j=jMin,jMax
367 DO i=iMin,iMax
368 guDiss(i,j) = guDiss(i,j)
369 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
370 & *recip_rAw(i,j,bi,bj)
371 & *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
372 ENDDO
373 ENDDO
374 ENDIF
375
376 C-- No-slip and drag BCs appear as body forces in cell abutting topography
377 IF ( no_slip_sides ) THEN
378 C- No-slip BCs impose a drag at walls...
379 CALL MOM_U_SIDEDRAG( bi, bj, k,
380 I uFld, del2u, h0FacZ,
381 I viscAh_Z, viscA4_Z,
382 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
383 O vF,
384 I myThid )
385 DO j=jMin,jMax
386 DO i=iMin,iMax
387 guDiss(i,j) = guDiss(i,j)+vF(i,j)
388 ENDDO
389 ENDDO
390 ENDIF
391
392 C- No-slip BCs impose a drag at bottom
393 IF ( bottomDragTerms ) THEN
394 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
395 DO j=jMin,jMax
396 DO i=iMin,iMax
397 guDiss(i,j) = guDiss(i,j)+vF(i,j)
398 ENDDO
399 ENDDO
400 ENDIF
401 #ifdef ALLOW_SHELFICE
402 IF ( useShelfIce.AND.bottomDragTerms ) THEN
403 CALL SHELFICE_U_DRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
404 DO j=jMin,jMax
405 DO i=iMin,iMax
406 guDiss(i,j) = guDiss(i,j) + vF(i,j)
407 ENDDO
408 ENDDO
409 ENDIF
410 #endif /* ALLOW_SHELFICE */
411
412 C--- Other dissipation terms in Meridional momentum equation
413
414 C-- Vertical flux (fVer is at upper face of "v" cell)
415 C Eddy component of vertical flux (interior component only) -> vrF
416 IF ( .NOT.implicitViscosity ) THEN
417 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
418 C Combine fluxes -> fVerV
419 DO j=jMin,jMax
420 DO i=iMin,iMax
421 fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
422 ENDDO
423 ENDDO
424 C-- Tendency is minus divergence of the fluxes
425 DO j=jMin,jMax
426 DO i=iMin,iMax
427 gvDiss(i,j) = gvDiss(i,j)
428 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
429 & *recip_rAs(i,j,bi,bj)
430 & *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
431 ENDDO
432 ENDDO
433 ENDIF
434
435 C-- No-slip and drag BCs appear as body forces in cell abutting topography
436 IF ( no_slip_sides ) THEN
437 C- No-slip BCs impose a drag at walls...
438 CALL MOM_V_SIDEDRAG( bi, bj, k,
439 I vFld, del2v, h0FacZ,
440 I viscAh_Z, viscA4_Z,
441 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
442 O vF,
443 I myThid )
444 DO j=jMin,jMax
445 DO i=iMin,iMax
446 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
447 ENDDO
448 ENDDO
449 ENDIF
450
451 C- No-slip BCs impose a drag at bottom
452 IF ( bottomDragTerms ) THEN
453 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
454 DO j=jMin,jMax
455 DO i=iMin,iMax
456 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
457 ENDDO
458 ENDDO
459 ENDIF
460 #ifdef ALLOW_SHELFICE
461 IF (useShelfIce.AND.bottomDragTerms ) THEN
462 CALL SHELFICE_V_DRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
463 DO j=jMin,jMax
464 DO i=iMin,iMax
465 gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
466 ENDDO
467 ENDDO
468 ENDIF
469 #endif /* ALLOW_SHELFICE */
470
471 C-- if (momViscosity) end of block.
472 ENDIF
473
474 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
475 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
476
477 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
478
479 C--- Prepare for Advection & Coriolis terms:
480 C- calculate absolute vorticity
481 IF (useAbsVorticity)
482 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
483
484 C-- Horizontal Coriolis terms
485 c IF (useCoriolis .AND. .NOT.useCDscheme
486 c & .AND. .NOT. useAbsVorticity) THEN
487 C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
488 IF ( useCoriolis .AND.
489 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
490 & ) THEN
491 IF (useAbsVorticity) THEN
492 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
493 & uCf,myThid)
494 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
495 & vCf,myThid)
496 ELSE
497 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
498 & uCf,vCf,myThid)
499 ENDIF
500 DO j=jMin,jMax
501 DO i=iMin,iMax
502 gU(i,j,k,bi,bj) = uCf(i,j)
503 gV(i,j,k,bi,bj) = vCf(i,j)
504 ENDDO
505 ENDDO
506 IF ( writeDiag ) THEN
507 IF (snapshot_mdsio) THEN
508 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
509 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
510 ENDIF
511 #ifdef ALLOW_MNC
512 IF (useMNC .AND. snapshot_mnc) THEN
513 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
514 & offsets, myThid)
515 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
516 & offsets, myThid)
517 ENDIF
518 #endif /* ALLOW_MNC */
519 ENDIF
520 #ifdef ALLOW_DIAGNOSTICS
521 IF ( useDiagnostics ) THEN
522 CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
523 CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
524 ENDIF
525 #endif /* ALLOW_DIAGNOSTICS */
526 ELSE
527 DO j=jMin,jMax
528 DO i=iMin,iMax
529 gU(i,j,k,bi,bj) = 0. _d 0
530 gV(i,j,k,bi,bj) = 0. _d 0
531 ENDDO
532 ENDDO
533 ENDIF
534
535 IF (momAdvection) THEN
536 C-- Horizontal advection of relative (or absolute) vorticity
537 IF ( (highOrderVorticity.OR.upwindVorticity)
538 & .AND.useAbsVorticity ) THEN
539 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
540 & uCf,myThid)
541 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
542 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
543 & uCf,myThid)
544 ELSEIF ( useAbsVorticity ) THEN
545 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
546 & uCf,myThid)
547 ELSE
548 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
549 & uCf,myThid)
550 ENDIF
551 DO j=jMin,jMax
552 DO i=iMin,iMax
553 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
554 ENDDO
555 ENDDO
556 IF ( (highOrderVorticity.OR.upwindVorticity)
557 & .AND.useAbsVorticity ) THEN
558 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
559 & vCf,myThid)
560 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
561 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
562 & vCf,myThid)
563 ELSEIF ( useAbsVorticity ) THEN
564 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
565 & vCf,myThid)
566 ELSE
567 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
568 & vCf,myThid)
569 ENDIF
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 ENDDO
575
576 IF ( writeDiag ) THEN
577 IF (snapshot_mdsio) THEN
578 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
579 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
580 ENDIF
581 #ifdef ALLOW_MNC
582 IF (useMNC .AND. snapshot_mnc) THEN
583 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
584 & offsets, myThid)
585 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
586 & offsets, myThid)
587 ENDIF
588 #endif /* ALLOW_MNC */
589 ENDIF
590
591 #ifdef ALLOW_TIMEAVE
592 IF (taveFreq.GT.0.) THEN
593 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
594 & Nr, k, bi, bj, myThid)
595 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
596 & Nr, k, bi, bj, myThid)
597 ENDIF
598 #endif /* ALLOW_TIMEAVE */
599 #ifdef ALLOW_DIAGNOSTICS
600 IF ( useDiagnostics ) THEN
601 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
602 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
603 ENDIF
604 #endif /* ALLOW_DIAGNOSTICS */
605
606 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
607 IF ( .NOT. momImplVertAdv ) THEN
608 CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
609 DO j=jMin,jMax
610 DO i=iMin,iMax
611 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
612 ENDDO
613 ENDDO
614 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
615 DO j=jMin,jMax
616 DO i=iMin,iMax
617 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
618 ENDDO
619 ENDDO
620 #ifdef ALLOW_DIAGNOSTICS
621 IF ( useDiagnostics ) THEN
622 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
623 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
624 ENDIF
625 #endif /* ALLOW_DIAGNOSTICS */
626 ENDIF
627
628 C-- Bernoulli term
629 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
630 DO j=jMin,jMax
631 DO i=iMin,iMax
632 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
633 ENDDO
634 ENDDO
635 CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
636 DO j=jMin,jMax
637 DO i=iMin,iMax
638 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
639 ENDDO
640 ENDDO
641 IF ( writeDiag ) THEN
642 IF (snapshot_mdsio) THEN
643 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
644 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
645 ENDIF
646 #ifdef ALLOW_MNC
647 IF (useMNC .AND. snapshot_mnc) THEN
648 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
649 & offsets, myThid)
650 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
651 & offsets, myThid)
652 ENDIF
653 #endif /* ALLOW_MNC */
654 ENDIF
655
656 C-- end if momAdvection
657 ENDIF
658
659 C-- 3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
660 IF ( use3dCoriolis ) THEN
661 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
662 DO j=jMin,jMax
663 DO i=iMin,iMax
664 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
665 ENDDO
666 ENDDO
667 IF ( usingCurvilinearGrid ) THEN
668 C- presently, non zero angleSinC array only supported with Curvilinear-Grid
669 CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
670 DO j=jMin,jMax
671 DO i=iMin,iMax
672 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
673 ENDDO
674 ENDDO
675 ENDIF
676 ENDIF
677
678 C-- Non-Hydrostatic (spherical) metric terms
679 IF ( useNHMTerms ) THEN
680 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
681 DO j=jMin,jMax
682 DO i=iMin,iMax
683 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
684 ENDDO
685 ENDDO
686 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
687 DO j=jMin,jMax
688 DO i=iMin,iMax
689 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
690 ENDDO
691 ENDDO
692 ENDIF
693
694 C-- Set du/dt & dv/dt on boundaries to zero
695 DO j=jMin,jMax
696 DO i=iMin,iMax
697 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
698 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
699 ENDDO
700 ENDDO
701
702 #ifdef ALLOW_DEBUG
703 IF ( debugLevel .GE. debLevC
704 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
705 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
706 & .AND. useCubedSphereExchange ) THEN
707 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
708 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
709 ENDIF
710 #endif /* ALLOW_DEBUG */
711
712 IF ( writeDiag ) THEN
713 IF (useBiharmonicVisc) THEN
714 CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
715 & bi,bj,k, myIter, myThid )
716 CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
717 & bi,bj,k, myIter, myThid )
718 CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
719 & bi,bj,k, myIter, myThid )
720 CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
721 & bi,bj,k, myIter, myThid )
722 ENDIF
723 IF (snapshot_mdsio) THEN
724 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
725 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
726 CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
727 CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
728 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
729 CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
730 & bi,bj,k, myIter, myThid )
731 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
732 CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
733 ENDIF
734 #ifdef ALLOW_MNC
735 IF (useMNC .AND. snapshot_mnc) THEN
736 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
737 & offsets, myThid)
738 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
739 & offsets, myThid)
740 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
741 & offsets, myThid)
742 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
743 & offsets, myThid)
744 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
745 & offsets, myThid)
746 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
747 & offsets, myThid)
748 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
749 & offsets, myThid)
750 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
751 & offsets, myThid)
752 ENDIF
753 #endif /* ALLOW_MNC */
754 ENDIF
755
756 #ifdef ALLOW_DIAGNOSTICS
757 IF ( useDiagnostics ) THEN
758 CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
759 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
760 IF (momViscosity) THEN
761 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
762 CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
763 CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
764 ENDIF
765 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
766 CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
767 CALL DIAGNOSTICS_FILL(strainBC,'Strain ',k,1,2,bi,bj,myThid)
768 ENDIF
769 CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
770 & 'Um_Advec',k,1,2,bi,bj,myThid)
771 CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
772 & 'Vm_Advec',k,1,2,bi,bj,myThid)
773 ENDIF
774 #endif /* ALLOW_DIAGNOSTICS */
775
776 #endif /* ALLOW_MOM_VECINV */
777
778 RETURN
779 END

  ViewVC Help
Powered by ViewVC 1.1.22