/[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.77 - (show annotations) (download)
Thu Sep 10 18:08:51 2015 UTC (9 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint66a, checkpoint65o
Changes since 1.76: +17 -13 lines
- add anelastic and deep-atmosphere geometry factor in pkg/mom_vecinv ; this
  allows to use Vector-Invariant form in deep atmos and anelastic formulation

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.76 2015/01/03 23:58: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 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
370 DO j=jMin,jMax
371 DO i=iMin,iMax
372 guDiss(i,j) = guDiss(i,j)
373 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
374 & *recip_rAw(i,j,bi,bj)
375 & *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
376 & *recip_deepFac2C(k)*recip_rhoFacC(k)
377 ENDDO
378 ENDDO
379 ENDIF
380
381 C-- No-slip and drag BCs appear as body forces in cell abutting topography
382 IF ( no_slip_sides ) THEN
383 C- No-slip BCs impose a drag at walls...
384 CALL MOM_U_SIDEDRAG( bi, bj, k,
385 I uFld, del2u, h0FacZ,
386 I viscAh_Z, viscA4_Z,
387 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
388 O vF,
389 I myThid )
390 DO j=jMin,jMax
391 DO i=iMin,iMax
392 guDiss(i,j) = guDiss(i,j)+vF(i,j)
393 ENDDO
394 ENDDO
395 ENDIF
396
397 C- No-slip BCs impose a drag at bottom
398 IF ( bottomDragTerms ) THEN
399 CALL MOM_U_BOTTOMDRAG( bi, bj, k,
400 I uFld, vFld, KE, kappaRU,
401 O vF,
402 I myThid )
403 DO j=jMin,jMax
404 DO i=iMin,iMax
405 guDiss(i,j) = guDiss(i,j)+vF(i,j)
406 ENDDO
407 ENDDO
408 ENDIF
409 #ifdef ALLOW_SHELFICE
410 IF ( useShelfIce ) THEN
411 CALL SHELFICE_U_DRAG( bi, bj, k,
412 I uFld, vFld, KE, kappaRU,
413 O vF,
414 I myThid )
415 DO j=jMin,jMax
416 DO i=iMin,iMax
417 guDiss(i,j) = guDiss(i,j) + vF(i,j)
418 ENDDO
419 ENDDO
420 ENDIF
421 #endif /* ALLOW_SHELFICE */
422
423 C--- Other dissipation terms in Meridional momentum equation
424
425 C-- Vertical flux (fVer is at upper face of "v" cell)
426 C Eddy component of vertical flux (interior component only) -> vrF
427 IF ( .NOT.implicitViscosity ) THEN
428 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,kappaRV,vrF,myThid)
429 C Combine fluxes -> fVerV
430 DO j=jMin,jMax
431 DO i=iMin,iMax
432 fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
433 ENDDO
434 ENDDO
435 C-- Tendency is minus divergence of the fluxes
436 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
437 DO j=jMin,jMax
438 DO i=iMin,iMax
439 gvDiss(i,j) = gvDiss(i,j)
440 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
441 & *recip_rAs(i,j,bi,bj)
442 & *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
443 & *recip_deepFac2C(k)*recip_rhoFacC(k)
444 ENDDO
445 ENDDO
446 ENDIF
447
448 C-- No-slip and drag BCs appear as body forces in cell abutting topography
449 IF ( no_slip_sides ) THEN
450 C- No-slip BCs impose a drag at walls...
451 CALL MOM_V_SIDEDRAG( bi, bj, k,
452 I vFld, del2v, h0FacZ,
453 I viscAh_Z, viscA4_Z,
454 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
455 O vF,
456 I myThid )
457 DO j=jMin,jMax
458 DO i=iMin,iMax
459 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
460 ENDDO
461 ENDDO
462 ENDIF
463
464 C- No-slip BCs impose a drag at bottom
465 IF ( bottomDragTerms ) THEN
466 CALL MOM_V_BOTTOMDRAG( bi, bj, k,
467 I uFld, vFld, KE, kappaRV,
468 O vF,
469 I myThid )
470 DO j=jMin,jMax
471 DO i=iMin,iMax
472 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
473 ENDDO
474 ENDDO
475 ENDIF
476 #ifdef ALLOW_SHELFICE
477 IF ( useShelfIce ) THEN
478 CALL SHELFICE_V_DRAG( bi, bj, k,
479 I uFld, vFld, KE, kappaRV,
480 O vF,
481 I myThid )
482 DO j=jMin,jMax
483 DO i=iMin,iMax
484 gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
485 ENDDO
486 ENDDO
487 ENDIF
488 #endif /* ALLOW_SHELFICE */
489
490 C-- if (momViscosity) end of block.
491 ENDIF
492
493 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
494 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
495
496 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
497
498 C--- Prepare for Advection & Coriolis terms:
499 C- calculate absolute vorticity
500 IF (useAbsVorticity)
501 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
502
503 C-- Horizontal Coriolis terms
504 c IF (useCoriolis .AND. .NOT.useCDscheme
505 c & .AND. .NOT. useAbsVorticity) THEN
506 C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
507 IF ( useCoriolis .AND.
508 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
509 & ) THEN
510 IF (useAbsVorticity) THEN
511 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,omega3,hFacZ,r_hFacZ,
512 & uCf,myThid)
513 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,omega3,hFacZ,r_hFacZ,
514 & vCf,myThid)
515 ELSE
516 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
517 & uCf,vCf,myThid)
518 ENDIF
519 DO j=jMin,jMax
520 DO i=iMin,iMax
521 gU(i,j,k,bi,bj) = uCf(i,j)
522 gV(i,j,k,bi,bj) = vCf(i,j)
523 ENDDO
524 ENDDO
525 IF ( writeDiag ) THEN
526 IF (snapshot_mdsio) THEN
527 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
528 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
529 ENDIF
530 #ifdef ALLOW_MNC
531 IF (useMNC .AND. snapshot_mnc) THEN
532 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
533 & offsets, myThid)
534 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
535 & offsets, myThid)
536 ENDIF
537 #endif /* ALLOW_MNC */
538 ENDIF
539 #ifdef ALLOW_DIAGNOSTICS
540 IF ( useDiagnostics ) THEN
541 CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
542 CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
543 ENDIF
544 #endif /* ALLOW_DIAGNOSTICS */
545 ELSE
546 DO j=jMin,jMax
547 DO i=iMin,iMax
548 gU(i,j,k,bi,bj) = 0. _d 0
549 gV(i,j,k,bi,bj) = 0. _d 0
550 ENDDO
551 ENDDO
552 ENDIF
553
554 IF (momAdvection) THEN
555 C-- Horizontal advection of relative (or absolute) vorticity
556 IF ( (highOrderVorticity.OR.upwindVorticity)
557 & .AND.useAbsVorticity ) THEN
558 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
559 & uCf,myThid)
560 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
561 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
562 & uCf,myThid)
563 ELSEIF ( useAbsVorticity ) THEN
564 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,omega3,hFacZ,r_hFacZ,
565 & uCf,myThid)
566 ELSE
567 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
568 & uCf,myThid)
569 ENDIF
570 DO j=jMin,jMax
571 DO i=iMin,iMax
572 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
573 ENDDO
574 ENDDO
575 IF ( (highOrderVorticity.OR.upwindVorticity)
576 & .AND.useAbsVorticity ) THEN
577 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,uFld,omega3,r_hFacZ,
578 & vCf,myThid)
579 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
580 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,uFld,vort3, r_hFacZ,
581 & vCf,myThid)
582 ELSEIF ( useAbsVorticity ) THEN
583 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,omega3,hFacZ,r_hFacZ,
584 & vCf,myThid)
585 ELSE
586 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
587 & vCf,myThid)
588 ENDIF
589 DO j=jMin,jMax
590 DO i=iMin,iMax
591 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
592 ENDDO
593 ENDDO
594
595 IF ( writeDiag ) THEN
596 IF (snapshot_mdsio) THEN
597 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
598 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
599 ENDIF
600 #ifdef ALLOW_MNC
601 IF (useMNC .AND. snapshot_mnc) THEN
602 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
603 & offsets, myThid)
604 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
605 & offsets, myThid)
606 ENDIF
607 #endif /* ALLOW_MNC */
608 ENDIF
609
610 #ifdef ALLOW_TIMEAVE
611 IF (taveFreq.GT.0.) THEN
612 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
613 & Nr, k, bi, bj, myThid)
614 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
615 & Nr, k, bi, bj, myThid)
616 ENDIF
617 #endif /* ALLOW_TIMEAVE */
618 #ifdef ALLOW_DIAGNOSTICS
619 IF ( useDiagnostics ) THEN
620 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
621 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
622 ENDIF
623 #endif /* ALLOW_DIAGNOSTICS */
624
625 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
626 IF ( .NOT. momImplVertAdv ) THEN
627 CALL MOM_VI_U_VERTSHEAR(bi,bj,k,uVel,wVel,uCf,myThid)
628 DO j=jMin,jMax
629 DO i=iMin,iMax
630 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
631 ENDDO
632 ENDDO
633 CALL MOM_VI_V_VERTSHEAR(bi,bj,k,vVel,wVel,vCf,myThid)
634 DO j=jMin,jMax
635 DO i=iMin,iMax
636 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
637 ENDDO
638 ENDDO
639 #ifdef ALLOW_DIAGNOSTICS
640 IF ( useDiagnostics ) THEN
641 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
642 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
643 ENDIF
644 #endif /* ALLOW_DIAGNOSTICS */
645 ENDIF
646
647 C-- Bernoulli term
648 CALL MOM_VI_U_GRAD_KE(bi,bj,k,KE,uCf,myThid)
649 DO j=jMin,jMax
650 DO i=iMin,iMax
651 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
652 ENDDO
653 ENDDO
654 CALL MOM_VI_V_GRAD_KE(bi,bj,k,KE,vCf,myThid)
655 DO j=jMin,jMax
656 DO i=iMin,iMax
657 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
658 ENDDO
659 ENDDO
660 IF ( writeDiag ) THEN
661 IF (snapshot_mdsio) THEN
662 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
663 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
664 ENDIF
665 #ifdef ALLOW_MNC
666 IF (useMNC .AND. snapshot_mnc) THEN
667 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
668 & offsets, myThid)
669 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
670 & offsets, myThid)
671 ENDIF
672 #endif /* ALLOW_MNC */
673 ENDIF
674
675 C-- end if momAdvection
676 ENDIF
677
678 C-- 3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
679 IF ( use3dCoriolis ) THEN
680 CALL MOM_U_CORIOLIS_NH(bi,bj,k,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 IF ( usingCurvilinearGrid ) THEN
687 C- presently, non zero angleSinC array only supported with Curvilinear-Grid
688 CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
689 DO j=jMin,jMax
690 DO i=iMin,iMax
691 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
692 ENDDO
693 ENDDO
694 ENDIF
695 ENDIF
696
697 C-- Non-Hydrostatic (spherical) metric terms
698 IF ( useNHMTerms ) THEN
699 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
700 DO j=jMin,jMax
701 DO i=iMin,iMax
702 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
703 ENDDO
704 ENDDO
705 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
706 DO j=jMin,jMax
707 DO i=iMin,iMax
708 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
709 ENDDO
710 ENDDO
711 ENDIF
712
713 C-- Set du/dt & dv/dt on boundaries to zero
714 DO j=jMin,jMax
715 DO i=iMin,iMax
716 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
717 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
718 ENDDO
719 ENDDO
720
721 #ifdef ALLOW_DEBUG
722 IF ( debugLevel .GE. debLevC
723 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
724 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
725 & .AND. useCubedSphereExchange ) THEN
726 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
727 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
728 ENDIF
729 #endif /* ALLOW_DEBUG */
730
731 IF ( writeDiag ) THEN
732 IF (useBiharmonicVisc) THEN
733 CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
734 & bi,bj,k, myIter, myThid )
735 CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
736 & bi,bj,k, myIter, myThid )
737 CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
738 & bi,bj,k, myIter, myThid )
739 CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
740 & bi,bj,k, myIter, myThid )
741 ENDIF
742 IF (snapshot_mdsio) THEN
743 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
744 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
745 CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
746 CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
747 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
748 CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
749 & bi,bj,k, myIter, myThid )
750 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
751 CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
752 ENDIF
753 #ifdef ALLOW_MNC
754 IF (useMNC .AND. snapshot_mnc) THEN
755 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
756 & offsets, myThid)
757 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
758 & offsets, myThid)
759 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
760 & offsets, myThid)
761 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
762 & offsets, myThid)
763 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
764 & offsets, myThid)
765 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
766 & offsets, myThid)
767 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
768 & offsets, myThid)
769 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
770 & offsets, myThid)
771 ENDIF
772 #endif /* ALLOW_MNC */
773 ENDIF
774
775 #ifdef ALLOW_DIAGNOSTICS
776 IF ( useDiagnostics ) THEN
777 CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
778 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
779 IF (momViscosity) THEN
780 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
781 ENDIF
782 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
783 CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
784 CALL DIAGNOSTICS_FILL(strainBC,'Strain ',k,1,2,bi,bj,myThid)
785 ENDIF
786 CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
787 & 'Um_Advec',k,1,2,bi,bj,myThid)
788 CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
789 & 'Vm_Advec',k,1,2,bi,bj,myThid)
790 ENDIF
791 #endif /* ALLOW_DIAGNOSTICS */
792
793 #endif /* ALLOW_MOM_VECINV */
794
795 RETURN
796 END

  ViewVC Help
Powered by ViewVC 1.1.22