/[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.53 - (hide annotations) (download)
Thu Sep 29 12:19:52 2005 UTC (18 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint57t_post, checkpint57u_post
Changes since 1.52: +21 -15 lines
 o make mnc honor the writeBinaryPrec flag for all the non-pickup and
   non-diagnostics output types

1 edhill 1.53 C $Header: /u/gcmpack/MITgcm/pkg/mom_vecinv/mom_vecinv.F,v 1.52 2005/09/28 15:53:20 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     SUBROUTINE MOM_VECINV(
7     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 adcroft 1.1
42     C == Routine arguments ==
43 jmc 1.31 C fVerU :: Flux of momentum in the vertical direction, out of the upper
44     C fVerV :: face of a cell K ( flux into the cell above ).
45     C guDiss :: dissipation tendency (all explicit terms), u component
46     C gvDiss :: dissipation tendency (all explicit terms), v component
47 adcroft 1.1 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
48     C results will be set.
49     C kUp, kDown - Index for upper and lower layers.
50     C myThid - Instance number for this innvocation of CALC_MOM_RHS
51     _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
52     _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
53     _RL fVerU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
54     _RL fVerV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
55 jmc 1.31 _RL guDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
56     _RL gvDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
57 adcroft 1.1 INTEGER kUp,kDown
58 jmc 1.15 _RL myTime
59 adcroft 1.2 INTEGER myIter
60 adcroft 1.1 INTEGER myThid
61     INTEGER bi,bj,iMin,iMax,jMin,jMax
62    
63 edhill 1.11 #ifdef ALLOW_MOM_VECINV
64 jmc 1.7
65 adcroft 1.2 C == Functions ==
66 jmc 1.38 LOGICAL DIFFERENT_MULTIPLE
67     EXTERNAL DIFFERENT_MULTIPLE
68 adcroft 1.2
69 adcroft 1.1 C == Local variables ==
70     _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71     _RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72     _RL uCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73     _RL vCf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 jmc 1.29 c _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 adcroft 1.1 _RL del2u(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76     _RL del2v(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 adcroft 1.3 _RL tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78     _RL strain(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 adcroft 1.1 _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL dStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL zStar(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     C I,J,K - Loop counters
86     INTEGER i,j,k
87     C xxxFac - On-off tracer parameters used for switching terms off.
88     _RL ArDudrFac
89 jmc 1.29 c _RL mtFacU
90 adcroft 1.1 _RL ArDvdrFac
91 jmc 1.29 c _RL mtFacV
92 adcroft 1.1 LOGICAL bottomDragTerms
93 jmc 1.15 LOGICAL writeDiag
94 adcroft 1.1 _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL omega3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96     _RL vort3(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97     _RL hDiv(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 baylor 1.50 _RL viscAh_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99     _RL viscAh_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100     _RL viscA4_Z(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101     _RL viscA4_D(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102     LOGICAL harmonic,biharmonic,useVariableViscosity
103 adcroft 1.1
104 edhill 1.25 #ifdef ALLOW_MNC
105     INTEGER offsets(9)
106 edhill 1.53 CHARACTER*(1) pf
107 edhill 1.25 #endif
108    
109 heimbach 1.9 #ifdef ALLOW_AUTODIFF_TAMC
110     C-- only the kDown part of fverU/V is set in this subroutine
111     C-- the kUp is still required
112     C-- In the case of mom_fluxform Kup is set as well
113     C-- (at least in part)
114     fVerU(1,1,kUp) = fVerU(1,1,kUp)
115     fVerV(1,1,kUp) = fVerV(1,1,kUp)
116     #endif
117    
118 jmc 1.38 writeDiag = DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock)
119 adcroft 1.1
120 edhill 1.24 #ifdef ALLOW_MNC
121     IF (useMNC .AND. snapshot_mnc .AND. writeDiag) THEN
122 edhill 1.53 IF ( writeBinaryPrec .EQ. precFloat64 ) THEN
123     pf(1:1) = 'D'
124     ELSE
125     pf(1:1) = 'R'
126     ENDIF
127 edhill 1.25 IF ((bi .EQ. 1).AND.(bj .EQ. 1).AND.(k .EQ. 1)) THEN
128     CALL MNC_CW_SET_UDIM('mom_vi', -1, myThid)
129 edhill 1.39 CALL MNC_CW_RL_W_S('D','mom_vi',0,0,'T',myTime,myThid)
130 edhill 1.25 CALL MNC_CW_SET_UDIM('mom_vi', 0, myThid)
131 edhill 1.39 CALL MNC_CW_I_W_S('I','mom_vi',0,0,'iter',myIter,myThid)
132 edhill 1.25 ENDIF
133     DO i = 1,9
134     offsets(i) = 0
135     ENDDO
136     offsets(3) = k
137     C write(*,*) 'offsets = ',(offsets(i),i=1,9)
138 edhill 1.24 ENDIF
139     #endif /* ALLOW_MNC */
140    
141 adcroft 1.1 C Initialise intermediate terms
142     DO J=1-OLy,sNy+OLy
143     DO I=1-OLx,sNx+OLx
144 jmc 1.31 vF(i,j) = 0.
145     vrF(i,j) = 0.
146 adcroft 1.1 uCf(i,j) = 0.
147     vCf(i,j) = 0.
148 jmc 1.31 c mT(i,j) = 0.
149 adcroft 1.1 del2u(i,j) = 0.
150     del2v(i,j) = 0.
151     dStar(i,j) = 0.
152     zStar(i,j) = 0.
153 jmc 1.31 guDiss(i,j)= 0.
154     gvDiss(i,j)= 0.
155 adcroft 1.1 vort3(i,j) = 0.
156 jmc 1.31 omega3(i,j)= 0.
157     ke(i,j) = 0.
158 baylor 1.50 viscAh_Z(i,j) = 0.
159     viscAh_D(i,j) = 0.
160     viscA4_Z(i,j) = 0.
161     viscA4_D(i,j) = 0.
162    
163 heimbach 1.8 #ifdef ALLOW_AUTODIFF_TAMC
164     strain(i,j) = 0. _d 0
165     tension(i,j) = 0. _d 0
166     #endif
167 adcroft 1.1 ENDDO
168     ENDDO
169    
170     C-- Term by term tracer parmeters
171     C o U momentum equation
172     ArDudrFac = vfFacMom*1.
173 jmc 1.29 c mTFacU = mtFacMom*1.
174 adcroft 1.1 C o V momentum equation
175     ArDvdrFac = vfFacMom*1.
176 jmc 1.29 c mTFacV = mtFacMom*1.
177 adcroft 1.1
178     IF ( no_slip_bottom
179     & .OR. bottomDragQuadratic.NE.0.
180     & .OR. bottomDragLinear.NE.0.) THEN
181     bottomDragTerms=.TRUE.
182     ELSE
183     bottomDragTerms=.FALSE.
184     ENDIF
185    
186     C-- Calculate open water fraction at vorticity points
187     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
188    
189     C Make local copies of horizontal flow field
190     DO j=1-OLy,sNy+OLy
191     DO i=1-OLx,sNx+OLx
192     uFld(i,j) = uVel(i,j,k,bi,bj)
193     vFld(i,j) = vVel(i,j,k,bi,bj)
194     ENDDO
195     ENDDO
196    
197 jmc 1.7 C note (jmc) : Dissipation and Vort3 advection do not necesary
198     C use the same maskZ (and hFacZ) => needs 2 call(s)
199     c CALL MOM_VI_HFACZ_DISS(bi,bj,k,hFacZ,r_hFacZ,myThid)
200    
201 jmc 1.52 CALL MOM_CALC_KE(bi,bj,k,selectKEscheme,uFld,vFld,KE,myThid)
202 adcroft 1.1
203 adcroft 1.17 CALL MOM_CALC_HDIV(bi,bj,k,2,uFld,vFld,hDiv,myThid)
204 adcroft 1.1
205 adcroft 1.18 CALL MOM_CALC_RELVORT3(bi,bj,k,uFld,vFld,hFacZ,vort3,myThid)
206 adcroft 1.1
207 adcroft 1.20 IF (useAbsVorticity)
208     & CALL MOM_CALC_ABSVORT3(bi,bj,k,vort3,omega3,myThid)
209 adcroft 1.1
210     IF (momViscosity) THEN
211 baylor 1.50 C Calculate Viscosities
212 jmc 1.52 CALL MOM_CALC_TENSION(bi,bj,k,uFld,vFld,tension,myThid)
213    
214     CALL MOM_CALC_STRAIN(bi,bj,k,uFld,vFld,hFacZ,strain,myThid)
215    
216     CALL MOM_CALC_VISC(
217 baylor 1.50 I bi,bj,k,
218     O viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
219     O harmonic,biharmonic,useVariableViscosity,
220     I hDiv,vort3,tension,strain,KE,hfacZ,
221     I myThid)
222    
223 adcroft 1.1 C Calculate del^2 u and del^2 v for bi-harmonic term
224 baylor 1.50 IF (biharmonic) THEN
225 adcroft 1.2 CALL MOM_VI_DEL2UV(bi,bj,k,hDiv,vort3,hFacZ,
226     O del2u,del2v,
227     & myThid)
228 jmc 1.48 CALL MOM_CALC_HDIV(bi,bj,k,2,del2u,del2v,dStar,myThid)
229     CALL MOM_CALC_RELVORT3(bi,bj,k,
230     & del2u,del2v,hFacZ,zStar,myThid)
231 adcroft 1.2 ENDIF
232 baylor 1.47
233 adcroft 1.1 C Calculate dissipation terms for U and V equations
234 baylor 1.47
235     C in terms of tension and strain
236     IF (useStrainTensionVisc) THEN
237     CALL MOM_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
238     I hFacZ,
239 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
240     I harmonic,biharmonic,useVariableViscosity,
241 baylor 1.47 O guDiss,gvDiss,
242     I myThid)
243     ELSE
244 adcroft 1.2 C in terms of vorticity and divergence
245 baylor 1.47 CALL MOM_VI_HDISSIP(bi,bj,k,hDiv,vort3,tension,strain,KE,
246     I hFacZ,dStar,zStar,
247 baylor 1.50 I viscAh_Z,viscAh_D,viscA4_Z,viscA4_D,
248     I harmonic,biharmonic,useVariableViscosity,
249 jmc 1.31 O guDiss,gvDiss,
250 baylor 1.47 & myThid)
251 adcroft 1.3 ENDIF
252 adcroft 1.1 ENDIF
253    
254 jmc 1.7 C- Return to standard hfacZ (min-4) and mask vort3 accordingly:
255     c CALL MOM_VI_MASK_VORT3(bi,bj,k,hFacZ,r_hFacZ,vort3,myThid)
256    
257 adcroft 1.1 C---- Zonal momentum equation starts here
258    
259     C-- Vertical flux (fVer is at upper face of "u" cell)
260    
261     C Eddy component of vertical flux (interior component only) -> vrF
262 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
263 jmc 1.44 CALL MOM_U_RVISCFLUX(bi,bj,k+1,uVel,KappaRU,vrF,myThid)
264 adcroft 1.1
265     C Combine fluxes
266 jmc 1.31 DO j=jMin,jMax
267     DO i=iMin,iMax
268     fVerU(i,j,kDown) = ArDudrFac*vrF(i,j)
269     ENDDO
270 adcroft 1.1 ENDDO
271    
272 jmc 1.31 C-- Tendency is minus divergence of the fluxes
273     DO j=2-Oly,sNy+Oly-1
274     DO i=2-Olx,sNx+Olx-1
275     guDiss(i,j) = guDiss(i,j)
276 adcroft 1.1 & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
277     & *recip_rAw(i,j,bi,bj)
278     & *(
279 jmc 1.42 & fVerU(i,j,kDown) - fVerU(i,j,kUp)
280     & )*rkSign
281 jmc 1.31 ENDDO
282 adcroft 1.1 ENDDO
283 jmc 1.31 ENDIF
284 adcroft 1.1
285     C-- No-slip and drag BCs appear as body forces in cell abutting topography
286     IF (momViscosity.AND.no_slip_sides) THEN
287     C- No-slip BCs impose a drag at walls...
288 baylor 1.50 CALL MOM_U_SIDEDRAG(
289     I bi,bj,k,
290     I uFld, del2u, hFacZ,
291     I viscAh_Z,viscA4_Z,
292     I harmonic,biharmonic,useVariableViscosity,
293     O vF,
294     I myThid)
295 adcroft 1.1 DO j=jMin,jMax
296     DO i=iMin,iMax
297 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
298 adcroft 1.1 ENDDO
299     ENDDO
300     ENDIF
301 heimbach 1.8
302 adcroft 1.1 C- No-slip BCs impose a drag at bottom
303     IF (momViscosity.AND.bottomDragTerms) THEN
304     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
305     DO j=jMin,jMax
306     DO i=iMin,iMax
307 jmc 1.31 guDiss(i,j) = guDiss(i,j)+vF(i,j)
308 adcroft 1.1 ENDDO
309     ENDDO
310     ENDIF
311    
312     C-- Metric terms for curvilinear grid systems
313     c IF (usingSphericalPolarMTerms) THEN
314     C o Spherical polar grid metric terms
315     c CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
316     c DO j=jMin,jMax
317     c DO i=iMin,iMax
318     c gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
319     c ENDDO
320     c ENDDO
321     c ENDIF
322    
323     C---- Meridional momentum equation starts here
324    
325     C-- Vertical flux (fVer is at upper face of "v" cell)
326    
327     C Eddy component of vertical flux (interior component only) -> vrF
328 jmc 1.31 IF (momViscosity.AND..NOT.implicitViscosity) THEN
329 jmc 1.44 CALL MOM_V_RVISCFLUX(bi,bj,k+1,vVel,KappaRV,vrF,myThid)
330 adcroft 1.1
331     C Combine fluxes -> fVerV
332 jmc 1.31 DO j=jMin,jMax
333     DO i=iMin,iMax
334     fVerV(i,j,kDown) = ArDvdrFac*vrF(i,j)
335     ENDDO
336 adcroft 1.1 ENDDO
337    
338 jmc 1.31 C-- Tendency is minus divergence of the fluxes
339     DO j=jMin,jMax
340     DO i=iMin,iMax
341     gvDiss(i,j) = gvDiss(i,j)
342 adcroft 1.1 & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
343     & *recip_rAs(i,j,bi,bj)
344     & *(
345 jmc 1.42 & fVerV(i,j,kDown) - fVerV(i,j,kUp)
346     & )*rkSign
347 jmc 1.31 ENDDO
348 adcroft 1.1 ENDDO
349 jmc 1.31 ENDIF
350 adcroft 1.1
351     C-- No-slip and drag BCs appear as body forces in cell abutting topography
352     IF (momViscosity.AND.no_slip_sides) THEN
353     C- No-slip BCs impose a drag at walls...
354 baylor 1.50 CALL MOM_V_SIDEDRAG(
355     I bi,bj,k,
356     I vFld, del2v, hFacZ,
357     I viscAh_Z,viscA4_Z,
358     I harmonic,biharmonic,useVariableViscosity,
359     O vF,
360     I myThid)
361 adcroft 1.1 DO j=jMin,jMax
362     DO i=iMin,iMax
363 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
364 adcroft 1.1 ENDDO
365     ENDDO
366     ENDIF
367     C- No-slip BCs impose a drag at bottom
368     IF (momViscosity.AND.bottomDragTerms) THEN
369     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
370     DO j=jMin,jMax
371     DO i=iMin,iMax
372 jmc 1.31 gvDiss(i,j) = gvDiss(i,j)+vF(i,j)
373 adcroft 1.1 ENDDO
374     ENDDO
375     ENDIF
376    
377     C-- Metric terms for curvilinear grid systems
378     c IF (usingSphericalPolarMTerms) THEN
379     C o Spherical polar grid metric terms
380     c CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
381     c DO j=jMin,jMax
382     c DO i=iMin,iMax
383     c gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
384     c ENDDO
385     c ENDDO
386     c ENDIF
387    
388 jmc 1.5 C-- Horizontal Coriolis terms
389 jmc 1.37 c IF (useCoriolis .AND. .NOT.useCDscheme
390     c & .AND. .NOT. useAbsVorticity) THEN
391     C- jmc: change it to keep the Coriolis terms when useAbsVorticity=T & momAdvection=F
392 jmc 1.46 IF ( useCoriolis .AND.
393 jmc 1.37 & .NOT.( useCDscheme .OR. useAbsVorticity.AND.momAdvection )
394     & ) THEN
395     IF (useAbsVorticity) THEN
396     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
397     & uCf,myThid)
398     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
399     & vCf,myThid)
400     ELSE
401     CALL MOM_VI_CORIOLIS(bi,bj,k,uFld,vFld,hFacZ,r_hFacZ,
402     & uCf,vCf,myThid)
403     ENDIF
404 jmc 1.5 DO j=jMin,jMax
405     DO i=iMin,iMax
406 jmc 1.43 gU(i,j,k,bi,bj) = uCf(i,j)
407     gV(i,j,k,bi,bj) = vCf(i,j)
408 jmc 1.5 ENDDO
409 adcroft 1.1 ENDDO
410 jmc 1.46
411 jmc 1.15 IF ( writeDiag ) THEN
412 edhill 1.24 IF (snapshot_mdsio) THEN
413     CALL WRITE_LOCAL_RL('fV','I10',1,uCf,bi,bj,k,myIter,myThid)
414     CALL WRITE_LOCAL_RL('fU','I10',1,vCf,bi,bj,k,myIter,myThid)
415     ENDIF
416     #ifdef ALLOW_MNC
417     IF (useMNC .AND. snapshot_mnc) THEN
418 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fV', uCf,
419 edhill 1.25 & offsets, myThid)
420 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'fU', vCf,
421 edhill 1.25 & offsets, myThid)
422 edhill 1.24 ENDIF
423     #endif /* ALLOW_MNC */
424 jmc 1.15 ENDIF
425 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
426     IF ( useDiagnostics ) THEN
427     CALL DIAGNOSTICS_FILL(uCf,'Um_Cori ',k,1,2,bi,bj,myThid)
428     CALL DIAGNOSTICS_FILL(vCf,'Vm_Cori ',k,1,2,bi,bj,myThid)
429     ENDIF
430     #endif /* ALLOW_DIAGNOSTICS */
431    
432 jmc 1.31 ELSE
433     DO j=jMin,jMax
434     DO i=iMin,iMax
435 jmc 1.43 gU(i,j,k,bi,bj) = 0. _d 0
436     gV(i,j,k,bi,bj) = 0. _d 0
437 jmc 1.31 ENDDO
438     ENDDO
439 jmc 1.5 ENDIF
440 adcroft 1.1
441 jmc 1.5 IF (momAdvection) THEN
442 jmc 1.41 C-- Horizontal advection of relative (or absolute) vorticity
443     IF (highOrderVorticity.AND.useAbsVorticity) THEN
444     CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,omega3,r_hFacZ,
445 adcroft 1.20 & uCf,myThid)
446 jmc 1.40 ELSEIF (highOrderVorticity) THEN
447 jmc 1.41 CALL MOM_VI_U_CORIOLIS_C4(bi,bj,k,vFld,vort3, r_hFacZ,
448     & uCf,myThid)
449     ELSEIF (useAbsVorticity) THEN
450     CALL MOM_VI_U_CORIOLIS(bi,bj,K,vFld,omega3,hFacZ,r_hFacZ,
451 jmc 1.40 & uCf,myThid)
452 adcroft 1.20 ELSE
453 jmc 1.41 CALL MOM_VI_U_CORIOLIS(bi,bj,k,vFld,vort3, hFacZ,r_hFacZ,
454 adcroft 1.20 & uCf,myThid)
455     ENDIF
456 jmc 1.5 DO j=jMin,jMax
457     DO i=iMin,iMax
458     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
459     ENDDO
460 adcroft 1.1 ENDDO
461 jmc 1.41 IF (highOrderVorticity.AND.useAbsVorticity) THEN
462     CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,omega3,r_hFacZ,
463 adcroft 1.20 & vCf,myThid)
464 jmc 1.40 ELSEIF (highOrderVorticity) THEN
465 jmc 1.41 CALL MOM_VI_V_CORIOLIS_C4(bi,bj,K,uFld,vort3, r_hFacZ,
466     & vCf,myThid)
467     ELSEIF (useAbsVorticity) THEN
468     CALL MOM_VI_V_CORIOLIS(bi,bj,K,uFld,omega3,hFacZ,r_hFacZ,
469 jmc 1.40 & vCf,myThid)
470 adcroft 1.20 ELSE
471 jmc 1.41 CALL MOM_VI_V_CORIOLIS(bi,bj,k,uFld,vort3, hFacZ,r_hFacZ,
472 adcroft 1.20 & vCf,myThid)
473     ENDIF
474 jmc 1.5 DO j=jMin,jMax
475     DO i=iMin,iMax
476     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
477     ENDDO
478 adcroft 1.1 ENDDO
479    
480 jmc 1.15 IF ( writeDiag ) THEN
481 edhill 1.24 IF (snapshot_mdsio) THEN
482     CALL WRITE_LOCAL_RL('zV','I10',1,uCf,bi,bj,k,myIter,myThid)
483     CALL WRITE_LOCAL_RL('zU','I10',1,vCf,bi,bj,k,myIter,myThid)
484     ENDIF
485     #ifdef ALLOW_MNC
486     IF (useMNC .AND. snapshot_mnc) THEN
487 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zV', uCf,
488 edhill 1.25 & offsets, myThid)
489 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'zU', vCf,
490 edhill 1.25 & offsets, myThid)
491 edhill 1.24 ENDIF
492     #endif /* ALLOW_MNC */
493 jmc 1.15 ENDIF
494 edhill 1.24
495 jmc 1.7 #ifdef ALLOW_TIMEAVE
496     IF (taveFreq.GT.0.) THEN
497     CALL TIMEAVE_CUMUL_1K1T(uZetatave,vCf,deltaTClock,
498     & Nr, k, bi, bj, myThid)
499     CALL TIMEAVE_CUMUL_1K1T(vZetatave,uCf,deltaTClock,
500     & Nr, k, bi, bj, myThid)
501     ENDIF
502 dimitri 1.13 #endif /* ALLOW_TIMEAVE */
503 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
504     IF ( useDiagnostics ) THEN
505     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvZ3',k,1,2,bi,bj,myThid)
506     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvZ3',k,1,2,bi,bj,myThid)
507     ENDIF
508     #endif /* ALLOW_DIAGNOSTICS */
509 jmc 1.7
510 jmc 1.5 C-- Vertical shear terms (-w*du/dr & -w*dv/dr)
511 jmc 1.12 IF ( .NOT. momImplVertAdv ) THEN
512     CALL MOM_VI_U_VERTSHEAR(bi,bj,K,uVel,wVel,uCf,myThid)
513     DO j=jMin,jMax
514     DO i=iMin,iMax
515     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
516     ENDDO
517 jmc 1.5 ENDDO
518 jmc 1.12 CALL MOM_VI_V_VERTSHEAR(bi,bj,K,vVel,wVel,vCf,myThid)
519     DO j=jMin,jMax
520     DO i=iMin,iMax
521     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
522     ENDDO
523 jmc 1.5 ENDDO
524 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
525     IF ( useDiagnostics ) THEN
526     CALL DIAGNOSTICS_FILL(uCf,'Um_AdvRe',k,1,2,bi,bj,myThid)
527     CALL DIAGNOSTICS_FILL(vCf,'Vm_AdvRe',k,1,2,bi,bj,myThid)
528     ENDIF
529     #endif /* ALLOW_DIAGNOSTICS */
530 jmc 1.12 ENDIF
531 adcroft 1.1
532     C-- Bernoulli term
533 jmc 1.5 CALL MOM_VI_U_GRAD_KE(bi,bj,K,KE,uCf,myThid)
534     DO j=jMin,jMax
535     DO i=iMin,iMax
536     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+uCf(i,j)
537     ENDDO
538     ENDDO
539     CALL MOM_VI_V_GRAD_KE(bi,bj,K,KE,vCf,myThid)
540     DO j=jMin,jMax
541     DO i=iMin,iMax
542     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vCf(i,j)
543     ENDDO
544 adcroft 1.1 ENDDO
545 jmc 1.15 IF ( writeDiag ) THEN
546 edhill 1.24 IF (snapshot_mdsio) THEN
547     CALL WRITE_LOCAL_RL('KEx','I10',1,uCf,bi,bj,k,myIter,myThid)
548     CALL WRITE_LOCAL_RL('KEy','I10',1,vCf,bi,bj,k,myIter,myThid)
549     ENDIF
550     #ifdef ALLOW_MNC
551     IF (useMNC .AND. snapshot_mnc) THEN
552 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEx', uCf,
553 edhill 1.25 & offsets, myThid)
554 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj, 'KEy', vCf,
555 edhill 1.25 & offsets, myThid)
556     ENDIF
557 edhill 1.24 #endif /* ALLOW_MNC */
558 jmc 1.15 ENDIF
559    
560 jmc 1.5 C-- end if momAdvection
561     ENDIF
562    
563     C-- Set du/dt & dv/dt on boundaries to zero
564 adcroft 1.1 DO j=jMin,jMax
565     DO i=iMin,iMax
566 jmc 1.5 gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
567     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
568 adcroft 1.1 ENDDO
569     ENDDO
570 jmc 1.5
571 jmc 1.22 #ifdef ALLOW_DEBUG
572     IF ( debugLevel .GE. debLevB
573     & .AND. k.EQ.4 .AND. myIter.EQ.nIter0
574     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
575     & .AND. useCubedSphereExchange ) THEN
576 jmc 1.23 CALL DEBUG_CS_CORNER_UV( ' uDiss,vDiss from MOM_VECINV',
577 jmc 1.31 & guDiss,gvDiss, k, standardMessageUnit,bi,bj,myThid )
578 jmc 1.22 ENDIF
579     #endif /* ALLOW_DEBUG */
580 adcroft 1.2
581 jmc 1.15 IF ( writeDiag ) THEN
582 edhill 1.24 IF (snapshot_mdsio) THEN
583     CALL WRITE_LOCAL_RL('Ds','I10',1,strain,bi,bj,k,myIter,myThid)
584     CALL WRITE_LOCAL_RL('Dt','I10',1,tension,bi,bj,k,myIter,
585     & myThid)
586 jmc 1.31 CALL WRITE_LOCAL_RL('Du','I10',1,guDiss,bi,bj,k,myIter,myThid)
587     CALL WRITE_LOCAL_RL('Dv','I10',1,gvDiss,bi,bj,k,myIter,myThid)
588 edhill 1.24 CALL WRITE_LOCAL_RL('Z3','I10',1,vort3,bi,bj,k,myIter,myThid)
589     CALL WRITE_LOCAL_RL('W3','I10',1,omega3,bi,bj,k,myIter,myThid)
590     CALL WRITE_LOCAL_RL('KE','I10',1,KE,bi,bj,k,myIter,myThid)
591 jmc 1.46 CALL WRITE_LOCAL_RL('D','I10',1,hDiv,bi,bj,k,myIter,myThid)
592 edhill 1.24 ENDIF
593     #ifdef ALLOW_MNC
594     IF (useMNC .AND. snapshot_mnc) THEN
595 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Ds',strain,
596 edhill 1.25 & offsets, myThid)
597 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dt',tension,
598 edhill 1.25 & offsets, myThid)
599 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Du',guDiss,
600 edhill 1.25 & offsets, myThid)
601 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Dv',gvDiss,
602 edhill 1.25 & offsets, myThid)
603 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'Z3',vort3,
604 edhill 1.25 & offsets, myThid)
605 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'W3',omega3,
606 edhill 1.25 & offsets, myThid)
607 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'KE',KE,
608 edhill 1.25 & offsets, myThid)
609 edhill 1.53 CALL MNC_CW_RL_W_OFFSET(pf,'mom_vi',bi,bj,'D', hDiv,
610 edhill 1.25 & offsets, myThid)
611 edhill 1.24 ENDIF
612     #endif /* ALLOW_MNC */
613 adcroft 1.1 ENDIF
614 jmc 1.41
615 jmc 1.46 #ifdef ALLOW_DIAGNOSTICS
616     IF ( useDiagnostics ) THEN
617     CALL DIAGNOSTICS_FILL(KE, 'momKE ',k,1,2,bi,bj,myThid)
618     CALL DIAGNOSTICS_FILL(hDiv, 'momHDiv ',k,1,2,bi,bj,myThid)
619     CALL DIAGNOSTICS_FILL(vort3, 'momVort3',k,1,2,bi,bj,myThid)
620     CALL DIAGNOSTICS_FILL(gU(1-Olx,1-Oly,k,bi,bj),
621     & 'Um_Advec',k,1,2,bi,bj,myThid)
622     CALL DIAGNOSTICS_FILL(gV(1-Olx,1-Oly,k,bi,bj),
623     & 'Vm_Advec',k,1,2,bi,bj,myThid)
624     IF (momViscosity) THEN
625 jmc 1.52 CALL DIAGNOSTICS_FILL(strain,'Strain ',k,1,2,bi,bj,myThid)
626     CALL DIAGNOSTICS_FILL(tension,'Tension ',k,1,2,bi,bj,myThid)
627 jmc 1.46 CALL DIAGNOSTICS_FILL(guDiss,'Um_Diss ',k,1,2,bi,bj,myThid)
628     CALL DIAGNOSTICS_FILL(gvDiss,'Vm_Diss ',k,1,2,bi,bj,myThid)
629     ENDIF
630     ENDIF
631     #endif /* ALLOW_DIAGNOSTICS */
632    
633 edhill 1.11 #endif /* ALLOW_MOM_VECINV */
634 adcroft 1.1
635     RETURN
636     END

  ViewVC Help
Powered by ViewVC 1.1.22