/[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.14 - (show annotations) (download)
Fri Oct 10 23:00:01 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51j_post
Changes since 1.13: +5 -1 lines
Added some AD-related initialisations in mom_vecinv/ mom_fluxform/

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_fluxform/mom_fluxform.F,v 1.13 2003/10/09 04:19:20 edhill 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 C I,J,K - Loop counters
126 C rVelMaskOverride - Factor for imposing special surface boundary conditions
127 C ( set according to free-surface condition ).
128 C hFacROpen - Lopped cell factos used tohold fraction of open
129 C hFacRClosed and closed cell wall.
130 _RL rVelMaskOverride
131 C xxxFac - On-off tracer parameters used for switching terms off.
132 _RL uDudxFac
133 _RL AhDudxFac
134 _RL A4DuxxdxFac
135 _RL vDudyFac
136 _RL AhDudyFac
137 _RL A4DuyydyFac
138 _RL rVelDudrFac
139 _RL ArDudrFac
140 _RL fuFac
141 _RL phxFac
142 _RL mtFacU
143 _RL uDvdxFac
144 _RL AhDvdxFac
145 _RL A4DvxxdxFac
146 _RL vDvdyFac
147 _RL AhDvdyFac
148 _RL A4DvyydyFac
149 _RL rVelDvdrFac
150 _RL ArDvdrFac
151 _RL fvFac
152 _RL phyFac
153 _RL vForcFac
154 _RL mtFacV
155 INTEGER km1,kp1
156 _RL wVelBottomOverride
157 LOGICAL bottomDragTerms
158 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
159 CEOP
160
161 km1=MAX(1,k-1)
162 kp1=MIN(Nr,k+1)
163 rVelMaskOverride=1.
164 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
165 wVelBottomOverride=1.
166 IF (k.EQ.Nr) wVelBottomOverride=0.
167
168 C Initialise intermediate terms
169 DO J=1-OLy,sNy+OLy
170 DO I=1-OLx,sNx+OLx
171 aF(i,j) = 0.
172 vF(i,j) = 0.
173 v4F(i,j) = 0.
174 vrF(i,j) = 0.
175 cF(i,j) = 0.
176 mT(i,j) = 0.
177 pF(i,j) = 0.
178 fZon(i,j) = 0.
179 fMer(i,j) = 0.
180 rTransU(i,j) = 0.
181 rTransV(i,j) = 0.
182 fVerU(i,j,1) = 0. _d 0
183 fVerU(i,j,2) = 0. _d 0
184 fVerV(i,j,1) = 0. _d 0
185 fVerV(i,j,2) = 0. _d 0
186 ENDDO
187 ENDDO
188
189 C-- Term by term tracer parmeters
190 C o U momentum equation
191 uDudxFac = afFacMom*1.
192 AhDudxFac = vfFacMom*1.
193 A4DuxxdxFac = vfFacMom*1.
194 vDudyFac = afFacMom*1.
195 AhDudyFac = vfFacMom*1.
196 A4DuyydyFac = vfFacMom*1.
197 rVelDudrFac = afFacMom*1.
198 ArDudrFac = vfFacMom*1.
199 mTFacU = mtFacMom*1.
200 fuFac = cfFacMom*1.
201 phxFac = pfFacMom*1.
202 C o V momentum equation
203 uDvdxFac = afFacMom*1.
204 AhDvdxFac = vfFacMom*1.
205 A4DvxxdxFac = vfFacMom*1.
206 vDvdyFac = afFacMom*1.
207 AhDvdyFac = vfFacMom*1.
208 A4DvyydyFac = vfFacMom*1.
209 rVelDvdrFac = afFacMom*1.
210 ArDvdrFac = vfFacMom*1.
211 mTFacV = mtFacMom*1.
212 fvFac = cfFacMom*1.
213 phyFac = pfFacMom*1.
214 vForcFac = foFacMom*1.
215
216 IF ( no_slip_bottom
217 & .OR. bottomDragQuadratic.NE.0.
218 & .OR. bottomDragLinear.NE.0.) THEN
219 bottomDragTerms=.TRUE.
220 ELSE
221 bottomDragTerms=.FALSE.
222 ENDIF
223
224 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
225 IF (staggerTimeStep) THEN
226 phxFac = 0.
227 phyFac = 0.
228 ENDIF
229
230 C-- Calculate open water fraction at vorticity points
231 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
232
233 C---- Calculate common quantities used in both U and V equations
234 C Calculate tracer cell face open areas
235 DO j=1-OLy,sNy+OLy
236 DO i=1-OLx,sNx+OLx
237 xA(i,j) = _dyG(i,j,bi,bj)
238 & *drF(k)*_hFacW(i,j,k,bi,bj)
239 yA(i,j) = _dxG(i,j,bi,bj)
240 & *drF(k)*_hFacS(i,j,k,bi,bj)
241 ENDDO
242 ENDDO
243
244 C Make local copies of horizontal flow field
245 DO j=1-OLy,sNy+OLy
246 DO i=1-OLx,sNx+OLx
247 uFld(i,j) = uVel(i,j,k,bi,bj)
248 vFld(i,j) = vVel(i,j,k,bi,bj)
249 ENDDO
250 ENDDO
251
252 C Calculate velocity field "volume transports" through tracer cell faces.
253 DO j=1-OLy,sNy+OLy
254 DO i=1-OLx,sNx+OLx
255 uTrans(i,j) = uFld(i,j)*xA(i,j)
256 vTrans(i,j) = vFld(i,j)*yA(i,j)
257 ENDDO
258 ENDDO
259
260 CALL MOM_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)
261
262 C--- First call (k=1): compute vertical adv. flux fVerU(kUp) & fVerV(kUp)
263 IF (momAdvection.AND.k.EQ.1) THEN
264
265 C- Calculate vertical transports above U & V points (West & South face):
266 CALL MOM_CALC_RTRANS( k, bi, bj,
267 O rTransU, rTransV,
268 I myTime, myIter, myThid)
269
270 C- Free surface correction term (flux at k=1)
271 CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,rTransU,af,myThid)
272 DO j=jMin,jMax
273 DO i=iMin,iMax
274 fVerU(i,j,kUp) = af(i,j)
275 ENDDO
276 ENDDO
277
278 CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,rTransV,af,myThid)
279 DO j=jMin,jMax
280 DO i=iMin,iMax
281 fVerV(i,j,kUp) = af(i,j)
282 ENDDO
283 ENDDO
284
285 C--- endif momAdvection & k=1
286 ENDIF
287
288
289 C--- Calculate vertical transports (at k+1) below U & V points :
290 IF (momAdvection) THEN
291 CALL MOM_CALC_RTRANS( k+1, bi, bj,
292 O rTransU, rTransV,
293 I myTime, myIter, myThid)
294 ENDIF
295
296
297 C---- Zonal momentum equation starts here
298
299 C Bi-harmonic term del^2 U -> v4F
300 IF (momViscosity .AND. viscA4.NE.0. )
301 & CALL MOM_U_DEL2U(bi,bj,k,uFld,hFacZ,v4f,myThid)
302
303 C--- Calculate mean and eddy fluxes between cells for zonal flow.
304
305 C-- Zonal flux (fZon is at east face of "u" cell)
306
307 C Mean flow component of zonal flux -> aF
308 IF (momAdvection)
309 & CALL MOM_U_ADV_UU(bi,bj,k,uTrans,uFld,aF,myThid)
310
311 C Laplacian and bi-harmonic terms -> vF
312 IF (momViscosity)
313 & CALL MOM_U_XVISCFLUX(bi,bj,k,uFld,v4F,vF,myThid)
314
315 C Combine fluxes -> fZon
316 DO j=jMin,jMax
317 DO i=iMin,iMax
318 fZon(i,j) = uDudxFac*aF(i,j) + AhDudxFac*vF(i,j)
319 ENDDO
320 ENDDO
321
322 C-- Meridional flux (fMer is at south face of "u" cell)
323
324 C Mean flow component of meridional flux
325 IF (momAdvection)
326 & CALL MOM_U_ADV_VU(bi,bj,k,vTrans,uFld,aF,myThid)
327
328 C Laplacian and bi-harmonic term
329 IF (momViscosity)
330 & CALL MOM_U_YVISCFLUX(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
331
332 C Combine fluxes -> fMer
333 DO j=jMin,jMax+1
334 DO i=iMin,iMax
335 fMer(i,j) = vDudyFac*aF(i,j) + AhDudyFac*vF(i,j)
336 ENDDO
337 ENDDO
338
339 C-- Vertical flux (fVer is at upper face of "u" cell)
340
341 C Mean flow component of vertical flux (at k+1) -> aF
342 IF (momAdvection)
343 & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,rTransU,af,myThid)
344
345 C Eddy component of vertical flux (interior component only) -> vrF
346 IF (momViscosity.AND..NOT.implicitViscosity)
347 & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
348
349 C Combine fluxes
350 DO j=jMin,jMax
351 DO i=iMin,iMax
352 fVerU(i,j,kDown) = rVelDudrFac*aF(i,j) + ArDudrFac*vrF(i,j)
353 ENDDO
354 ENDDO
355
356 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
357 DO j=jMin,jMax
358 DO i=iMin,iMax
359 gU(i,j,k,bi,bj) =
360 #ifdef OLD_UV_GEOM
361 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/
362 & ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) )
363 #else
364 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
365 & *recip_rAw(i,j,bi,bj)
366 #endif
367 & *(fZon(i,j ) - fZon(i-1,j)
368 & +fMer(i,j+1) - fMer(i ,j)
369 & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
370 & )
371 & - phxFac*dPhiHydX(i,j)
372 ENDDO
373 ENDDO
374
375 #ifdef NONLIN_FRSURF
376 C-- account for 3.D divergence of the flow in rStar coordinate:
377 IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
378 DO j=jMin,jMax
379 DO i=iMin,iMax
380 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
381 & - (rStarExpW(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
382 & *uVel(i,j,k,bi,bj)
383 ENDDO
384 ENDDO
385 ENDIF
386 IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
387 DO j=jMin,jMax
388 DO i=iMin,iMax
389 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)
390 & - rStarDhWDt(i,j,bi,bj)*uVel(i,j,k,bi,bj)
391 ENDDO
392 ENDDO
393 ENDIF
394 #endif /* NONLIN_FRSURF */
395
396 C-- No-slip and drag BCs appear as body forces in cell abutting topography
397 IF (momViscosity.AND.no_slip_sides) THEN
398 C- No-slip BCs impose a drag at walls...
399 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
400 DO j=jMin,jMax
401 DO i=iMin,iMax
402 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
403 ENDDO
404 ENDDO
405 ENDIF
406 C- No-slip BCs impose a drag at bottom
407 IF (momViscosity.AND.bottomDragTerms) THEN
408 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
409 DO j=jMin,jMax
410 DO i=iMin,iMax
411 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
412 ENDDO
413 ENDDO
414 ENDIF
415
416 C-- Forcing term (moved to timestep.F)
417 c IF (momForcing)
418 c & CALL EXTERNAL_FORCING_U(
419 c I iMin,iMax,jMin,jMax,bi,bj,k,
420 c I myTime,myThid)
421
422 C-- Metric terms for curvilinear grid systems
423 IF (useNHMTerms) THEN
424 C o Non-hydrosatic metric terms
425 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
426 DO j=jMin,jMax
427 DO i=iMin,iMax
428 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
429 ENDDO
430 ENDDO
431 ENDIF
432 IF (usingSphericalPolarMTerms) THEN
433 CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)
434 DO j=jMin,jMax
435 DO i=iMin,iMax
436 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
437 ENDDO
438 ENDDO
439 ENDIF
440
441 C-- Set du/dt on boundaries to zero
442 DO j=jMin,jMax
443 DO i=iMin,iMax
444 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
445 ENDDO
446 ENDDO
447
448
449 C---- Meridional momentum equation starts here
450
451 C Bi-harmonic term del^2 V -> v4F
452 IF (momViscosity .AND. viscA4.NE.0. )
453 & CALL MOM_V_DEL2V(bi,bj,k,vFld,hFacZ,v4f,myThid)
454
455 C--- Calculate mean and eddy fluxes between cells for meridional flow.
456
457 C-- Zonal flux (fZon is at west face of "v" cell)
458
459 C Mean flow component of zonal flux -> aF
460 IF (momAdvection)
461 & CALL MOM_V_ADV_UV(bi,bj,k,uTrans,vFld,af,myThid)
462
463 C Laplacian and bi-harmonic terms -> vF
464 IF (momViscosity)
465 & CALL MOM_V_XVISCFLUX(bi,bj,k,vFld,v4f,hFacZ,vf,myThid)
466
467 C Combine fluxes -> fZon
468 DO j=jMin,jMax
469 DO i=iMin,iMax+1
470 fZon(i,j) = uDvdxFac*aF(i,j) + AhDvdxFac*vF(i,j)
471 ENDDO
472 ENDDO
473
474 C-- Meridional flux (fMer is at north face of "v" cell)
475
476 C Mean flow component of meridional flux
477 IF (momAdvection)
478 & CALL MOM_V_ADV_VV(bi,bj,k,vTrans,vFld,af,myThid)
479
480 C Laplacian and bi-harmonic term
481 IF (momViscosity)
482 & CALL MOM_V_YVISCFLUX(bi,bj,k,vFld,v4f,vf,myThid)
483
484 C Combine fluxes -> fMer
485 DO j=jMin,jMax
486 DO i=iMin,iMax
487 fMer(i,j) = vDvdyFac*aF(i,j) + AhDvdyFac*vF(i,j)
488 ENDDO
489 ENDDO
490
491 C-- Vertical flux (fVer is at upper face of "v" cell)
492
493 C o Mean flow component of vertical flux
494 IF (momAdvection)
495 & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,rTransV,af,myThid)
496
497 C Eddy component of vertical flux (interior component only) -> vrF
498 IF (momViscosity.AND..NOT.implicitViscosity)
499 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
500
501 C Combine fluxes -> fVerV
502 DO j=jMin,jMax
503 DO i=iMin,iMax
504 fVerV(i,j,kDown) = rVelDvdrFac*aF(i,j) + ArDvdrFac*vrF(i,j)
505 ENDDO
506 ENDDO
507
508 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
509 DO j=jMin,jMax
510 DO i=iMin,iMax
511 gV(i,j,k,bi,bj) =
512 #ifdef OLD_UV_GEOM
513 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/
514 & ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) )
515 #else
516 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
517 & *recip_rAs(i,j,bi,bj)
518 #endif
519 & *(fZon(i+1,j) - fZon(i,j )
520 & +fMer(i,j ) - fMer(i,j-1)
521 & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
522 & )
523 & - phyFac*dPhiHydY(i,j)
524 ENDDO
525 ENDDO
526
527 #ifdef NONLIN_FRSURF
528 C-- account for 3.D divergence of the flow in rStar coordinate:
529 IF ( momAdvection .AND. select_rStar.GT.0 ) THEN
530 DO j=jMin,jMax
531 DO i=iMin,iMax
532 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
533 & - (rStarExpS(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf
534 & *vVel(i,j,k,bi,bj)
535 ENDDO
536 ENDDO
537 ENDIF
538 IF ( momAdvection .AND. select_rStar.LT.0 ) THEN
539 DO j=jMin,jMax
540 DO i=iMin,iMax
541 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)
542 & - rStarDhSDt(i,j,bi,bj)*vVel(i,j,k,bi,bj)
543 ENDDO
544 ENDDO
545 ENDIF
546 #endif /* NONLIN_FRSURF */
547
548 C-- No-slip and drag BCs appear as body forces in cell abutting topography
549 IF (momViscosity.AND.no_slip_sides) THEN
550 C- No-slip BCs impose a drag at walls...
551 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,v4F,hFacZ,vF,myThid)
552 DO j=jMin,jMax
553 DO i=iMin,iMax
554 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
555 ENDDO
556 ENDDO
557 ENDIF
558 C- No-slip BCs impose a drag at bottom
559 IF (momViscosity.AND.bottomDragTerms) THEN
560 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
561 DO j=jMin,jMax
562 DO i=iMin,iMax
563 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
564 ENDDO
565 ENDDO
566 ENDIF
567
568 C-- Forcing term (moved to timestep.F)
569 c IF (momForcing)
570 c & CALL EXTERNAL_FORCING_V(
571 c I iMin,iMax,jMin,jMax,bi,bj,k,
572 c I myTime,myThid)
573
574 C-- Metric terms for curvilinear grid systems
575 IF (useNHMTerms) THEN
576 C o Spherical polar grid metric terms
577 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
578 DO j=jMin,jMax
579 DO i=iMin,iMax
580 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
581 ENDDO
582 ENDDO
583 ENDIF
584 IF (usingSphericalPolarMTerms) THEN
585 CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)
586 DO j=jMin,jMax
587 DO i=iMin,iMax
588 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
589 ENDDO
590 ENDDO
591 ENDIF
592
593 C-- Set dv/dt on boundaries to zero
594 DO j=jMin,jMax
595 DO i=iMin,iMax
596 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
597 ENDDO
598 ENDDO
599
600 C-- Coriolis term
601 C Note. As coded here, coriolis will not work with "thin walls"
602 c IF (useCDscheme) THEN
603 c CALL MOM_CDSCHEME(bi,bj,k,dPhiHydX,dPhiHydY,myThid)
604 c ELSE
605 IF (.NOT.useCDscheme) THEN
606 CALL MOM_U_CORIOLIS(bi,bj,k,vFld,cf,myThid)
607 DO j=jMin,jMax
608 DO i=iMin,iMax
609 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
610 ENDDO
611 ENDDO
612 CALL MOM_V_CORIOLIS(bi,bj,k,uFld,cf,myThid)
613 DO j=jMin,jMax
614 DO i=iMin,iMax
615 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j)
616 ENDDO
617 ENDDO
618 ENDIF
619
620 IF (nonHydrostatic.OR.quasiHydrostatic) THEN
621 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,cf,myThid)
622 DO j=jMin,jMax
623 DO i=iMin,iMax
624 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
625 ENDDO
626 ENDDO
627 ENDIF
628
629 RETURN
630 END

  ViewVC Help
Powered by ViewVC 1.1.22