/[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.24 - (show annotations) (download)
Thu Oct 7 21:52:29 2004 UTC (19 years, 7 months ago) by edhill
Branch: MAIN
Changes since 1.23: +67 -16 lines
 o mnc-ify the MOM_VECINV() diagFreq output

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.23 2004/09/24 17:02:34 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 dPhiHydX,dPhiHydY,KappaRU,KappaRV,
9 U fVerU, fVerV,
10 I myTime, myIter, myThid)
11 C /==========================================================\
12 C | S/R MOM_VECINV |
13 C | o Form the right hand-side of the momentum equation. |
14 C |==========================================================|
15 C | Terms are evaluated one layer at a time working from |
16 C | the bottom to the top. The vertically integrated |
17 C | barotropic flow tendency term is evluated by summing the |
18 C | tendencies. |
19 C | Notes: |
20 C | We have not sorted out an entirely satisfactory formula |
21 C | for the diffusion equation bc with lopping. The present |
22 C | form produces a diffusive flux that does not scale with |
23 C | open-area. Need to do something to solidfy this and to |
24 C | deal "properly" with thin walls. |
25 C \==========================================================/
26 IMPLICIT NONE
27
28 C == Global variables ==
29 #include "SIZE.h"
30 #include "DYNVARS.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #ifdef ALLOW_MNC
34 #include "MNC_PARAMS.h"
35 #endif
36 #include "GRID.h"
37 #ifdef ALLOW_TIMEAVE
38 #include "TIMEAVE_STATV.h"
39 #endif
40
41 C == Routine arguments ==
42 C fVerU - Flux of momentum in the vertical
43 C fVerV direction out of the upper face of a cell K
44 C ( flux into the cell above ).
45 C dPhiHydX,Y :: Gradient (X & Y dir.) of Hydrostatic Potential
46 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
47 C results will be set.
48 C kUp, kDown - Index for upper and lower layers.
49 C myThid - Instance number for this innvocation of CALC_MOM_RHS
50 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
51 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
52 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
53 _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
54 _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
55 _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
56 INTEGER kUp,kDown
57 _RL myTime
58 INTEGER myIter
59 INTEGER myThid
60 INTEGER bi,bj,iMin,iMax,jMin,jMax
61
62 #ifdef ALLOW_MOM_VECINV
63
64 C == Functions ==
65 LOGICAL DIFFERENT_MULTIPLE
66 EXTERNAL DIFFERENT_MULTIPLE
67
68 C == Local variables ==
69 _RL aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
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 _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 _RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 C I,J,K - Loop counters
91 INTEGER i,j,k
92 C rVelMaskOverride - Factor for imposing special surface boundary conditions
93 C ( set according to free-surface condition ).
94 C hFacROpen - Lopped cell factos used tohold fraction of open
95 C hFacRClosed and closed cell wall.
96 _RL rVelMaskOverride
97 C xxxFac - On-off tracer parameters used for switching terms off.
98 _RL uDudxFac
99 _RL AhDudxFac
100 _RL A4DuxxdxFac
101 _RL vDudyFac
102 _RL AhDudyFac
103 _RL A4DuyydyFac
104 _RL rVelDudrFac
105 _RL ArDudrFac
106 _RL fuFac
107 _RL phxFac
108 _RL mtFacU
109 _RL uDvdxFac
110 _RL AhDvdxFac
111 _RL A4DvxxdxFac
112 _RL vDvdyFac
113 _RL AhDvdyFac
114 _RL A4DvyydyFac
115 _RL rVelDvdrFac
116 _RL ArDvdrFac
117 _RL fvFac
118 _RL phyFac
119 _RL vForcFac
120 _RL mtFacV
121 _RL wVelBottomOverride
122 LOGICAL bottomDragTerms
123 LOGICAL writeDiag
124 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125 _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126 _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127 _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128
129 #ifdef ALLOW_AUTODIFF_TAMC
130 C-- only the kDown part of fverU/V is set in this subroutine
131 C-- the kUp is still required
132 C-- In the case of mom_fluxform Kup is set as well
133 C-- (at least in part)
134 fVerU(1,1,kUp) = fVerU(1,1,kUp)
135 fVerV(1,1,kUp) = fVerV(1,1,kUp)
136 #endif
137
138 rVelMaskOverride=1.
139 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
140 wVelBottomOverride=1.
141 IF (k.EQ.Nr) wVelBottomOverride=0.
142 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime,
143 & myTime-deltaTClock)
144
145 #ifdef ALLOW_MNC
146 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
147 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
148 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
149 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
150 ENDIF
151 #endif /* ALLOW_MNC */
152
153 C Initialise intermediate terms
154 DO J=1-OLy,sNy+OLy
155 DO I=1-OLx,sNx+OLx
156 aF(i,j) = 0.
157 vF(i,j) = 0.
158 vrF(i,j) = 0.
159 uCf(i,j) = 0.
160 vCf(i,j) = 0.
161 mT(i,j) = 0.
162 pF(i,j) = 0.
163 del2u(i,j) = 0.
164 del2v(i,j) = 0.
165 dStar(i,j) = 0.
166 zStar(i,j) = 0.
167 uDiss(i,j) = 0.
168 vDiss(i,j) = 0.
169 vort3(i,j) = 0.
170 omega3(i,j) = 0.
171 ke(i,j) = 0.
172 #ifdef ALLOW_AUTODIFF_TAMC
173 strain(i,j) = 0. _d 0
174 tension(i,j) = 0. _d 0
175 #endif
176 ENDDO
177 ENDDO
178
179 C-- Term by term tracer parmeters
180 C o U momentum equation
181 uDudxFac = afFacMom*1.
182 AhDudxFac = vfFacMom*1.
183 A4DuxxdxFac = vfFacMom*1.
184 vDudyFac = afFacMom*1.
185 AhDudyFac = vfFacMom*1.
186 A4DuyydyFac = vfFacMom*1.
187 rVelDudrFac = afFacMom*1.
188 ArDudrFac = vfFacMom*1.
189 mTFacU = mtFacMom*1.
190 fuFac = cfFacMom*1.
191 phxFac = pfFacMom*1.
192 C o V momentum equation
193 uDvdxFac = afFacMom*1.
194 AhDvdxFac = vfFacMom*1.
195 A4DvxxdxFac = vfFacMom*1.
196 vDvdyFac = afFacMom*1.
197 AhDvdyFac = vfFacMom*1.
198 A4DvyydyFac = vfFacMom*1.
199 rVelDvdrFac = afFacMom*1.
200 ArDvdrFac = vfFacMom*1.
201 mTFacV = mtFacMom*1.
202 fvFac = cfFacMom*1.
203 phyFac = pfFacMom*1.
204 vForcFac = foFacMom*1.
205
206 IF ( no_slip_bottom
207 & .OR. bottomDragQuadratic.NE.0.
208 & .OR. bottomDragLinear.NE.0.) THEN
209 bottomDragTerms=.TRUE.
210 ELSE
211 bottomDragTerms=.FALSE.
212 ENDIF
213
214 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
215 IF (staggerTimeStep) THEN
216 phxFac = 0.
217 phyFac = 0.
218 ENDIF
219
220 C-- Calculate open water fraction at vorticity points
221 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
222
223 C---- Calculate common quantities used in both U and V equations
224 C Calculate tracer cell face open areas
225 DO j=1-OLy,sNy+OLy
226 DO i=1-OLx,sNx+OLx
227 xA(i,j) = _dyG(i,j,bi,bj)
228 & *drF(k)*_hFacW(i,j,k,bi,bj)
229 yA(i,j) = _dxG(i,j,bi,bj)
230 & *drF(k)*_hFacS(i,j,k,bi,bj)
231 ENDDO
232 ENDDO
233
234 C Make local copies of horizontal flow field
235 DO j=1-OLy,sNy+OLy
236 DO i=1-OLx,sNx+OLx
237 uFld(i,j) = uVel(i,j,k,bi,bj)
238 vFld(i,j) = vVel(i,j,k,bi,bj)
239 ENDDO
240 ENDDO
241
242 C note (jmc) : Dissipation and Vort3 advection do not necesary
243 C use the same maskZ (and hFacZ) => needs 2 call(s)
244 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
245
246 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
247
248 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
249
250 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
251
252 IF (useAbsVorticity)
253 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
254
255 IF (momViscosity) THEN
256 C Calculate del^2 u and del^2 v for bi-harmonic term
257 IF (viscA4.NE.0.
258 & .OR. viscA4Grid.NE.0.
259 & .OR. viscC4leith.NE.0.
260 & ) THEN
261 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
262 O del2u,del2v,
263 & myThid)
264 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
265 CALL MOM_CALC_RELVORT3(
266 & bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)
267 ENDIF
268 C Calculate dissipation terms for U and V equations
269 C in terms of vorticity and divergence
270 IF (viscAh.NE.0. .OR. viscA4.NE.0.
271 & .OR. viscAhGrid.NE.0. .OR. viscA4Grid.NE.0.
272 & .OR. viscC2leith.NE.0. .OR. viscC4leith.NE.0.
273 & ) THEN
274 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,
275 O uDiss,vDiss,
276 & myThid)
277 ENDIF
278 C or in terms of tension and strain
279 IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
280 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
281 O tension,
282 I myThid)
283 CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
284 O strain,
285 I myThid)
286 CALL MOM_HDISSIP(bi,bj,k,
287 I tension,strain,hFacZ,viscAtension,viscAstrain,
288 O uDiss,vDiss,
289 I myThid)
290 ENDIF
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---- Zonal momentum equation starts here
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)
302 & CALL MOM_U_RVISCFLUX(bi,bj,k,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 + coriolis + pressure term
312 DO j=2-Oly,sNy+Oly-1
313 DO i=2-Olx,sNx+Olx-1
314 gU(i,j,k,bi,bj) = uDiss(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,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
319 & )
320 & - phxFac*dPhiHydX(i,j)
321 ENDDO
322 ENDDO
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(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
328 DO j=jMin,jMax
329 DO i=iMin,iMax
330 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
331 ENDDO
332 ENDDO
333 ENDIF
334
335 C- No-slip BCs impose a drag at bottom
336 IF (momViscosity.AND.bottomDragTerms) THEN
337 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
338 DO j=jMin,jMax
339 DO i=iMin,iMax
340 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
341 ENDDO
342 ENDDO
343 ENDIF
344
345 C-- Metric terms for curvilinear grid systems
346 c IF (usingSphericalPolarMTerms) THEN
347 C o Spherical polar grid metric terms
348 c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
349 c DO j=jMin,jMax
350 c DO i=iMin,iMax
351 c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
352 c ENDDO
353 c ENDDO
354 c ENDIF
355
356 C---- Meridional momentum equation starts here
357
358 C-- Vertical flux (fVer is at upper face of "v" cell)
359
360 C Eddy component of vertical flux (interior component only) -> vrF
361 IF (momViscosity.AND..NOT.implicitViscosity)
362 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
363
364 C Combine fluxes -> fVerV
365 DO j=jMin,jMax
366 DO i=iMin,iMax
367 fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
368 ENDDO
369 ENDDO
370
371 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
372 DO j=jMin,jMax
373 DO i=iMin,iMax
374 gV(i,j,k,bi,bj) = vDiss(i,j)
375 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
376 & *recip_rAs(i,j,bi,bj)
377 & *(
378 & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
379 & )
380 & - phyFac*dPhiHydY(i,j)
381 ENDDO
382 ENDDO
383
384 C-- No-slip and drag BCs appear as body forces in cell abutting topography
385 IF (momViscosity.AND.no_slip_sides) THEN
386 C- No-slip BCs impose a drag at walls...
387 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
388 DO j=jMin,jMax
389 DO i=iMin,iMax
390 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
391 ENDDO
392 ENDDO
393 ENDIF
394 C- No-slip BCs impose a drag at bottom
395 IF (momViscosity.AND.bottomDragTerms) THEN
396 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
397 DO j=jMin,jMax
398 DO i=iMin,iMax
399 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
400 ENDDO
401 ENDDO
402 ENDIF
403
404 C-- Metric terms for curvilinear grid systems
405 c IF (usingSphericalPolarMTerms) THEN
406 C o Spherical polar grid metric terms
407 c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
408 c DO j=jMin,jMax
409 c DO i=iMin,iMax
410 c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
411 c ENDDO
412 c ENDDO
413 c ENDIF
414
415 C-- Horizontal Coriolis terms
416 IF (useCoriolis .AND. .NOT.useCDscheme
417 & .AND. .NOT. useAbsVorticity) THEN
418 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
419 & uCf,vCf,myThid)
420 DO j=jMin,jMax
421 DO i=iMin,iMax
422 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
423 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
424 ENDDO
425 ENDDO
426 IF ( writeDiag ) THEN
427 IF (snapshot_mdsio) THEN
428 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
429 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
430 ENDIF
431 #ifdef ALLOW_MNC
432 IF (useMNC .AND. snapshot_mnc) THEN
433 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'fV', uCf, myThid)
434 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'fU', vCf, myThid)
435 ENDIF
436 #endif /* ALLOW_MNC */
437 ENDIF
438 ENDIF
439
440 IF (momAdvection) THEN
441 C-- Horizontal advection of relative vorticity
442 IF (useAbsVorticity) THEN
443 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
444 & uCf,myThid)
445 ELSE
446 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3,hFacZ,r_hFacZ,
447 & uCf,myThid)
448 ENDIF
449 c CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid)
450 DO j=jMin,jMax
451 DO i=iMin,iMax
452 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
453 ENDDO
454 ENDDO
455 IF (useAbsVorticity) THEN
456 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
457 & vCf,myThid)
458 ELSE
459 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3,hFacZ,r_hFacZ,
460 & vCf,myThid)
461 ENDIF
462 c CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid)
463 DO j=jMin,jMax
464 DO i=iMin,iMax
465 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
466 ENDDO
467 ENDDO
468
469 IF ( writeDiag ) THEN
470 IF (snapshot_mdsio) THEN
471 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
472 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
473 ENDIF
474 #ifdef ALLOW_MNC
475 IF (useMNC .AND. snapshot_mnc) THEN
476 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'zV', uCf, myThid)
477 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'zU', vCf, myThid)
478 ENDIF
479 #endif /* ALLOW_MNC */
480 ENDIF
481
482 #ifdef ALLOW_TIMEAVE
483 #ifndef HRCUBE
484 IF (taveFreq.GT.0.) THEN
485 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
486 & Nr, k, bi, bj, myThid)
487 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
488 & Nr, k, bi, bj, myThid)
489 ENDIF
490 #endif /* ndef HRCUBE */
491 #endif /* ALLOW_TIMEAVE */
492
493 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
494 IF ( .NOT. momImplVertAdv ) THEN
495 CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
496 DO j=jMin,jMax
497 DO i=iMin,iMax
498 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
499 ENDDO
500 ENDDO
501 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
502 DO j=jMin,jMax
503 DO i=iMin,iMax
504 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
505 ENDDO
506 ENDDO
507 ENDIF
508
509 C-- Bernoulli term
510 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
511 DO j=jMin,jMax
512 DO i=iMin,iMax
513 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
514 ENDDO
515 ENDDO
516 CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
517 DO j=jMin,jMax
518 DO i=iMin,iMax
519 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
520 ENDDO
521 ENDDO
522 IF ( writeDiag ) THEN
523 IF (snapshot_mdsio) THEN
524 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
525 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
526 ENDIF
527 #ifdef ALLOW_MNC
528 IF (useMNC .AND. snapshot_mnc) THEN
529 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'KEx', uCf, myThid)
530 CALL MNC_CW_RL_W('D','mom_vi',0,0, 'KEy', vCf, myThid)
531 ENDIF
532 #endif /* ALLOW_MNC */
533 ENDIF
534
535 C-- end if momAdvection
536 ENDIF
537
538 C-- Set du/dt & dv/dt on boundaries to zero
539 DO j=jMin,jMax
540 DO i=iMin,iMax
541 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
542 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
543 ENDDO
544 ENDDO
545
546 #ifdef ALLOW_DEBUG
547 IF ( debugLevel .GE. debLevB
548 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
549 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
550 & .AND. useCubedSphereExchange ) THEN
551 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
552 & uDiss,vDiss, k, standardMessageUnit,bi,bj,myThid )
553 ENDIF
554 #endif /* ALLOW_DEBUG */
555
556 IF ( writeDiag ) THEN
557 IF (snapshot_mdsio) THEN
558 CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
559 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,
560 & myThid)
561 CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid)
562 CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid)
563 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
564 CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
565 CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
566 CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid)
567 ENDIF
568 #ifdef ALLOW_MNC
569 IF (useMNC .AND. snapshot_mnc) THEN
570 CALL MNC_CW_RL_W('D','mom_vi',0,0,'Ds',strain, myThid)
571 CALL MNC_CW_RL_W('D','mom_vi',0,0,'Dt',tension, myThid)
572 CALL MNC_CW_RL_W('D','mom_vi',0,0,'Du',uDiss, myThid)
573 CALL MNC_CW_RL_W('D','mom_vi',0,0,'Dv',vDiss, myThid)
574 CALL MNC_CW_RL_W('D','mom_vi',0,0,'Z3',vort3, myThid)
575 CALL MNC_CW_RL_W('D','mom_vi',0,0,'W3',omega3, myThid)
576 CALL MNC_CW_RL_W('D','mom_vi',0,0,'KE',KE, myThid)
577 CALL MNC_CW_RL_W('D','mom_vi',0,0,'D', hdiv, myThid)
578 ENDIF
579 #endif /* ALLOW_MNC */
580 ENDIF
581
582 #endif /* ALLOW_MOM_VECINV */
583
584 RETURN
585 END

  ViewVC Help
Powered by ViewVC 1.1.22