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

Annotation of /MITgcm/pkg/mom_fluxform/mom_fluxform.F

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


Revision 1.2 - (hide annotations) (download)
Fri Aug 17 18:40:30 2001 UTC (22 years, 9 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre9, checkpoint40pre8, checkpoint40
Changes since 1.1: +5 -4 lines
Added method for dumping intermediate local arrays:
 mdsio_writetile - same as mdsio_writefield except works from inside bi,bj loop
 mdsio_writelocal - same as mdsio_writetile except works for local arrays
 write_local_r? - higher-level wrapper for mdsio_writelocal

Controlled by diagFreq. Defaults to zero (ie. no dumps)

Example given at end of mom_vecinv.F that dumps some local arrays.

1 adcroft 1.2 C $Header: /u/gcmpack/models/MITgcmUV/pkg/mom_fluxform/mom_fluxform.F,v 1.1 2001/08/16 17:16:03 adcroft Exp $
2     C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6     SUBROUTINE MOM_FLUXFORM(
7     I bi,bj,iMin,iMax,jMin,jMax,k,kUp,kDown,
8     I phi_hyd,KappaRU,KappaRV,
9     U fVerU, fVerV,
10 adcroft 1.2 I myCurrentTime, myIter, myThid)
11 adcroft 1.1 C /==========================================================\
12     C | S/R MOM_FLUXFORM |
13     C | o Form the right hand-side of the momentum equation. |
14     C |==========================================================|
15     C | Terms are evaluated one layer at a time working from |
16     C | the bottom to the top. The vertically integrated |
17     C | barotropic flow tendency term is evluated by summing the |
18     C | tendencies. |
19     C | Notes: |
20     C | We have not sorted out an entirely satisfactory formula |
21     C | for the diffusion equation bc with lopping. The present |
22     C | form produces a diffusive flux that does not scale with |
23     C | open-area. Need to do something to solidfy this and to |
24     C | deal "properly" with thin walls. |
25     C \==========================================================/
26     IMPLICIT NONE
27    
28     C == Global variables ==
29     #include "SIZE.h"
30     #include "DYNVARS.h"
31     #include "FFIELDS.h"
32     #include "EEPARAMS.h"
33     #include "PARAMS.h"
34     #include "GRID.h"
35     #include "SURFACE.h"
36    
37     C == Routine arguments ==
38     C fZon - Work array for flux of momentum in the east-west
39     C direction at the west face of a cell.
40     C fMer - Work array for flux of momentum in the north-south
41     C direction at the south face of a cell.
42     C fVerU - Flux of momentum in the vertical
43     C fVerV direction out of the upper face of a cell K
44     C ( flux into the cell above ).
45     C phi_hyd - Hydrostatic pressure
46     C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
47     C results will be set.
48     C kUp, kDown - Index for upper and lower layers.
49     C myThid - Instance number for this innvocation of CALC_MOM_RHS
50     _RL phi_hyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
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     INTEGER kUp,kDown
56 adcroft 1.2 _RL myCurrentTime
57     INTEGER myIter
58 adcroft 1.1 INTEGER myThid
59     INTEGER bi,bj,iMin,iMax,jMin,jMax
60    
61     C == Local variables ==
62     C ab15, ab05 - Weights for Adams-Bashforth time stepping scheme.
63     C i,j,k - Loop counters
64     C wMaskOverride - Land sea flag override for top layer.
65     C afFacMom - Tracer parameters for turning terms
66     C vfFacMom on and off.
67     C pfFacMom afFacMom - Advective terms
68     C cfFacMom vfFacMom - Eddy viscosity terms
69     C mTFacMom pfFacMom - Pressure terms
70     C cfFacMom - Coriolis terms
71     C foFacMom - Forcing
72     C mTFacMom - Metric term
73     C vF - Temporary holding viscous term (Laplacian)
74     C v4F - Temporary holding viscous term (Biharmonic)
75     C cF - Temporary holding coriolis term.
76     C mT - Temporary holding metric terms(s).
77     C pF - Temporary holding pressure|potential gradient terms.
78     C uDudxFac, AhDudxFac, etc ... individual term tracer parameters
79     _RL aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81     _RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL vrF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL cF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88     _RS hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89     _RS r_hFacZ(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90     _RS xA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91     _RS yA(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92     _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL uFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL vFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96     C I,J,K - Loop counters
97     INTEGER i,j,k
98     C rVelMaskOverride - Factor for imposing special surface boundary conditions
99     C ( set according to free-surface condition ).
100     C hFacROpen - Lopped cell factos used tohold fraction of open
101     C hFacRClosed and closed cell wall.
102     _RL rVelMaskOverride
103     C xxxFac - On-off tracer parameters used for switching terms off.
104     _RL uDudxFac
105     _RL AhDudxFac
106     _RL A4DuxxdxFac
107     _RL vDudyFac
108     _RL AhDudyFac
109     _RL A4DuyydyFac
110     _RL rVelDudrFac
111     _RL ArDudrFac
112     _RL fuFac
113     _RL phxFac
114     _RL mtFacU
115     _RL uDvdxFac
116     _RL AhDvdxFac
117     _RL A4DvxxdxFac
118     _RL vDvdyFac
119     _RL AhDvdyFac
120     _RL A4DvyydyFac
121     _RL rVelDvdrFac
122     _RL ArDvdrFac
123     _RL fvFac
124     _RL phyFac
125     _RL vForcFac
126     _RL mtFacV
127     C ab05, ab15 - Adams-Bashforth time-stepping weights.
128     _RL ab05, ab15
129     INTEGER km1,kp1
130     _RL wVelBottomOverride
131     LOGICAL bottomDragTerms
132     _RL KE(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133    
134     km1=MAX(1,k-1)
135     kp1=MIN(Nr,k+1)
136     rVelMaskOverride=1.
137     IF ( k .EQ. 1 ) rVelMaskOverride=freeSurfFac
138     wVelBottomOverride=1.
139     IF (k.EQ.Nr) wVelBottomOverride=0.
140    
141     C Initialise intermediate terms
142     DO J=1-OLy,sNy+OLy
143     DO I=1-OLx,sNx+OLx
144     aF(i,j) = 0.
145     vF(i,j) = 0.
146     v4F(i,j) = 0.
147     vrF(i,j) = 0.
148     cF(i,j) = 0.
149     mT(i,j) = 0.
150     pF(i,j) = 0.
151     fZon(i,j) = 0.
152     fMer(i,j) = 0.
153     ENDDO
154     ENDDO
155    
156     C-- Term by term tracer parmeters
157     C o U momentum equation
158     uDudxFac = afFacMom*1.
159     AhDudxFac = vfFacMom*1.
160     A4DuxxdxFac = vfFacMom*1.
161     vDudyFac = afFacMom*1.
162     AhDudyFac = vfFacMom*1.
163     A4DuyydyFac = vfFacMom*1.
164     rVelDudrFac = afFacMom*1.
165     ArDudrFac = vfFacMom*1.
166     mTFacU = mtFacMom*1.
167     fuFac = cfFacMom*1.
168     phxFac = pfFacMom*1.
169     C o V momentum equation
170     uDvdxFac = afFacMom*1.
171     AhDvdxFac = vfFacMom*1.
172     A4DvxxdxFac = vfFacMom*1.
173     vDvdyFac = afFacMom*1.
174     AhDvdyFac = vfFacMom*1.
175     A4DvyydyFac = vfFacMom*1.
176     rVelDvdrFac = afFacMom*1.
177     ArDvdrFac = vfFacMom*1.
178     mTFacV = mtFacMom*1.
179     fvFac = cfFacMom*1.
180     phyFac = pfFacMom*1.
181     vForcFac = foFacMom*1.
182    
183     IF ( no_slip_bottom
184     & .OR. bottomDragQuadratic.NE.0.
185     & .OR. bottomDragLinear.NE.0.) THEN
186     bottomDragTerms=.TRUE.
187     ELSE
188     bottomDragTerms=.FALSE.
189     ENDIF
190    
191     C-- with stagger time stepping, grad Phi_Hyp is directly incoporated in TIMESTEP
192     IF (staggerTimeStep) THEN
193     phxFac = 0.
194     phyFac = 0.
195     ENDIF
196    
197     C-- Adams-Bashforth weighting factors
198     ab15 = 1.5 _d 0 + abEps
199     ab05 = -0.5 _d 0 - abEps
200    
201     C-- Calculate open water fraction at vorticity points
202     CALL MOM_CALC_HFACZ(bi,bj,k,hFacZ,r_hFacZ,myThid)
203    
204     C---- Calculate common quantities used in both U and V equations
205     C Calculate tracer cell face open areas
206     DO j=1-OLy,sNy+OLy
207     DO i=1-OLx,sNx+OLx
208     xA(i,j) = _dyG(i,j,bi,bj)
209     & *drF(k)*_hFacW(i,j,k,bi,bj)
210     yA(i,j) = _dxG(i,j,bi,bj)
211     & *drF(k)*_hFacS(i,j,k,bi,bj)
212     ENDDO
213     ENDDO
214    
215     C Make local copies of horizontal flow field
216     DO j=1-OLy,sNy+OLy
217     DO i=1-OLx,sNx+OLx
218     uFld(i,j) = uVel(i,j,k,bi,bj)
219     vFld(i,j) = vVel(i,j,k,bi,bj)
220     ENDDO
221     ENDDO
222    
223     C Calculate velocity field "volume transports" through tracer cell faces.
224     DO j=1-OLy,sNy+OLy
225     DO i=1-OLx,sNx+OLx
226     uTrans(i,j) = uFld(i,j)*xA(i,j)
227     vTrans(i,j) = vFld(i,j)*yA(i,j)
228     ENDDO
229     ENDDO
230    
231     CALL MOM_CALC_KE(bi,bj,k,uFld,vFld,KE,myThid)
232    
233     C---- Zonal momentum equation starts here
234    
235     C Bi-harmonic term del^2 U -> v4F
236     IF (momViscosity)
237     & CALL MOM_U_DEL2U(bi,bj,k,uFld,hFacZ,v4f,myThid)
238    
239     C--- Calculate mean and eddy fluxes between cells for zonal flow.
240    
241     C-- Zonal flux (fZon is at east face of "u" cell)
242    
243     C Mean flow component of zonal flux -> aF
244     IF (momAdvection)
245     & CALL MOM_U_ADV_UU(bi,bj,k,uTrans,uFld,aF,myThid)
246    
247     C Laplacian and bi-harmonic terms -> vF
248     IF (momViscosity)
249     & CALL MOM_U_XVISCFLUX(bi,bj,k,uFld,v4F,vF,myThid)
250    
251     C Combine fluxes -> fZon
252     DO j=jMin,jMax
253     DO i=iMin,iMax
254     fZon(i,j) = uDudxFac*aF(i,j) + AhDudxFac*vF(i,j)
255     ENDDO
256     ENDDO
257    
258     C-- Meridional flux (fMer is at south face of "u" cell)
259    
260     C Mean flow component of meridional flux
261     IF (momAdvection)
262     & CALL MOM_U_ADV_VU(bi,bj,k,vTrans,uFld,aF,myThid)
263    
264     C Laplacian and bi-harmonic term
265     IF (momViscosity)
266     & CALL MOM_U_YVISCFLUX(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
267    
268     C Combine fluxes -> fMer
269     DO j=jMin,jMax
270     DO i=iMin,iMax
271     fMer(i,j) = vDudyFac*aF(i,j) + AhDudyFac*vF(i,j)
272     ENDDO
273     ENDDO
274    
275     C-- Vertical flux (fVer is at upper face of "u" cell)
276    
277     C-- Free surface correction term (flux at k=1)
278     IF (momAdvection.AND.k.EQ.1) THEN
279     CALL MOM_U_ADV_WU(bi,bj,k,uVel,wVel,af,myThid)
280     DO j=jMin,jMax
281     DO i=iMin,iMax
282     fVerU(i,j,kUp) = af(i,j)
283     ENDDO
284     ENDDO
285     ENDIF
286     C Mean flow component of vertical flux (at k+1) -> aF
287     IF (momAdvection)
288     & CALL MOM_U_ADV_WU(bi,bj,k+1,uVel,wVel,af,myThid)
289    
290     C Eddy component of vertical flux (interior component only) -> vrF
291     IF (momViscosity.AND..NOT.implicitViscosity)
292     & CALL MOM_U_RVISCFLUX(bi,bj,k,uVel,KappaRU,vrF,myThid)
293    
294     C Combine fluxes
295     DO j=jMin,jMax
296     DO i=iMin,iMax
297     fVerU(i,j,kDown) = rVelDudrFac*aF(i,j) + ArDudrFac*vrF(i,j)
298     ENDDO
299     ENDDO
300    
301     C--- Hydrostatic term ( -1/rhoConst . dphi/dx )
302     IF (momPressureForcing) THEN
303     DO j=jMin,jMax
304     DO i=iMin,iMax
305     pf(i,j) = - _recip_dxC(i,j,bi,bj)
306     & *(phi_hyd(i,j,k)-phi_hyd(i-1,j,k))
307     ENDDO
308     ENDDO
309     ENDIF
310    
311     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
312     DO j=jMin,jMax
313     DO i=iMin,iMax
314     gU(i,j,k,bi,bj) =
315     #ifdef OLD_UV_GEOM
316     & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)/
317     & ( 0.5 _d 0*(rA(i,j,bi,bj)+rA(i-1,j,bi,bj)) )
318     #else
319     & -_recip_hFacW(i,j,k,bi,bj)*recip_drF(k)
320     & *recip_rAw(i,j,bi,bj)
321     #endif
322     & *(fZon(i,j ) - fZon(i-1,j)
323     & +fMer(i,j+1) - fMer(i ,j)
324     & +fVerU(i,j,kUp)*rkFac - fVerU(i,j,kDown)*rkFac
325     & )
326     & _PHM( +phxFac * pf(i,j) )
327     ENDDO
328     ENDDO
329    
330     C-- No-slip and drag BCs appear as body forces in cell abutting topography
331     IF (momViscosity.AND.no_slip_sides) THEN
332     C- No-slip BCs impose a drag at walls...
333     CALL MOM_U_SIDEDRAG(bi,bj,k,uFld,v4F,hFacZ,vF,myThid)
334     DO j=jMin,jMax
335     DO i=iMin,iMax
336     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
337     ENDDO
338     ENDDO
339     ENDIF
340     C- No-slip BCs impose a drag at bottom
341     IF (momViscosity.AND.bottomDragTerms) THEN
342     CALL MOM_U_BOTTOMDRAG(bi,bj,k,uFld,KE,KappaRU,vF,myThid)
343     DO j=jMin,jMax
344     DO i=iMin,iMax
345     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+vF(i,j)
346     ENDDO
347     ENDDO
348     ENDIF
349    
350     C-- Forcing term
351     IF (momForcing)
352     & CALL EXTERNAL_FORCING_U(
353     I iMin,iMax,jMin,jMax,bi,bj,k,
354     I myCurrentTime,myThid)
355    
356     C-- Metric terms for curvilinear grid systems
357     IF (usingSphericalPolarMTerms) THEN
358     C o Spherical polar grid metric terms
359     CALL MOM_U_METRIC_NH(bi,bj,k,uFld,wVel,mT,myThid)
360     DO j=jMin,jMax
361     DO i=iMin,iMax
362     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
363     ENDDO
364     ENDDO
365     CALL MOM_U_METRIC_SPHERE(bi,bj,k,uFld,vFld,mT,myThid)
366     DO j=jMin,jMax
367     DO i=iMin,iMax
368     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+mTFacU*mT(i,j)
369     ENDDO
370     ENDDO
371     ENDIF
372    
373     C-- Set du/dt on boundaries to zero
374     DO j=jMin,jMax
375     DO i=iMin,iMax
376     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)*_maskW(i,j,k,bi,bj)
377     ENDDO
378     ENDDO
379    
380    
381     C---- Meridional momentum equation starts here
382    
383     C Bi-harmonic term del^2 V -> v4F
384     IF (momViscosity)
385     & CALL MOM_V_DEL2V(bi,bj,k,vFld,hFacZ,v4f,myThid)
386    
387     C--- Calculate mean and eddy fluxes between cells for meridional flow.
388    
389     C-- Zonal flux (fZon is at west face of "v" cell)
390    
391     C Mean flow component of zonal flux -> aF
392     IF (momAdvection)
393     & CALL MOM_V_ADV_UV(bi,bj,k,uTrans,vFld,af,myThid)
394    
395     C Laplacian and bi-harmonic terms -> vF
396     IF (momViscosity)
397     & CALL MOM_V_XVISCFLUX(bi,bj,k,vFld,v4f,hFacZ,vf,myThid)
398    
399     C Combine fluxes -> fZon
400     DO j=jMin,jMax
401     DO i=iMin,iMax
402     fZon(i,j) = uDvdxFac*aF(i,j) + AhDvdxFac*vF(i,j)
403     ENDDO
404     ENDDO
405    
406     C-- Meridional flux (fMer is at north face of "v" cell)
407    
408     C Mean flow component of meridional flux
409     IF (momAdvection)
410     & CALL MOM_V_ADV_VV(bi,bj,k,vTrans,vFld,af,myThid)
411    
412     C Laplacian and bi-harmonic term
413     IF (momViscosity)
414     & CALL MOM_V_YVISCFLUX(bi,bj,k,vFld,v4f,vf,myThid)
415    
416     C Combine fluxes -> fMer
417     DO j=jMin,jMax
418     DO i=iMin,iMax
419     fMer(i,j) = vDvdyFac*aF(i,j) + AhDvdyFac*vF(i,j)
420     ENDDO
421     ENDDO
422    
423     C-- Vertical flux (fVer is at upper face of "v" cell)
424    
425     C-- Free surface correction term (flux at k=1)
426     IF (momAdvection.AND.k.EQ.1) THEN
427     CALL MOM_V_ADV_WV(bi,bj,k,vVel,wVel,af,myThid)
428     DO j=jMin,jMax
429     DO i=iMin,iMax
430     fVerV(i,j,kUp) = af(i,j)
431     ENDDO
432     ENDDO
433     ENDIF
434     C o Mean flow component of vertical flux
435     IF (momAdvection)
436     & CALL MOM_V_ADV_WV(bi,bj,k+1,vVel,wVel,af,myThid)
437    
438     C Eddy component of vertical flux (interior component only) -> vrF
439     IF (momViscosity.AND..NOT.implicitViscosity)
440     & CALL MOM_V_RVISCFLUX(bi,bj,k,vVel,KappaRV,vrf,myThid)
441    
442     C Combine fluxes -> fVerV
443     DO j=jMin,jMax
444     DO i=iMin,iMax
445     fVerV(i,j,kDown) = rVelDvdrFac*aF(i,j) + ArDvdrFac*vrF(i,j)
446     ENDDO
447     ENDDO
448    
449     C--- Hydorstatic term (-1/rhoConst . dphi/dy )
450     IF (momPressureForcing) THEN
451     DO j=jMin,jMax
452     DO i=iMin,iMax
453     pF(i,j) = -_recip_dyC(i,j,bi,bj)
454     & *(phi_hyd(i,j,k)-phi_hyd(i,j-1,k))
455     ENDDO
456     ENDDO
457     ENDIF
458    
459     C-- Tendency is minus divergence of the fluxes + coriolis + pressure term
460     DO j=jMin,jMax
461     DO i=iMin,iMax
462     gV(i,j,k,bi,bj) =
463     #ifdef OLD_UV_GEOM
464     & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)/
465     & ( 0.5 _d 0*(_rA(i,j,bi,bj)+_rA(i,j-1,bi,bj)) )
466     #else
467     & -_recip_hFacS(i,j,k,bi,bj)*recip_drF(k)
468     & *recip_rAs(i,j,bi,bj)
469     #endif
470     & *(fZon(i+1,j) - fZon(i,j )
471     & +fMer(i,j ) - fMer(i,j-1)
472     & +fVerV(i,j,kUp)*rkFac - fVerV(i,j,kDown)*rkFac
473     & )
474     & _PHM( +phyFac*pf(i,j) )
475     ENDDO
476     ENDDO
477    
478     C-- No-slip and drag BCs appear as body forces in cell abutting topography
479     IF (momViscosity.AND.no_slip_sides) THEN
480     C- No-slip BCs impose a drag at walls...
481     CALL MOM_V_SIDEDRAG(bi,bj,k,vFld,v4F,hFacZ,vF,myThid)
482     DO j=jMin,jMax
483     DO i=iMin,iMax
484     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
485     ENDDO
486     ENDDO
487     ENDIF
488     C- No-slip BCs impose a drag at bottom
489     IF (momViscosity.AND.bottomDragTerms) THEN
490     CALL MOM_V_BOTTOMDRAG(bi,bj,k,vFld,KE,KappaRV,vF,myThid)
491     DO j=jMin,jMax
492     DO i=iMin,iMax
493     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+vF(i,j)
494     ENDDO
495     ENDDO
496     ENDIF
497    
498     C-- Forcing term
499     IF (momForcing)
500     & CALL EXTERNAL_FORCING_V(
501     I iMin,iMax,jMin,jMax,bi,bj,k,
502     I myCurrentTime,myThid)
503    
504     C-- Metric terms for curvilinear grid systems
505     IF (usingSphericalPolarMTerms) THEN
506     C o Spherical polar grid metric terms
507     CALL MOM_V_METRIC_NH(bi,bj,k,vFld,wVel,mT,myThid)
508     DO j=jMin,jMax
509     DO i=iMin,iMax
510     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
511     ENDDO
512     ENDDO
513     CALL MOM_V_METRIC_SPHERE(bi,bj,k,uFld,mT,myThid)
514     DO j=jMin,jMax
515     DO i=iMin,iMax
516     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+mTFacV*mT(i,j)
517     ENDDO
518     ENDDO
519     ENDIF
520    
521     C-- Set dv/dt on boundaries to zero
522     DO j=jMin,jMax
523     DO i=iMin,iMax
524     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)*_maskS(i,j,k,bi,bj)
525     ENDDO
526     ENDDO
527    
528     C-- Coriolis term
529     C Note. As coded here, coriolis will not work with "thin walls"
530     #ifdef INCLUDE_CD_CODE
531     CALL MOM_CDSCHEME(bi,bj,k,phi_hyd,myThid)
532     #else
533     CALL MOM_U_CORIOLIS(bi,bj,k,vFld,cf,myThid)
534     DO j=jMin,jMax
535     DO i=iMin,iMax
536     gU(i,j,k,bi,bj) = gU(i,j,k,bi,bj)+fuFac*cf(i,j)
537     ENDDO
538     ENDDO
539     CALL MOM_V_CORIOLIS(bi,bj,k,uFld,cf,myThid)
540     DO j=jMin,jMax
541     DO i=iMin,iMax
542     gV(i,j,k,bi,bj) = gV(i,j,k,bi,bj)+fvFac*cf(i,j)
543     ENDDO
544     ENDDO
545     #endif /* INCLUDE_CD_CODE */
546    
547     RETURN
548     END

  ViewVC Help
Powered by ViewVC 1.1.22