/[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.29 - (show annotations) (download)
Fri Nov 5 18:39:15 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
Changes since 1.28: +7 -58 lines
remove unused variables

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.28 2004/11/02 01:04:08 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 vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
70 _RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71 _RL uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 c _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84 _RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 _RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 C I,J,K - Loop counters
87 INTEGER i,j,k
88 C rVelMaskOverride - Factor for imposing special surface boundary conditions
89 C ( set according to free-surface condition ).
90 C hFacROpen - Lopped cell factos used tohold fraction of open
91 C hFacRClosed and closed cell wall.
92 _RL rVelMaskOverride
93 C xxxFac - On-off tracer parameters used for switching terms off.
94 _RL ArDudrFac
95 _RL phxFac
96 c _RL mtFacU
97 _RL ArDvdrFac
98 _RL phyFac
99 c _RL mtFacV
100 _RL wVelBottomOverride
101 LOGICAL bottomDragTerms
102 LOGICAL writeDiag
103 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107
108 #ifdef ALLOW_MNC
109 INTEGER offsets(9)
110 #endif
111
112 #ifdef ALLOW_AUTODIFF_TAMC
113 C-- only the kDown part of fverU/V is set in this subroutine
114 C-- the kUp is still required
115 C-- In the case of mom_fluxform Kup is set as well
116 C-- (at least in part)
117 fVerU(1,1,kUp) = fVerU(1,1,kUp)
118 fVerV(1,1,kUp) = fVerV(1,1,kUp)
119 #endif
120
121 rVelMaskOverride=1.
122 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
123 wVelBottomOverride=1.
124 IF (k.EQ.Nr) wVelBottomOverride=0.
125 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime,
126 & myTime-deltaTClock)
127
128 #ifdef ALLOW_MNC
129 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
130 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
131 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
132 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
133 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
134 ENDIF
135 DO i = 1,9
136 offsets(i) = 0
137 ENDDO
138 offsets(3) = k
139 C write(*,*) 'offsets = ',(offsets(i),i=1,9)
140 ENDIF
141 #endif /* ALLOW_MNC */
142
143 C Initialise intermediate terms
144 DO J=1-OLy,sNy+OLy
145 DO I=1-OLx,sNx+OLx
146 vF(i,j) = 0.
147 vrF(i,j) = 0.
148 uCf(i,j) = 0.
149 vCf(i,j) = 0.
150 c mT(i,j) = 0.
151 del2u(i,j) = 0.
152 del2v(i,j) = 0.
153 dStar(i,j) = 0.
154 zStar(i,j) = 0.
155 uDiss(i,j) = 0.
156 vDiss(i,j) = 0.
157 vort3(i,j) = 0.
158 omega3(i,j) = 0.
159 ke(i,j) = 0.
160 #ifdef ALLOW_AUTODIFF_TAMC
161 strain(i,j) = 0. _d 0
162 tension(i,j) = 0. _d 0
163 #endif
164 ENDDO
165 ENDDO
166
167 C-- Term by term tracer parmeters
168 C o U momentum equation
169 ArDudrFac = vfFacMom*1.
170 c mTFacU = mtFacMom*1.
171 phxFac = pfFacMom*1.
172 C o V momentum equation
173 ArDvdrFac = vfFacMom*1.
174 c mTFacV = mtFacMom*1.
175 phyFac = pfFacMom*1.
176
177 IF ( no_slip_bottom
178 & .OR. bottomDragQuadratic.NE.0.
179 & .OR. bottomDragLinear.NE.0.) THEN
180 bottomDragTerms=.TRUE.
181 ELSE
182 bottomDragTerms=.FALSE.
183 ENDIF
184
185 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
186 IF (staggerTimeStep) THEN
187 phxFac = 0.
188 phyFac = 0.
189 ENDIF
190
191 C-- Calculate open water fraction at vorticity points
192 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
193
194 C Make local copies of horizontal flow field
195 DO j=1-OLy,sNy+OLy
196 DO i=1-OLx,sNx+OLx
197 uFld(i,j) = uVel(i,j,k,bi,bj)
198 vFld(i,j) = vVel(i,j,k,bi,bj)
199 ENDDO
200 ENDDO
201
202 C note (jmc) : Dissipation and Vort3 advection do not necesary
203 C use the same maskZ (and hFacZ) => needs 2 call(s)
204 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
205
206 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
207
208 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
209
210 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
211
212 IF (useAbsVorticity)
213 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
214
215 IF (momViscosity) THEN
216 C Calculate del^2 u and del^2 v for bi-harmonic term
217 IF (viscA4.NE.0.
218 & .OR. viscA4Grid.NE.0.
219 & .OR. viscC4leith.NE.0.
220 & ) THEN
221 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
222 O del2u,del2v,
223 & myThid)
224 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
225 CALL MOM_CALC_RELVORT3(
226 & bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)
227 ENDIF
228 C Calculate dissipation terms for U and V equations
229 C in terms of vorticity and divergence
230 IF ( viscAhD.NE.0. .OR. viscAhZ.NE.0.
231 & .OR. viscA4D.NE.0. .OR. viscA4Z.NE.0.
232 & .OR. viscAhGrid.NE.0. .OR. viscA4Grid.NE.0.
233 & .OR. viscC2leith.NE.0. .OR. viscC4leith.NE.0.
234 & ) THEN
235 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,
236 O uDiss,vDiss,
237 & myThid)
238 ENDIF
239 C or in terms of tension and strain
240 IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
241 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
242 O tension,
243 I myThid)
244 CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
245 O strain,
246 I myThid)
247 CALL MOM_HDISSIP(bi,bj,k,
248 I tension,strain,hFacZ,viscAtension,viscAstrain,
249 O uDiss,vDiss,
250 I myThid)
251 ENDIF
252 ENDIF
253
254 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
255 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
256
257 C---- Zonal momentum equation starts here
258
259 C-- Vertical flux (fVer is at upper face of "u" cell)
260
261 C Eddy component of vertical flux (interior component only) -> vrF
262 IF (momViscosity.AND..NOT.implicitViscosity)
263 & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
264
265 C Combine fluxes
266 DO j=jMin,jMax
267 DO i=iMin,iMax
268 fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
269 ENDDO
270 ENDDO
271
272 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
273 DO j=2-Oly,sNy+Oly-1
274 DO i=2-Olx,sNx+Olx-1
275 gU(i,j,k,bi,bj) = uDiss(i,j)
276 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
277 & *recip_rAw(i,j,bi,bj)
278 & *(
279 & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
280 & )
281 & - phxFac*dPhiHydX(i,j)
282 ENDDO
283 ENDDO
284
285 C-- No-slip and drag BCs appear as body forces in cell abutting topography
286 IF (momViscosity.AND.no_slip_sides) THEN
287 C- No-slip BCs impose a drag at walls...
288 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
289 DO j=jMin,jMax
290 DO i=iMin,iMax
291 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
292 ENDDO
293 ENDDO
294 ENDIF
295
296 C- No-slip BCs impose a drag at bottom
297 IF (momViscosity.AND.bottomDragTerms) THEN
298 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
299 DO j=jMin,jMax
300 DO i=iMin,iMax
301 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
302 ENDDO
303 ENDDO
304 ENDIF
305
306 C-- Metric terms for curvilinear grid systems
307 c IF (usingSphericalPolarMTerms) THEN
308 C o Spherical polar grid metric terms
309 c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
310 c DO j=jMin,jMax
311 c DO i=iMin,iMax
312 c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
313 c ENDDO
314 c ENDDO
315 c ENDIF
316
317 C---- Meridional momentum equation starts here
318
319 C-- Vertical flux (fVer is at upper face of "v" cell)
320
321 C Eddy component of vertical flux (interior component only) -> vrF
322 IF (momViscosity.AND..NOT.implicitViscosity)
323 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
324
325 C Combine fluxes -> fVerV
326 DO j=jMin,jMax
327 DO i=iMin,iMax
328 fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
329 ENDDO
330 ENDDO
331
332 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
333 DO j=jMin,jMax
334 DO i=iMin,iMax
335 gV(i,j,k,bi,bj) = vDiss(i,j)
336 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
337 & *recip_rAs(i,j,bi,bj)
338 & *(
339 & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
340 & )
341 & - phyFac*dPhiHydY(i,j)
342 ENDDO
343 ENDDO
344
345 C-- No-slip and drag BCs appear as body forces in cell abutting topography
346 IF (momViscosity.AND.no_slip_sides) THEN
347 C- No-slip BCs impose a drag at walls...
348 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
349 DO j=jMin,jMax
350 DO i=iMin,iMax
351 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
352 ENDDO
353 ENDDO
354 ENDIF
355 C- No-slip BCs impose a drag at bottom
356 IF (momViscosity.AND.bottomDragTerms) THEN
357 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
358 DO j=jMin,jMax
359 DO i=iMin,iMax
360 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
361 ENDDO
362 ENDDO
363 ENDIF
364
365 C-- Metric terms for curvilinear grid systems
366 c IF (usingSphericalPolarMTerms) THEN
367 C o Spherical polar grid metric terms
368 c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
369 c DO j=jMin,jMax
370 c DO i=iMin,iMax
371 c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
372 c ENDDO
373 c ENDDO
374 c ENDIF
375
376 C-- Horizontal Coriolis terms
377 IF (useCoriolis .AND. .NOT.useCDscheme
378 & .AND. .NOT. useAbsVorticity) THEN
379 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
380 & uCf,vCf,myThid)
381 DO j=jMin,jMax
382 DO i=iMin,iMax
383 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
384 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
385 ENDDO
386 ENDDO
387 IF ( writeDiag ) THEN
388 IF (snapshot_mdsio) THEN
389 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
390 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
391 ENDIF
392 #ifdef ALLOW_MNC
393 IF (useMNC .AND. snapshot_mnc) THEN
394 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fV', uCf,
395 & offsets, myThid)
396 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fU', vCf,
397 & offsets, myThid)
398 ENDIF
399 #endif /* ALLOW_MNC */
400 ENDIF
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 HRCUBE
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 HRCUBE */
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 & uDiss,vDiss, 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,uDiss,bi,bj,k,myIter,myThid)
529 CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,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',uDiss,
542 & offsets, myThid)
543 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dv',vDiss,
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