/[MITgcm]/MITgcm/pkg/mom_fluxform/mom_fluxform.F
ViewVC logotype

Contents of /MITgcm/pkg/mom_fluxform/mom_fluxform.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.22 - (show annotations) (download)
Wed Jun 22 00:32:15 2005 UTC (18 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57n_post, checkpoint57l_post, checkpoint57j_post, checkpoint57o_post, checkpoint57k_post
Changes since 1.21: +3 -3 lines
"rkSign" replaces "-rkFac" (<- removed)

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_fluxform/mom_fluxform.F,v 1.21 2004/10/29 16:25:37 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
30 CBOP
31 C !ROUTINE: MOM_FLUXFORM
32
33 C !INTERFACE: ==========================================================
34 SUBROUTINE MOM_FLUXFORM(
35 I bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
36 I dPhihydX,dPhiHydY,KappaRU,KappaRV,
37 U fVerU, fVerV,
38 I myTime,myIter,myThid)
39
40 C !DESCRIPTION:
41 C Calculates all the horizontal accelerations except for the implicit surface
42 C pressure gradient and implciit vertical viscosity.
43
44 C !USES: ===============================================================
45 C == Global variables ==
46 IMPLICIT NONE
47 #include "SIZE.h"
48 #include "DYNVARS.h"
49 #include "FFIELDS.h"
50 #include "EEPARAMS.h"
51 #include "PARAMS.h"
52 #include "GRID.h"
53 #include "SURFACE.h"
54
55 C !INPUT PARAMETERS: ===================================================
56 C bi,bj :: tile indices
57 C iMin,iMax,jMin,jMAx :: loop ranges
58 C k :: vertical level
59 C kUp :: =1 or 2 for consecutive k
60 C kDown :: =2 or 1 for consecutive k
61 C dPhiHydX,Y :: Gradient (X & Y dir.) of Hydrostatic Potential
62 C KappaRU :: vertical viscosity
63 C KappaRV :: vertical viscosity
64 C fVerU :: vertical flux of U, 2 1/2 dim for pipe-lining
65 C fVerV :: vertical flux of V, 2 1/2 dim for pipe-lining
66 C myTime :: current time
67 C myIter :: current time-step number
68 C myThid :: thread number
69 INTEGER bi,bj,iMin,iMax,jMin,jMax
70 INTEGER k,kUp,kDown
71 _RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
74 _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
75 _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
76 _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
77 _RL myTime
78 INTEGER myIter
79 INTEGER myThid
80
81 C !OUTPUT PARAMETERS: ==================================================
82 C None - updates gU() and gV() in common blocks
83
84 C !LOCAL VARIABLES: ====================================================
85 C i,j :: loop indices
86 C aF :: advective flux
87 C vF :: viscous flux
88 C v4F :: bi-harmonic viscous flux
89 C vrF :: vertical viscous flux
90 C cF :: Coriolis acceleration
91 C mT :: Metric terms
92 C pF :: Pressure gradient
93 C fZon :: zonal fluxes
94 C fMer :: meridional fluxes
95 INTEGER i,j
96 _RL aF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL vF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RL cF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RL mT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL pF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 C wMaskOverride - Land sea flag override for top layer.
106 C afFacMom - Tracer parameters for turning terms
107 C vfFacMom on and off.
108 C pfFacMom afFacMom - Advective terms
109 C cfFacMom vfFacMom - Eddy viscosity terms
110 C mTFacMom pfFacMom - Pressure terms
111 C cfFacMom - Coriolis terms
112 C foFacMom - Forcing
113 C mTFacMom - Metric term
114 C uDudxFac, AhDudxFac, etc ... individual term tracer parameters
115 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123 _RL rTransU(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124 _RL rTransV(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126 c _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127 c _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128 c _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 c _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 c _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 c _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134 C I,J,K - Loop counters
135 C rVelMaskOverride - Factor for imposing special surface boundary conditions
136 C ( set according to free-surface condition ).
137 C hFacROpen - Lopped cell factos used tohold fraction of open
138 C hFacRClosed and closed cell wall.
139 _RL rVelMaskOverride
140 C xxxFac - On-off tracer parameters used for switching terms off.
141 _RL uDudxFac
142 _RL AhDudxFac
143 _RL A4DuxxdxFac
144 _RL vDudyFac
145 _RL AhDudyFac
146 _RL A4DuyydyFac
147 _RL rVelDudrFac
148 _RL ArDudrFac
149 _RL fuFac
150 _RL phxFac
151 _RL mtFacU
152 _RL uDvdxFac
153 _RL AhDvdxFac
154 _RL A4DvxxdxFac
155 _RL vDvdyFac
156 _RL AhDvdyFac
157 _RL A4DvyydyFac
158 _RL rVelDvdrFac
159 _RL ArDvdrFac
160 _RL fvFac
161 _RL phyFac
162 _RL vForcFac
163 _RL mtFacV
164 INTEGER km1,kp1
165 _RL wVelBottomOverride
166 LOGICAL bottomDragTerms
167 CEOP
168
169 km1=MAX(1,k-1)
170 kp1=MIN(Nr,k+1)
171 rVelMaskOverride=1.
172 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
173 wVelBottomOverride=1.
174 IF (k.EQ.Nr) wVelBottomOverride=0.
175
176 C Initialise intermediate terms
177 DO J=1-OLy,sNy+OLy
178 DO I=1-OLx,sNx+OLx
179 aF(i,j) = 0.
180 vF(i,j) = 0.
181 v4F(i,j) = 0.
182 vrF(i,j) = 0.
183 cF(i,j) = 0.
184 mT(i,j) = 0.
185 pF(i,j) = 0.
186 fZon(i,j) = 0.
187 fMer(i,j) = 0.
188 rTransU(i,j) = 0.
189 rTransV(i,j) = 0.
190 strain(i,j) = 0.
191 tension(i,j) = 0.
192 ENDDO
193 ENDDO
194
195 C-- Term by term tracer parmeters
196 C o U momentum equation
197 uDudxFac = afFacMom*1.
198 AhDudxFac = vfFacMom*1.
199 A4DuxxdxFac = vfFacMom*1.
200 vDudyFac = afFacMom*1.
201 AhDudyFac = vfFacMom*1.
202 A4DuyydyFac = vfFacMom*1.
203 rVelDudrFac = afFacMom*1.
204 ArDudrFac = vfFacMom*1.
205 mTFacU = mtFacMom*1.
206 fuFac = cfFacMom*1.
207 phxFac = pfFacMom*1.
208 C o V momentum equation
209 uDvdxFac = afFacMom*1.
210 AhDvdxFac = vfFacMom*1.
211 A4DvxxdxFac = vfFacMom*1.
212 vDvdyFac = afFacMom*1.
213 AhDvdyFac = vfFacMom*1.
214 A4DvyydyFac = vfFacMom*1.
215 rVelDvdrFac = afFacMom*1.
216 ArDvdrFac = vfFacMom*1.
217 mTFacV = mtFacMom*1.
218 fvFac = cfFacMom*1.
219 phyFac = pfFacMom*1.
220 vForcFac = foFacMom*1.
221
222 IF ( no_slip_bottom
223 & .OR. bottomDragQuadratic.NE.0.
224 & .OR. bottomDragLinear.NE.0.) THEN
225 bottomDragTerms=.TRUE.
226 ELSE
227 bottomDragTerms=.FALSE.
228 ENDIF
229
230 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
231 IF (staggerTimeStep) THEN
232 phxFac = 0.
233 phyFac = 0.
234 ENDIF
235
236 C-- Calculate open water fraction at vorticity points
237 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
238
239 C---- Calculate common quantities used in both U and V equations
240 C Calculate tracer cell face open areas
241 DO j=1-OLy,sNy+OLy
242 DO i=1-OLx,sNx+OLx
243 xA(i,j) = _dyG(i,j,bi,bj)
244 & *drF(k)*_hFacW(i,j,k,bi,bj)
245 yA(i,j) = _dxG(i,j,bi,bj)
246 & *drF(k)*_hFacS(i,j,k,bi,bj)
247 ENDDO
248 ENDDO
249
250 C Make local copies of horizontal flow field
251 DO j=1-OLy,sNy+OLy
252 DO i=1-OLx,sNx+OLx
253 uFld(i,j) = uVel(i,j,k,bi,bj)
254 vFld(i,j) = vVel(i,j,k,bi,bj)
255 ENDDO
256 ENDDO
257
258 C Calculate velocity field "volume transports" through tracer cell faces.
259 DO j=1-OLy,sNy+OLy
260 DO i=1-OLx,sNx+OLx
261 uTrans(i,j) = uFld(i,j)*xA(i,j)
262 vTrans(i,j) = vFld(i,j)*yA(i,j)
263 ENDDO
264 ENDDO
265
266 CALL MOM_CALC_KE(bi,bj,k,3,uFld,vFld,KE,myThid)
267
268 IF (viscAstrain.NE.0. .OR. viscAtension.NE.0.) THEN
269 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,
270 O tension,
271 I myThid)
272 CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,
273 O strain,
274 I myThid)
275 ENDIF
276
277 C--- First call (k=1): compute vertical adv. flux fVerU(kUp) & fVerV(kUp)
278 IF (momAdvection.AND.k.EQ.1) THEN
279
280 C- Calculate vertical transports above U & V points (West & South face):
281 CALL MOM_CALC_RTRANS( k, bi, bj,
282 O rTransU, rTransV,
283 I myTime, myIter, myThid)
284
285 C- Free surface correction term (flux at k=1)
286 CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,rTransU,af,myThid)
287 DO j=jMin,jMax
288 DO i=iMin,iMax
289 fVerU(i,j,kUp) = af(i,j)
290 ENDDO
291 ENDDO
292
293 CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,rTransV,af,myThid)
294 DO j=jMin,jMax
295 DO i=iMin,iMax
296 fVerV(i,j,kUp) = af(i,j)
297 ENDDO
298 ENDDO
299
300 C--- endif momAdvection & k=1
301 ENDIF
302
303
304 C--- Calculate vertical transports (at k+1) below U & V points :
305 IF (momAdvection) THEN
306 CALL MOM_CALC_RTRANS( k+1, bi, bj,
307 O rTransU, rTransV,
308 I myTime, myIter, myThid)
309 ENDIF
310
311 c IF (momViscosity) THEN
312 c & CALL MOM_CALC_VISCOSITY(bi,bj,k,
313 c I uFld,vFld,
314 c O viscAh_D,viscAh_Z,myThid)
315
316 C---- Zonal momentum equation starts here
317
318 C Bi-harmonic term del^2 U -> v4F
319 IF (momViscosity .AND. viscA4.NE.0. )
320 & CALL MOM_U_DEL2U(bi,bj,k,uFld,hFacZ,v4f,myThid)
321
322 C--- Calculate mean and eddy fluxes between cells for zonal flow.
323
324 C-- Zonal flux (fZon is at east face of "u" cell)
325
326 C Mean flow component of zonal flux -> aF
327 IF (momAdvection)
328 & CALL MOM_U_ADV_UU(bi,bj,k,uTrans,uFld,aF,myThid)
329
330 C Laplacian and bi-harmonic terms -> vF
331 IF (momViscosity)
332 & CALL MOM_U_XVISCFLUX(bi,bj,k,uFld,v4F,vF,myThid)
333
334 C Combine fluxes -> fZon
335 DO j=jMin,jMax
336 DO i=iMin,iMax
337 fZon(i,j) = uDudxFac*aF(i,j) + AhDudxFac*vF(i,j)
338 ENDDO
339 ENDDO
340
341 C-- Meridional flux (fMer is at south face of "u" cell)
342
343 C Mean flow component of meridional flux
344 IF (momAdvection)
345 & CALL MOM_U_ADV_VU(bi,bj,k,vTrans,uFld,aF,myThid)
346
347 C Laplacian and bi-harmonic term
348 IF (momViscosity)
349 & CALL MOM_U_YVISCFLUX(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
350
351 C Combine fluxes -> fMer
352 DO j=jMin,jMax+1
353 DO i=iMin,iMax
354 fMer(i,j) = vDudyFac*aF(i,j) + AhDudyFac*vF(i,j)
355 ENDDO
356 ENDDO
357
358 C-- Vertical flux (fVer is at upper face of "u" cell)
359
360 C Mean flow component of vertical flux (at k+1) -> aF
361 IF (momAdvection)
362 & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,rTransU,af,myThid)
363
364 C Eddy component of vertical flux (interior component only) -> vrF
365 IF (momViscosity.AND..NOT.implicitViscosity)
366 & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
367
368 C Combine fluxes
369 DO j=jMin,jMax
370 DO i=iMin,iMax
371 fVerU(i,j,kDown) = rVelDudrFac*aF(i,j) + ArDudrFac*vrF(i,j)
372 ENDDO
373 ENDDO
374
375 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
376 DO j=jMin,jMax
377 DO i=iMin,iMax
378 gU(i,j,k,bi,bj) =
379 #ifdef OLD_UV_GEOM
380 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/
381 & ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) )
382 #else
383 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
384 & *recip_rAw(i,j,bi,bj)
385 #endif
386 & *(fZon(i,j ) - fZon(i-1,j)
387 & +fMer(i,j+1) - fMer(i ,j)
388 & -fVerU(i,j,kUp)*rkSign + fVerU(i,j,kDown)*rkSign
389 & )
390 & - phxFac*dPhiHydX(i,j)
391 ENDDO
392 ENDDO
393
394 #ifdef NONLIN_FRSURF
395 C-- account for 3.D divergence of the flow in rStar coordinate:
396 IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
397 DO j=jMin,jMax
398 DO i=iMin,iMax
399 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
400 & - (rStarExpW(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
401 & *uVel(i,j,k,bi,bj)
402 ENDDO
403 ENDDO
404 ENDIF
405 IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
406 DO j=jMin,jMax
407 DO i=iMin,iMax
408 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
409 & - rStarDhWDt(i,j,bi,bj)*uVel(i,j,k,bi,bj)
410 ENDDO
411 ENDDO
412 ENDIF
413 #endif /* NONLIN_FRSURF */
414
415 C-- No-slip and drag BCs appear as body forces in cell abutting topography
416 IF (momViscosity.AND.no_slip_sides) THEN
417 C- No-slip BCs impose a drag at walls...
418 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
419 DO j=jMin,jMax
420 DO i=iMin,iMax
421 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
422 ENDDO
423 ENDDO
424 ENDIF
425 C- No-slip BCs impose a drag at bottom
426 IF (momViscosity.AND.bottomDragTerms) THEN
427 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
428 DO j=jMin,jMax
429 DO i=iMin,iMax
430 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
431 ENDDO
432 ENDDO
433 ENDIF
434
435 C-- Forcing term (moved to timestep.F)
436 c IF (momForcing)
437 c & CALL EXTERNAL_FORCING_U(
438 c I iMin,iMax,jMin,jMax,bi,bj,k,
439 c I myTime,myThid)
440
441 C-- Metric terms for curvilinear grid systems
442 IF (useNHMTerms) THEN
443 C o Non-hydrosatic metric terms
444 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
445 DO j=jMin,jMax
446 DO i=iMin,iMax
447 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
448 ENDDO
449 ENDDO
450 ENDIF
451 IF (usingSphericalPolarMTerms) THEN
452 CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)
453 DO j=jMin,jMax
454 DO i=iMin,iMax
455 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
456 ENDDO
457 ENDDO
458
459 ENDIF
460 IF (usingCylindricalGrid) THEN
461 CALL MOM_U_METRIC_CYLINDER(bi,bj,k,uFld,vFld,mT,myThid)
462 DO j=jMin,jMax
463 DO i=iMin,iMax
464 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
465 ENDDO
466 ENDDO
467
468 ENDIF
469 C-- Set du/dt on boundaries to zero
470 DO j=jMin,jMax
471 DO i=iMin,iMax
472 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
473 ENDDO
474 ENDDO
475
476
477 C---- Meridional momentum equation starts here
478
479 C Bi-harmonic term del^2 V -> v4F
480 IF (momViscosity .AND. viscA4.NE.0. )
481 & CALL MOM_V_DEL2V(bi,bj,k,vFld,hFacZ,v4f,myThid)
482
483 C--- Calculate mean and eddy fluxes between cells for meridional flow.
484
485 C-- Zonal flux (fZon is at west face of "v" cell)
486
487 C Mean flow component of zonal flux -> aF
488 IF (momAdvection)
489 & CALL MOM_V_ADV_UV(bi,bj,k,uTrans,vFld,af,myThid)
490
491 C Laplacian and bi-harmonic terms -> vF
492 IF (momViscosity)
493 & CALL MOM_V_XVISCFLUX(bi,bj,k,vFld,v4f,hFacZ,vf,myThid)
494
495 C Combine fluxes -> fZon
496 DO j=jMin,jMax
497 DO i=iMin,iMax+1
498 fZon(i,j) = uDvdxFac*aF(i,j) + AhDvdxFac*vF(i,j)
499 ENDDO
500 ENDDO
501
502 C-- Meridional flux (fMer is at north face of "v" cell)
503
504 C Mean flow component of meridional flux
505 IF (momAdvection)
506 & CALL MOM_V_ADV_VV(bi,bj,k,vTrans,vFld,af,myThid)
507
508 C Laplacian and bi-harmonic term
509 IF (momViscosity)
510 & CALL MOM_V_YVISCFLUX(bi,bj,k,vFld,v4f,vf,myThid)
511
512 C Combine fluxes -> fMer
513 DO j=jMin,jMax
514 DO i=iMin,iMax
515 fMer(i,j) = vDvdyFac*aF(i,j) + AhDvdyFac*vF(i,j)
516 ENDDO
517 ENDDO
518
519 C-- Vertical flux (fVer is at upper face of "v" cell)
520
521 C o Mean flow component of vertical flux
522 IF (momAdvection)
523 & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,rTransV,af,myThid)
524
525 C Eddy component of vertical flux (interior component only) -> vrF
526 IF (momViscosity.AND..NOT.implicitViscosity)
527 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
528
529 C Combine fluxes -> fVerV
530 DO j=jMin,jMax
531 DO i=iMin,iMax
532 fVerV(i,j,kDown) = rVelDvdrFac*aF(i,j) + ArDvdrFac*vrF(i,j)
533 ENDDO
534 ENDDO
535
536 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
537 DO j=jMin,jMax
538 DO i=iMin,iMax
539 gV(i,j,k,bi,bj) =
540 #ifdef OLD_UV_GEOM
541 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/
542 & ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) )
543 #else
544 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
545 & *recip_rAs(i,j,bi,bj)
546 #endif
547 & *(fZon(i+1,j) - fZon(i,j )
548 & +fMer(i,j ) - fMer(i,j-1)
549 & -fVerV(i,j,kUp)*rkSign + fVerV(i,j,kDown)*rkSign
550 & )
551 & - phyFac*dPhiHydY(i,j)
552 ENDDO
553 ENDDO
554
555 #ifdef NONLIN_FRSURF
556 C-- account for 3.D divergence of the flow in rStar coordinate:
557 IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
558 DO j=jMin,jMax
559 DO i=iMin,iMax
560 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
561 & - (rStarExpS(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
562 & *vVel(i,j,k,bi,bj)
563 ENDDO
564 ENDDO
565 ENDIF
566 IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
567 DO j=jMin,jMax
568 DO i=iMin,iMax
569 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
570 & - rStarDhSDt(i,j,bi,bj)*vVel(i,j,k,bi,bj)
571 ENDDO
572 ENDDO
573 ENDIF
574 #endif /* NONLIN_FRSURF */
575
576 C-- No-slip and drag BCs appear as body forces in cell abutting topography
577 IF (momViscosity.AND.no_slip_sides) THEN
578 C- No-slip BCs impose a drag at walls...
579 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,v4F,hFacZ,vF,myThid)
580 DO j=jMin,jMax
581 DO i=iMin,iMax
582 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
583 ENDDO
584 ENDDO
585 ENDIF
586 C- No-slip BCs impose a drag at bottom
587 IF (momViscosity.AND.bottomDragTerms) THEN
588 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
589 DO j=jMin,jMax
590 DO i=iMin,iMax
591 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
592 ENDDO
593 ENDDO
594 ENDIF
595
596 C-- Forcing term (moved to timestep.F)
597 c IF (momForcing)
598 c & CALL EXTERNAL_FORCING_V(
599 c I iMin,iMax,jMin,jMax,bi,bj,k,
600 c I myTime,myThid)
601
602 C-- Metric terms for curvilinear grid systems
603 IF (useNHMTerms) THEN
604 C o Spherical polar grid metric terms
605 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
606 DO j=jMin,jMax
607 DO i=iMin,iMax
608 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
609 ENDDO
610 ENDDO
611 ENDIF
612 IF (usingSphericalPolarMTerms) THEN
613 CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)
614 DO j=jMin,jMax
615 DO i=iMin,iMax
616 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
617 ENDDO
618 ENDDO
619 ENDIF
620 IF (usingCylindricalGrid) THEN
621 CALL MOM_V_METRIC_CYLINDER(bi,bj,k,uFld,vFld,mT,myThid)
622 DO j=jMin,jMax
623 DO i=iMin,iMax
624 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
625 ENDDO
626 ENDDO
627 ENDIF
628
629 C-- Set dv/dt on boundaries to zero
630 DO j=jMin,jMax
631 DO i=iMin,iMax
632 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
633 ENDDO
634 ENDDO
635
636 C-- Coriolis term
637 C Note. As coded here, coriolis will not work with "thin walls"
638 c IF (useCDscheme) THEN
639 c CALL MOM_CDSCHEME(bi,bj,k,dPhiHydX,dPhiHydY,myThid)
640 c ELSE
641 IF (.NOT.useCDscheme) THEN
642 CALL MOM_U_CORIOLIS(bi,bj,k,vFld,cf,myThid)
643 DO j=jMin,jMax
644 DO i=iMin,iMax
645 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
646 ENDDO
647 ENDDO
648 CALL MOM_V_CORIOLIS(bi,bj,k,uFld,cf,myThid)
649 DO j=jMin,jMax
650 DO i=iMin,iMax
651 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j)
652 ENDDO
653 ENDDO
654 ENDIF
655
656 IF (nonHydrostatic.OR.quasiHydrostatic) THEN
657 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,cf,myThid)
658 DO j=jMin,jMax
659 DO i=iMin,iMax
660 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
661 ENDDO
662 ENDDO
663 ENDIF
664
665 RETURN
666 END

  ViewVC Help
Powered by ViewVC 1.1.22