/[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.76 - (show annotations) (download)
Sat Jan 3 23:58:53 2015 UTC (9 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65j, checkpoint65k, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m
Changes since 1.75: +26 -14 lines
- add one argument (the other velocity component) to S/R MOM_U/V_BOTTOMDRAG
  and S/R SHELFICE_U/V_DRAG
- remove condition on bottomDragTerms when calling SHELFICE_U/V_DRAG
  (similar to mom_fluxform calls).

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

  ViewVC Help
Powered by ViewVC 1.1.22