/[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.48 - (show annotations) (download)
Mon Sep 19 19:58:05 2005 UTC (18 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57s_post
Changes since 1.47: +4 -1 lines
fix recent modification (calculation of viscosities in mom_calc_visc.F)

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.47 2005/09/16 19:32:20 baylor 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 - Instance number for this innvocation of CALC_MOM_RHS
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 c _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 C I,J,K - Loop counters
86 INTEGER i,j,k
87 C xxxFac - On-off tracer parameters used for switching terms off.
88 _RL ArDudrFac
89 c _RL mtFacU
90 _RL ArDvdrFac
91 c _RL mtFacV
92 LOGICAL bottomDragTerms
93 LOGICAL writeDiag
94 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98
99 #ifdef ALLOW_MNC
100 INTEGER offsets(9)
101 #endif
102
103 #ifdef ALLOW_AUTODIFF_TAMC
104 C-- only the kDown part of fverU/V is set in this subroutine
105 C-- the kUp is still required
106 C-- In the case of mom_fluxform Kup is set as well
107 C-- (at least in part)
108 fVerU(1,1,kUp) = fVerU(1,1,kUp)
109 fVerV(1,1,kUp) = fVerV(1,1,kUp)
110 #endif
111
112 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
113
114 #ifdef ALLOW_MNC
115 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
116 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
117 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
118 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
119 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
120 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
121 ENDIF
122 DO i = 1,9
123 offsets(i) = 0
124 ENDDO
125 offsets(3) = k
126 C write(*,*) 'offsets = ',(offsets(i),i=1,9)
127 ENDIF
128 #endif /* ALLOW_MNC */
129
130 C Initialise intermediate terms
131 DO J=1-OLy,sNy+OLy
132 DO I=1-OLx,sNx+OLx
133 vF(i,j) = 0.
134 vrF(i,j) = 0.
135 uCf(i,j) = 0.
136 vCf(i,j) = 0.
137 c mT(i,j) = 0.
138 del2u(i,j) = 0.
139 del2v(i,j) = 0.
140 dStar(i,j) = 0.
141 zStar(i,j) = 0.
142 guDiss(i,j)= 0.
143 gvDiss(i,j)= 0.
144 vort3(i,j) = 0.
145 omega3(i,j)= 0.
146 ke(i,j) = 0.
147 #ifdef ALLOW_AUTODIFF_TAMC
148 strain(i,j) = 0. _d 0
149 tension(i,j) = 0. _d 0
150 #endif
151 ENDDO
152 ENDDO
153
154 C-- Term by term tracer parmeters
155 C o U momentum equation
156 ArDudrFac = vfFacMom*1.
157 c mTFacU = mtFacMom*1.
158 C o V momentum equation
159 ArDvdrFac = vfFacMom*1.
160 c mTFacV = mtFacMom*1.
161
162 IF ( no_slip_bottom
163 & .OR. bottomDragQuadratic.NE.0.
164 & .OR. bottomDragLinear.NE.0.) THEN
165 bottomDragTerms=.TRUE.
166 ELSE
167 bottomDragTerms=.FALSE.
168 ENDIF
169
170 C-- Calculate open water fraction at vorticity points
171 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
172
173 C Make local copies of horizontal flow field
174 DO j=1-OLy,sNy+OLy
175 DO i=1-OLx,sNx+OLx
176 uFld(i,j) = uVel(i,j,k,bi,bj)
177 vFld(i,j) = vVel(i,j,k,bi,bj)
178 ENDDO
179 ENDDO
180
181 C note (jmc) : Dissipation and Vort3 advection do not necesary
182 C use the same maskZ (and hFacZ) => needs 2 call(s)
183 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
184
185 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
186
187 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
188
189 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
190
191 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
192
193 CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
194
195 IF (useAbsVorticity)
196 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
197
198 IF (momViscosity) THEN
199 C Calculate del^2 u and del^2 v for bi-harmonic term
200 IF ( (viscA4.NE.0. .AND. no_slip_sides)
201 & .OR. viscA4D.NE.0. .OR. viscA4Z.NE.0.
202 & .OR. viscA4Grid.NE.0.
203 & .OR. viscC4leith.NE.0.
204 & .OR. viscC4leithD.NE.0.
205 & .OR. viscC4smag.NE.0.
206 & ) THEN
207 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
208 O del2u,del2v,
209 & myThid)
210 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
211 CALL MOM_CALC_RELVORT3(bi,bj,k,
212 & del2u,del2v,hFacZ,zStar,myThid)
213 ENDIF
214
215 C Calculate dissipation terms for U and V equations
216
217 C in terms of tension and strain
218 IF (useStrainTensionVisc) THEN
219 CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
220 I hFacZ,
221 O guDiss,gvDiss,
222 I myThid)
223 ELSE
224 C in terms of vorticity and divergence
225 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
226 I hFacZ,dStar,zStar,
227 O guDiss,gvDiss,
228 & myThid)
229 ENDIF
230 ENDIF
231
232 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
233 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
234
235 C---- Zonal momentum equation starts here
236
237 C-- Vertical flux (fVer is at upper face of "u" cell)
238
239 C Eddy component of vertical flux (interior component only) -> vrF
240 IF (momViscosity.AND..NOT.implicitViscosity) THEN
241 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
242
243 C Combine fluxes
244 DO j=jMin,jMax
245 DO i=iMin,iMax
246 fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
247 ENDDO
248 ENDDO
249
250 C-- Tendency is minus divergence of the fluxes
251 DO j=2-Oly,sNy+Oly-1
252 DO i=2-Olx,sNx+Olx-1
253 guDiss(i,j) = guDiss(i,j)
254 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
255 & *recip_rAw(i,j,bi,bj)
256 & *(
257 & fVerU(i,j,kDown) - fVerU(i,j,kUp)
258 & )*rkSign
259 ENDDO
260 ENDDO
261 ENDIF
262
263 C-- No-slip and drag BCs appear as body forces in cell abutting topography
264 IF (momViscosity.AND.no_slip_sides) THEN
265 C- No-slip BCs impose a drag at walls...
266 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
267 DO j=jMin,jMax
268 DO i=iMin,iMax
269 guDiss(i,j) = guDiss(i,j)+vF(i,j)
270 ENDDO
271 ENDDO
272 ENDIF
273
274 C- No-slip BCs impose a drag at bottom
275 IF (momViscosity.AND.bottomDragTerms) THEN
276 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
277 DO j=jMin,jMax
278 DO i=iMin,iMax
279 guDiss(i,j) = guDiss(i,j)+vF(i,j)
280 ENDDO
281 ENDDO
282 ENDIF
283
284 C-- Metric terms for curvilinear grid systems
285 c IF (usingSphericalPolarMTerms) THEN
286 C o Spherical polar grid metric terms
287 c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
288 c DO j=jMin,jMax
289 c DO i=iMin,iMax
290 c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
291 c ENDDO
292 c ENDDO
293 c ENDIF
294
295 C---- Meridional momentum equation starts here
296
297 C-- Vertical flux (fVer is at upper face of "v" cell)
298
299 C Eddy component of vertical flux (interior component only) -> vrF
300 IF (momViscosity.AND..NOT.implicitViscosity) THEN
301 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
302
303 C Combine fluxes -> fVerV
304 DO j=jMin,jMax
305 DO i=iMin,iMax
306 fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
307 ENDDO
308 ENDDO
309
310 C-- Tendency is minus divergence of the fluxes
311 DO j=jMin,jMax
312 DO i=iMin,iMax
313 gvDiss(i,j) = gvDiss(i,j)
314 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
315 & *recip_rAs(i,j,bi,bj)
316 & *(
317 & fVerV(i,j,kDown) - fVerV(i,j,kUp)
318 & )*rkSign
319 ENDDO
320 ENDDO
321 ENDIF
322
323 C-- No-slip and drag BCs appear as body forces in cell abutting topography
324 IF (momViscosity.AND.no_slip_sides) THEN
325 C- No-slip BCs impose a drag at walls...
326 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
327 DO j=jMin,jMax
328 DO i=iMin,iMax
329 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
330 ENDDO
331 ENDDO
332 ENDIF
333 C- No-slip BCs impose a drag at bottom
334 IF (momViscosity.AND.bottomDragTerms) THEN
335 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
336 DO j=jMin,jMax
337 DO i=iMin,iMax
338 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
339 ENDDO
340 ENDDO
341 ENDIF
342
343 C-- Metric terms for curvilinear grid systems
344 c IF (usingSphericalPolarMTerms) THEN
345 C o Spherical polar grid metric terms
346 c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
347 c DO j=jMin,jMax
348 c DO i=iMin,iMax
349 c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
350 c ENDDO
351 c ENDDO
352 c ENDIF
353
354 C-- Horizontal Coriolis terms
355 c IF (useCoriolis .AND. .NOT.useCDscheme
356 c & .AND. .NOT. useAbsVorticity) THEN
357 C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
358 IF ( useCoriolis .AND.
359 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
360 & ) THEN
361 IF (useAbsVorticity) THEN
362 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
363 & uCf,myThid)
364 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
365 & vCf,myThid)
366 ELSE
367 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
368 & uCf,vCf,myThid)
369 ENDIF
370 DO j=jMin,jMax
371 DO i=iMin,iMax
372 gU(i,j,k,bi,bj) = uCf(i,j)
373 gV(i,j,k,bi,bj) = vCf(i,j)
374 ENDDO
375 ENDDO
376
377 IF ( writeDiag ) THEN
378 IF (snapshot_mdsio) THEN
379 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
380 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
381 ENDIF
382 #ifdef ALLOW_MNC
383 IF (useMNC .AND. snapshot_mnc) THEN
384 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fV', uCf,
385 & offsets, myThid)
386 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fU', vCf,
387 & offsets, myThid)
388 ENDIF
389 #endif /* ALLOW_MNC */
390 ENDIF
391 #ifdef ALLOW_DIAGNOSTICS
392 IF ( useDiagnostics ) THEN
393 CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
394 CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
395 ENDIF
396 #endif /* ALLOW_DIAGNOSTICS */
397
398 ELSE
399 DO j=jMin,jMax
400 DO i=iMin,iMax
401 gU(i,j,k,bi,bj) = 0. _d 0
402 gV(i,j,k,bi,bj) = 0. _d 0
403 ENDDO
404 ENDDO
405 ENDIF
406
407 IF (momAdvection) THEN
408 C-- Horizontal advection of relative (or absolute) vorticity
409 IF (highOrderVorticity.AND.useAbsVorticity) THEN
410 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
411 & uCf,myThid)
412 ELSEIF (highOrderVorticity) THEN
413 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
414 & uCf,myThid)
415 ELSEIF (useAbsVorticity) THEN
416 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
417 & uCf,myThid)
418 ELSE
419 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
420 & uCf,myThid)
421 ENDIF
422 DO j=jMin,jMax
423 DO i=iMin,iMax
424 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
425 ENDDO
426 ENDDO
427 IF (highOrderVorticity.AND.useAbsVorticity) THEN
428 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
429 & vCf,myThid)
430 ELSEIF (highOrderVorticity) THEN
431 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
432 & vCf,myThid)
433 ELSEIF (useAbsVorticity) THEN
434 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
435 & vCf,myThid)
436 ELSE
437 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
438 & vCf,myThid)
439 ENDIF
440 DO j=jMin,jMax
441 DO i=iMin,iMax
442 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
443 ENDDO
444 ENDDO
445
446 IF ( writeDiag ) THEN
447 IF (snapshot_mdsio) THEN
448 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
449 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
450 ENDIF
451 #ifdef ALLOW_MNC
452 IF (useMNC .AND. snapshot_mnc) THEN
453 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zV', uCf,
454 & offsets, myThid)
455 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zU', vCf,
456 & offsets, myThid)
457 ENDIF
458 #endif /* ALLOW_MNC */
459 ENDIF
460
461 #ifdef ALLOW_TIMEAVE
462 IF (taveFreq.GT.0.) THEN
463 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
464 & Nr, k, bi, bj, myThid)
465 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
466 & Nr, k, bi, bj, myThid)
467 ENDIF
468 #endif /* ALLOW_TIMEAVE */
469 #ifdef ALLOW_DIAGNOSTICS
470 IF ( useDiagnostics ) THEN
471 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
472 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
473 ENDIF
474 #endif /* ALLOW_DIAGNOSTICS */
475
476 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
477 IF ( .NOT. momImplVertAdv ) THEN
478 CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
479 DO j=jMin,jMax
480 DO i=iMin,iMax
481 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
482 ENDDO
483 ENDDO
484 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
485 DO j=jMin,jMax
486 DO i=iMin,iMax
487 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
488 ENDDO
489 ENDDO
490 #ifdef ALLOW_DIAGNOSTICS
491 IF ( useDiagnostics ) THEN
492 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
493 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
494 ENDIF
495 #endif /* ALLOW_DIAGNOSTICS */
496 ENDIF
497
498 C-- Bernoulli term
499 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
500 DO j=jMin,jMax
501 DO i=iMin,iMax
502 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
503 ENDDO
504 ENDDO
505 CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
506 DO j=jMin,jMax
507 DO i=iMin,iMax
508 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
509 ENDDO
510 ENDDO
511 IF ( writeDiag ) THEN
512 IF (snapshot_mdsio) THEN
513 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
514 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
515 ENDIF
516 #ifdef ALLOW_MNC
517 IF (useMNC .AND. snapshot_mnc) THEN
518 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEx', uCf,
519 & offsets, myThid)
520 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEy', vCf,
521 & offsets, myThid)
522 ENDIF
523 #endif /* ALLOW_MNC */
524 ENDIF
525
526 C-- end if momAdvection
527 ENDIF
528
529 C-- Set du/dt & dv/dt on boundaries to zero
530 DO j=jMin,jMax
531 DO i=iMin,iMax
532 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
533 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
534 ENDDO
535 ENDDO
536
537 #ifdef ALLOW_DEBUG
538 IF ( debugLevel .GE. debLevB
539 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
540 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
541 & .AND. useCubedSphereExchange ) THEN
542 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
543 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
544 ENDIF
545 #endif /* ALLOW_DEBUG */
546
547 IF ( writeDiag ) THEN
548 IF (snapshot_mdsio) THEN
549 CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
550 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,
551 & myThid)
552 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss,bi,bj,k,myIter,myThid)
553 CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss,bi,bj,k,myIter,myThid)
554 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
555 CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
556 CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
557 CALL WRITE_LOCAL_RL('D','I10',1,hDiv,bi,bj,k,myIter,myThid)
558 ENDIF
559 #ifdef ALLOW_MNC
560 IF (useMNC .AND. snapshot_mnc) THEN
561 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Ds',strain,
562 & offsets, myThid)
563 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dt',tension,
564 & offsets, myThid)
565 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Du',guDiss,
566 & offsets, myThid)
567 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dv',gvDiss,
568 & offsets, myThid)
569 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Z3',vort3,
570 & offsets, myThid)
571 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'W3',omega3,
572 & offsets, myThid)
573 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'KE',KE,
574 & offsets, myThid)
575 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'D', hDiv,
576 & offsets, myThid)
577 ENDIF
578 #endif /* ALLOW_MNC */
579 ENDIF
580
581 #ifdef ALLOW_DIAGNOSTICS
582 IF ( useDiagnostics ) THEN
583 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
584 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
585 CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid)
586 CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),
587 & 'Um_Advec',k,1,2,bi,bj,myThid)
588 CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),
589 & 'Vm_Advec',k,1,2,bi,bj,myThid)
590 IF (momViscosity) THEN
591 CALL DIAGNOSTICS_FILL(guDiss,'Um_Diss ',k,1,2,bi,bj,myThid)
592 CALL DIAGNOSTICS_FILL(gvDiss,'Vm_Diss ',k,1,2,bi,bj,myThid)
593 ENDIF
594 ENDIF
595 #endif /* ALLOW_DIAGNOSTICS */
596
597 #endif /* ALLOW_MOM_VECINV */
598
599 RETURN
600 END

  ViewVC Help
Powered by ViewVC 1.1.22