/[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.59 - (hide annotations) (download)
Tue Jul 18 03:23:30 2006 UTC (17 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58r_post, checkpoint58n_post
Changes since 1.58: +30 -1 lines
vort3 partial recomputation warning treated.

1 heimbach 1.59 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.58 2006/07/13 03:02:48 jmc Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4 adcroft 1.21 #include "MOM_VECINV_OPTIONS.h"
5 adcroft 1.1
6 jmc 1.57 SUBROUTINE MOM_VECINV(
7 adcroft 1.1 I bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
8 jmc 1.43 I KappaRU, KappaRV,
9 adcroft 1.1 U fVerU, fVerV,
10 jmc 1.31 O guDiss, gvDiss,
11 jmc 1.15 I myTime, myIter, myThid)
12 adcroft 1.1 C /==========================================================\
13     C | S/R MOM_VECINV |
14     C | o Form the right hand-side of the momentum equation. |
15     C |==========================================================|
16     C | Terms are evaluated one layer at a time working from |
17     C | the bottom to the top. The vertically integrated |
18     C | barotropic flow tendency term is evluated by summing the |
19     C | tendencies. |
20     C | Notes: |
21     C | We have not sorted out an entirely satisfactory formula |
22     C | for the diffusion equation bc with lopping. The present |
23     C | form produces a diffusive flux that does not scale with |
24     C | open-area. Need to do something to solidfy this and to |
25     C | deal "properly" with thin walls. |
26     C \==========================================================/
27     IMPLICIT NONE
28    
29     C == Global variables ==
30     #include "SIZE.h"
31     #include "DYNVARS.h"
32     #include "EEPARAMS.h"
33     #include "PARAMS.h"
34 edhill 1.27 #ifdef ALLOW_MNC
35     #include "MNC_PARAMS.h"
36     #endif
37 adcroft 1.1 #include "GRID.h"
38 jmc 1.7 #ifdef ALLOW_TIMEAVE
39     #include "TIMEAVE_STATV.h"
40     #endif
41 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
42     # include "tamc.h"
43     # include "tamc_keys.h"
44     #endif
45 adcroft 1.1
46     C == Routine arguments ==
47 jmc 1.31 C fVerU :: Flux of momentum in the vertical direction, out of the upper
48     C fVerV :: face of a cell K ( flux into the cell above ).
49     C guDiss :: dissipation tendency (all explicit terms), u component
50     C gvDiss :: dissipation tendency (all explicit terms), v component
51 adcroft 1.1 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
52     C results will be set.
53     C kUp, kDown - Index for upper and lower layers.
54 jmc 1.54 C myThid :: my Thread Id number
55 adcroft 1.1 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
56     _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
57     _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
58     _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
59 jmc 1.31 _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
60     _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61 adcroft 1.1 INTEGER kUp,kDown
62 jmc 1.15 _RL myTime
63 adcroft 1.2 INTEGER myIter
64 adcroft 1.1 INTEGER myThid
65     INTEGER bi,bj,iMin,iMax,jMin,jMax
66    
67 edhill 1.11 #ifdef ALLOW_MOM_VECINV
68 jmc 1.7
69 adcroft 1.2 C == Functions ==
70 jmc 1.38 LOGICAL DIFFERENT_MULTIPLE
71     EXTERNAL DIFFERENT_MULTIPLE
72 adcroft 1.2
73 adcroft 1.1 C == Local variables ==
74     _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 jmc 1.54 _RL vrF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76     _RL uCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77     _RL vCf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78     _RS hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79     _RS r_hFacZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL del2u (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL del2v (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL dStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL zStar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RL tension (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL strain (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     _RL KE (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89     _RL omega3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90     _RL vort3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91     _RL hDiv (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92     _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96     C i,j,k :: Loop counters
97 adcroft 1.1 INTEGER i,j,k
98     C xxxFac - On-off tracer parameters used for switching terms off.
99     _RL ArDudrFac
100     _RL ArDvdrFac
101 jmc 1.54 _RL sideMaskFac
102 adcroft 1.1 LOGICAL bottomDragTerms
103 jmc 1.15 LOGICAL writeDiag
104 baylor 1.50 LOGICAL harmonic,biharmonic,useVariableViscosity
105 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
106     INTEGER imomkey
107     #endif
108 adcroft 1.1
109 edhill 1.25 #ifdef ALLOW_MNC
110     INTEGER offsets(9)
111 edhill 1.53 CHARACTER*(1) pf
112 edhill 1.25 #endif
113    
114 heimbach 1.9 #ifdef ALLOW_AUTODIFF_TAMC
115     C-- only the kDown part of fverU/V is set in this subroutine
116     C-- the kUp is still required
117     C-- In the case of mom_fluxform Kup is set as well
118     C-- (at least in part)
119     fVerU(1,1,kUp) = fVerU(1,1,kUp)
120     fVerV(1,1,kUp) = fVerV(1,1,kUp)
121     #endif
122    
123 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
124     act0 = k - 1
125     max0 = Nr
126     act1 = bi - myBxLo(myThid)
127     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
128     act2 = bj - myByLo(myThid)
129     max2 = myByHi(myThid) - myByLo(myThid) + 1
130     act3 = myThid - 1
131     max3 = nTx*nTy
132     act4 = ikey_dynamics - 1
133     imomkey = (act0 + 1)
134     & + act1*max0
135     & + act2*max0*max1
136     & + act3*max0*max1*max2
137     & + act4*max0*max1*max2*max3
138     #endif /* ALLOW_AUTODIFF_TAMC */
139    
140 jmc 1.38 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
141 adcroft 1.1
142 edhill 1.24 #ifdef ALLOW_MNC
143     IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
144 edhill 1.53 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
145     pf(1:1) = 'D'
146     ELSE
147     pf(1:1) = 'R'
148     ENDIF
149 edhill 1.25 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
150     CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
151 edhill 1.39 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
152 edhill 1.25 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
153 edhill 1.39 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
154 edhill 1.25 ENDIF
155     DO i = 1,9
156     offsets(i) = 0
157     ENDDO
158     offsets(3) = k
159     C write(*,*) 'offsets = ',(offsets(i),i=1,9)
160 edhill 1.24 ENDIF
161     #endif /* ALLOW_MNC */
162    
163 adcroft 1.1 C Initialise intermediate terms
164     DO J=1-OLy,sNy+OLy
165     DO I=1-OLx,sNx+OLx
166 jmc 1.31 vF(i,j) = 0.
167     vrF(i,j) = 0.
168 adcroft 1.1 uCf(i,j) = 0.
169     vCf(i,j) = 0.
170     del2u(i,j) = 0.
171     del2v(i,j) = 0.
172     dStar(i,j) = 0.
173     zStar(i,j) = 0.
174 jmc 1.31 guDiss(i,j)= 0.
175     gvDiss(i,j)= 0.
176 adcroft 1.1 vort3(i,j) = 0.
177 jmc 1.31 omega3(i,j)= 0.
178 jmc 1.54 KE(i,j) = 0.
179 baylor 1.50 viscAh_Z(i,j) = 0.
180     viscAh_D(i,j) = 0.
181     viscA4_Z(i,j) = 0.
182     viscA4_D(i,j) = 0.
183    
184 heimbach 1.8 #ifdef ALLOW_AUTODIFF_TAMC
185     strain(i,j) = 0. _d 0
186     tension(i,j) = 0. _d 0
187 heimbach 1.55 hFacZ(i,j) = 0. _d 0
188 heimbach 1.8 #endif
189 adcroft 1.1 ENDDO
190     ENDDO
191    
192     C-- Term by term tracer parmeters
193     C o U momentum equation
194     ArDudrFac = vfFacMom*1.
195     C o V momentum equation
196     ArDvdrFac = vfFacMom*1.
197    
198 jmc 1.54 C note: using standard stencil (no mask) results in under-estimating
199     C vorticity at a no-slip boundary by a factor of 2 = sideDragFactor
200     IF ( no_slip_sides ) THEN
201     sideMaskFac = sideDragFactor
202     ELSE
203     sideMaskFac = 0. _d 0
204     ENDIF
205    
206 adcroft 1.1 IF ( no_slip_bottom
207     & .OR. bottomDragQuadratic.NE.0.
208     & .OR. bottomDragLinear.NE.0.) THEN
209     bottomDragTerms=.TRUE.
210     ELSE
211     bottomDragTerms=.FALSE.
212     ENDIF
213    
214     C-- Calculate open water fraction at vorticity points
215     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
216    
217     C Make local copies of horizontal flow field
218     DO j=1-OLy,sNy+OLy
219     DO i=1-OLx,sNx+OLx
220     uFld(i,j) = uVel(i,j,k,bi,bj)
221     vFld(i,j) = vVel(i,j,k,bi,bj)
222     ENDDO
223     ENDDO
224    
225 jmc 1.7 C note (jmc) : Dissipation and Vort3 advection do not necesary
226     C use the same maskZ (and hFacZ) => needs 2 call(s)
227     c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
228    
229 jmc 1.52 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
230 adcroft 1.1
231 jmc 1.54 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
232 adcroft 1.1
233 jmc 1.54 IF (momViscosity) THEN
234 jmc 1.57 C-- For viscous term, compute horizontal divergence, tension & strain
235 jmc 1.54 C and mask relative vorticity (free-slip case):
236 adcroft 1.1
237 heimbach 1.59 #ifdef ALLOW_AUTODIFF_TAMC
238     CADJ STORE vort3(:,:) =
239     CADJ & comlev1_bibj_k, key = imomkey, byte = isbyte
240     #endif
241    
242 jmc 1.54 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
243 adcroft 1.1
244 jmc 1.52 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
245    
246     CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
247    
248 jmc 1.54 C- account for no-slip / free-slip BC:
249     DO j=1-Oly,sNy+Oly
250     DO i=1-Olx,sNx+Olx
251     IF ( hFacZ(i,j).EQ.0. ) THEN
252     vort3(i,j) = sideMaskFac*vort3(i,j)
253     strain(i,j) = sideMaskFac*strain(i,j)
254     ENDIF
255     ENDDO
256     ENDDO
257    
258     C-- Calculate Viscosities
259 jmc 1.52 CALL MOM_CALC_VISC(
260 baylor 1.50 I bi,bj,k,
261     O viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
262     O harmonic,biharmonic,useVariableViscosity,
263     I hDiv,vort3,tension,strain,KE,hfacZ,
264     I myThid)
265    
266 adcroft 1.1 C Calculate del^2 u and del^2 v for bi-harmonic term
267 baylor 1.50 IF (biharmonic) THEN
268 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
269     O del2u,del2v,
270     & myThid)
271 jmc 1.48 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
272     CALL MOM_CALC_RELVORT3(bi,bj,k,
273     & del2u,del2v,hFacZ,zStar,myThid)
274 adcroft 1.2 ENDIF
275 baylor 1.47
276 jmc 1.54 C- Strain diagnostics:
277     IF ( writeDiag ) THEN
278     IF (snapshot_mdsio) THEN
279     CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
280     ENDIF
281     #ifdef ALLOW_MNC
282     IF (useMNC .AND. snapshot_mnc) THEN
283     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
284     & offsets, myThid)
285     ENDIF
286     #endif /* ALLOW_MNC */
287     ENDIF
288     #ifdef ALLOW_DIAGNOSTICS
289     IF ( useDiagnostics ) THEN
290     CALL DIAGNOSTICS_FILL(strain, 'Strain ',k,1,2,bi,bj,myThid)
291     ENDIF
292     #endif /* ALLOW_DIAGNOSTICS */
293    
294     C--- Calculate dissipation terms for U and V equations
295 baylor 1.47
296     C in terms of tension and strain
297     IF (useStrainTensionVisc) THEN
298 jmc 1.54 C mask strain as if free-slip since side-drag is computed separately
299     DO j=1-Oly,sNy+Oly
300     DO i=1-Olx,sNx+Olx
301     IF ( hFacZ(i,j).EQ.0. ) strain(i,j) = 0. _d 0
302     ENDDO
303     ENDDO
304 baylor 1.47 CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
305     I hFacZ,
306 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
307     I harmonic,biharmonic,useVariableViscosity,
308 baylor 1.47 O guDiss,gvDiss,
309     I myThid)
310     ELSE
311 adcroft 1.2 C in terms of vorticity and divergence
312 baylor 1.47 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
313     I hFacZ,dStar,zStar,
314 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
315     I harmonic,biharmonic,useVariableViscosity,
316 jmc 1.31 O guDiss,gvDiss,
317 jmc 1.57 & myThid)
318 adcroft 1.3 ENDIF
319 jmc 1.54 C-- if (momViscosity) end of block.
320 adcroft 1.1 ENDIF
321    
322 jmc 1.7 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
323     c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
324    
325 jmc 1.54 C--- Other dissipation terms in Zonal momentum equation
326 adcroft 1.1
327     C-- Vertical flux (fVer is at upper face of "u" cell)
328    
329     C Eddy component of vertical flux (interior component only) -> vrF
330 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
331 jmc 1.44 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
332 adcroft 1.1
333     C Combine fluxes
334 jmc 1.31 DO j=jMin,jMax
335     DO i=iMin,iMax
336     fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
337     ENDDO
338 adcroft 1.1 ENDDO
339    
340 jmc 1.31 C-- Tendency is minus divergence of the fluxes
341     DO j=2-Oly,sNy+Oly-1
342     DO i=2-Olx,sNx+Olx-1
343     guDiss(i,j) = guDiss(i,j)
344 adcroft 1.1 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
345     & *recip_rAw(i,j,bi,bj)
346     & *(
347 jmc 1.42 & fVerU(i,j,kDown) - fVerU(i,j,kUp)
348     & )*rkSign
349 jmc 1.31 ENDDO
350 adcroft 1.1 ENDDO
351 jmc 1.31 ENDIF
352 adcroft 1.1
353 jmc 1.57 C-- No-slip and drag BCs appear as body forces in cell abutting topography
354 adcroft 1.1 IF (momViscosity.AND.no_slip_sides) THEN
355     C- No-slip BCs impose a drag at walls...
356 baylor 1.50 CALL MOM_U_SIDEDRAG(
357     I bi,bj,k,
358     I uFld, del2u, hFacZ,
359     I viscAh_Z,viscA4_Z,
360     I harmonic,biharmonic,useVariableViscosity,
361     O vF,
362     I myThid)
363 adcroft 1.1 DO j=jMin,jMax
364     DO i=iMin,iMax
365 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
366 adcroft 1.1 ENDDO
367     ENDDO
368     ENDIF
369     C- No-slip BCs impose a drag at bottom
370     IF (momViscosity.AND.bottomDragTerms) THEN
371     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
372     DO j=jMin,jMax
373     DO i=iMin,iMax
374 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
375 adcroft 1.1 ENDDO
376     ENDDO
377     ENDIF
378 mlosch 1.56 #ifdef ALLOW_SHELFICE
379     IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
380     CALL SHELFICE_U_DRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
381     DO j=jMin,jMax
382     DO i=iMin,iMax
383     guDiss(i,j) = guDiss(i,j) + vF(i,j)
384     ENDDO
385     ENDDO
386     ENDIF
387     #endif /* ALLOW_SHELFICE */
388    
389 adcroft 1.1
390 jmc 1.54 C--- Other dissipation terms in Meridional momentum equation
391 adcroft 1.1
392     C-- Vertical flux (fVer is at upper face of "v" cell)
393    
394     C Eddy component of vertical flux (interior component only) -> vrF
395 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
396 jmc 1.44 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
397 adcroft 1.1
398     C Combine fluxes -> fVerV
399 jmc 1.31 DO j=jMin,jMax
400     DO i=iMin,iMax
401     fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
402     ENDDO
403 adcroft 1.1 ENDDO
404    
405 jmc 1.31 C-- Tendency is minus divergence of the fluxes
406     DO j=jMin,jMax
407     DO i=iMin,iMax
408     gvDiss(i,j) = gvDiss(i,j)
409 adcroft 1.1 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
410     & *recip_rAs(i,j,bi,bj)
411     & *(
412 jmc 1.42 & fVerV(i,j,kDown) - fVerV(i,j,kUp)
413     & )*rkSign
414 jmc 1.31 ENDDO
415 adcroft 1.1 ENDDO
416 jmc 1.31 ENDIF
417 adcroft 1.1
418 jmc 1.57 C-- No-slip and drag BCs appear as body forces in cell abutting topography
419 adcroft 1.1 IF (momViscosity.AND.no_slip_sides) THEN
420     C- No-slip BCs impose a drag at walls...
421 baylor 1.50 CALL MOM_V_SIDEDRAG(
422     I bi,bj,k,
423     I vFld, del2v, hFacZ,
424     I viscAh_Z,viscA4_Z,
425     I harmonic,biharmonic,useVariableViscosity,
426     O vF,
427     I myThid)
428 adcroft 1.1 DO j=jMin,jMax
429     DO i=iMin,iMax
430 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
431 adcroft 1.1 ENDDO
432     ENDDO
433     ENDIF
434     C- No-slip BCs impose a drag at bottom
435     IF (momViscosity.AND.bottomDragTerms) THEN
436     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
437     DO j=jMin,jMax
438     DO i=iMin,iMax
439 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
440 adcroft 1.1 ENDDO
441     ENDDO
442     ENDIF
443 mlosch 1.56 #ifdef ALLOW_SHELFICE
444     IF (useShelfIce.AND.momViscosity.AND.bottomDragTerms) THEN
445     CALL SHELFICE_V_DRAG(bi,bj,k,vFld,KE,KappaRU,vF,myThid)
446     DO j=jMin,jMax
447     DO i=iMin,iMax
448     gvDiss(i,j) = gvDiss(i,j) + vF(i,j)
449     ENDDO
450     ENDDO
451     ENDIF
452     #endif /* ALLOW_SHELFICE */
453    
454 adcroft 1.1
455 jmc 1.54 C- Vorticity diagnostics:
456     IF ( writeDiag ) THEN
457     IF (snapshot_mdsio) THEN
458     CALL WRITE_LOCAL_RL('Z3','I10',1,vort3, bi,bj,k,myIter,myThid)
459     ENDIF
460     #ifdef ALLOW_MNC
461     IF (useMNC .AND. snapshot_mnc) THEN
462     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
463     & offsets, myThid)
464     ENDIF
465     #endif /* ALLOW_MNC */
466     ENDIF
467     #ifdef ALLOW_DIAGNOSTICS
468     IF ( useDiagnostics ) THEN
469     CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid)
470     ENDIF
471     #endif /* ALLOW_DIAGNOSTICS */
472    
473     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
474    
475     C--- Prepare for Advection & Coriolis terms:
476     C- Mask relative vorticity and calculate absolute vorticity
477     DO j=1-Oly,sNy+Oly
478     DO i=1-Olx,sNx+Olx
479     IF ( hFacZ(i,j).EQ.0. ) vort3(i,j) = 0.
480     ENDDO
481     ENDDO
482     IF (useAbsVorticity)
483     & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
484 adcroft 1.1
485 jmc 1.5 C-- Horizontal Coriolis terms
486 jmc 1.37 c IF (useCoriolis .AND. .NOT.useCDscheme
487     c & .AND. .NOT. useAbsVorticity) THEN
488     C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
489 jmc 1.46 IF ( useCoriolis .AND.
490 jmc 1.37 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
491     & ) THEN
492     IF (useAbsVorticity) THEN
493     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
494     & uCf,myThid)
495     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
496     & vCf,myThid)
497     ELSE
498     CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
499     & uCf,vCf,myThid)
500     ENDIF
501 jmc 1.5 DO j=jMin,jMax
502     DO i=iMin,iMax
503 jmc 1.43 gU(i,j,k,bi,bj) = uCf(i,j)
504     gV(i,j,k,bi,bj) = vCf(i,j)
505 jmc 1.5 ENDDO
506 adcroft 1.1 ENDDO
507 jmc 1.15 IF ( writeDiag ) THEN
508 edhill 1.24 IF (snapshot_mdsio) THEN
509     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
510     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
511     ENDIF
512     #ifdef ALLOW_MNC
513     IF (useMNC .AND. snapshot_mnc) THEN
514 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
515 edhill 1.25 & offsets, myThid)
516 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
517 edhill 1.25 & offsets, myThid)
518 edhill 1.24 ENDIF
519     #endif /* ALLOW_MNC */
520 jmc 1.15 ENDIF
521 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
522     IF ( useDiagnostics ) THEN
523     CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
524     CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
525     ENDIF
526     #endif /* ALLOW_DIAGNOSTICS */
527 jmc 1.31 ELSE
528     DO j=jMin,jMax
529     DO i=iMin,iMax
530 jmc 1.43 gU(i,j,k,bi,bj) = 0. _d 0
531     gV(i,j,k,bi,bj) = 0. _d 0
532 jmc 1.31 ENDDO
533     ENDDO
534 jmc 1.5 ENDIF
535 adcroft 1.1
536 jmc 1.5 IF (momAdvection) THEN
537 jmc 1.41 C-- Horizontal advection of relative (or absolute) vorticity
538     IF (highOrderVorticity.AND.useAbsVorticity) THEN
539     CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
540 adcroft 1.20 & uCf,myThid)
541 jmc 1.40 ELSEIF (highOrderVorticity) THEN
542 jmc 1.41 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
543     & uCf,myThid)
544     ELSEIF (useAbsVorticity) THEN
545     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
546 jmc 1.40 & uCf,myThid)
547 adcroft 1.20 ELSE
548 jmc 1.41 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
549 adcroft 1.20 & uCf,myThid)
550     ENDIF
551 jmc 1.5 DO j=jMin,jMax
552     DO i=iMin,iMax
553     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
554     ENDDO
555 adcroft 1.1 ENDDO
556 jmc 1.41 IF (highOrderVorticity.AND.useAbsVorticity) THEN
557     CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
558 adcroft 1.20 & vCf,myThid)
559 jmc 1.40 ELSEIF (highOrderVorticity) THEN
560 jmc 1.41 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
561     & vCf,myThid)
562     ELSEIF (useAbsVorticity) THEN
563     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
564 jmc 1.40 & vCf,myThid)
565 adcroft 1.20 ELSE
566 jmc 1.41 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
567 adcroft 1.20 & vCf,myThid)
568     ENDIF
569 jmc 1.5 DO j=jMin,jMax
570     DO i=iMin,iMax
571     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
572     ENDDO
573 adcroft 1.1 ENDDO
574    
575 jmc 1.15 IF ( writeDiag ) THEN
576 edhill 1.24 IF (snapshot_mdsio) THEN
577     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
578     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
579     ENDIF
580     #ifdef ALLOW_MNC
581     IF (useMNC .AND. snapshot_mnc) THEN
582 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
583 edhill 1.25 & offsets, myThid)
584 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
585 edhill 1.25 & offsets, myThid)
586 edhill 1.24 ENDIF
587     #endif /* ALLOW_MNC */
588 jmc 1.15 ENDIF
589 edhill 1.24
590 jmc 1.7 #ifdef ALLOW_TIMEAVE
591     IF (taveFreq.GT.0.) THEN
592     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
593     & Nr, k, bi, bj, myThid)
594     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
595     & Nr, k, bi, bj, myThid)
596     ENDIF
597 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
598 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
599     IF ( useDiagnostics ) THEN
600     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
601     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
602     ENDIF
603     #endif /* ALLOW_DIAGNOSTICS */
604 jmc 1.7
605 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
606 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
607     CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
608     DO j=jMin,jMax
609     DO i=iMin,iMax
610     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
611     ENDDO
612 jmc 1.5 ENDDO
613 jmc 1.12 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
614     DO j=jMin,jMax
615     DO i=iMin,iMax
616     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
617     ENDDO
618 jmc 1.5 ENDDO
619 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
620     IF ( useDiagnostics ) THEN
621     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
622     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
623     ENDIF
624     #endif /* ALLOW_DIAGNOSTICS */
625 jmc 1.12 ENDIF
626 adcroft 1.1
627     C-- Bernoulli term
628 jmc 1.5 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
629     DO j=jMin,jMax
630     DO i=iMin,iMax
631     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
632     ENDDO
633     ENDDO
634     CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
635     DO j=jMin,jMax
636     DO i=iMin,iMax
637     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
638     ENDDO
639 adcroft 1.1 ENDDO
640 jmc 1.15 IF ( writeDiag ) THEN
641 edhill 1.24 IF (snapshot_mdsio) THEN
642     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
643     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
644     ENDIF
645     #ifdef ALLOW_MNC
646     IF (useMNC .AND. snapshot_mnc) THEN
647 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
648 edhill 1.25 & offsets, myThid)
649 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
650 edhill 1.25 & offsets, myThid)
651 jmc 1.54 ENDIF
652 edhill 1.24 #endif /* ALLOW_MNC */
653 jmc 1.15 ENDIF
654    
655 jmc 1.5 C-- end if momAdvection
656     ENDIF
657    
658 jmc 1.57 C-- 3.D Coriolis term (horizontal momentum, Eastward component: -f'*w)
659 jmc 1.58 IF ( use3dCoriolis ) THEN
660 jmc 1.57 CALL MOM_U_CORIOLIS_NH(bi,bj,k,wVel,uCf,myThid)
661     DO j=jMin,jMax
662     DO i=iMin,iMax
663     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
664     ENDDO
665     ENDDO
666     IF ( usingCurvilinearGrid ) THEN
667     C- presently, non zero angleSinC array only supported with Curvilinear-Grid
668     CALL MOM_V_CORIOLIS_NH(bi,bj,k,wVel,vCf,myThid)
669     DO j=jMin,jMax
670     DO i=iMin,iMax
671     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
672     ENDDO
673     ENDDO
674     ENDIF
675     ENDIF
676    
677     C-- Non-Hydrostatic (spherical) metric terms
678     IF ( useNHMTerms ) THEN
679     CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,uCf,myThid)
680     DO j=jMin,jMax
681     DO i=iMin,iMax
682     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
683     ENDDO
684     ENDDO
685     CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,vCf,myThid)
686     DO j=jMin,jMax
687     DO i=iMin,iMax
688     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
689     ENDDO
690     ENDDO
691     ENDIF
692 jmc 1.54
693 jmc 1.5 C-- Set du/dt & dv/dt on boundaries to zero
694 adcroft 1.1 DO j=jMin,jMax
695     DO i=iMin,iMax
696 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
697     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
698 adcroft 1.1 ENDDO
699     ENDDO
700 jmc 1.5
701 jmc 1.22 #ifdef ALLOW_DEBUG
702     IF ( debugLevel .GE. debLevB
703     & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
704     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
705     & .AND. useCubedSphereExchange ) THEN
706 jmc 1.23 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
707 jmc 1.31 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
708 jmc 1.22 ENDIF
709     #endif /* ALLOW_DEBUG */
710 adcroft 1.2
711 jmc 1.15 IF ( writeDiag ) THEN
712 edhill 1.24 IF (snapshot_mdsio) THEN
713 jmc 1.54 CALL WRITE_LOCAL_RL('W3','I10',1,omega3, bi,bj,k,myIter,myThid)
714     CALL WRITE_LOCAL_RL('KE','I10',1,KE, bi,bj,k,myIter,myThid)
715     CALL WRITE_LOCAL_RL('D', 'I10',1,hDiv, bi,bj,k,myIter,myThid)
716     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,myThid)
717     CALL WRITE_LOCAL_RL('Du','I10',1,guDiss, bi,bj,k,myIter,myThid)
718     CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss, bi,bj,k,myIter,myThid)
719 edhill 1.24 ENDIF
720     #ifdef ALLOW_MNC
721     IF (useMNC .AND. snapshot_mnc) THEN
722 jmc 1.54 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
723     & offsets, myThid)
724     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
725     & offsets, myThid)
726     CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
727 edhill 1.25 & offsets, myThid)
728 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
729 edhill 1.25 & offsets, myThid)
730 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
731 edhill 1.25 & offsets, myThid)
732 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
733 edhill 1.25 & offsets, myThid)
734 edhill 1.24 ENDIF
735     #endif /* ALLOW_MNC */
736 adcroft 1.1 ENDIF
737 jmc 1.41
738 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
739     IF ( useDiagnostics ) THEN
740 jmc 1.54 CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
741 jmc 1.46 IF (momViscosity) THEN
742 jmc 1.54 CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
743 jmc 1.52 CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
744 jmc 1.54 CALL DIAGNOSTICS_FILL(guDiss, 'Um_Diss ',k,1,2,bi,bj,myThid)
745     CALL DIAGNOSTICS_FILL(gvDiss, 'Vm_Diss ',k,1,2,bi,bj,myThid)
746 jmc 1.46 ENDIF
747 jmc 1.54 CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),
748     & 'Um_Advec',k,1,2,bi,bj,myThid)
749     CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),
750     & 'Vm_Advec',k,1,2,bi,bj,myThid)
751 jmc 1.46 ENDIF
752     #endif /* ALLOW_DIAGNOSTICS */
753    
754 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
755 adcroft 1.1
756     RETURN
757     END

  ViewVC Help
Powered by ViewVC 1.1.22