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

Annotation of /MITgcm/pkg/mom_vecinv/mom_vecinv.F

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


Revision 1.80 - (hide annotations) (download)
Fri Apr 28 17:17:14 2017 UTC (7 years, 1 month ago) by mlosch
Branch: MAIN
Changes since 1.79: +21 -9 lines
 pass these runtime flags as formal parameters to
  s/r mom_vi_u/v_coriolis, mom_vi_u/v_coriolis_c4, so that these routines
  can also be used in pkg/seaice:
  selectVortScheme, highOrderVorticity, upwindVorticity, useJamartMomAdv

1 mlosch 1.80 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.79 2017/03/06 20:05:36 gforget Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4 adcroft 1.21 #include "MOM_VECINV_OPTIONS.h"
5 jmc 1.73 #ifdef ALLOW_AUTODIFF
6     # include "AUTODIFF_OPTIONS.h"
7     #endif
8 jmc 1.69 #ifdef ALLOW_MOM_COMMON
9     # include "MOM_COMMON_OPTIONS.h"
10     #endif
11 adcroft 1.1
12 jmc 1.57 SUBROUTINE MOM_VECINV(
13 jmc 1.66 I bi,bj,k,iMin,iMax,jMin,jMax,
14 jmc 1.75 I kappaRU, kappaRV,
15 jmc 1.66 I fVerUkm, fVerVkm,
16     O fVerUkp, fVerVkp,
17 jmc 1.31 O guDiss, gvDiss,
18 jmc 1.66 I myTime, myIter, myThid )
19     C *==========================================================*
20 adcroft 1.1 C | S/R MOM_VECINV |
21     C | o Form the right hand-side of the momentum equation. |
22 jmc 1.66 C *==========================================================*
23 adcroft 1.1 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 jmc 1.66 C *==========================================================*
34 adcroft 1.1 IMPLICIT NONE
35    
36     C == Global variables ==
37     #include "SIZE.h"
38     #include "EEPARAMS.h"
39     #include "PARAMS.h"
40 jmc 1.69 #include "GRID.h"
41 jmc 1.71 #include "SURFACE.h"
42 jmc 1.69 #include "DYNVARS.h"
43     #ifdef ALLOW_MOM_COMMON
44     # include "MOM_VISC.h"
45 edhill 1.27 #endif
46 jmc 1.7 #ifdef ALLOW_TIMEAVE
47 jmc 1.69 # include "TIMEAVE_STATV.h"
48     #endif
49     #ifdef ALLOW_MNC
50     # include "MNC_PARAMS.h"
51 jmc 1.7 #endif
52 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
53     # include "tamc.h"
54     # include "tamc_keys.h"
55     #endif
56 adcroft 1.1
57     C == Routine arguments ==
58 jmc 1.66 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 jmc 1.77 C fVerV :: face of a cell k ( flux into the cell above ).
63 jmc 1.66 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 jmc 1.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 jmc 1.66 _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 jmc 1.31 _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83 jmc 1.15 _RL myTime
84 adcroft 1.2 INTEGER myIter
85 adcroft 1.1 INTEGER myThid
86    
87 edhill 1.11 #ifdef ALLOW_MOM_VECINV
88 jmc 1.7
89 adcroft 1.2 C == Functions ==
90 jmc 1.38 LOGICAL DIFFERENT_MULTIPLE
91     EXTERNAL DIFFERENT_MULTIPLE
92 adcroft 1.2
93 adcroft 1.1 C == Local variables ==
94 jmc 1.72 C strainBC :: same as strain but account for no-slip BC
95     C vort3BC :: same as vort3 but account for no-slip BC
96 adcroft 1.1 _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 jmc 1.54 _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 jmc 1.71 _RS h0FacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 jmc 1.54 _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 jmc 1.72 _RL strainBC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112 jmc 1.54 _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 jmc 1.72 _RL vort3BC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 jmc 1.54 _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 jmc 1.66 C i,j :: Loop counters
122     INTEGER i,j
123     C xxxFac :: On-off tracer parameters used for switching terms off.
124 adcroft 1.1 _RL ArDudrFac
125     _RL ArDvdrFac
126 jmc 1.54 _RL sideMaskFac
127 adcroft 1.1 LOGICAL bottomDragTerms
128 jmc 1.15 LOGICAL writeDiag
129 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
130     INTEGER imomkey
131     #endif
132 adcroft 1.1
133 edhill 1.25 #ifdef ALLOW_MNC
134     INTEGER offsets(9)
135 edhill 1.53 CHARACTER*(1) pf
136 edhill 1.25 #endif
137    
138 jmc 1.73 #ifdef ALLOW_AUTODIFF
139 heimbach 1.9 C-- only the kDown part of fverU/V is set in this subroutine
140     C-- the kUp is still required
141 jmc 1.77 C-- In the case of mom_fluxform kUp is set as well
142 heimbach 1.9 C-- (at least in part)
143 jmc 1.66 fVerUkm(1,1) = fVerUkm(1,1)
144     fVerVkm(1,1) = fVerVkm(1,1)
145 heimbach 1.9 #endif
146    
147 heimbach 1.59 #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 jmc 1.38 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
165 adcroft 1.1
166 edhill 1.24 #ifdef ALLOW_MNC
167     IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
168 edhill 1.53 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
169     pf(1:1) = 'D'
170     ELSE
171     pf(1:1) = 'R'
172     ENDIF
173 edhill 1.25 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
174     CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
175 edhill 1.39 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
176 edhill 1.25 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
177 edhill 1.39 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
178 edhill 1.25 ENDIF
179     DO i = 1,9
180     offsets(i) = 0
181     ENDDO
182     offsets(3) = k
183 jmc 1.61 c write(*,*) 'offsets = ',(offsets(i),i=1,9)
184 edhill 1.24 ENDIF
185     #endif /* ALLOW_MNC */
186    
187 jmc 1.61 C-- Initialise intermediate terms
188     DO j=1-OLy,sNy+OLy
189     DO i=1-OLx,sNx+OLx
190 jmc 1.31 vF(i,j) = 0.
191     vrF(i,j) = 0.
192 adcroft 1.1 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 jmc 1.31 guDiss(i,j)= 0.
199     gvDiss(i,j)= 0.
200 adcroft 1.1 vort3(i,j) = 0.
201 jmc 1.31 omega3(i,j)= 0.
202 jmc 1.54 KE(i,j) = 0.
203 jmc 1.61 C- need to initialise hDiv for MOM_VI_DEL2UV(call FILL_CS_CORNER_TR_RL)
204     hDiv(i,j) = 0.
205 jmc 1.70 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 heimbach 1.8 strain(i,j) = 0. _d 0
210 jmc 1.72 strainBC(i,j)= 0. _d 0
211 heimbach 1.8 tension(i,j) = 0. _d 0
212 jmc 1.73 #ifdef ALLOW_AUTODIFF
213 heimbach 1.55 hFacZ(i,j) = 0. _d 0
214 heimbach 1.8 #endif
215 adcroft 1.1 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 jmc 1.54 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 jmc 1.78 IF ( selectImplicitDrag.EQ.0 .AND.
233     & ( no_slip_bottom
234 jmc 1.76 & .OR. selectBotDragQuadr.GE.0
235 jmc 1.78 & .OR. bottomDragLinear.NE.0. ) ) THEN
236 adcroft 1.1 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 gforget 1.79 #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 jmc 1.7 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 jmc 1.52 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
272 adcroft 1.1
273 jmc 1.54 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
274 adcroft 1.1
275 gforget 1.79 #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 jmc 1.72 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 gforget 1.79 #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 jmc 1.54 IF (momViscosity) THEN
303 jmc 1.57 C-- For viscous term, compute horizontal divergence, tension & strain
304 jmc 1.54 C and mask relative vorticity (free-slip case):
305 adcroft 1.1
306 jmc 1.71 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 gforget 1.79 #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 jmc 1.54 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
331 adcroft 1.1
332 jmc 1.72 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 jmc 1.54 ENDDO
345 jmc 1.72 ENDIF
346 jmc 1.54
347 gforget 1.79 #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 jmc 1.70 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 jmc 1.72 C- uses vort3BC & strainBC which account for no-slip / free-slip BC
369 jmc 1.70 CALL MOM_CALC_VISC( bi, bj, k,
370     O viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
371 jmc 1.72 I hDiv, vort3BC, tension, strainBC, KE, hfacZ,
372 jmc 1.70 I myThid )
373     ENDIF
374 baylor 1.50
375 gforget 1.79 #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 adcroft 1.1 C Calculate del^2 u and del^2 v for bi-harmonic term
396 jmc 1.69 IF (useBiharmonicVisc) THEN
397 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
398     O del2u,del2v,
399 jmc 1.72 I myThid)
400 gforget 1.79 #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 jmc 1.48 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 jmc 1.54 ENDIF
410    
411 gforget 1.79 #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 jmc 1.54 C--- Calculate dissipation terms for U and V equations
423 baylor 1.47
424 jmc 1.72 C- in terms of tension and strain
425 baylor 1.47 IF (useStrainTensionVisc) THEN
426 jmc 1.72 C use masked strain as if free-slip since side-drag is computed separately
427 jmc 1.69 CALL MOM_HDISSIP( bi, bj, k,
428 jmc 1.72 I tension, strain, hFacZ,
429 jmc 1.69 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
430     I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
431     O guDiss, gvDiss,
432     I myThid )
433 baylor 1.47 ELSE
434 jmc 1.72 C- in terms of vorticity and divergence
435 jmc 1.69 CALL MOM_VI_HDISSIP( bi, bj, k,
436 jmc 1.72 I hDiv, vort3, dStar, zStar, hFacZ,
437 jmc 1.69 I viscAh_Z, viscAh_D, viscA4_Z, viscA4_D,
438     I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
439     O guDiss, gvDiss,
440 jmc 1.72 I myThid )
441 adcroft 1.3 ENDIF
442 jmc 1.7
443 jmc 1.54 C--- Other dissipation terms in Zonal momentum equation
444 adcroft 1.1
445     C-- Vertical flux (fVer is at upper face of "u" cell)
446     C Eddy component of vertical flux (interior component only) -> vrF
447 jmc 1.70 IF ( .NOT.implicitViscosity ) THEN
448 jmc 1.75 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,kappaRU,vrF,myThid)
449 adcroft 1.1 C Combine fluxes
450 jmc 1.31 DO j=jMin,jMax
451     DO i=iMin,iMax
452 jmc 1.66 fVerUkp(i,j) = ArDudrFac*vrF(i,j)
453 jmc 1.31 ENDDO
454 adcroft 1.1 ENDDO
455 gforget 1.79
456     #ifdef ALLOW_AUTODIFF_TAMC
457     CADJ STORE fVerUkp(:,:) =
458     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
459     #endif
460    
461 jmc 1.31 C-- Tendency is minus divergence of the fluxes
462 jmc 1.77 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
463 jmc 1.67 DO j=jMin,jMax
464     DO i=iMin,iMax
465 jmc 1.31 guDiss(i,j) = guDiss(i,j)
466 adcroft 1.1 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
467     & *recip_rAw(i,j,bi,bj)
468 jmc 1.66 & *( fVerUkp(i,j) - fVerUkm(i,j) )*rkSign
469 jmc 1.77 & *recip_deepFac2C(k)*recip_rhoFacC(k)
470 jmc 1.31 ENDDO
471 adcroft 1.1 ENDDO
472 jmc 1.31 ENDIF
473 adcroft 1.1
474 jmc 1.57 C-- No-slip and drag BCs appear as body forces in cell abutting topography
475 jmc 1.70 IF ( no_slip_sides ) THEN
476 adcroft 1.1 C- No-slip BCs impose a drag at walls...
477 jmc 1.69 CALL MOM_U_SIDEDRAG( bi, bj, k,
478 jmc 1.71 I uFld, del2u, h0FacZ,
479 jmc 1.69 I viscAh_Z, viscA4_Z,
480     I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
481     O vF,
482     I myThid )
483 adcroft 1.1 DO j=jMin,jMax
484     DO i=iMin,iMax
485 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
486 adcroft 1.1 ENDDO
487     ENDDO
488     ENDIF
489 jmc 1.70
490 adcroft 1.1 C- No-slip BCs impose a drag at bottom
491 jmc 1.70 IF ( bottomDragTerms ) THEN
492 jmc 1.76 CALL MOM_U_BOTTOMDRAG( bi, bj, k,
493     I uFld, vFld, KE, kappaRU,
494     O vF,
495     I myThid )
496 adcroft 1.1 DO j=jMin,jMax
497     DO i=iMin,iMax
498 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
499 adcroft 1.1 ENDDO
500     ENDDO
501     ENDIF
502 mlosch 1.56 #ifdef ALLOW_SHELFICE
503 jmc 1.76 IF ( useShelfIce ) THEN
504     CALL SHELFICE_U_DRAG( bi, bj, k,
505     I uFld, vFld, KE, kappaRU,
506     O vF,
507     I myThid )
508 mlosch 1.56 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 jmc 1.54 C--- Other dissipation terms in Meridional momentum equation
517 adcroft 1.1
518     C-- Vertical flux (fVer is at upper face of "v" cell)
519     C Eddy component of vertical flux (interior component only) -> vrF
520 jmc 1.70 IF ( .NOT.implicitViscosity ) THEN
521 jmc 1.75 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,kappaRV,vrF,myThid)
522 adcroft 1.1 C Combine fluxes -> fVerV
523 jmc 1.31 DO j=jMin,jMax
524     DO i=iMin,iMax
525 jmc 1.66 fVerVkp(i,j) = ArDvdrFac*vrF(i,j)
526 jmc 1.31 ENDDO
527 adcroft 1.1 ENDDO
528 gforget 1.79 #ifdef ALLOW_AUTODIFF_TAMC
529     CADJ STORE fVerVkp(:,:) =
530     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
531     #endif
532 jmc 1.31 C-- Tendency is minus divergence of the fluxes
533 jmc 1.77 C vert.visc.flx is scaled by deepFac2F (deep-atmos) and rhoFac (anelastic)
534 jmc 1.31 DO j=jMin,jMax
535     DO i=iMin,iMax
536     gvDiss(i,j) = gvDiss(i,j)
537 adcroft 1.1 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
538 jmc 1.66 & *recip_rAs(i,j,bi,bj)
539     & *( fVerVkp(i,j) - fVerVkm(i,j) )*rkSign
540 jmc 1.77 & *recip_deepFac2C(k)*recip_rhoFacC(k)
541 jmc 1.31 ENDDO
542 adcroft 1.1 ENDDO
543 jmc 1.31 ENDIF
544 adcroft 1.1
545 jmc 1.57 C-- No-slip and drag BCs appear as body forces in cell abutting topography
546 jmc 1.70 IF ( no_slip_sides ) THEN
547 adcroft 1.1 C- No-slip BCs impose a drag at walls...
548 jmc 1.69 CALL MOM_V_SIDEDRAG( bi, bj, k,
549 jmc 1.71 I vFld, del2v, h0FacZ,
550 jmc 1.69 I viscAh_Z, viscA4_Z,
551     I useHarmonicVisc, useBiharmonicVisc, useVariableVisc,
552     O vF,
553     I myThid )
554 adcroft 1.1 DO j=jMin,jMax
555     DO i=iMin,iMax
556 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
557 adcroft 1.1 ENDDO
558     ENDDO
559     ENDIF
560 jmc 1.70
561 adcroft 1.1 C- No-slip BCs impose a drag at bottom
562 jmc 1.70 IF ( bottomDragTerms ) THEN
563 jmc 1.76 CALL MOM_V_BOTTOMDRAG( bi, bj, k,
564     I uFld, vFld, KE, kappaRV,
565     O vF,
566     I myThid )
567 adcroft 1.1 DO j=jMin,jMax
568     DO i=iMin,iMax
569 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
570 adcroft 1.1 ENDDO
571     ENDDO
572     ENDIF
573 mlosch 1.56 #ifdef ALLOW_SHELFICE
574 jmc 1.76 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 mlosch 1.56 #endif /* ALLOW_SHELFICE */
586    
587 jmc 1.70 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 adcroft 1.1
593 jmc 1.54 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
594    
595     C--- Prepare for Advection & Coriolis terms:
596 jmc 1.72 C- calculate absolute vorticity
597 jmc 1.54 IF (useAbsVorticity)
598     & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
599 adcroft 1.1
600 gforget 1.79 #ifdef ALLOW_AUTODIFF_TAMC
601     CADJ STORE omega3(:,:) =
602     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
603     #endif
604    
605 jmc 1.5 C-- Horizontal Coriolis terms
606 jmc 1.37 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 jmc 1.46 IF ( useCoriolis .AND.
610 jmc 1.37 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
611     & ) THEN
612     IF (useAbsVorticity) THEN
613 jmc 1.77 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,omega3,hFacZ,r_hFacZ,
614 jmc 1.37 & uCf,myThid)
615 jmc 1.77 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,omega3,hFacZ,r_hFacZ,
616 jmc 1.37 & vCf,myThid)
617     ELSE
618     CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
619     & uCf,vCf,myThid)
620     ENDIF
621 jmc 1.5 DO j=jMin,jMax
622     DO i=iMin,iMax
623 jmc 1.43 gU(i,j,k,bi,bj) = uCf(i,j)
624     gV(i,j,k,bi,bj) = vCf(i,j)
625 jmc 1.5 ENDDO
626 adcroft 1.1 ENDDO
627 jmc 1.15 IF ( writeDiag ) THEN
628 edhill 1.24 IF (snapshot_mdsio) THEN
629     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
630     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
631     ENDIF
632     #ifdef ALLOW_MNC
633     IF (useMNC .AND. snapshot_mnc) THEN
634 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
635 edhill 1.25 & offsets, myThid)
636 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
637 edhill 1.25 & offsets, myThid)
638 edhill 1.24 ENDIF
639     #endif /* ALLOW_MNC */
640 jmc 1.15 ENDIF
641 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
642     IF ( useDiagnostics ) THEN
643     CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
644     CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
645     ENDIF
646     #endif /* ALLOW_DIAGNOSTICS */
647 jmc 1.31 ELSE
648     DO j=jMin,jMax
649     DO i=iMin,iMax
650 jmc 1.43 gU(i,j,k,bi,bj) = 0. _d 0
651     gV(i,j,k,bi,bj) = 0. _d 0
652 jmc 1.31 ENDDO
653     ENDDO
654 jmc 1.5 ENDIF
655 adcroft 1.1
656 gforget 1.79 #ifdef ALLOW_AUTODIFF_TAMC
657     CADJ STORE ucf(:,:) =
658     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
659     CADJ STORE vcf(:,:) =
660     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
661     #endif
662    
663 jmc 1.5 IF (momAdvection) THEN
664 jmc 1.41 C-- Horizontal advection of relative (or absolute) vorticity
665 jmc 1.62 IF ( (highOrderVorticity.OR.upwindVorticity)
666     & .AND.useAbsVorticity ) THEN
667 mlosch 1.80 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
668     & highOrderVorticity,upwindVorticity,
669     & vFld,omega3,r_hFacZ,
670 adcroft 1.20 & uCf,myThid)
671 jmc 1.62 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
672 mlosch 1.80 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,selectVortScheme,
673     & highOrderVorticity, upwindVorticity,
674     & vFld,vort3, r_hFacZ,
675 jmc 1.41 & uCf,myThid)
676 jmc 1.62 ELSEIF ( useAbsVorticity ) THEN
677 mlosch 1.80 CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
678     & vFld,omega3,hFacZ,r_hFacZ,
679 jmc 1.40 & uCf,myThid)
680 adcroft 1.20 ELSE
681 mlosch 1.80 CALL MOM_VI_U_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
682     & vFld,vort3, hFacZ,r_hFacZ,
683 adcroft 1.20 & uCf,myThid)
684     ENDIF
685 jmc 1.5 DO j=jMin,jMax
686     DO i=iMin,iMax
687     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
688     ENDDO
689 adcroft 1.1 ENDDO
690 jmc 1.62 IF ( (highOrderVorticity.OR.upwindVorticity)
691     & .AND.useAbsVorticity ) THEN
692 mlosch 1.80 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
693     & highOrderVorticity, upwindVorticity,
694     & uFld,omega3,r_hFacZ,
695 adcroft 1.20 & vCf,myThid)
696 jmc 1.62 ELSEIF ( (highOrderVorticity.OR.upwindVorticity) ) THEN
697 mlosch 1.80 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,k,selectVortScheme,
698     & highOrderVorticity, upwindVorticity,
699     & uFld,vort3, r_hFacZ,
700 jmc 1.41 & vCf,myThid)
701 jmc 1.62 ELSEIF ( useAbsVorticity ) THEN
702 mlosch 1.80 CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
703     & uFld,omega3,hFacZ,r_hFacZ,
704 jmc 1.40 & vCf,myThid)
705 adcroft 1.20 ELSE
706 mlosch 1.80 CALL MOM_VI_V_CORIOLIS(bi,bj,k,selectVortScheme,useJamartMomAdv,
707     & uFld,vort3, hFacZ,r_hFacZ,
708 adcroft 1.20 & vCf,myThid)
709     ENDIF
710 jmc 1.5 DO j=jMin,jMax
711     DO i=iMin,iMax
712     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
713     ENDDO
714 adcroft 1.1 ENDDO
715    
716 gforget 1.79 #ifdef ALLOW_AUTODIFF_TAMC
717     CADJ STORE ucf(:,:) =
718     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
719     CADJ STORE vcf(:,:) =
720     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
721     #endif
722    
723 jmc 1.15 IF ( writeDiag ) THEN
724 edhill 1.24 IF (snapshot_mdsio) THEN
725     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
726     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
727     ENDIF
728     #ifdef ALLOW_MNC
729     IF (useMNC .AND. snapshot_mnc) THEN
730 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
731 edhill 1.25 & offsets, myThid)
732 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
733 edhill 1.25 & offsets, myThid)
734 edhill 1.24 ENDIF
735     #endif /* ALLOW_MNC */
736 jmc 1.15 ENDIF
737 edhill 1.24
738 jmc 1.7 #ifdef ALLOW_TIMEAVE
739     IF (taveFreq.GT.0.) THEN
740     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
741     & Nr, k, bi, bj, myThid)
742     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
743     & Nr, k, bi, bj, myThid)
744     ENDIF
745 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
746 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
747     IF ( useDiagnostics ) THEN
748     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
749     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
750     ENDIF
751     #endif /* ALLOW_DIAGNOSTICS */
752 jmc 1.7
753 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
754 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
755 jmc 1.77 CALL MOM_VI_U_VERTSHEAR(bi,bj,k,uVel,wVel,uCf,myThid)
756 jmc 1.12 DO j=jMin,jMax
757     DO i=iMin,iMax
758     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
759     ENDDO
760 jmc 1.5 ENDDO
761 jmc 1.77 CALL MOM_VI_V_VERTSHEAR(bi,bj,k,vVel,wVel,vCf,myThid)
762 jmc 1.12 DO j=jMin,jMax
763     DO i=iMin,iMax
764     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
765     ENDDO
766 jmc 1.5 ENDDO
767 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
768     IF ( useDiagnostics ) THEN
769     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
770     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
771     ENDIF
772     #endif /* ALLOW_DIAGNOSTICS */
773 jmc 1.12 ENDIF
774 adcroft 1.1
775     C-- Bernoulli term
776 jmc 1.77 CALL MOM_VI_U_GRAD_KE(bi,bj,k,KE,uCf,myThid)
777 jmc 1.5 DO j=jMin,jMax
778     DO i=iMin,iMax
779     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
780     ENDDO
781     ENDDO
782 jmc 1.77 CALL MOM_VI_V_GRAD_KE(bi,bj,k,KE,vCf,myThid)
783 jmc 1.5 DO j=jMin,jMax
784     DO i=iMin,iMax
785     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
786     ENDDO
787 adcroft 1.1 ENDDO
788 jmc 1.15 IF ( writeDiag ) THEN
789 edhill 1.24 IF (snapshot_mdsio) THEN
790     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
791     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
792     ENDIF
793     #ifdef ALLOW_MNC
794     IF (useMNC .AND. snapshot_mnc) THEN
795 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
796 edhill 1.25 & offsets, myThid)
797 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
798 edhill 1.25 & offsets, myThid)
799 jmc 1.54 ENDIF
800 edhill 1.24 #endif /* ALLOW_MNC */
801 jmc 1.15 ENDIF
802    
803 jmc 1.5 C-- end if momAdvection
804     ENDIF
805    
806 jmc 1.63 C-- 3.D Coriolis term (horizontal momentum, Eastward component: -fprime*w)
807 jmc 1.58 IF ( use3dCoriolis ) THEN
808 jmc 1.57 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
809     DO j=jMin,jMax
810     DO i=iMin,iMax
811     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
812     ENDDO
813     ENDDO
814     IF ( usingCurvilinearGrid ) THEN
815     C- presently, non zero angleSinC array only supported with Curvilinear-Grid
816     CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
817     DO j=jMin,jMax
818     DO i=iMin,iMax
819     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
820     ENDDO
821     ENDDO
822     ENDIF
823     ENDIF
824    
825     C-- Non-Hydrostatic (spherical) metric terms
826     IF ( useNHMTerms ) THEN
827     CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
828     DO j=jMin,jMax
829     DO i=iMin,iMax
830     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
831     ENDDO
832     ENDDO
833     CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
834     DO j=jMin,jMax
835     DO i=iMin,iMax
836     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
837     ENDDO
838     ENDDO
839     ENDIF
840 jmc 1.54
841 jmc 1.5 C-- Set du/dt & dv/dt on boundaries to zero
842 adcroft 1.1 DO j=jMin,jMax
843     DO i=iMin,iMax
844 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
845     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
846 adcroft 1.1 ENDDO
847     ENDDO
848 jmc 1.5
849 jmc 1.22 #ifdef ALLOW_DEBUG
850 jmc 1.65 IF ( debugLevel .GE. debLevC
851 jmc 1.22 & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
852     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
853     & .AND. useCubedSphereExchange ) THEN
854 jmc 1.23 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
855 jmc 1.31 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
856 jmc 1.22 ENDIF
857     #endif /* ALLOW_DEBUG */
858 adcroft 1.2
859 jmc 1.15 IF ( writeDiag ) THEN
860 jmc 1.72 IF (useBiharmonicVisc) THEN
861     CALL WRITE_LOCAL_RL( 'del2u', 'I10', 1, del2u,
862     & bi,bj,k, myIter, myThid )
863     CALL WRITE_LOCAL_RL( 'del2v', 'I10', 1, del2v,
864     & bi,bj,k, myIter, myThid )
865     CALL WRITE_LOCAL_RL( 'dStar', 'I10', 1, dStar,
866     & bi,bj,k, myIter, myThid )
867     CALL WRITE_LOCAL_RL( 'zStar', 'I10', 1, zStar,
868     & bi,bj,k, myIter, myThid )
869     ENDIF
870 edhill 1.24 IF (snapshot_mdsio) THEN
871 jmc 1.54 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
872 jmc 1.72 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3BC,bi,bj,k,myIter,myThid)
873 jmc 1.54 CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
874     CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
875     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
876 jmc 1.72 CALL WRITE_LOCAL_RL( 'Ds', 'I10', 1, strainBC,
877     & bi,bj,k, myIter, myThid )
878 jmc 1.54 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
879     CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
880 edhill 1.24 ENDIF
881     #ifdef ALLOW_MNC
882     IF (useMNC .AND. snapshot_mnc) THEN
883 jmc 1.54 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
884     & offsets, myThid)
885 jmc 1.72 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3BC,
886     & offsets, myThid)
887 jmc 1.54 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
888     & offsets, myThid)
889     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
890 edhill 1.25 & offsets, myThid)
891 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
892 edhill 1.25 & offsets, myThid)
893 jmc 1.72 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strainBC,
894     & offsets, myThid)
895 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
896 edhill 1.25 & offsets, myThid)
897 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
898 edhill 1.25 & offsets, myThid)
899 edhill 1.24 ENDIF
900     #endif /* ALLOW_MNC */
901 adcroft 1.1 ENDIF
902 jmc 1.41
903 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
904     IF ( useDiagnostics ) THEN
905 jmc 1.72 CALL DIAGNOSTICS_FILL(vort3BC,'momVort3',k,1,2,bi,bj,myThid)
906 jmc 1.54 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
907 jmc 1.46 IF (momViscosity) THEN
908 jmc 1.54 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
909 jmc 1.46 ENDIF
910 jmc 1.72 IF ( useVariableVisc .OR. useStrainTensionVisc ) THEN
911     CALL DIAGNOSTICS_FILL(tension, 'Tension ',k,1,2,bi,bj,myThid)
912     CALL DIAGNOSTICS_FILL(strainBC,'Strain ',k,1,2,bi,bj,myThid)
913     ENDIF
914 jmc 1.66 CALL DIAGNOSTICS_FILL(gU(1-OLx,1-OLy,k,bi,bj),
915 jmc 1.54 & 'Um_Advec',k,1,2,bi,bj,myThid)
916 jmc 1.66 CALL DIAGNOSTICS_FILL(gV(1-OLx,1-OLy,k,bi,bj),
917 jmc 1.54 & 'Vm_Advec',k,1,2,bi,bj,myThid)
918 jmc 1.46 ENDIF
919     #endif /* ALLOW_DIAGNOSTICS */
920    
921 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
922 adcroft 1.1
923     RETURN
924     END

  ViewVC Help
Powered by ViewVC 1.1.22