1 |
C $Header: /u/gcmpack/MITgcm/pkg/mom_fluxform/mom_fluxform.F,v 1.53 2014/12/24 19:12:56 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
CBOI |
5 |
C !TITLE: pkg/mom\_advdiff |
6 |
C !AUTHORS: adcroft@mit.edu |
7 |
C !INTRODUCTION: Flux-form Momentum Equations Package |
8 |
C |
9 |
C Package "mom\_fluxform" provides methods for calculating explicit terms |
10 |
C in the momentum equation cast in flux-form: |
11 |
C \begin{eqnarray*} |
12 |
C G^u & = & -\frac{1}{\rho} \partial_x \phi_h |
13 |
C -\nabla \cdot {\bf v} u |
14 |
C -fv |
15 |
C +\frac{1}{\rho} \nabla \cdot {\bf \tau}^x |
16 |
C + \mbox{metrics} |
17 |
C \\ |
18 |
C G^v & = & -\frac{1}{\rho} \partial_y \phi_h |
19 |
C -\nabla \cdot {\bf v} v |
20 |
C +fu |
21 |
C +\frac{1}{\rho} \nabla \cdot {\bf \tau}^y |
22 |
C + \mbox{metrics} |
23 |
C \end{eqnarray*} |
24 |
C where ${\bf v}=(u,v,w)$ and $\tau$, the stress tensor, includes surface |
25 |
C stresses as well as internal viscous stresses. |
26 |
CEOI |
27 |
|
28 |
#include "MOM_FLUXFORM_OPTIONS.h" |
29 |
#ifdef ALLOW_AUTODIFF |
30 |
# include "AUTODIFF_OPTIONS.h" |
31 |
#endif |
32 |
#ifdef ALLOW_MOM_COMMON |
33 |
# include "MOM_COMMON_OPTIONS.h" |
34 |
#endif |
35 |
|
36 |
CBOP |
37 |
C !ROUTINE: MOM_FLUXFORM |
38 |
|
39 |
C !INTERFACE: ========================================================== |
40 |
SUBROUTINE MOM_FLUXFORM( |
41 |
I bi,bj,k,iMin,iMax,jMin,jMax, |
42 |
I kappaRU, kappaRV, |
43 |
U fVerUkm, fVerVkm, |
44 |
O fVerUkp, fVerVkp, |
45 |
O guDiss, gvDiss, |
46 |
I myTime, myIter, myThid ) |
47 |
|
48 |
C !DESCRIPTION: |
49 |
C Calculates all the horizontal accelerations except for the implicit surface |
50 |
C pressure gradient and implicit vertical viscosity. |
51 |
|
52 |
C !USES: =============================================================== |
53 |
C == Global variables == |
54 |
IMPLICIT NONE |
55 |
#include "SIZE.h" |
56 |
#include "EEPARAMS.h" |
57 |
#include "PARAMS.h" |
58 |
#include "GRID.h" |
59 |
#include "DYNVARS.h" |
60 |
#include "FFIELDS.h" |
61 |
#include "SURFACE.h" |
62 |
#ifdef ALLOW_MOM_COMMON |
63 |
# include "MOM_VISC.h" |
64 |
#endif |
65 |
#ifdef ALLOW_AUTODIFF |
66 |
# include "tamc.h" |
67 |
# include "tamc_keys.h" |
68 |
# include "MOM_FLUXFORM.h" |
69 |
#endif |
70 |
|
71 |
C !INPUT PARAMETERS: =================================================== |
72 |
C bi,bj :: current tile indices |
73 |
C k :: current vertical level |
74 |
C iMin,iMax,jMin,jMax :: loop ranges |
75 |
C kappaRU :: vertical viscosity |
76 |
C kappaRV :: vertical viscosity |
77 |
C fVerUkm :: vertical advective flux of U, interface above (k-1/2) |
78 |
C fVerVkm :: vertical advective flux of V, interface above (k-1/2) |
79 |
C fVerUkp :: vertical advective flux of U, interface below (k+1/2) |
80 |
C fVerVkp :: vertical advective flux of V, interface below (k+1/2) |
81 |
C guDiss :: dissipation tendency (all explicit terms), u component |
82 |
C gvDiss :: dissipation tendency (all explicit terms), v component |
83 |
C myTime :: current time |
84 |
C myIter :: current time-step number |
85 |
C myThid :: my Thread Id number |
86 |
INTEGER bi,bj,k |
87 |
INTEGER iMin,iMax,jMin,jMax |
88 |
_RL kappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1) |
89 |
_RL kappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1) |
90 |
_RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
_RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
_RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
93 |
_RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
94 |
_RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
95 |
_RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
96 |
_RL myTime |
97 |
INTEGER myIter |
98 |
INTEGER myThid |
99 |
|
100 |
C !OUTPUT PARAMETERS: ================================================== |
101 |
C None - updates gU() and gV() in common blocks |
102 |
|
103 |
C !LOCAL VARIABLES: ==================================================== |
104 |
C i,j :: loop indices |
105 |
C vF :: viscous flux |
106 |
C v4F :: bi-harmonic viscous flux |
107 |
C cF :: Coriolis acceleration |
108 |
C mT :: Metric terms |
109 |
C fZon :: zonal fluxes |
110 |
C fMer :: meridional fluxes |
111 |
C fVrUp,fVrDw :: vertical viscous fluxes at interface k & k+1 |
112 |
INTEGER i,j |
113 |
#ifdef ALLOW_AUTODIFF_TAMC |
114 |
INTEGER imomkey |
115 |
#endif |
116 |
_RL vF(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
117 |
_RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
118 |
_RL cF(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
119 |
_RL mT(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
120 |
_RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
121 |
_RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
122 |
_RL fVrUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
123 |
_RL fVrDw(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
124 |
C afFacMom :: Tracer parameters for turning terms on and off. |
125 |
C vfFacMom |
126 |
C pfFacMom afFacMom - Advective terms |
127 |
C cfFacMom vfFacMom - Eddy viscosity terms |
128 |
C mtFacMom pfFacMom - Pressure terms |
129 |
C cfFacMom - Coriolis terms |
130 |
C foFacMom - Forcing |
131 |
C mtFacMom - Metric term |
132 |
C uDudxFac, AhDudxFac, etc ... individual term parameters for switching terms off |
133 |
_RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
134 |
_RS h0FacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
135 |
_RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
136 |
_RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
137 |
_RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
138 |
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
139 |
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
140 |
_RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
141 |
_RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
142 |
_RL rTransU(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
143 |
_RL rTransV(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
144 |
_RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
145 |
_RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
146 |
_RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
147 |
_RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
148 |
_RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
149 |
_RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
150 |
_RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
151 |
_RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
152 |
_RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
153 |
_RL uDudxFac |
154 |
_RL AhDudxFac |
155 |
_RL vDudyFac |
156 |
_RL AhDudyFac |
157 |
_RL rVelDudrFac |
158 |
_RL ArDudrFac |
159 |
_RL fuFac |
160 |
_RL mtFacU |
161 |
_RL mtNHFacU |
162 |
_RL uDvdxFac |
163 |
_RL AhDvdxFac |
164 |
_RL vDvdyFac |
165 |
_RL AhDvdyFac |
166 |
_RL rVelDvdrFac |
167 |
_RL ArDvdrFac |
168 |
_RL fvFac |
169 |
_RL mtFacV |
170 |
_RL mtNHFacV |
171 |
_RL sideMaskFac |
172 |
LOGICAL bottomDragTerms |
173 |
CEOP |
174 |
#ifdef MOM_BOUNDARY_CONSERVE |
175 |
COMMON / MOM_FLUXFORM_LOCAL / uBnd, vBnd |
176 |
_RL uBnd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
177 |
_RL vBnd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
178 |
#endif /* MOM_BOUNDARY_CONSERVE */ |
179 |
|
180 |
#ifdef ALLOW_AUTODIFF_TAMC |
181 |
act0 = k - 1 |
182 |
max0 = Nr |
183 |
act1 = bi - myBxLo(myThid) |
184 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
185 |
act2 = bj - myByLo(myThid) |
186 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
187 |
act3 = myThid - 1 |
188 |
max3 = nTx*nTy |
189 |
act4 = ikey_dynamics - 1 |
190 |
imomkey = (act0 + 1) |
191 |
& + act1*max0 |
192 |
& + act2*max0*max1 |
193 |
& + act3*max0*max1*max2 |
194 |
& + act4*max0*max1*max2*max3 |
195 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
196 |
|
197 |
C Initialise intermediate terms |
198 |
DO j=1-OLy,sNy+OLy |
199 |
DO i=1-OLx,sNx+OLx |
200 |
vF(i,j) = 0. |
201 |
v4F(i,j) = 0. |
202 |
cF(i,j) = 0. |
203 |
mT(i,j) = 0. |
204 |
fZon(i,j) = 0. |
205 |
fMer(i,j) = 0. |
206 |
fVrUp(i,j)= 0. |
207 |
fVrDw(i,j)= 0. |
208 |
rTransU(i,j)= 0. |
209 |
rTransV(i,j)= 0. |
210 |
c KE(i,j) = 0. |
211 |
hDiv(i,j) = 0. |
212 |
vort3(i,j) = 0. |
213 |
strain(i,j) = 0. |
214 |
tension(i,j)= 0. |
215 |
guDiss(i,j) = 0. |
216 |
gvDiss(i,j) = 0. |
217 |
ENDDO |
218 |
ENDDO |
219 |
|
220 |
C-- Term by term tracer parmeters |
221 |
C o U momentum equation |
222 |
uDudxFac = afFacMom*1. |
223 |
AhDudxFac = vfFacMom*1. |
224 |
vDudyFac = afFacMom*1. |
225 |
AhDudyFac = vfFacMom*1. |
226 |
rVelDudrFac = afFacMom*1. |
227 |
ArDudrFac = vfFacMom*1. |
228 |
mtFacU = mtFacMom*1. |
229 |
mtNHFacU = 1. |
230 |
fuFac = cfFacMom*1. |
231 |
C o V momentum equation |
232 |
uDvdxFac = afFacMom*1. |
233 |
AhDvdxFac = vfFacMom*1. |
234 |
vDvdyFac = afFacMom*1. |
235 |
AhDvdyFac = vfFacMom*1. |
236 |
rVelDvdrFac = afFacMom*1. |
237 |
ArDvdrFac = vfFacMom*1. |
238 |
mtFacV = mtFacMom*1. |
239 |
mtNHFacV = 1. |
240 |
fvFac = cfFacMom*1. |
241 |
|
242 |
IF (implicitViscosity) THEN |
243 |
ArDudrFac = 0. |
244 |
ArDvdrFac = 0. |
245 |
ENDIF |
246 |
|
247 |
C note: using standard stencil (no mask) results in under-estimating |
248 |
C vorticity at a no-slip boundary by a factor of 2 = sideDragFactor |
249 |
IF ( no_slip_sides ) THEN |
250 |
sideMaskFac = sideDragFactor |
251 |
ELSE |
252 |
sideMaskFac = 0. _d 0 |
253 |
ENDIF |
254 |
|
255 |
IF ( no_slip_bottom |
256 |
& .OR. selectBotDragQuadr.GE.0 |
257 |
& .OR. bottomDragLinear.NE.0.) THEN |
258 |
bottomDragTerms=.TRUE. |
259 |
ELSE |
260 |
bottomDragTerms=.FALSE. |
261 |
ENDIF |
262 |
|
263 |
C-- Calculate open water fraction at vorticity points |
264 |
CALL MOM_CALC_HFACZ( bi,bj,k,hFacZ,r_hFacZ,myThid ) |
265 |
|
266 |
C---- Calculate common quantities used in both U and V equations |
267 |
C Calculate tracer cell face open areas |
268 |
DO j=1-OLy,sNy+OLy |
269 |
DO i=1-OLx,sNx+OLx |
270 |
xA(i,j) = _dyG(i,j,bi,bj)*deepFacC(k) |
271 |
& *drF(k)*_hFacW(i,j,k,bi,bj) |
272 |
yA(i,j) = _dxG(i,j,bi,bj)*deepFacC(k) |
273 |
& *drF(k)*_hFacS(i,j,k,bi,bj) |
274 |
h0FacZ(i,j) = hFacZ(i,j) |
275 |
ENDDO |
276 |
ENDDO |
277 |
#ifdef NONLIN_FRSURF |
278 |
IF ( momViscosity .AND. no_slip_sides |
279 |
& .AND. nonlinFreeSurf.GT.0 ) THEN |
280 |
DO j=2-OLy,sNy+OLy |
281 |
DO i=2-OLx,sNx+OLx |
282 |
h0FacZ(i,j) = MIN( |
283 |
& MIN( h0FacW(i,j,k,bi,bj), h0FacW(i,j-1,k,bi,bj) ), |
284 |
& MIN( h0FacS(i,j,k,bi,bj), h0FacS(i-1,j,k,bi,bj) ) ) |
285 |
ENDDO |
286 |
ENDDO |
287 |
ENDIF |
288 |
#endif /* NONLIN_FRSURF */ |
289 |
|
290 |
C Make local copies of horizontal flow field |
291 |
DO j=1-OLy,sNy+OLy |
292 |
DO i=1-OLx,sNx+OLx |
293 |
uFld(i,j) = uVel(i,j,k,bi,bj) |
294 |
vFld(i,j) = vVel(i,j,k,bi,bj) |
295 |
ENDDO |
296 |
ENDDO |
297 |
|
298 |
C Calculate velocity field "volume transports" through tracer cell faces. |
299 |
C anelastic: transports are scaled by rhoFacC (~ mass transport) |
300 |
DO j=1-OLy,sNy+OLy |
301 |
DO i=1-OLx,sNx+OLx |
302 |
uTrans(i,j) = uFld(i,j)*xA(i,j)*rhoFacC(k) |
303 |
vTrans(i,j) = vFld(i,j)*yA(i,j)*rhoFacC(k) |
304 |
ENDDO |
305 |
ENDDO |
306 |
|
307 |
CALL MOM_CALC_KE( bi,bj,k,2,uFld,vFld,KE,myThid ) |
308 |
IF ( useVariableVisc ) THEN |
309 |
CALL MOM_CALC_HDIV( bi,bj,k,2,uFld,vFld,hDiv,myThid ) |
310 |
CALL MOM_CALC_RELVORT3( bi,bj,k,uFld,vFld,hFacZ,vort3,myThid ) |
311 |
CALL MOM_CALC_TENSION( bi,bj,k,uFld,vFld,tension,myThid ) |
312 |
CALL MOM_CALC_STRAIN( bi,bj,k,uFld,vFld,hFacZ,strain,myThid ) |
313 |
DO j=1-OLy,sNy+OLy |
314 |
DO i=1-OLx,sNx+OLx |
315 |
IF ( hFacZ(i,j).EQ.0. ) THEN |
316 |
vort3(i,j) = sideMaskFac*vort3(i,j) |
317 |
strain(i,j) = sideMaskFac*strain(i,j) |
318 |
ENDIF |
319 |
ENDDO |
320 |
ENDDO |
321 |
#ifdef ALLOW_DIAGNOSTICS |
322 |
IF ( useDiagnostics ) THEN |
323 |
CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid) |
324 |
CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid) |
325 |
CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid) |
326 |
CALL DIAGNOSTICS_FILL(strain, 'Strain ',k,1,2,bi,bj,myThid) |
327 |
ENDIF |
328 |
#endif |
329 |
ENDIF |
330 |
|
331 |
C--- First call (k=1): compute vertical adv. flux fVerUkm & fVerVkm |
332 |
IF (momAdvection.AND.k.EQ.1) THEN |
333 |
|
334 |
#ifdef MOM_BOUNDARY_CONSERVE |
335 |
CALL MOM_UV_BOUNDARY( bi, bj, k, |
336 |
I uVel, vVel, |
337 |
O uBnd(1-OLx,1-OLy,k,bi,bj), |
338 |
O vBnd(1-OLx,1-OLy,k,bi,bj), |
339 |
I myTime, myIter, myThid ) |
340 |
#endif /* MOM_BOUNDARY_CONSERVE */ |
341 |
|
342 |
C- Calculate vertical transports above U & V points (West & South face): |
343 |
|
344 |
#ifdef ALLOW_AUTODIFF_TAMC |
345 |
# ifdef NONLIN_FRSURF |
346 |
# ifndef DISABLE_RSTAR_CODE |
347 |
CADJ STORE dwtransc(:,:,bi,bj) = |
348 |
CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte |
349 |
CADJ STORE dwtransu(:,:,bi,bj) = |
350 |
CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte |
351 |
CADJ STORE dwtransv(:,:,bi,bj) = |
352 |
CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte |
353 |
# endif |
354 |
# endif /* NONLIN_FRSURF */ |
355 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
356 |
CALL MOM_CALC_RTRANS( k, bi, bj, |
357 |
O rTransU, rTransV, |
358 |
I myTime, myIter, myThid ) |
359 |
|
360 |
C- Free surface correction term (flux at k=1) |
361 |
CALL MOM_U_ADV_WU( bi,bj,k,uVel,wVel,rTransU, |
362 |
O fVerUkm, myThid ) |
363 |
|
364 |
CALL MOM_V_ADV_WV( bi,bj,k,vVel,wVel,rTransV, |
365 |
O fVerVkm, myThid ) |
366 |
|
367 |
C--- endif momAdvection & k=1 |
368 |
ENDIF |
369 |
|
370 |
C--- Calculate vertical transports (at k+1) below U & V points : |
371 |
IF (momAdvection) THEN |
372 |
CALL MOM_CALC_RTRANS( k+1, bi, bj, |
373 |
O rTransU, rTransV, |
374 |
I myTime, myIter, myThid ) |
375 |
ENDIF |
376 |
|
377 |
#ifdef MOM_BOUNDARY_CONSERVE |
378 |
IF ( momAdvection .AND. k.LT.Nr ) THEN |
379 |
CALL MOM_UV_BOUNDARY( bi, bj, k+1, |
380 |
I uVel, vVel, |
381 |
O uBnd(1-OLx,1-OLy,k+1,bi,bj), |
382 |
O vBnd(1-OLx,1-OLy,k+1,bi,bj), |
383 |
I myTime, myIter, myThid ) |
384 |
ENDIF |
385 |
#endif /* MOM_BOUNDARY_CONSERVE */ |
386 |
|
387 |
IF (momViscosity) THEN |
388 |
DO j=1-OLy,sNy+OLy |
389 |
DO i=1-OLx,sNx+OLx |
390 |
viscAh_D(i,j) = viscAhD |
391 |
viscAh_Z(i,j) = viscAhZ |
392 |
viscA4_D(i,j) = viscA4D |
393 |
viscA4_Z(i,j) = viscA4Z |
394 |
ENDDO |
395 |
ENDDO |
396 |
IF ( useVariableVisc ) THEN |
397 |
CALL MOM_CALC_VISC( bi, bj, k, |
398 |
O viscAh_Z, viscAh_D, viscA4_Z, viscA4_D, |
399 |
I hDiv, vort3, tension, strain, KE, hFacZ, |
400 |
I myThid ) |
401 |
ENDIF |
402 |
ENDIF |
403 |
|
404 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
405 |
|
406 |
C---- Zonal momentum equation starts here |
407 |
|
408 |
IF (momAdvection) THEN |
409 |
C--- Calculate mean fluxes (advection) between cells for zonal flow. |
410 |
|
411 |
#ifdef MOM_BOUNDARY_CONSERVE |
412 |
CALL MOM_U_ADV_UU( bi,bj,k,uTrans,uBnd(1-OLx,1-OLy,k,bi,bj), |
413 |
O fZon,myThid ) |
414 |
CALL MOM_U_ADV_VU( bi,bj,k,vTrans,uBnd(1-OLx,1-OLy,k,bi,bj), |
415 |
O fMer,myThid ) |
416 |
CALL MOM_U_ADV_WU( |
417 |
I bi,bj,k+1,uBnd,wVel,rTransU, |
418 |
O fVerUkp, myThid ) |
419 |
#else /* MOM_BOUNDARY_CONSERVE */ |
420 |
C-- Zonal flux (fZon is at east face of "u" cell) |
421 |
C Mean flow component of zonal flux -> fZon |
422 |
CALL MOM_U_ADV_UU( bi,bj,k,uTrans,uFld,fZon,myThid ) |
423 |
|
424 |
C-- Meridional flux (fMer is at south face of "u" cell) |
425 |
C Mean flow component of meridional flux -> fMer |
426 |
CALL MOM_U_ADV_VU( bi,bj,k,vTrans,uFld,fMer,myThid ) |
427 |
|
428 |
C-- Vertical flux (fVer is at upper face of "u" cell) |
429 |
C Mean flow component of vertical flux (at k+1) -> fVer |
430 |
CALL MOM_U_ADV_WU( |
431 |
I bi,bj,k+1,uVel,wVel,rTransU, |
432 |
O fVerUkp, myThid ) |
433 |
#endif /* MOM_BOUNDARY_CONSERVE */ |
434 |
|
435 |
C-- Tendency is minus divergence of the fluxes + coriolis + pressure term |
436 |
DO j=jMin,jMax |
437 |
DO i=iMin,iMax |
438 |
gU(i,j,k,bi,bj) = |
439 |
#ifdef OLD_UV_GEOM |
440 |
& -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/ |
441 |
& ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) ) |
442 |
#else |
443 |
& -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k) |
444 |
& *recip_rAw(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k) |
445 |
#endif |
446 |
& *( ( fZon(i,j ) - fZon(i-1,j) )*uDudxFac |
447 |
& +( fMer(i,j+1) - fMer(i, j) )*vDudyFac |
448 |
& +( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign*rVelDudrFac |
449 |
& ) |
450 |
ENDDO |
451 |
ENDDO |
452 |
|
453 |
#ifdef ALLOW_DIAGNOSTICS |
454 |
IF ( useDiagnostics ) THEN |
455 |
CALL DIAGNOSTICS_FILL( fZon, 'ADVx_Um ',k,1,2,bi,bj,myThid) |
456 |
CALL DIAGNOSTICS_FILL( fMer, 'ADVy_Um ',k,1,2,bi,bj,myThid) |
457 |
CALL DIAGNOSTICS_FILL(fVerUkm,'ADVrE_Um',k,1,2,bi,bj,myThid) |
458 |
ENDIF |
459 |
#endif |
460 |
|
461 |
#ifdef NONLIN_FRSURF |
462 |
C-- account for 3.D divergence of the flow in rStar coordinate: |
463 |
# ifndef DISABLE_RSTAR_CODE |
464 |
IF ( select_rStar.GT.0 ) THEN |
465 |
DO j=jMin,jMax |
466 |
DO i=iMin,iMax |
467 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) |
468 |
& - (rStarExpW(i,j,bi,bj) - 1. _d 0)/deltaTFreeSurf |
469 |
& *uVel(i,j,k,bi,bj) |
470 |
ENDDO |
471 |
ENDDO |
472 |
ENDIF |
473 |
IF ( select_rStar.LT.0 ) THEN |
474 |
DO j=jMin,jMax |
475 |
DO i=iMin,iMax |
476 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) |
477 |
& - rStarDhWDt(i,j,bi,bj)*uVel(i,j,k,bi,bj) |
478 |
ENDDO |
479 |
ENDDO |
480 |
ENDIF |
481 |
# endif /* DISABLE_RSTAR_CODE */ |
482 |
#endif /* NONLIN_FRSURF */ |
483 |
|
484 |
#ifdef ALLOW_ADDFLUID |
485 |
IF ( selectAddFluid.GE.1 ) THEN |
486 |
DO j=jMin,jMax |
487 |
DO i=iMin,iMax |
488 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj) |
489 |
& + uVel(i,j,k,bi,bj)*mass2rUnit*0.5 _d 0 |
490 |
& *( addMass(i-1,j,k,bi,bj) + addMass(i,j,k,bi,bj) ) |
491 |
& *_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)*recip_rhoFacC(k) |
492 |
& * recip_rAw(i,j,bi,bj)*recip_deepFac2C(k) |
493 |
ENDDO |
494 |
ENDDO |
495 |
ENDIF |
496 |
#endif /* ALLOW_ADDFLUID */ |
497 |
|
498 |
ELSE |
499 |
C- if momAdvection / else |
500 |
DO j=1-OLy,sNy+OLy |
501 |
DO i=1-OLx,sNx+OLx |
502 |
gU(i,j,k,bi,bj) = 0. _d 0 |
503 |
ENDDO |
504 |
ENDDO |
505 |
|
506 |
C- endif momAdvection. |
507 |
ENDIF |
508 |
|
509 |
IF (momViscosity) THEN |
510 |
C--- Calculate eddy fluxes (dissipation) between cells for zonal flow. |
511 |
|
512 |
C Bi-harmonic term del^2 U -> v4F |
513 |
IF ( useBiharmonicVisc ) |
514 |
& CALL MOM_U_DEL2U( bi, bj, k, uFld, hFacZ, h0FacZ, |
515 |
O v4f, myThid ) |
516 |
|
517 |
C Laplacian and bi-harmonic terms, Zonal Fluxes -> fZon |
518 |
CALL MOM_U_XVISCFLUX( bi,bj,k,uFld,v4F,fZon, |
519 |
I viscAh_D,viscA4_D,myThid ) |
520 |
|
521 |
C Laplacian and bi-harmonic termis, Merid Fluxes -> fMer |
522 |
CALL MOM_U_YVISCFLUX( bi,bj,k,uFld,v4F,hFacZ,fMer, |
523 |
I viscAh_Z,viscA4_Z,myThid ) |
524 |
|
525 |
C Eddy component of vertical flux (interior component only) -> fVrUp & fVrDw |
526 |
IF (.NOT.implicitViscosity) THEN |
527 |
CALL MOM_U_RVISCFLUX( bi,bj, k, uVel,kappaRU,fVrUp,myThid ) |
528 |
CALL MOM_U_RVISCFLUX( bi,bj,k+1,uVel,kappaRU,fVrDw,myThid ) |
529 |
ENDIF |
530 |
|
531 |
C-- Tendency is minus divergence of the fluxes |
532 |
C anelastic: hor.visc.fluxes are not scaled by rhoFac (by vert.visc.flx is) |
533 |
DO j=jMin,jMax |
534 |
DO i=iMin,iMax |
535 |
guDiss(i,j) = |
536 |
#ifdef OLD_UV_GEOM |
537 |
& -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/ |
538 |
& ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) ) |
539 |
#else |
540 |
& -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k) |
541 |
& *recip_rAw(i,j,bi,bj)*recip_deepFac2C(k) |
542 |
#endif |
543 |
& *( ( fZon(i,j ) - fZon(i-1,j) )*AhDudxFac |
544 |
& +( fMer(i,j+1) - fMer(i, j) )*AhDudyFac |
545 |
& +( fVrDw(i,j) - fVrUp(i,j) )*rkSign*ArDudrFac |
546 |
& *recip_rhoFacC(k) |
547 |
& ) |
548 |
ENDDO |
549 |
ENDDO |
550 |
|
551 |
#ifdef ALLOW_DIAGNOSTICS |
552 |
IF ( useDiagnostics ) THEN |
553 |
CALL DIAGNOSTICS_FILL(fZon, 'VISCx_Um',k,1,2,bi,bj,myThid) |
554 |
CALL DIAGNOSTICS_FILL(fMer, 'VISCy_Um',k,1,2,bi,bj,myThid) |
555 |
IF (.NOT.implicitViscosity) |
556 |
& CALL DIAGNOSTICS_FILL(fVrUp,'VISrE_Um',k,1,2,bi,bj,myThid) |
557 |
ENDIF |
558 |
#endif |
559 |
|
560 |
C-- No-slip and drag BCs appear as body forces in cell abutting topography |
561 |
IF (no_slip_sides) THEN |
562 |
C- No-slip BCs impose a drag at walls... |
563 |
CALL MOM_U_SIDEDRAG( bi, bj, k, |
564 |
I uFld, v4f, h0FacZ, |
565 |
I viscAh_Z, viscA4_Z, |
566 |
I useHarmonicVisc, useBiharmonicVisc, useVariableVisc, |
567 |
O vF, |
568 |
I myThid ) |
569 |
DO j=jMin,jMax |
570 |
DO i=iMin,iMax |
571 |
gUdiss(i,j) = gUdiss(i,j) + vF(i,j) |
572 |
ENDDO |
573 |
ENDDO |
574 |
ENDIF |
575 |
C- No-slip BCs impose a drag at bottom |
576 |
IF (bottomDragTerms) THEN |
577 |
CALL MOM_U_BOTTOMDRAG( bi, bj, k, |
578 |
I uFld, vFld, KE, kappaRU, |
579 |
O vF, |
580 |
I myThid ) |
581 |
DO j=jMin,jMax |
582 |
DO i=iMin,iMax |
583 |
gUdiss(i,j) = gUdiss(i,j) + vF(i,j) |
584 |
ENDDO |
585 |
ENDDO |
586 |
ENDIF |
587 |
|
588 |
#ifdef ALLOW_SHELFICE |
589 |
IF (useShelfIce) THEN |
590 |
CALL SHELFICE_U_DRAG( bi, bj, k, |
591 |
I uFld, vFld, KE, kappaRU, |
592 |
O vF, |
593 |
I myThid ) |
594 |
DO j=jMin,jMax |
595 |
DO i=iMin,iMax |
596 |
gUdiss(i,j) = gUdiss(i,j) + vF(i,j) |
597 |
ENDDO |
598 |
ENDDO |
599 |
ENDIF |
600 |
#endif /* ALLOW_SHELFICE */ |
601 |
|
602 |
C- endif momViscosity |
603 |
ENDIF |
604 |
|
605 |
C-- Forcing term (moved to timestep.F) |
606 |
c IF (momForcing) |
607 |
c & CALL EXTERNAL_FORCING_U( |
608 |
c I iMin,iMax,jMin,jMax,bi,bj,k, |
609 |
c I myTime,myThid) |
610 |
|
611 |
C-- Metric terms for curvilinear grid systems |
612 |
IF (useNHMTerms) THEN |
613 |
C o Non-Hydrostatic (spherical) metric terms |
614 |
CALL MOM_U_METRIC_NH( bi,bj,k,uFld,wVel,mT,myThid ) |
615 |
DO j=jMin,jMax |
616 |
DO i=iMin,iMax |
617 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mtNHFacU*mT(i,j) |
618 |
ENDDO |
619 |
ENDDO |
620 |
ENDIF |
621 |
IF ( usingSphericalPolarGrid .AND. metricTerms ) THEN |
622 |
C o Spherical polar grid metric terms |
623 |
CALL MOM_U_METRIC_SPHERE( bi,bj,k,uFld,vFld,mT,myThid ) |
624 |
DO j=jMin,jMax |
625 |
DO i=iMin,iMax |
626 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mtFacU*mT(i,j) |
627 |
ENDDO |
628 |
ENDDO |
629 |
ENDIF |
630 |
IF ( usingCylindricalGrid .AND. metricTerms ) THEN |
631 |
C o Cylindrical grid metric terms |
632 |
CALL MOM_U_METRIC_CYLINDER( bi,bj,k,uFld,vFld,mT,myThid ) |
633 |
DO j=jMin,jMax |
634 |
DO i=iMin,iMax |
635 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mtFacU*mT(i,j) |
636 |
ENDDO |
637 |
ENDDO |
638 |
ENDIF |
639 |
|
640 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
641 |
|
642 |
C---- Meridional momentum equation starts here |
643 |
|
644 |
IF (momAdvection) THEN |
645 |
|
646 |
#ifdef MOM_BOUNDARY_CONSERVE |
647 |
CALL MOM_V_ADV_UV( bi,bj,k,uTrans,vBnd(1-OLx,1-OLy,k,bi,bj), |
648 |
O fZon,myThid ) |
649 |
CALL MOM_V_ADV_VV( bi,bj,k,vTrans,vBnd(1-OLx,1-OLy,k,bi,bj), |
650 |
O fMer,myThid ) |
651 |
CALL MOM_V_ADV_WV( bi,bj,k+1,vBnd,wVel,rTransV, |
652 |
O fVerVkp, myThid ) |
653 |
#else /* MOM_BOUNDARY_CONSERVE */ |
654 |
C--- Calculate mean fluxes (advection) between cells for meridional flow. |
655 |
C Mean flow component of zonal flux -> fZon |
656 |
CALL MOM_V_ADV_UV( bi,bj,k,uTrans,vFld,fZon,myThid ) |
657 |
|
658 |
C-- Meridional flux (fMer is at north face of "v" cell) |
659 |
C Mean flow component of meridional flux -> fMer |
660 |
CALL MOM_V_ADV_VV( bi,bj,k,vTrans,vFld,fMer,myThid ) |
661 |
|
662 |
C-- Vertical flux (fVer is at upper face of "v" cell) |
663 |
C Mean flow component of vertical flux (at k+1) -> fVerV |
664 |
CALL MOM_V_ADV_WV( bi,bj,k+1,vVel,wVel,rTransV, |
665 |
O fVerVkp, myThid ) |
666 |
#endif /* MOM_BOUNDARY_CONSERVE */ |
667 |
|
668 |
C-- Tendency is minus divergence of the fluxes + coriolis + pressure term |
669 |
DO j=jMin,jMax |
670 |
DO i=iMin,iMax |
671 |
gV(i,j,k,bi,bj) = |
672 |
#ifdef OLD_UV_GEOM |
673 |
& -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/ |
674 |
& ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) ) |
675 |
#else |
676 |
& -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k) |
677 |
& *recip_rAs(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k) |
678 |
#endif |
679 |
& *( ( fZon(i+1,j) - fZon(i,j ) )*uDvdxFac |
680 |
& +( fMer(i, j) - fMer(i,j-1) )*vDvdyFac |
681 |
& +( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign*rVelDvdrFac |
682 |
& ) |
683 |
ENDDO |
684 |
ENDDO |
685 |
|
686 |
#ifdef ALLOW_DIAGNOSTICS |
687 |
IF ( useDiagnostics ) THEN |
688 |
CALL DIAGNOSTICS_FILL( fZon, 'ADVx_Vm ',k,1,2,bi,bj,myThid) |
689 |
CALL DIAGNOSTICS_FILL( fMer, 'ADVy_Vm ',k,1,2,bi,bj,myThid) |
690 |
CALL DIAGNOSTICS_FILL(fVerVkm,'ADVrE_Vm',k,1,2,bi,bj,myThid) |
691 |
ENDIF |
692 |
#endif |
693 |
|
694 |
#ifdef NONLIN_FRSURF |
695 |
C-- account for 3.D divergence of the flow in rStar coordinate: |
696 |
# ifndef DISABLE_RSTAR_CODE |
697 |
IF ( select_rStar.GT.0 ) THEN |
698 |
DO j=jMin,jMax |
699 |
DO i=iMin,iMax |
700 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) |
701 |
& - (rStarExpS(i,j,bi,bj) - 1. _d 0)/deltaTFreeSurf |
702 |
& *vVel(i,j,k,bi,bj) |
703 |
ENDDO |
704 |
ENDDO |
705 |
ENDIF |
706 |
IF ( select_rStar.LT.0 ) THEN |
707 |
DO j=jMin,jMax |
708 |
DO i=iMin,iMax |
709 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) |
710 |
& - rStarDhSDt(i,j,bi,bj)*vVel(i,j,k,bi,bj) |
711 |
ENDDO |
712 |
ENDDO |
713 |
ENDIF |
714 |
# endif /* DISABLE_RSTAR_CODE */ |
715 |
#endif /* NONLIN_FRSURF */ |
716 |
|
717 |
#ifdef ALLOW_ADDFLUID |
718 |
IF ( selectAddFluid.GE.1 ) THEN |
719 |
DO j=jMin,jMax |
720 |
DO i=iMin,iMax |
721 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj) |
722 |
& + vVel(i,j,k,bi,bj)*mass2rUnit*0.5 _d 0 |
723 |
& *( addMass(i,j-1,k,bi,bj) + addMass(i,j,k,bi,bj) ) |
724 |
& *_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)*recip_rhoFacC(k) |
725 |
& * recip_rAs(i,j,bi,bj)*recip_deepFac2C(k) |
726 |
ENDDO |
727 |
ENDDO |
728 |
ENDIF |
729 |
#endif /* ALLOW_ADDFLUID */ |
730 |
|
731 |
ELSE |
732 |
C- if momAdvection / else |
733 |
DO j=1-OLy,sNy+OLy |
734 |
DO i=1-OLx,sNx+OLx |
735 |
gV(i,j,k,bi,bj) = 0. _d 0 |
736 |
ENDDO |
737 |
ENDDO |
738 |
|
739 |
C- endif momAdvection. |
740 |
ENDIF |
741 |
|
742 |
IF (momViscosity) THEN |
743 |
C--- Calculate eddy fluxes (dissipation) between cells for meridional flow. |
744 |
C Bi-harmonic term del^2 V -> v4F |
745 |
IF ( useBiharmonicVisc ) |
746 |
& CALL MOM_V_DEL2V( bi, bj, k, vFld, hFacZ, h0FacZ, |
747 |
O v4f, myThid ) |
748 |
|
749 |
C Laplacian and bi-harmonic terms, Zonal Fluxes -> fZon |
750 |
CALL MOM_V_XVISCFLUX( bi,bj,k,vFld,v4f,hFacZ,fZon, |
751 |
I viscAh_Z,viscA4_Z,myThid ) |
752 |
|
753 |
C Laplacian and bi-harmonic termis, Merid Fluxes -> fMer |
754 |
CALL MOM_V_YVISCFLUX( bi,bj,k,vFld,v4f,fMer, |
755 |
I viscAh_D,viscA4_D,myThid ) |
756 |
|
757 |
C Eddy component of vertical flux (interior component only) -> fVrUp & fVrDw |
758 |
IF (.NOT.implicitViscosity) THEN |
759 |
CALL MOM_V_RVISCFLUX( bi,bj, k, vVel,KappaRV,fVrUp,myThid ) |
760 |
CALL MOM_V_RVISCFLUX( bi,bj,k+1,vVel,KappaRV,fVrDw,myThid ) |
761 |
ENDIF |
762 |
|
763 |
C-- Tendency is minus divergence of the fluxes + coriolis + pressure term |
764 |
C anelastic: hor.visc.fluxes are not scaled by rhoFac (by vert.visc.flx is) |
765 |
DO j=jMin,jMax |
766 |
DO i=iMin,iMax |
767 |
gvDiss(i,j) = |
768 |
#ifdef OLD_UV_GEOM |
769 |
& -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/ |
770 |
& ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) ) |
771 |
#else |
772 |
& -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k) |
773 |
& *recip_rAs(i,j,bi,bj)*recip_deepFac2C(k) |
774 |
#endif |
775 |
& *( ( fZon(i+1,j) - fZon(i,j ) )*AhDvdxFac |
776 |
& +( fMer(i, j) - fMer(i,j-1) )*AhDvdyFac |
777 |
& +( fVrDw(i,j) - fVrUp(i,j) )*rkSign*ArDvdrFac |
778 |
& *recip_rhoFacC(k) |
779 |
& ) |
780 |
ENDDO |
781 |
ENDDO |
782 |
|
783 |
#ifdef ALLOW_DIAGNOSTICS |
784 |
IF ( useDiagnostics ) THEN |
785 |
CALL DIAGNOSTICS_FILL(fZon, 'VISCx_Vm',k,1,2,bi,bj,myThid) |
786 |
CALL DIAGNOSTICS_FILL(fMer, 'VISCy_Vm',k,1,2,bi,bj,myThid) |
787 |
IF (.NOT.implicitViscosity) |
788 |
& CALL DIAGNOSTICS_FILL(fVrUp,'VISrE_Vm',k,1,2,bi,bj,myThid) |
789 |
ENDIF |
790 |
#endif |
791 |
|
792 |
C-- No-slip and drag BCs appear as body forces in cell abutting topography |
793 |
IF (no_slip_sides) THEN |
794 |
C- No-slip BCs impose a drag at walls... |
795 |
CALL MOM_V_SIDEDRAG( bi, bj, k, |
796 |
I vFld, v4f, h0FacZ, |
797 |
I viscAh_Z, viscA4_Z, |
798 |
I useHarmonicVisc, useBiharmonicVisc, useVariableVisc, |
799 |
O vF, |
800 |
I myThid ) |
801 |
DO j=jMin,jMax |
802 |
DO i=iMin,iMax |
803 |
gvDiss(i,j) = gvDiss(i,j) + vF(i,j) |
804 |
ENDDO |
805 |
ENDDO |
806 |
ENDIF |
807 |
C- No-slip BCs impose a drag at bottom |
808 |
IF (bottomDragTerms) THEN |
809 |
CALL MOM_V_BOTTOMDRAG( bi, bj, k, |
810 |
I uFld, vFld, KE, kappaRV, |
811 |
O vF, |
812 |
I myThid ) |
813 |
DO j=jMin,jMax |
814 |
DO i=iMin,iMax |
815 |
gvDiss(i,j) = gvDiss(i,j) + vF(i,j) |
816 |
ENDDO |
817 |
ENDDO |
818 |
ENDIF |
819 |
|
820 |
#ifdef ALLOW_SHELFICE |
821 |
IF (useShelfIce) THEN |
822 |
CALL SHELFICE_V_DRAG( bi, bj, k, |
823 |
I uFld, vFld, KE, kappaRV, |
824 |
O vF, |
825 |
I myThid ) |
826 |
DO j=jMin,jMax |
827 |
DO i=iMin,iMax |
828 |
gvDiss(i,j) = gvDiss(i,j) + vF(i,j) |
829 |
ENDDO |
830 |
ENDDO |
831 |
ENDIF |
832 |
#endif /* ALLOW_SHELFICE */ |
833 |
|
834 |
C- endif momViscosity |
835 |
ENDIF |
836 |
|
837 |
C-- Forcing term (moved to timestep.F) |
838 |
c IF (momForcing) |
839 |
c & CALL EXTERNAL_FORCING_V( |
840 |
c I iMin,iMax,jMin,jMax,bi,bj,k, |
841 |
c I myTime,myThid) |
842 |
|
843 |
C-- Metric terms for curvilinear grid systems |
844 |
IF (useNHMTerms) THEN |
845 |
C o Non-Hydrostatic (spherical) metric terms |
846 |
CALL MOM_V_METRIC_NH( bi,bj,k,vFld,wVel,mT,myThid ) |
847 |
DO j=jMin,jMax |
848 |
DO i=iMin,iMax |
849 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mtNHFacV*mT(i,j) |
850 |
ENDDO |
851 |
ENDDO |
852 |
ENDIF |
853 |
IF ( usingSphericalPolarGrid .AND. metricTerms ) THEN |
854 |
C o Spherical polar grid metric terms |
855 |
CALL MOM_V_METRIC_SPHERE( bi,bj,k,uFld,mT,myThid ) |
856 |
DO j=jMin,jMax |
857 |
DO i=iMin,iMax |
858 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mtFacV*mT(i,j) |
859 |
ENDDO |
860 |
ENDDO |
861 |
ENDIF |
862 |
IF ( usingCylindricalGrid .AND. metricTerms ) THEN |
863 |
C o Cylindrical grid metric terms |
864 |
CALL MOM_V_METRIC_CYLINDER( bi,bj,k,uFld,vFld,mT,myThid ) |
865 |
DO j=jMin,jMax |
866 |
DO i=iMin,iMax |
867 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mtFacV*mT(i,j) |
868 |
ENDDO |
869 |
ENDDO |
870 |
ENDIF |
871 |
|
872 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
873 |
|
874 |
C-- Coriolis term (call to CD_CODE_SCHEME has been moved to timestep.F) |
875 |
IF (.NOT.useCDscheme) THEN |
876 |
CALL MOM_U_CORIOLIS( bi,bj,k,vFld,cf,myThid ) |
877 |
DO j=jMin,jMax |
878 |
DO i=iMin,iMax |
879 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j) |
880 |
ENDDO |
881 |
ENDDO |
882 |
#ifdef ALLOW_DIAGNOSTICS |
883 |
IF ( useDiagnostics ) |
884 |
& CALL DIAGNOSTICS_FILL(cf,'Um_Cori ',k,1,2,bi,bj,myThid) |
885 |
#endif |
886 |
CALL MOM_V_CORIOLIS( bi,bj,k,uFld,cf,myThid ) |
887 |
DO j=jMin,jMax |
888 |
DO i=iMin,iMax |
889 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j) |
890 |
ENDDO |
891 |
ENDDO |
892 |
#ifdef ALLOW_DIAGNOSTICS |
893 |
IF ( useDiagnostics ) |
894 |
& CALL DIAGNOSTICS_FILL(cf,'Vm_Cori ',k,1,2,bi,bj,myThid) |
895 |
#endif |
896 |
ENDIF |
897 |
|
898 |
C-- 3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w) |
899 |
IF ( use3dCoriolis ) THEN |
900 |
CALL MOM_U_CORIOLIS_NH( bi,bj,k,wVel,cf,myThid ) |
901 |
DO j=jMin,jMax |
902 |
DO i=iMin,iMax |
903 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j) |
904 |
ENDDO |
905 |
ENDDO |
906 |
IF ( usingCurvilinearGrid ) THEN |
907 |
C- presently, non zero angleSinC array only supported with Curvilinear-Grid |
908 |
CALL MOM_V_CORIOLIS_NH( bi,bj,k,wVel,cf,myThid ) |
909 |
DO j=jMin,jMax |
910 |
DO i=iMin,iMax |
911 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j) |
912 |
ENDDO |
913 |
ENDDO |
914 |
ENDIF |
915 |
ENDIF |
916 |
|
917 |
C-- Set du/dt & dv/dt on boundaries to zero |
918 |
DO j=jMin,jMax |
919 |
DO i=iMin,iMax |
920 |
gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj) |
921 |
guDiss(i,j) = guDiss(i,j) *_maskW(i,j,k,bi,bj) |
922 |
gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj) |
923 |
gvDiss(i,j) = gvDiss(i,j) *_maskS(i,j,k,bi,bj) |
924 |
ENDDO |
925 |
ENDDO |
926 |
|
927 |
#ifdef ALLOW_DIAGNOSTICS |
928 |
IF ( useDiagnostics ) THEN |
929 |
CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid) |
930 |
CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj), |
931 |
& 'Um_Advec',k,1,2,bi,bj,myThid) |
932 |
CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj), |
933 |
& 'Vm_Advec',k,1,2,bi,bj,myThid) |
934 |
ENDIF |
935 |
#endif /* ALLOW_DIAGNOSTICS */ |
936 |
|
937 |
RETURN |
938 |
END |