/[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.36 - (hide annotations) (download)
Wed Apr 6 18:42:12 2005 UTC (19 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57g_pre, checkpoint57g_post, checkpoint57f_post
Changes since 1.35: +5 -5 lines
use baseTime as time origin ; DIFF_BASE_MULTIPLE replaces DIFFERENT_MULTIPLE

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

  ViewVC Help
Powered by ViewVC 1.1.22