/[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.32 - (show annotations) (download)
Sat Dec 4 05:59:50 2004 UTC (19 years, 5 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint57, checkpoint57a_post, checkpoint57a_pre
Changes since 1.31: +3 -3 lines
Added CPP option MINIMAL_TAVE_OUTPUT for minimal time-averaged output:
S, T, U, V, W, ETA, and phiHydLow.

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

  ViewVC Help
Powered by ViewVC 1.1.22