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