/[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.4 - (show annotations) (download)
Fri Sep 28 16:49:54 2001 UTC (22 years, 7 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint44e_post, checkpoint46l_post, checkpoint46g_pre, release1_p13_pre, checkpoint46f_post, checkpoint44f_post, checkpoint46b_post, checkpoint43a-release1mods, release1_p13, checkpoint46l_pre, chkpt44d_post, release1_p8, release1_p9, release1_p1, release1_p2, release1_p3, release1_p4, release1_p5, release1_p6, release1_p7, checkpoint44e_pre, release1_b1, checkpoint43, release1_chkpt44d_post, release1_p11, icebear5, icebear4, icebear3, icebear2, checkpoint46d_pre, release1-branch_tutorials, checkpoint45d_post, checkpoint46j_pre, chkpt44a_post, checkpoint44h_pre, checkpoint46a_post, checkpoint46j_post, checkpoint46k_post, chkpt44c_pre, checkpoint45a_post, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, checkpoint44g_post, checkpoint46e_pre, checkpoint45b_post, checkpoint46b_pre, release1-branch-end, release1_final_v1, checkpoint46c_pre, checkpoint46, checkpoint44b_post, checkpoint46h_pre, checkpoint46a_pre, checkpoint45c_post, ecco_ice2, ecco_ice1, checkpoint44h_post, checkpoint46g_post, release1_p12_pre, ecco_c44_e22, ecco_c44_e25, chkpt44a_pre, checkpoint46i_post, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, checkpoint46c_post, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, checkpoint46e_post, release1_beta1, checkpoint44b_pre, checkpoint44, checkpoint45, checkpoint46h_post, chkpt44c_post, checkpoint44f_pre, checkpoint46d_post, release1-branch_branchpoint
Branch point for: c24_e25_ice, release1_final, release1-branch, release1, ecco-branch, release1_50yr, icebear, release1_coupled
Changes since 1.3: +2 -3 lines
Changes for structuing protex document.

1 C $Header: /u/gcmpack/models/MITgcmUV/pkg/mom_fluxform/mom_fluxform.F,v 1.3 2001/09/26 19:05:21 adcroft 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 "CPP_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 phi_hyd,KappaRU,KappaRV,
37 U fVerU, fVerV,
38 I myCurrentTime,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 phi_hyd :: hydrostatic pressure (perturbation)
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 myCurrentTime :: 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 phi_hyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
72 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
73 _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
74 _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
75 _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
76 _RL myCurrentTime
77 INTEGER myIter
78 INTEGER myThid
79
80 C !OUTPUT PARAMETERS: ==================================================
81 C None - updates gU() and gV() in common blocks
82
83 C !LOCAL VARIABLES: ====================================================
84 C i,j :: loop indices
85 C aF :: advective flux
86 C vF :: viscous flux
87 C v4F :: bi-harmonic viscous flux
88 C vrF :: vertical viscous flux
89 C cF :: Coriolis acceleration
90 C mT :: Metric terms
91 C pF :: Pressure gradient
92 C fZon :: zonal fluxes
93 C fMer :: meridional fluxes
94 INTEGER i,j
95 _RL aF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL vF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL cF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RL mT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RL pF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 C wMaskOverride - Land sea flag override for top layer.
105 C afFacMom - Tracer parameters for turning terms
106 C vfFacMom on and off.
107 C pfFacMom afFacMom - Advective terms
108 C cfFacMom vfFacMom - Eddy viscosity terms
109 C mTFacMom pfFacMom - Pressure terms
110 C cfFacMom - Coriolis terms
111 C foFacMom - Forcing
112 C mTFacMom - Metric term
113 C uDudxFac, AhDudxFac, etc ... individual term tracer parameters
114 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115 _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 C I,J,K - Loop counters
123 C rVelMaskOverride - Factor for imposing special surface boundary conditions
124 C ( set according to free-surface condition ).
125 C hFacROpen - Lopped cell factos used tohold fraction of open
126 C hFacRClosed and closed cell wall.
127 _RL rVelMaskOverride
128 C xxxFac - On-off tracer parameters used for switching terms off.
129 _RL uDudxFac
130 _RL AhDudxFac
131 _RL A4DuxxdxFac
132 _RL vDudyFac
133 _RL AhDudyFac
134 _RL A4DuyydyFac
135 _RL rVelDudrFac
136 _RL ArDudrFac
137 _RL fuFac
138 _RL phxFac
139 _RL mtFacU
140 _RL uDvdxFac
141 _RL AhDvdxFac
142 _RL A4DvxxdxFac
143 _RL vDvdyFac
144 _RL AhDvdyFac
145 _RL A4DvyydyFac
146 _RL rVelDvdrFac
147 _RL ArDvdrFac
148 _RL fvFac
149 _RL phyFac
150 _RL vForcFac
151 _RL mtFacV
152 INTEGER km1,kp1
153 _RL wVelBottomOverride
154 LOGICAL bottomDragTerms
155 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156 CEOP
157
158 km1=MAX(1,k-1)
159 kp1=MIN(Nr,k+1)
160 rVelMaskOverride=1.
161 IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
162 wVelBottomOverride=1.
163 IF (k.EQ.Nr) wVelBottomOverride=0.
164
165 C Initialise intermediate terms
166 DO J=1-OLy,sNy+OLy
167 DO I=1-OLx,sNx+OLx
168 aF(i,j) = 0.
169 vF(i,j) = 0.
170 v4F(i,j) = 0.
171 vrF(i,j) = 0.
172 cF(i,j) = 0.
173 mT(i,j) = 0.
174 pF(i,j) = 0.
175 fZon(i,j) = 0.
176 fMer(i,j) = 0.
177 ENDDO
178 ENDDO
179
180 C-- Term by term tracer parmeters
181 C o U momentum equation
182 uDudxFac = afFacMom*1.
183 AhDudxFac = vfFacMom*1.
184 A4DuxxdxFac = vfFacMom*1.
185 vDudyFac = afFacMom*1.
186 AhDudyFac = vfFacMom*1.
187 A4DuyydyFac = vfFacMom*1.
188 rVelDudrFac = afFacMom*1.
189 ArDudrFac = vfFacMom*1.
190 mTFacU = mtFacMom*1.
191 fuFac = cfFacMom*1.
192 phxFac = pfFacMom*1.
193 C o V momentum equation
194 uDvdxFac = afFacMom*1.
195 AhDvdxFac = vfFacMom*1.
196 A4DvxxdxFac = vfFacMom*1.
197 vDvdyFac = afFacMom*1.
198 AhDvdyFac = vfFacMom*1.
199 A4DvyydyFac = vfFacMom*1.
200 rVelDvdrFac = afFacMom*1.
201 ArDvdrFac = vfFacMom*1.
202 mTFacV = mtFacMom*1.
203 fvFac = cfFacMom*1.
204 phyFac = pfFacMom*1.
205 vForcFac = foFacMom*1.
206
207 IF ( no_slip_bottom
208 & .OR. bottomDragQuadratic.NE.0.
209 & .OR. bottomDragLinear.NE.0.) THEN
210 bottomDragTerms=.TRUE.
211 ELSE
212 bottomDragTerms=.FALSE.
213 ENDIF
214
215 C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
216 IF (staggerTimeStep) THEN
217 phxFac = 0.
218 phyFac = 0.
219 ENDIF
220
221 C-- Calculate open water fraction at vorticity points
222 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
223
224 C---- Calculate common quantities used in both U and V equations
225 C Calculate tracer cell face open areas
226 DO j=1-OLy,sNy+OLy
227 DO i=1-OLx,sNx+OLx
228 xA(i,j) = _dyG(i,j,bi,bj)
229 & *drF(k)*_hFacW(i,j,k,bi,bj)
230 yA(i,j) = _dxG(i,j,bi,bj)
231 & *drF(k)*_hFacS(i,j,k,bi,bj)
232 ENDDO
233 ENDDO
234
235 C Make local copies of horizontal flow field
236 DO j=1-OLy,sNy+OLy
237 DO i=1-OLx,sNx+OLx
238 uFld(i,j) = uVel(i,j,k,bi,bj)
239 vFld(i,j) = vVel(i,j,k,bi,bj)
240 ENDDO
241 ENDDO
242
243 C Calculate velocity field "volume transports" through tracer cell faces.
244 DO j=1-OLy,sNy+OLy
245 DO i=1-OLx,sNx+OLx
246 uTrans(i,j) = uFld(i,j)*xA(i,j)
247 vTrans(i,j) = vFld(i,j)*yA(i,j)
248 ENDDO
249 ENDDO
250
251 CALL MOM_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)
252
253 C---- Zonal momentum equation starts here
254
255 C Bi-harmonic term del^2 U -> v4F
256 IF (momViscosity)
257 & CALL MOM_U_DEL2U(bi,bj,k,uFld,hFacZ,v4f,myThid)
258
259 C--- Calculate mean and eddy fluxes between cells for zonal flow.
260
261 C-- Zonal flux (fZon is at east face of "u" cell)
262
263 C Mean flow component of zonal flux -> aF
264 IF (momAdvection)
265 & CALL MOM_U_ADV_UU(bi,bj,k,uTrans,uFld,aF,myThid)
266
267 C Laplacian and bi-harmonic terms -> vF
268 IF (momViscosity)
269 & CALL MOM_U_XVISCFLUX(bi,bj,k,uFld,v4F,vF,myThid)
270
271 C Combine fluxes -> fZon
272 DO j=jMin,jMax
273 DO i=iMin,iMax
274 fZon(i,j) = uDudxFac*aF(i,j) + AhDudxFac*vF(i,j)
275 ENDDO
276 ENDDO
277
278 C-- Meridional flux (fMer is at south face of "u" cell)
279
280 C Mean flow component of meridional flux
281 IF (momAdvection)
282 & CALL MOM_U_ADV_VU(bi,bj,k,vTrans,uFld,aF,myThid)
283
284 C Laplacian and bi-harmonic term
285 IF (momViscosity)
286 & CALL MOM_U_YVISCFLUX(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
287
288 C Combine fluxes -> fMer
289 DO j=jMin,jMax
290 DO i=iMin,iMax
291 fMer(i,j) = vDudyFac*aF(i,j) + AhDudyFac*vF(i,j)
292 ENDDO
293 ENDDO
294
295 C-- Vertical flux (fVer is at upper face of "u" cell)
296
297 C-- Free surface correction term (flux at k=1)
298 IF (momAdvection.AND.k.EQ.1) THEN
299 CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,af,myThid)
300 DO j=jMin,jMax
301 DO i=iMin,iMax
302 fVerU(i,j,kUp) = af(i,j)
303 ENDDO
304 ENDDO
305 ENDIF
306 C Mean flow component of vertical flux (at k+1) -> aF
307 IF (momAdvection)
308 & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,af,myThid)
309
310 C Eddy component of vertical flux (interior component only) -> vrF
311 IF (momViscosity.AND..NOT.implicitViscosity)
312 & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
313
314 C Combine fluxes
315 DO j=jMin,jMax
316 DO i=iMin,iMax
317 fVerU(i,j,kDown) = rVelDudrFac*aF(i,j) + ArDudrFac*vrF(i,j)
318 ENDDO
319 ENDDO
320
321 C--- Hydrostatic term ( -1/rhoConst . dphi/dx )
322 IF (momPressureForcing) THEN
323 DO j=jMin,jMax
324 DO i=iMin,iMax
325 pf(i,j) = - _recip_dxC(i,j,bi,bj)
326 & *(phi_hyd(i,j,k)-phi_hyd(i-1,j,k))
327 ENDDO
328 ENDDO
329 ENDIF
330
331 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
332 DO j=jMin,jMax
333 DO i=iMin,iMax
334 gU(i,j,k,bi,bj) =
335 #ifdef OLD_UV_GEOM
336 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/
337 & ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) )
338 #else
339 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
340 & *recip_rAw(i,j,bi,bj)
341 #endif
342 & *(fZon(i,j ) - fZon(i-1,j)
343 & +fMer(i,j+1) - fMer(i ,j)
344 & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
345 & )
346 & _PHM( +phxFac * pf(i,j) )
347 ENDDO
348 ENDDO
349
350 C-- No-slip and drag BCs appear as body forces in cell abutting topography
351 IF (momViscosity.AND.no_slip_sides) THEN
352 C- No-slip BCs impose a drag at walls...
353 CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
354 DO j=jMin,jMax
355 DO i=iMin,iMax
356 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
357 ENDDO
358 ENDDO
359 ENDIF
360 C- No-slip BCs impose a drag at bottom
361 IF (momViscosity.AND.bottomDragTerms) THEN
362 CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
363 DO j=jMin,jMax
364 DO i=iMin,iMax
365 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
366 ENDDO
367 ENDDO
368 ENDIF
369
370 C-- Forcing term
371 IF (momForcing)
372 & CALL EXTERNAL_FORCING_U(
373 I iMin,iMax,jMin,jMax,bi,bj,k,
374 I myCurrentTime,myThid)
375
376 C-- Metric terms for curvilinear grid systems
377 IF (usingSphericalPolarMTerms) THEN
378 C o Spherical polar grid metric terms
379 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
380 DO j=jMin,jMax
381 DO i=iMin,iMax
382 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
383 ENDDO
384 ENDDO
385 CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)
386 DO j=jMin,jMax
387 DO i=iMin,iMax
388 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
389 ENDDO
390 ENDDO
391 ENDIF
392
393 C-- Set du/dt on boundaries to zero
394 DO j=jMin,jMax
395 DO i=iMin,iMax
396 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
397 ENDDO
398 ENDDO
399
400
401 C---- Meridional momentum equation starts here
402
403 C Bi-harmonic term del^2 V -> v4F
404 IF (momViscosity)
405 & CALL MOM_V_DEL2V(bi,bj,k,vFld,hFacZ,v4f,myThid)
406
407 C--- Calculate mean and eddy fluxes between cells for meridional flow.
408
409 C-- Zonal flux (fZon is at west face of "v" cell)
410
411 C Mean flow component of zonal flux -> aF
412 IF (momAdvection)
413 & CALL MOM_V_ADV_UV(bi,bj,k,uTrans,vFld,af,myThid)
414
415 C Laplacian and bi-harmonic terms -> vF
416 IF (momViscosity)
417 & CALL MOM_V_XVISCFLUX(bi,bj,k,vFld,v4f,hFacZ,vf,myThid)
418
419 C Combine fluxes -> fZon
420 DO j=jMin,jMax
421 DO i=iMin,iMax
422 fZon(i,j) = uDvdxFac*aF(i,j) + AhDvdxFac*vF(i,j)
423 ENDDO
424 ENDDO
425
426 C-- Meridional flux (fMer is at north face of "v" cell)
427
428 C Mean flow component of meridional flux
429 IF (momAdvection)
430 & CALL MOM_V_ADV_VV(bi,bj,k,vTrans,vFld,af,myThid)
431
432 C Laplacian and bi-harmonic term
433 IF (momViscosity)
434 & CALL MOM_V_YVISCFLUX(bi,bj,k,vFld,v4f,vf,myThid)
435
436 C Combine fluxes -> fMer
437 DO j=jMin,jMax
438 DO i=iMin,iMax
439 fMer(i,j) = vDvdyFac*aF(i,j) + AhDvdyFac*vF(i,j)
440 ENDDO
441 ENDDO
442
443 C-- Vertical flux (fVer is at upper face of "v" cell)
444
445 C-- Free surface correction term (flux at k=1)
446 IF (momAdvection.AND.k.EQ.1) THEN
447 CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,af,myThid)
448 DO j=jMin,jMax
449 DO i=iMin,iMax
450 fVerV(i,j,kUp) = af(i,j)
451 ENDDO
452 ENDDO
453 ENDIF
454 C o Mean flow component of vertical flux
455 IF (momAdvection)
456 & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,af,myThid)
457
458 C Eddy component of vertical flux (interior component only) -> vrF
459 IF (momViscosity.AND..NOT.implicitViscosity)
460 & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
461
462 C Combine fluxes -> fVerV
463 DO j=jMin,jMax
464 DO i=iMin,iMax
465 fVerV(i,j,kDown) = rVelDvdrFac*aF(i,j) + ArDvdrFac*vrF(i,j)
466 ENDDO
467 ENDDO
468
469 C--- Hydorstatic term (-1/rhoConst . dphi/dy )
470 IF (momPressureForcing) THEN
471 DO j=jMin,jMax
472 DO i=iMin,iMax
473 pF(i,j) = -_recip_dyC(i,j,bi,bj)
474 & *(phi_hyd(i,j,k)-phi_hyd(i,j-1,k))
475 ENDDO
476 ENDDO
477 ENDIF
478
479 C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
480 DO j=jMin,jMax
481 DO i=iMin,iMax
482 gV(i,j,k,bi,bj) =
483 #ifdef OLD_UV_GEOM
484 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/
485 & ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) )
486 #else
487 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
488 & *recip_rAs(i,j,bi,bj)
489 #endif
490 & *(fZon(i+1,j) - fZon(i,j )
491 & +fMer(i,j ) - fMer(i,j-1)
492 & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
493 & )
494 & _PHM( +phyFac*pf(i,j) )
495 ENDDO
496 ENDDO
497
498 C-- No-slip and drag BCs appear as body forces in cell abutting topography
499 IF (momViscosity.AND.no_slip_sides) THEN
500 C- No-slip BCs impose a drag at walls...
501 CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,v4F,hFacZ,vF,myThid)
502 DO j=jMin,jMax
503 DO i=iMin,iMax
504 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
505 ENDDO
506 ENDDO
507 ENDIF
508 C- No-slip BCs impose a drag at bottom
509 IF (momViscosity.AND.bottomDragTerms) THEN
510 CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
511 DO j=jMin,jMax
512 DO i=iMin,iMax
513 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
514 ENDDO
515 ENDDO
516 ENDIF
517
518 C-- Forcing term
519 IF (momForcing)
520 & CALL EXTERNAL_FORCING_V(
521 I iMin,iMax,jMin,jMax,bi,bj,k,
522 I myCurrentTime,myThid)
523
524 C-- Metric terms for curvilinear grid systems
525 IF (usingSphericalPolarMTerms) THEN
526 C o Spherical polar grid metric terms
527 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
528 DO j=jMin,jMax
529 DO i=iMin,iMax
530 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
531 ENDDO
532 ENDDO
533 CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)
534 DO j=jMin,jMax
535 DO i=iMin,iMax
536 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
537 ENDDO
538 ENDDO
539 ENDIF
540
541 C-- Set dv/dt on boundaries to zero
542 DO j=jMin,jMax
543 DO i=iMin,iMax
544 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
545 ENDDO
546 ENDDO
547
548 C-- Coriolis term
549 C Note. As coded here, coriolis will not work with "thin walls"
550 #ifdef INCLUDE_CD_CODE
551 CALL MOM_CDSCHEME(bi,bj,k,phi_hyd,myThid)
552 #else
553 CALL MOM_U_CORIOLIS(bi,bj,k,vFld,cf,myThid)
554 DO j=jMin,jMax
555 DO i=iMin,iMax
556 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
557 ENDDO
558 ENDDO
559 CALL MOM_V_CORIOLIS(bi,bj,k,uFld,cf,myThid)
560 DO j=jMin,jMax
561 DO i=iMin,iMax
562 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j)
563 ENDDO
564 ENDDO
565 #endif /* INCLUDE_CD_CODE */
566
567 RETURN
568 END

  ViewVC Help
Powered by ViewVC 1.1.22