/[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.58 - (show annotations) (download)
Thu Jul 13 03:02:48 2006 UTC (17 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58m_post
Changes since 1.57: +2 -2 lines
test use3dCoriolis (where needed)

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

  ViewVC Help
Powered by ViewVC 1.1.22