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

Contents of /MITgcm/pkg/mom_vecinv/mom_vecinv.F

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


Revision 1.81 - (show annotations) (download)
Sat Apr 29 16:11:38 2017 UTC (7 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.80: +5 -3 lines
fix previous modif

1 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.80 2017/04/28 17:17:14 mlosch Exp $
2 C $Name: $
3
4 #include "MOM_VECINV_OPTIONS.h"
5 #ifdef ALLOW_AUTODIFF
6 # include "AUTODIFF_OPTIONS.h"
7 #endif
8 #ifdef ALLOW_MOM_COMMON
9 # include "MOM_COMMON_OPTIONS.h"
10 #endif
11
12 SUBROUTINE MOM_VECINV(
13 I bi,bj,k,iMin,iMax,jMin,jMax,
14 I kappaRU, kappaRV,
15 I fVerUkm, fVerVkm,
16 O fVerUkp, fVerVkp,
17 O guDiss, gvDiss,
18 I myTime, myIter, myThid )
19 C *==========================================================*
20 C | S/R MOM_VECINV |
21 C | o Form the right hand-side of the momentum equation. |
22 C *==========================================================*
23 C | Terms are evaluated one layer at a time working from |
24 C | the bottom to the top. The vertically integrated |
25 C | barotropic flow tendency term is evluated by summing the |
26 C | tendencies. |
27 C | Notes: |
28 C | We have not sorted out an entirely satisfactory formula |
29 C | for the diffusion equation bc with lopping. The present |
30 C | form produces a diffusive flux that does not scale with |
31 C | open-area. Need to do something to solidfy this and to |
32 C | deal "properly" with thin walls. |
33 C *==========================================================*
34 IMPLICIT NONE
35
36 C == Global variables ==
37 #include "SIZE.h"
38 #include "EEPARAMS.h"
39 #include "PARAMS.h"
40 #include "GRID.h"
41 #include "SURFACE.h"
42 #include "DYNVARS.h"
43 #ifdef ALLOW_MOM_COMMON
44 # include "MOM_VISC.h"
45 #endif
46 #ifdef ALLOW_TIMEAVE
47 # include "TIMEAVE_STATV.h"
48 #endif
49 #ifdef ALLOW_MNC
50 # include "MNC_PARAMS.h"
51 #endif
52 #ifdef ALLOW_AUTODIFF_TAMC
53 # include "tamc.h"
54 # include "tamc_keys.h"
55 #endif
56
57 C == Routine arguments ==
58 C bi,bj :: current tile indices
59 C k :: current vertical level
60 C iMin,iMax,jMin,jMax :: loop ranges
61 C fVerU :: Flux of momentum in the vertical direction, out of the upper
62 C fVerV :: face of a cell k ( flux into the cell above ).
63 C fVerUkm :: vertical viscous flux of U, interface above (k-1/2)
64 C fVerVkm :: vertical viscous flux of V, interface above (k-1/2)
65 C fVerUkp :: vertical viscous flux of U, interface below (k+1/2)
66 C fVerVkp :: vertical viscous flux of V, interface below (k+1/2)
67
68 C guDiss :: dissipation tendency (all explicit terms), u component
69 C gvDiss :: dissipation tendency (all explicit terms), v component
70 C myTime :: current time
71 C myIter :: current time-step number
72 C myThid :: my Thread Id number
73 INTEGER bi,bj,k
74 INTEGER iMin,iMax,jMin,jMax
75 _RL kappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
76 _RL kappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr+1)
77 _RL fVerUkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL fVerVkm(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 _RL fVerUkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80 _RL fVerVkp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 _RL myTime
84 INTEGER myIter
85 INTEGER myThid
86
87 #ifdef ALLOW_MOM_VECINV
88
89 C == Functions ==
90 LOGICAL DIFFERENT_MULTIPLE
91 EXTERNAL DIFFERENT_MULTIPLE
92
93 C == Local variables ==
94 C strainBC :: same as strain but account for no-slip BC
95 C vort3BC :: same as vort3 but account for no-slip BC
96 _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RS hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RS h0FacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL del2u (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL del2v (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL dStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 _RL zStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110 _RL strain (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111 _RL strainBC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112 _RL KE (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113 _RL omega3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114 _RL vort3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115 _RL vort3BC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RL hDiv (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 C i,j :: Loop counters
122 INTEGER i,j
123 C xxxFac :: On-off tracer parameters used for switching terms off.
124 _RL ArDudrFac
125 _RL ArDvdrFac
126 _RL sideMaskFac
127 LOGICAL bottomDragTerms
128 LOGICAL writeDiag
129 #ifdef ALLOW_AUTODIFF_TAMC
130 INTEGER imomkey
131 #endif
132
133 #ifdef ALLOW_MNC
134 INTEGER offsets(9)
135 CHARACTER*(1) pf
136 #endif
137
138 #ifdef ALLOW_AUTODIFF
139 C-- only the kDown part of fverU/V is set in this subroutine
140 C-- the kUp is still required
141 C-- In the case of mom_fluxform kUp is set as well
142 C-- (at least in part)
143 fVerUkm(1,1) = fVerUkm(1,1)
144 fVerVkm(1,1) = fVerVkm(1,1)
145 #endif
146
147 #ifdef ALLOW_AUTODIFF_TAMC
148 act0 = k - 1
149 max0 = Nr
150 act1 = bi - myBxLo(myThid)
151 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
152 act2 = bj - myByLo(myThid)
153 max2 = myByHi(myThid) - myByLo(myThid) + 1
154 act3 = myThid - 1
155 max3 = nTx*nTy
156 act4 = ikey_dynamics - 1
157 imomkey = (act0 + 1)
158 & + act1*max0
159 & + act2*max0*max1
160 & + act3*max0*max1*max2
161 & + act4*max0*max1*max2*max3
162 #endif /* ALLOW_AUTODIFF_TAMC */
163
164 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
165
166 #ifdef ALLOW_MNC
167 IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
168 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
169 pf(1:1) = 'D'
170 ELSE
171 pf(1:1) = 'R'
172 ENDIF
173 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
174 CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
175 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
176 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
177 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
178 ENDIF
179 DO i = 1,9
180 offsets(i) = 0
181 ENDDO
182 offsets(3) = k
183 c write(*,*) 'offsets = ',(offsets(i),i=1,9)
184 ENDIF
185 #endif /* ALLOW_MNC */
186
187 C-- Initialise intermediate terms
188 DO j=1-OLy,sNy+OLy
189 DO i=1-OLx,sNx+OLx
190 vF(i,j) = 0.
191 vrF(i,j) = 0.
192 uCf(i,j) = 0.
193 vCf(i,j) = 0.
194 del2u(i,j) = 0.
195 del2v(i,j) = 0.
196 dStar(i,j) = 0.
197 zStar(i,j) = 0.
198 guDiss(i,j)= 0.
199 gvDiss(i,j)= 0.
200 vort3(i,j) = 0.
201 omega3(i,j)= 0.
202 KE(i,j) = 0.
203 C- need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
204 hDiv(i,j) = 0.
205 c viscAh_Z(i,j) = 0.
206 c viscAh_D(i,j) = 0.
207 c viscA4_Z(i,j) = 0.
208 c viscA4_D(i,j) = 0.
209 strain(i,j) = 0. _d 0
210 strainBC(i,j)= 0. _d 0
211 tension(i,j) = 0. _d 0
212 #ifdef ALLOW_AUTODIFF
213 hFacZ(i,j) = 0. _d 0
214 #endif
215 ENDDO
216 ENDDO
217
218 C-- Term by term tracer parmeters
219 C o U momentum equation
220 ArDudrFac = vfFacMom*1.
221 C o V momentum equation
222 ArDvdrFac = vfFacMom*1.
223
224 C note: using standard stencil (no mask) results in under-estimating
225 C vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
226 IF ( no_slip_sides ) THEN
227 sideMaskFac = sideDragFactor
228 ELSE
229 sideMaskFac = 0. _d 0
230 ENDIF
231
232 IF ( selectImplicitDrag.EQ.0 .AND.
233 & ( no_slip_bottom
234 & .OR. selectBotDragQuadr.GE.0
235 & .OR. bottomDragLinear.NE.0. ) ) THEN
236 bottomDragTerms=.TRUE.
237 ELSE
238 bottomDragTerms=.FALSE.
239 ENDIF
240
241 C-- Calculate open water fraction at vorticity points
242 CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
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 #ifdef ALLOW_AUTODIFF_TAMC
253 CADJ STORE ufld(:,:) =
254 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
255 CADJ STORE vfld(:,:) =
256 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
257 CADJ STORE hFacZ(:,:) =
258 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
259 CADJ STORE r_hFacZ(:,:) =
260 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
261 CADJ STORE fverukm(:,:) =
262 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
263 CADJ STORE fvervkm(:,:) =
264 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
265 #endif
266
267 C note (jmc) : Dissipation and Vort3 advection do not necesary
268 C use the same maskZ (and hFacZ) => needs 2 call(s)
269 c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
270
271 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
272
273 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
274
275 #ifdef ALLOW_AUTODIFF_TAMC
276 CADJ STORE ke(:,:) =
277 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
278 CADJ STORE vort3(:,:) =
279 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
280 CADJ STORE vort3bc(:,:) =
281 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
282 #endif
283
284 C- mask vort3 and account for no-slip / free-slip BC in vort3BC:
285 DO j=1-OLy,sNy+OLy
286 DO i=1-OLx,sNx+OLx
287 vort3BC(i,j) = vort3(i,j)
288 IF ( hFacZ(i,j).EQ.zeroRS ) THEN
289 vort3BC(i,j) = sideMaskFac*vort3BC(i,j)
290 vort3(i,j) = 0.
291 ENDIF
292 ENDDO
293 ENDDO
294
295 #ifdef ALLOW_AUTODIFF_TAMC
296 CADJ STORE vort3(:,:) =
297 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
298 CADJ STORE vort3bc(:,:) =
299 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
300 #endif
301
302 IF (momViscosity) THEN
303 C-- For viscous term, compute horizontal divergence, tension & strain
304 C and mask relative vorticity (free-slip case):
305
306 DO j=1-OLy,sNy+OLy
307 DO i=1-OLx,sNx+OLx
308 h0FacZ(i,j) = hFacZ(i,j)
309 ENDDO
310 ENDDO
311 #ifdef NONLIN_FRSURF
312 IF ( no_slip_sides .AND. nonlinFreeSurf.GT.0 ) THEN
313 DO j=2-OLy,sNy+OLy
314 DO i=2-OLx,sNx+OLx
315 h0FacZ(i,j) = MIN(
316 & MIN( h0FacW(i,j,k,bi,bj), h0FacW(i,j-1,k,bi,bj) ),
317 & MIN( h0FacS(i,j,k,bi,bj), h0FacS(i-1,j,k,bi,bj) ) )
318 ENDDO
319 ENDDO
320 ENDIF
321 #endif /* NONLIN_FRSURF */
322
323 #ifdef ALLOW_AUTODIFF_TAMC
324 CADJ STORE h0FacZ(:,:) =
325 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
326 CADJ STORE hFacZ(:,:) =
327 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
328 #endif
329
330 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
331
332 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
333 CALL MOM_CALC_TENSION( bi,bj,k,uFld,vFld,tension,myThid )
334 CALL MOM_CALC_STRAIN( bi,bj,k,uFld,vFld,hFacZ,strain,myThid )
335 C- mask strain and account for no-slip / free-slip BC in strainBC:
336 DO j=1-OLy,sNy+OLy
337 DO i=1-OLx,sNx+OLx
338 strainBC(i,j) = strain(i,j)
339 IF ( hFacZ(i,j).EQ.zeroRS ) THEN
340 strainBC(i,j) = sideMaskFac*strainBC(i,j)
341 strain(i,j) = 0.
342 ENDIF
343 ENDDO
344 ENDDO
345 ENDIF
346
347 #ifdef ALLOW_AUTODIFF_TAMC
348 CADJ STORE hdiv(:,:) =
349 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
350 CADJ STORE tension(:,:) =
351 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
352 CADJ STORE strain(:,:) =
353 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
354 CADJ STORE strainbc(:,:) =
355 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
356 #endif
357
358 C-- Calculate Lateral Viscosities
359 DO j=1-OLy,sNy+OLy
360 DO i=1-OLx,sNx+OLx
361 viscAh_D(i,j) = viscAhD
362 viscAh_Z(i,j) = viscAhZ
363 viscA4_D(i,j) = viscA4D
364 viscA4_Z(i,j) = viscA4Z
365 ENDDO
366 ENDDO
367 IF ( useVariableVisc ) THEN
368 C- uses vort3BC & strainBC which account for no-slip / free-slip BC
369 CALL MOM_CALC_VISC( bi, bj, k,
370 O viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
371 I hDiv, vort3BC, tension, strainBC, KE, hfacZ,
372 I myThid )
373 ENDIF
374
375 #ifdef ALLOW_AUTODIFF_TAMC
376 CADJ STORE viscAh_Z(:,:) =
377 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
378 CADJ STORE viscAh_D(:,:) =
379 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
380 CADJ STORE viscA4_Z(:,:) =
381 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
382 CADJ STORE viscA4_D(:,:) =
383 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
384 #endif
385
386 #ifdef ALLOW_AUTODIFF_TAMC
387 CADJ STORE hDiv(:,:) =
388 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
389 CADJ STORE vort3(:,:) =
390 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
391 CADJ STORE hFacZ(:,:) =
392 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
393 #endif
394
395 C Calculate del^2 u and del^2 v for bi-harmonic term
396 IF (useBiharmonicVisc) THEN
397 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
398 O del2u,del2v,
399 I myThid)
400 #ifdef ALLOW_AUTODIFF_TAMC
401 CADJ STORE del2u(:,:) =
402 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
403 CADJ STORE del2v(:,:) =
404 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
405 #endif
406 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
407 CALL MOM_CALC_RELVORT3(bi,bj,k,
408 & del2u,del2v,hFacZ,zStar,myThid)
409 ENDIF
410
411 #ifdef ALLOW_AUTODIFF_TAMC
412 CADJ STORE del2u(:,:) =
413 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
414 CADJ STORE del2v(:,:) =
415 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
416 CADJ STORE dStar(:,:) =
417 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
418 CADJ STORE zStar(:,:) =
419 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
420 #endif
421
422 C--- Calculate dissipation terms for U and V equations
423
424 C- in terms of tension and strain
425 IF (useStrainTensionVisc) THEN
426 C use masked strain as if free-slip since side-drag is computed separately
427 CALL MOM_HDISSIP( bi, bj, k,
428 I tension, strain, hFacZ,
429 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
430 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
431 O guDiss, gvDiss,
432 I myThid )
433 ELSE
434 C- in terms of vorticity and divergence
435 CALL MOM_VI_HDISSIP( bi, bj, k,
436 I hDiv, vort3, dStar, zStar, hFacZ,
437 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
438 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
439 O guDiss, gvDiss,
440 I myThid )
441 ENDIF
442
443 C--- Other dissipation terms in Zonal momentum equation
444
445 C-- Vertical flux (fVer is at upper face of "u" cell)
446 C Eddy component of vertical flux (interior component only) -> vrF
447 IF ( .NOT.implicitViscosity ) THEN
448 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,kappaRU,vrF,myThid)
449 C Combine fluxes
450 DO j=jMin,jMax
451 DO i=iMin,iMax
452 fVerUkp(i,j) = ArDudrFac*vrF(i,j)
453 ENDDO
454 ENDDO
455
456 #ifdef ALLOW_AUTODIFF_TAMC
457 CADJ STORE fVerUkp(:,:) =
458 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
459 #endif
460
461 C-- Tendency is minus divergence of the fluxes
462 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
463 DO j=jMin,jMax
464 DO i=iMin,iMax
465 guDiss(i,j) = guDiss(i,j)
466 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
467 & *recip_rAw(i,j,bi,bj)
468 & *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
469 & *recip_deepFac2C(k)*recip_rhoFacC(k)
470 ENDDO
471 ENDDO
472 ENDIF
473
474 C-- No-slip and drag BCs appear as body forces in cell abutting topography
475 IF ( no_slip_sides ) THEN
476 C- No-slip BCs impose a drag at walls...
477 CALL MOM_U_SIDEDRAG( bi, bj, k,
478 I uFld, del2u, h0FacZ,
479 I viscAh_Z, viscA4_Z,
480 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
481 O vF,
482 I myThid )
483 DO j=jMin,jMax
484 DO i=iMin,iMax
485 guDiss(i,j) = guDiss(i,j)+vF(i,j)
486 ENDDO
487 ENDDO
488 ENDIF
489
490 C- No-slip BCs impose a drag at bottom
491 IF ( bottomDragTerms ) THEN
492 CALL MOM_U_BOTTOMDRAG( bi, bj, k,
493 I uFld, vFld, KE, kappaRU,
494 O vF,
495 I myThid )
496 DO j=jMin,jMax
497 DO i=iMin,iMax
498 guDiss(i,j) = guDiss(i,j)+vF(i,j)
499 ENDDO
500 ENDDO
501 ENDIF
502 #ifdef ALLOW_SHELFICE
503 IF ( useShelfIce ) THEN
504 CALL SHELFICE_U_DRAG( bi, bj, k,
505 I uFld, vFld, KE, kappaRU,
506 O vF,
507 I myThid )
508 DO j=jMin,jMax
509 DO i=iMin,iMax
510 guDiss(i,j) = guDiss(i,j) + vF(i,j)
511 ENDDO
512 ENDDO
513 ENDIF
514 #endif /* ALLOW_SHELFICE */
515
516 C--- Other dissipation terms in Meridional momentum equation
517
518 C-- Vertical flux (fVer is at upper face of "v" cell)
519 C Eddy component of vertical flux (interior component only) -> vrF
520 IF ( .NOT.implicitViscosity ) THEN
521 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,kappaRV,vrF,myThid)
522 C Combine fluxes -> fVerV
523 DO j=jMin,jMax
524 DO i=iMin,iMax
525 fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
526 ENDDO
527 ENDDO
528 #ifdef ALLOW_AUTODIFF_TAMC
529 CADJ STORE fVerVkp(:,:) =
530 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
531 #endif
532 C-- Tendency is minus divergence of the fluxes
533 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
534 DO j=jMin,jMax
535 DO i=iMin,iMax
536 gvDiss(i,j) = gvDiss(i,j)
537 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
538 & *recip_rAs(i,j,bi,bj)
539 & *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
540 & *recip_deepFac2C(k)*recip_rhoFacC(k)
541 ENDDO
542 ENDDO
543 ENDIF
544
545 C-- No-slip and drag BCs appear as body forces in cell abutting topography
546 IF ( no_slip_sides ) THEN
547 C- No-slip BCs impose a drag at walls...
548 CALL MOM_V_SIDEDRAG( bi, bj, k,
549 I vFld, del2v, h0FacZ,
550 I viscAh_Z, viscA4_Z,
551 I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
552 O vF,
553 I myThid )
554 DO j=jMin,jMax
555 DO i=iMin,iMax
556 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
557 ENDDO
558 ENDDO
559 ENDIF
560
561 C- No-slip BCs impose a drag at bottom
562 IF ( bottomDragTerms ) THEN
563 CALL MOM_V_BOTTOMDRAG( bi, bj, k,
564 I uFld, vFld, KE, kappaRV,
565 O vF,
566 I myThid )
567 DO j=jMin,jMax
568 DO i=iMin,iMax
569 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
570 ENDDO
571 ENDDO
572 ENDIF
573 #ifdef ALLOW_SHELFICE
574 IF ( useShelfIce ) THEN
575 CALL SHELFICE_V_DRAG( bi, bj, k,
576 I uFld, vFld, KE, kappaRV,
577 O vF,
578 I myThid )
579 DO j=jMin,jMax
580 DO i=iMin,iMax
581 gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
582 ENDDO
583 ENDDO
584 ENDIF
585 #endif /* ALLOW_SHELFICE */
586
587 C-- if (momViscosity) end of block.
588 ENDIF
589
590 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
591 c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
592
593 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
594
595 C--- Prepare for Advection & Coriolis terms:
596 C- calculate absolute vorticity
597 IF (useAbsVorticity)
598 & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
599
600 #ifdef ALLOW_AUTODIFF_TAMC
601 CADJ STORE omega3(:,:) =
602 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
603 #endif
604
605 C-- Horizontal Coriolis terms
606 c IF (useCoriolis .AND. .NOT.useCDscheme
607 c & .AND. .NOT. useAbsVorticity) THEN
608 C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
609 IF ( useCoriolis .AND.
610 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
611 & ) THEN
612 IF (useAbsVorticity) THEN
613 CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
614 & vFld,omega3,hFacZ,r_hFacZ,
615 & uCf,myThid)
616 CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
617 & uFld,omega3,hFacZ,r_hFacZ,
618 & vCf,myThid)
619 ELSE
620 CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
621 & uCf,vCf,myThid)
622 ENDIF
623 DO j=jMin,jMax
624 DO i=iMin,iMax
625 gU(i,j,k,bi,bj) = uCf(i,j)
626 gV(i,j,k,bi,bj) = vCf(i,j)
627 ENDDO
628 ENDDO
629 IF ( writeDiag ) THEN
630 IF (snapshot_mdsio) THEN
631 CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
632 CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
633 ENDIF
634 #ifdef ALLOW_MNC
635 IF (useMNC .AND. snapshot_mnc) THEN
636 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
637 & offsets, myThid)
638 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
639 & offsets, myThid)
640 ENDIF
641 #endif /* ALLOW_MNC */
642 ENDIF
643 #ifdef ALLOW_DIAGNOSTICS
644 IF ( useDiagnostics ) THEN
645 CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
646 CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
647 ENDIF
648 #endif /* ALLOW_DIAGNOSTICS */
649 ELSE
650 DO j=jMin,jMax
651 DO i=iMin,iMax
652 gU(i,j,k,bi,bj) = 0. _d 0
653 gV(i,j,k,bi,bj) = 0. _d 0
654 ENDDO
655 ENDDO
656 ENDIF
657
658 #ifdef ALLOW_AUTODIFF_TAMC
659 CADJ STORE ucf(:,:) =
660 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
661 CADJ STORE vcf(:,:) =
662 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
663 #endif
664
665 IF (momAdvection) THEN
666 C-- Horizontal advection of relative (or absolute) vorticity
667 IF ( (highOrderVorticity.OR.upwindVorticity)
668 & .AND.useAbsVorticity ) THEN
669 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
670 & highOrderVorticity,upwindVorticity,
671 & vFld,omega3,r_hFacZ,
672 & uCf,myThid)
673 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
674 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
675 & highOrderVorticity, upwindVorticity,
676 & vFld,vort3, r_hFacZ,
677 & uCf,myThid)
678 ELSEIF ( useAbsVorticity ) THEN
679 CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
680 & vFld,omega3,hFacZ,r_hFacZ,
681 & uCf,myThid)
682 ELSE
683 CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
684 & vFld,vort3, hFacZ,r_hFacZ,
685 & uCf,myThid)
686 ENDIF
687 DO j=jMin,jMax
688 DO i=iMin,iMax
689 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
690 ENDDO
691 ENDDO
692 IF ( (highOrderVorticity.OR.upwindVorticity)
693 & .AND.useAbsVorticity ) THEN
694 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
695 & highOrderVorticity, upwindVorticity,
696 & uFld,omega3,r_hFacZ,
697 & vCf,myThid)
698 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
699 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
700 & highOrderVorticity, upwindVorticity,
701 & uFld,vort3, r_hFacZ,
702 & vCf,myThid)
703 ELSEIF ( useAbsVorticity ) THEN
704 CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
705 & uFld,omega3,hFacZ,r_hFacZ,
706 & vCf,myThid)
707 ELSE
708 CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
709 & uFld,vort3, hFacZ,r_hFacZ,
710 & vCf,myThid)
711 ENDIF
712 DO j=jMin,jMax
713 DO i=iMin,iMax
714 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
715 ENDDO
716 ENDDO
717
718 #ifdef ALLOW_AUTODIFF_TAMC
719 CADJ STORE ucf(:,:) =
720 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
721 CADJ STORE vcf(:,:) =
722 CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
723 #endif
724
725 IF ( writeDiag ) THEN
726 IF (snapshot_mdsio) THEN
727 CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
728 CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
729 ENDIF
730 #ifdef ALLOW_MNC
731 IF (useMNC .AND. snapshot_mnc) THEN
732 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
733 & offsets, myThid)
734 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
735 & offsets, myThid)
736 ENDIF
737 #endif /* ALLOW_MNC */
738 ENDIF
739
740 #ifdef ALLOW_TIMEAVE
741 IF (taveFreq.GT.0.) THEN
742 CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
743 & Nr, k, bi, bj, myThid)
744 CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
745 & Nr, k, bi, bj, myThid)
746 ENDIF
747 #endif /* ALLOW_TIMEAVE */
748 #ifdef ALLOW_DIAGNOSTICS
749 IF ( useDiagnostics ) THEN
750 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
751 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
752 ENDIF
753 #endif /* ALLOW_DIAGNOSTICS */
754
755 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
756 IF ( .NOT. momImplVertAdv ) THEN
757 CALL MOM_VI_U_VERTSHEAR(bi,bj,k,uVel,wVel,uCf,myThid)
758 DO j=jMin,jMax
759 DO i=iMin,iMax
760 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
761 ENDDO
762 ENDDO
763 CALL MOM_VI_V_VERTSHEAR(bi,bj,k,vVel,wVel,vCf,myThid)
764 DO j=jMin,jMax
765 DO i=iMin,iMax
766 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
767 ENDDO
768 ENDDO
769 #ifdef ALLOW_DIAGNOSTICS
770 IF ( useDiagnostics ) THEN
771 CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
772 CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
773 ENDIF
774 #endif /* ALLOW_DIAGNOSTICS */
775 ENDIF
776
777 C-- Bernoulli term
778 CALL MOM_VI_U_GRAD_KE(bi,bj,k,KE,uCf,myThid)
779 DO j=jMin,jMax
780 DO i=iMin,iMax
781 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
782 ENDDO
783 ENDDO
784 CALL MOM_VI_V_GRAD_KE(bi,bj,k,KE,vCf,myThid)
785 DO j=jMin,jMax
786 DO i=iMin,iMax
787 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
788 ENDDO
789 ENDDO
790 IF ( writeDiag ) THEN
791 IF (snapshot_mdsio) THEN
792 CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
793 CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
794 ENDIF
795 #ifdef ALLOW_MNC
796 IF (useMNC .AND. snapshot_mnc) THEN
797 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
798 & offsets, myThid)
799 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
800 & offsets, myThid)
801 ENDIF
802 #endif /* ALLOW_MNC */
803 ENDIF
804
805 C-- end if momAdvection
806 ENDIF
807
808 C-- 3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
809 IF ( use3dCoriolis ) THEN
810 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
811 DO j=jMin,jMax
812 DO i=iMin,iMax
813 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
814 ENDDO
815 ENDDO
816 IF ( usingCurvilinearGrid ) THEN
817 C- presently, non zero angleSinC array only supported with Curvilinear-Grid
818 CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
819 DO j=jMin,jMax
820 DO i=iMin,iMax
821 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
822 ENDDO
823 ENDDO
824 ENDIF
825 ENDIF
826
827 C-- Non-Hydrostatic (spherical) metric terms
828 IF ( useNHMTerms ) THEN
829 CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
830 DO j=jMin,jMax
831 DO i=iMin,iMax
832 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
833 ENDDO
834 ENDDO
835 CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
836 DO j=jMin,jMax
837 DO i=iMin,iMax
838 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
839 ENDDO
840 ENDDO
841 ENDIF
842
843 C-- Set du/dt & dv/dt on boundaries to zero
844 DO j=jMin,jMax
845 DO i=iMin,iMax
846 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
847 gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
848 ENDDO
849 ENDDO
850
851 #ifdef ALLOW_DEBUG
852 IF ( debugLevel .GE. debLevC
853 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
854 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
855 & .AND. useCubedSphereExchange ) THEN
856 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
857 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
858 ENDIF
859 #endif /* ALLOW_DEBUG */
860
861 IF ( writeDiag ) THEN
862 IF (useBiharmonicVisc) THEN
863 CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
864 & bi,bj,k, myIter, myThid )
865 CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
866 & bi,bj,k, myIter, myThid )
867 CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
868 & bi,bj,k, myIter, myThid )
869 CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
870 & bi,bj,k, myIter, myThid )
871 ENDIF
872 IF (snapshot_mdsio) THEN
873 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
874 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
875 CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
876 CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
877 CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
878 CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
879 & bi,bj,k, myIter, myThid )
880 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
881 CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
882 ENDIF
883 #ifdef ALLOW_MNC
884 IF (useMNC .AND. snapshot_mnc) THEN
885 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
886 & offsets, myThid)
887 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
888 & offsets, myThid)
889 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
890 & offsets, myThid)
891 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
892 & offsets, myThid)
893 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
894 & offsets, myThid)
895 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
896 & offsets, myThid)
897 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
898 & offsets, myThid)
899 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
900 & offsets, myThid)
901 ENDIF
902 #endif /* ALLOW_MNC */
903 ENDIF
904
905 #ifdef ALLOW_DIAGNOSTICS
906 IF ( useDiagnostics ) THEN
907 CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
908 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
909 IF (momViscosity) THEN
910 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
911 ENDIF
912 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
913 CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
914 CALL DIAGNOSTICS_FILL(strainBC,'Strain ',k,1,2,bi,bj,myThid)
915 ENDIF
916 CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
917 & 'Um_Advec',k,1,2,bi,bj,myThid)
918 CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
919 & 'Vm_Advec',k,1,2,bi,bj,myThid)
920 ENDIF
921 #endif /* ALLOW_DIAGNOSTICS */
922
923 #endif /* ALLOW_MOM_VECINV */
924
925 RETURN
926 END

  ViewVC Help
Powered by ViewVC 1.1.22