/[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.25 - (show annotations) (download)
Fri Oct 8 17:03:21 2004 UTC (19 years, 7 months ago) by edhill
Branch: MAIN
Changes since 1.24: +44 -19 lines
 o add ability of MNC to write local and "partial" (eg. 2D slices where
   the full 3D field is never actually stored) arrays to NetCDF files
   with the correct (that is, the complete multi-dimensional) set of
   array indicies
   - used in mom_vecinv() to write the diagFreq output
   - tested (demonstrated) in verification/aim.5l_cs

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.24 2004/10/07 21:52:29 edhill 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_MNC
130 INTEGER offsets(9)
131 #endif
132
133 #ifdef ALLOW_AUTODIFF_TAMC
134 C-- only the kDown part of fverU/V is set in this subroutine
135 C-- the kUp is still required
136 C-- In the case of mom_fluxform Kup is set as well
137 C-- (at least in part)
138 fVerU(1,1,kUp) = fVerU(1,1,kUp)
139 fVerV(1,1,kUp) = fVerV(1,1,kUp)
140 #endif
141
142 rVelMaskOverride=1.
143 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
144 wVelBottomOverride=1.
145 IF (k.EQ.Nr) wVelBottomOverride=0.
146 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime,
147 & myTime-deltaTClock)
148
149 #ifdef ALLOW_MNC
150 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
151 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
152 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
153 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
154 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
155 ENDIF
156 DO i = 1,9
157 offsets(i) = 0
158 ENDDO
159 offsets(3) = k
160 C write(*,*) 'offsets = ',(offsets(i),i=1,9)
161 ENDIF
162 #endif /* ALLOW_MNC */
163
164 C Initialise intermediate terms
165 DO J=1-OLy,sNy+OLy
166 DO I=1-OLx,sNx+OLx
167 aF(i,j) = 0.
168 vF(i,j) = 0.
169 vrF(i,j) = 0.
170 uCf(i,j) = 0.
171 vCf(i,j) = 0.
172 mT(i,j) = 0.
173 pF(i,j) = 0.
174 del2u(i,j) = 0.
175 del2v(i,j) = 0.
176 dStar(i,j) = 0.
177 zStar(i,j) = 0.
178 uDiss(i,j) = 0.
179 vDiss(i,j) = 0.
180 vort3(i,j) = 0.
181 omega3(i,j) = 0.
182 ke(i,j) = 0.
183 #ifdef ALLOW_AUTODIFF_TAMC
184 strain(i,j) = 0. _d 0
185 tension(i,j) = 0. _d 0
186 #endif
187 ENDDO
188 ENDDO
189
190 C-- Term by term tracer parmeters
191 C o U momentum equation
192 uDudxFac = afFacMom*1.
193 AhDudxFac = vfFacMom*1.
194 A4DuxxdxFac = vfFacMom*1.
195 vDudyFac = afFacMom*1.
196 AhDudyFac = vfFacMom*1.
197 A4DuyydyFac = vfFacMom*1.
198 rVelDudrFac = afFacMom*1.
199 ArDudrFac = vfFacMom*1.
200 mTFacU = mtFacMom*1.
201 fuFac = cfFacMom*1.
202 phxFac = pfFacMom*1.
203 C o V momentum equation
204 uDvdxFac = afFacMom*1.
205 AhDvdxFac = vfFacMom*1.
206 A4DvxxdxFac = vfFacMom*1.
207 vDvdyFac = afFacMom*1.
208 AhDvdyFac = vfFacMom*1.
209 A4DvyydyFac = vfFacMom*1.
210 rVelDvdrFac = afFacMom*1.
211 ArDvdrFac = vfFacMom*1.
212 mTFacV = mtFacMom*1.
213 fvFac = cfFacMom*1.
214 phyFac = pfFacMom*1.
215 vForcFac = foFacMom*1.
216
217 IF ( no_slip_bottom
218 & .OR. bottomDragQuadratic.NE.0.
219 & .OR. bottomDragLinear.NE.0.) THEN
220 bottomDragTerms=.TRUE.
221 ELSE
222 bottomDragTerms=.FALSE.
223 ENDIF
224
225 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
226 IF (staggerTimeStep) THEN
227 phxFac = 0.
228 phyFac = 0.
229 ENDIF
230
231 C-- Calculate open water fraction at vorticity points
232 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
233
234 C---- Calculate common quantities used in both U and V equations
235 C Calculate tracer cell face open areas
236 DO j=1-OLy,sNy+OLy
237 DO i=1-OLx,sNx+OLx
238 xA(i,j) = _dyG(i,j,bi,bj)
239 & *drF(k)*_hFacW(i,j,k,bi,bj)
240 yA(i,j) = _dxG(i,j,bi,bj)
241 & *drF(k)*_hFacS(i,j,k,bi,bj)
242 ENDDO
243 ENDDO
244
245 C Make local copies of horizontal flow field
246 DO j=1-OLy,sNy+OLy
247 DO i=1-OLx,sNx+OLx
248 uFld(i,j) = uVel(i,j,k,bi,bj)
249 vFld(i,j) = vVel(i,j,k,bi,bj)
250 ENDDO
251 ENDDO
252
253 C note (jmc) : Dissipation and Vort3 advection do not necesary
254 C use the same maskZ (and hFacZ) => needs 2 call(s)
255 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
256
257 CALL MOM_CALC_KE(bi,bj,k,2,uFld,vFld,KE,myThid)
258
259 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
260
261 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
262
263 IF (useAbsVorticity)
264 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
265
266 IF (momViscosity) THEN
267 C Calculate del^2 u and del^2 v for bi-harmonic term
268 IF (viscA4.NE.0.
269 & .OR. viscA4Grid.NE.0.
270 & .OR. viscC4leith.NE.0.
271 & ) THEN
272 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
273 O del2u,del2v,
274 & myThid)
275 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
276 CALL MOM_CALC_RELVORT3(
277 & bi,bj,k,del2u,del2v,hFacZ,zStar,myThid)
278 ENDIF
279 C Calculate dissipation terms for U and V equations
280 C in terms of vorticity and divergence
281 IF (viscAh.NE.0. .OR. viscA4.NE.0.
282 & .OR. viscAhGrid.NE.0. .OR. viscA4Grid.NE.0.
283 & .OR. viscC2leith.NE.0. .OR. viscC4leith.NE.0.
284 & ) THEN
285 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar,
286 O uDiss,vDiss,
287 & myThid)
288 ENDIF
289 C or in terms of tension and strain
290 IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
291 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
292 O tension,
293 I myThid)
294 CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
295 O strain,
296 I myThid)
297 CALL MOM_HDISSIP(bi,bj,k,
298 I tension,strain,hFacZ,viscAtension,viscAstrain,
299 O uDiss,vDiss,
300 I myThid)
301 ENDIF
302 ENDIF
303
304 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
305 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
306
307 C---- Zonal momentum equation starts here
308
309 C-- Vertical flux (fVer is at upper face of "u" cell)
310
311 C Eddy component of vertical flux (interior component only) -> vrF
312 IF (momViscosity.AND..NOT.implicitViscosity)
313 & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
314
315 C Combine fluxes
316 DO j=jMin,jMax
317 DO i=iMin,iMax
318 fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
319 ENDDO
320 ENDDO
321
322 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
323 DO j=2-Oly,sNy+Oly-1
324 DO i=2-Olx,sNx+Olx-1
325 gU(i,j,k,bi,bj) = uDiss(i,j)
326 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
327 & *recip_rAw(i,j,bi,bj)
328 & *(
329 & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
330 & )
331 & - phxFac*dPhiHydX(i,j)
332 ENDDO
333 ENDDO
334
335 C-- No-slip and drag BCs appear as body forces in cell abutting topography
336 IF (momViscosity.AND.no_slip_sides) THEN
337 C- No-slip BCs impose a drag at walls...
338 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid)
339 DO j=jMin,jMax
340 DO i=iMin,iMax
341 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
342 ENDDO
343 ENDDO
344 ENDIF
345
346 C- No-slip BCs impose a drag at bottom
347 IF (momViscosity.AND.bottomDragTerms) THEN
348 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
349 DO j=jMin,jMax
350 DO i=iMin,iMax
351 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
352 ENDDO
353 ENDDO
354 ENDIF
355
356 C-- Metric terms for curvilinear grid systems
357 c IF (usingSphericalPolarMTerms) THEN
358 C o Spherical polar grid metric terms
359 c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
360 c DO j=jMin,jMax
361 c DO i=iMin,iMax
362 c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
363 c ENDDO
364 c ENDDO
365 c ENDIF
366
367 C---- Meridional momentum equation starts here
368
369 C-- Vertical flux (fVer is at upper face of "v" cell)
370
371 C Eddy component of vertical flux (interior component only) -> vrF
372 IF (momViscosity.AND..NOT.implicitViscosity)
373 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
374
375 C Combine fluxes -> fVerV
376 DO j=jMin,jMax
377 DO i=iMin,iMax
378 fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
379 ENDDO
380 ENDDO
381
382 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
383 DO j=jMin,jMax
384 DO i=iMin,iMax
385 gV(i,j,k,bi,bj) = vDiss(i,j)
386 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
387 & *recip_rAs(i,j,bi,bj)
388 & *(
389 & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
390 & )
391 & - phyFac*dPhiHydY(i,j)
392 ENDDO
393 ENDDO
394
395 C-- No-slip and drag BCs appear as body forces in cell abutting topography
396 IF (momViscosity.AND.no_slip_sides) THEN
397 C- No-slip BCs impose a drag at walls...
398 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid)
399 DO j=jMin,jMax
400 DO i=iMin,iMax
401 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+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 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
411 ENDDO
412 ENDDO
413 ENDIF
414
415 C-- Metric terms for curvilinear grid systems
416 c IF (usingSphericalPolarMTerms) THEN
417 C o Spherical polar grid metric terms
418 c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
419 c DO j=jMin,jMax
420 c DO i=iMin,iMax
421 c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
422 c ENDDO
423 c ENDDO
424 c ENDIF
425
426 C-- Horizontal Coriolis terms
427 IF (useCoriolis .AND. .NOT.useCDscheme
428 & .AND. .NOT. useAbsVorticity) THEN
429 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
430 & uCf,vCf,myThid)
431 DO j=jMin,jMax
432 DO i=iMin,iMax
433 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
434 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
435 ENDDO
436 ENDDO
437 IF ( writeDiag ) THEN
438 IF (snapshot_mdsio) THEN
439 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
440 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
441 ENDIF
442 #ifdef ALLOW_MNC
443 IF (useMNC .AND. snapshot_mnc) THEN
444 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fV', uCf,
445 & offsets, myThid)
446 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'fU', vCf,
447 & offsets, myThid)
448 ENDIF
449 #endif /* ALLOW_MNC */
450 ENDIF
451 ENDIF
452
453 IF (momAdvection) THEN
454 C-- Horizontal advection of relative vorticity
455 IF (useAbsVorticity) THEN
456 CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
457 & uCf,myThid)
458 ELSE
459 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3,hFacZ,r_hFacZ,
460 & uCf,myThid)
461 ENDIF
462 c CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid)
463 DO j=jMin,jMax
464 DO i=iMin,iMax
465 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
466 ENDDO
467 ENDDO
468 IF (useAbsVorticity) THEN
469 CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
470 & vCf,myThid)
471 ELSE
472 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3,hFacZ,r_hFacZ,
473 & vCf,myThid)
474 ENDIF
475 c CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid)
476 DO j=jMin,jMax
477 DO i=iMin,iMax
478 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
479 ENDDO
480 ENDDO
481
482 IF ( writeDiag ) THEN
483 IF (snapshot_mdsio) THEN
484 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
485 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
486 ENDIF
487 #ifdef ALLOW_MNC
488 IF (useMNC .AND. snapshot_mnc) THEN
489 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zV', uCf,
490 & offsets, myThid)
491 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'zU', vCf,
492 & offsets, myThid)
493 ENDIF
494 #endif /* ALLOW_MNC */
495 ENDIF
496
497 #ifdef ALLOW_TIMEAVE
498 #ifndef HRCUBE
499 IF (taveFreq.GT.0.) THEN
500 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
501 & Nr, k, bi, bj, myThid)
502 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
503 & Nr, k, bi, bj, myThid)
504 ENDIF
505 #endif /* ndef HRCUBE */
506 #endif /* ALLOW_TIMEAVE */
507
508 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
509 IF ( .NOT. momImplVertAdv ) THEN
510 CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,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_VERTSHEAR(bi,bj,K,vVel,wVel,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 ENDIF
523
524 C-- Bernoulli term
525 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
526 DO j=jMin,jMax
527 DO i=iMin,iMax
528 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
529 ENDDO
530 ENDDO
531 CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
532 DO j=jMin,jMax
533 DO i=iMin,iMax
534 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
535 ENDDO
536 ENDDO
537 IF ( writeDiag ) THEN
538 IF (snapshot_mdsio) THEN
539 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
540 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
541 ENDIF
542 #ifdef ALLOW_MNC
543 IF (useMNC .AND. snapshot_mnc) THEN
544 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEx', uCf,
545 & offsets, myThid)
546 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj, 'KEy', vCf,
547 & offsets, myThid)
548 ENDIF
549 #endif /* ALLOW_MNC */
550 ENDIF
551
552 C-- end if momAdvection
553 ENDIF
554
555 C-- Set du/dt & dv/dt on boundaries to zero
556 DO j=jMin,jMax
557 DO i=iMin,iMax
558 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
559 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
560 ENDDO
561 ENDDO
562
563 #ifdef ALLOW_DEBUG
564 IF ( debugLevel .GE. debLevB
565 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
566 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
567 & .AND. useCubedSphereExchange ) THEN
568 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
569 & uDiss,vDiss, k, standardMessageUnit,bi,bj,myThid )
570 ENDIF
571 #endif /* ALLOW_DEBUG */
572
573 IF ( writeDiag ) THEN
574 IF (snapshot_mdsio) THEN
575 CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
576 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,
577 & myThid)
578 CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid)
579 CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid)
580 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
581 CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
582 CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
583 CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid)
584 ENDIF
585 #ifdef ALLOW_MNC
586 IF (useMNC .AND. snapshot_mnc) THEN
587 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Ds',strain,
588 & offsets, myThid)
589 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dt',tension,
590 & offsets, myThid)
591 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Du',uDiss,
592 & offsets, myThid)
593 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Dv',vDiss,
594 & offsets, myThid)
595 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'Z3',vort3,
596 & offsets, myThid)
597 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'W3',omega3,
598 & offsets, myThid)
599 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'KE',KE,
600 & offsets, myThid)
601 CALL MNC_CW_RL_W_OFFSET('D','mom_vi',bi,bj,'D', hdiv,
602 & offsets, myThid)
603 ENDIF
604 #endif /* ALLOW_MNC */
605 ENDIF
606
607 #endif /* ALLOW_MOM_VECINV */
608
609 RETURN
610 END

  ViewVC Help
Powered by ViewVC 1.1.22