1 |
C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.5 2003/04/11 13:35:03 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_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 myCurrentTime, 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 |
#include "GRID.h" |
34 |
|
35 |
C == Routine arguments == |
36 |
C fVerU - Flux of momentum in the vertical |
37 |
C fVerV direction out of the upper face of a cell K |
38 |
C ( flux into the cell above ). |
39 |
C dPhiHydX,Y :: Gradient (X & Y dir.) of Hydrostatic Potential |
40 |
C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation |
41 |
C results will be set. |
42 |
C kUp, kDown - Index for upper and lower layers. |
43 |
C myThid - Instance number for this innvocation of CALC_MOM_RHS |
44 |
_RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
45 |
_RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
46 |
_RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
47 |
_RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
48 |
_RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
49 |
_RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
50 |
INTEGER kUp,kDown |
51 |
_RL myCurrentTime |
52 |
INTEGER myIter |
53 |
INTEGER myThid |
54 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
55 |
|
56 |
C == Functions == |
57 |
LOGICAL DIFFERENT_MULTIPLE |
58 |
EXTERNAL DIFFERENT_MULTIPLE |
59 |
|
60 |
C == Local variables == |
61 |
_RL aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
62 |
_RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
63 |
_RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
_RL uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
_RL vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
66 |
_RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
67 |
_RL pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
68 |
_RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
69 |
_RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
70 |
_RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
71 |
_RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
72 |
_RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
73 |
_RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
74 |
_RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
75 |
_RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
76 |
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
77 |
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
78 |
_RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
79 |
_RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
80 |
_RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
81 |
_RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
82 |
_RL uDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
83 |
_RL vDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
84 |
C I,J,K - Loop counters |
85 |
INTEGER i,j,k |
86 |
C rVelMaskOverride - Factor for imposing special surface boundary conditions |
87 |
C ( set according to free-surface condition ). |
88 |
C hFacROpen - Lopped cell factos used tohold fraction of open |
89 |
C hFacRClosed and closed cell wall. |
90 |
_RL rVelMaskOverride |
91 |
C xxxFac - On-off tracer parameters used for switching terms off. |
92 |
_RL uDudxFac |
93 |
_RL AhDudxFac |
94 |
_RL A4DuxxdxFac |
95 |
_RL vDudyFac |
96 |
_RL AhDudyFac |
97 |
_RL A4DuyydyFac |
98 |
_RL rVelDudrFac |
99 |
_RL ArDudrFac |
100 |
_RL fuFac |
101 |
_RL phxFac |
102 |
_RL mtFacU |
103 |
_RL uDvdxFac |
104 |
_RL AhDvdxFac |
105 |
_RL A4DvxxdxFac |
106 |
_RL vDvdyFac |
107 |
_RL AhDvdyFac |
108 |
_RL A4DvyydyFac |
109 |
_RL rVelDvdrFac |
110 |
_RL ArDvdrFac |
111 |
_RL fvFac |
112 |
_RL phyFac |
113 |
_RL vForcFac |
114 |
_RL mtFacV |
115 |
INTEGER km1,kp1 |
116 |
_RL wVelBottomOverride |
117 |
LOGICAL bottomDragTerms |
118 |
_RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
119 |
_RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
120 |
_RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
121 |
_RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
122 |
|
123 |
km1=MAX(1,k-1) |
124 |
kp1=MIN(Nr,k+1) |
125 |
rVelMaskOverride=1. |
126 |
IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac |
127 |
wVelBottomOverride=1. |
128 |
IF (k.EQ.Nr) wVelBottomOverride=0. |
129 |
|
130 |
C Initialise intermediate terms |
131 |
DO J=1-OLy,sNy+OLy |
132 |
DO I=1-OLx,sNx+OLx |
133 |
aF(i,j) = 0. |
134 |
vF(i,j) = 0. |
135 |
vrF(i,j) = 0. |
136 |
uCf(i,j) = 0. |
137 |
vCf(i,j) = 0. |
138 |
mT(i,j) = 0. |
139 |
pF(i,j) = 0. |
140 |
del2u(i,j) = 0. |
141 |
del2v(i,j) = 0. |
142 |
dStar(i,j) = 0. |
143 |
zStar(i,j) = 0. |
144 |
uDiss(i,j) = 0. |
145 |
vDiss(i,j) = 0. |
146 |
vort3(i,j) = 0. |
147 |
omega3(i,j) = 0. |
148 |
ke(i,j) = 0. |
149 |
ENDDO |
150 |
ENDDO |
151 |
|
152 |
C-- Term by term tracer parmeters |
153 |
C o U momentum equation |
154 |
uDudxFac = afFacMom*1. |
155 |
AhDudxFac = vfFacMom*1. |
156 |
A4DuxxdxFac = vfFacMom*1. |
157 |
vDudyFac = afFacMom*1. |
158 |
AhDudyFac = vfFacMom*1. |
159 |
A4DuyydyFac = vfFacMom*1. |
160 |
rVelDudrFac = afFacMom*1. |
161 |
ArDudrFac = vfFacMom*1. |
162 |
mTFacU = mtFacMom*1. |
163 |
fuFac = cfFacMom*1. |
164 |
phxFac = pfFacMom*1. |
165 |
C o V momentum equation |
166 |
uDvdxFac = afFacMom*1. |
167 |
AhDvdxFac = vfFacMom*1. |
168 |
A4DvxxdxFac = vfFacMom*1. |
169 |
vDvdyFac = afFacMom*1. |
170 |
AhDvdyFac = vfFacMom*1. |
171 |
A4DvyydyFac = vfFacMom*1. |
172 |
rVelDvdrFac = afFacMom*1. |
173 |
ArDvdrFac = vfFacMom*1. |
174 |
mTFacV = mtFacMom*1. |
175 |
fvFac = cfFacMom*1. |
176 |
phyFac = pfFacMom*1. |
177 |
vForcFac = foFacMom*1. |
178 |
|
179 |
IF ( no_slip_bottom |
180 |
& .OR. bottomDragQuadratic.NE.0. |
181 |
& .OR. bottomDragLinear.NE.0.) THEN |
182 |
bottomDragTerms=.TRUE. |
183 |
ELSE |
184 |
bottomDragTerms=.FALSE. |
185 |
ENDIF |
186 |
|
187 |
C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP |
188 |
IF (staggerTimeStep) THEN |
189 |
phxFac = 0. |
190 |
phyFac = 0. |
191 |
ENDIF |
192 |
|
193 |
C-- Calculate open water fraction at vorticity points |
194 |
CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid) |
195 |
|
196 |
C---- Calculate common quantities used in both U and V equations |
197 |
C Calculate tracer cell face open areas |
198 |
DO j=1-OLy,sNy+OLy |
199 |
DO i=1-OLx,sNx+OLx |
200 |
xA(i,j) = _dyG(i,j,bi,bj) |
201 |
& *drF(k)*_hFacW(i,j,k,bi,bj) |
202 |
yA(i,j) = _dxG(i,j,bi,bj) |
203 |
& *drF(k)*_hFacS(i,j,k,bi,bj) |
204 |
ENDDO |
205 |
ENDDO |
206 |
|
207 |
C Make local copies of horizontal flow field |
208 |
DO j=1-OLy,sNy+OLy |
209 |
DO i=1-OLx,sNx+OLx |
210 |
uFld(i,j) = uVel(i,j,k,bi,bj) |
211 |
vFld(i,j) = vVel(i,j,k,bi,bj) |
212 |
ENDDO |
213 |
ENDDO |
214 |
|
215 |
C Calculate velocity field "volume transports" through tracer cell faces. |
216 |
DO j=1-OLy,sNy+OLy |
217 |
DO i=1-OLx,sNx+OLx |
218 |
uTrans(i,j) = uFld(i,j)*xA(i,j) |
219 |
vTrans(i,j) = vFld(i,j)*yA(i,j) |
220 |
ENDDO |
221 |
ENDDO |
222 |
|
223 |
CALL MOM_VI_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid) |
224 |
|
225 |
CALL MOM_VI_CALC_HDIV(bi,bj,k,uFld,vFld,hDiv,myThid) |
226 |
|
227 |
CALL MOM_VI_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid) |
228 |
|
229 |
c CALL MOM_VI_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid) |
230 |
|
231 |
IF (momViscosity) THEN |
232 |
C Calculate del^2 u and del^2 v for bi-harmonic term |
233 |
IF (viscA4.NE.0.) THEN |
234 |
CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ, |
235 |
O del2u,del2v, |
236 |
& myThid) |
237 |
CALL MOM_VI_CALC_HDIV(bi,bj,k,del2u,del2v,dStar,myThid) |
238 |
CALL MOM_VI_CALC_RELVORT3( |
239 |
& bi,bj,k,del2u,del2v,hFacZ,zStar,myThid) |
240 |
ENDIF |
241 |
C Calculate dissipation terms for U and V equations |
242 |
C in terms of vorticity and divergence |
243 |
IF (viscAh.NE.0. .OR. viscA4.NE.0.) THEN |
244 |
CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,hFacZ,dStar,zStar, |
245 |
O uDiss,vDiss, |
246 |
& myThid) |
247 |
ENDIF |
248 |
C or in terms of tension and strain |
249 |
IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN |
250 |
CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld, |
251 |
O tension, |
252 |
I myThid) |
253 |
CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ, |
254 |
O strain, |
255 |
I myThid) |
256 |
CALL MOM_HDISSIP(bi,bj,k, |
257 |
I tension,strain,hFacZ,viscAtension,viscAstrain, |
258 |
O uDiss,vDiss, |
259 |
I myThid) |
260 |
ENDIF |
261 |
ENDIF |
262 |
|
263 |
C---- Zonal momentum equation starts here |
264 |
|
265 |
C-- Vertical flux (fVer is at upper face of "u" cell) |
266 |
|
267 |
C Eddy component of vertical flux (interior component only) -> vrF |
268 |
IF (momViscosity.AND..NOT.implicitViscosity) |
269 |
& CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid) |
270 |
|
271 |
C Combine fluxes |
272 |
DO j=jMin,jMax |
273 |
DO i=iMin,iMax |
274 |
fVerU(i,j,kDown) = ArDudrFac*vrF(i,j) |
275 |
ENDDO |
276 |
ENDDO |
277 |
|
278 |
C-- Tendency is minus divergence of the fluxes + coriolis + pressure term |
279 |
DO j=2-Oly,sNy+Oly-1 |
280 |
DO i=2-Olx,sNx+Olx-1 |
281 |
gU(i,j,k,bi,bj) = uDiss(i,j) |
282 |
& -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k) |
283 |
& *recip_rAw(i,j,bi,bj) |
284 |
& *( |
285 |
& +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac |
286 |
& ) |
287 |
& - phxFac*dPhiHydX(i,j) |
288 |
ENDDO |
289 |
ENDDO |
290 |
|
291 |
C-- No-slip and drag BCs appear as body forces in cell abutting topography |
292 |
IF (momViscosity.AND.no_slip_sides) THEN |
293 |
C- No-slip BCs impose a drag at walls... |
294 |
CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,del2u,hFacZ,vF,myThid) |
295 |
DO j=jMin,jMax |
296 |
DO i=iMin,iMax |
297 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j) |
298 |
ENDDO |
299 |
ENDDO |
300 |
ENDIF |
301 |
C- No-slip BCs impose a drag at bottom |
302 |
IF (momViscosity.AND.bottomDragTerms) THEN |
303 |
CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid) |
304 |
DO j=jMin,jMax |
305 |
DO i=iMin,iMax |
306 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j) |
307 |
ENDDO |
308 |
ENDDO |
309 |
ENDIF |
310 |
|
311 |
C-- Forcing term (moved to timestep.F) |
312 |
c IF (momForcing) |
313 |
c & CALL EXTERNAL_FORCING_U( |
314 |
c I iMin,iMax,jMin,jMax,bi,bj,k, |
315 |
c I myCurrentTime,myThid) |
316 |
|
317 |
C-- Metric terms for curvilinear grid systems |
318 |
c IF (usingSphericalPolarMTerms) THEN |
319 |
C o Spherical polar grid metric terms |
320 |
c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid) |
321 |
c DO j=jMin,jMax |
322 |
c DO i=iMin,iMax |
323 |
c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j) |
324 |
c ENDDO |
325 |
c ENDDO |
326 |
c ENDIF |
327 |
|
328 |
|
329 |
C---- Meridional momentum equation starts here |
330 |
|
331 |
C-- Vertical flux (fVer is at upper face of "v" cell) |
332 |
|
333 |
C Eddy component of vertical flux (interior component only) -> vrF |
334 |
IF (momViscosity.AND..NOT.implicitViscosity) |
335 |
& CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid) |
336 |
|
337 |
C Combine fluxes -> fVerV |
338 |
DO j=jMin,jMax |
339 |
DO i=iMin,iMax |
340 |
fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j) |
341 |
ENDDO |
342 |
ENDDO |
343 |
|
344 |
C-- Tendency is minus divergence of the fluxes + coriolis + pressure term |
345 |
DO j=jMin,jMax |
346 |
DO i=iMin,iMax |
347 |
gV(i,j,k,bi,bj) = vDiss(i,j) |
348 |
& -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k) |
349 |
& *recip_rAs(i,j,bi,bj) |
350 |
& *( |
351 |
& +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac |
352 |
& ) |
353 |
& - phyFac*dPhiHydY(i,j) |
354 |
ENDDO |
355 |
ENDDO |
356 |
|
357 |
C-- No-slip and drag BCs appear as body forces in cell abutting topography |
358 |
IF (momViscosity.AND.no_slip_sides) THEN |
359 |
C- No-slip BCs impose a drag at walls... |
360 |
CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,del2v,hFacZ,vF,myThid) |
361 |
DO j=jMin,jMax |
362 |
DO i=iMin,iMax |
363 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j) |
364 |
ENDDO |
365 |
ENDDO |
366 |
ENDIF |
367 |
C- No-slip BCs impose a drag at bottom |
368 |
IF (momViscosity.AND.bottomDragTerms) THEN |
369 |
CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid) |
370 |
DO j=jMin,jMax |
371 |
DO i=iMin,iMax |
372 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j) |
373 |
ENDDO |
374 |
ENDDO |
375 |
ENDIF |
376 |
|
377 |
C-- Forcing term (moved to timestep.F) |
378 |
c IF (momForcing) |
379 |
c & CALL EXTERNAL_FORCING_V( |
380 |
c I iMin,iMax,jMin,jMax,bi,bj,k, |
381 |
c I myCurrentTime,myThid) |
382 |
|
383 |
C-- Metric terms for curvilinear grid systems |
384 |
c IF (usingSphericalPolarMTerms) THEN |
385 |
C o Spherical polar grid metric terms |
386 |
c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid) |
387 |
c DO j=jMin,jMax |
388 |
c DO i=iMin,iMax |
389 |
c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j) |
390 |
c ENDDO |
391 |
c ENDDO |
392 |
c ENDIF |
393 |
|
394 |
C-- Horizontal Coriolis terms |
395 |
IF (useCoriolis .AND. .NOT.useCDscheme) THEN |
396 |
CALL MOM_VI_CORIOLIS(bi,bj,K,uFld,vFld,omega3,r_hFacZ, |
397 |
& uCf,vCf,myThid) |
398 |
DO j=jMin,jMax |
399 |
DO i=iMin,iMax |
400 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j) |
401 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j) |
402 |
ENDDO |
403 |
ENDDO |
404 |
ENDIF |
405 |
|
406 |
IF (momAdvection) THEN |
407 |
C-- Horizontal advection of relative vorticity |
408 |
c CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,r_hFacZ,uCf,myThid) |
409 |
CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid) |
410 |
c CALL MOM_VI_U_CORIOLIS_C4(bi,bj,K,vFld,vort3,r_hFacZ,uCf,myThid) |
411 |
DO j=jMin,jMax |
412 |
DO i=iMin,iMax |
413 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j) |
414 |
ENDDO |
415 |
ENDDO |
416 |
c CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,r_hFacZ,vCf,myThid) |
417 |
CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid) |
418 |
c CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3,r_hFacZ,vCf,myThid) |
419 |
DO j=jMin,jMax |
420 |
DO i=iMin,iMax |
421 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j) |
422 |
ENDDO |
423 |
ENDDO |
424 |
|
425 |
C-- Vertical shear terms (-w*du/dr & -w*dv/dr) |
426 |
CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid) |
427 |
DO j=jMin,jMax |
428 |
DO i=iMin,iMax |
429 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j) |
430 |
ENDDO |
431 |
ENDDO |
432 |
CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid) |
433 |
DO j=jMin,jMax |
434 |
DO i=iMin,iMax |
435 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j) |
436 |
ENDDO |
437 |
ENDDO |
438 |
|
439 |
C-- Bernoulli term |
440 |
CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid) |
441 |
DO j=jMin,jMax |
442 |
DO i=iMin,iMax |
443 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j) |
444 |
ENDDO |
445 |
ENDDO |
446 |
CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid) |
447 |
DO j=jMin,jMax |
448 |
DO i=iMin,iMax |
449 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j) |
450 |
ENDDO |
451 |
ENDDO |
452 |
C-- end if momAdvection |
453 |
ENDIF |
454 |
|
455 |
C-- Set du/dt & dv/dt on boundaries to zero |
456 |
DO j=jMin,jMax |
457 |
DO i=iMin,iMax |
458 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj) |
459 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj) |
460 |
ENDDO |
461 |
ENDDO |
462 |
|
463 |
|
464 |
IF ( |
465 |
& DIFFERENT_MULTIPLE(diagFreq,myCurrentTime, |
466 |
& myCurrentTime-deltaTClock) |
467 |
& ) THEN |
468 |
CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid) |
469 |
CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid) |
470 |
CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid) |
471 |
CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid) |
472 |
CALL WRITE_LOCAL_RL('Du','I10',1,uDiss,bi,bj,k,myIter,myThid) |
473 |
CALL WRITE_LOCAL_RL('Dv','I10',1,vDiss,bi,bj,k,myIter,myThid) |
474 |
CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid) |
475 |
c CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid) |
476 |
CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid) |
477 |
CALL WRITE_LOCAL_RL('D','I10',1,hdiv,bi,bj,k,myIter,myThid) |
478 |
ENDIF |
479 |
|
480 |
RETURN |
481 |
END |