/[MITgcm]/MITgcm/model/src/calc_gw.F
ViewVC logotype

Annotation of /MITgcm/model/src/calc_gw.F

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


Revision 1.34 - (hide annotations) (download)
Thu Jul 13 21:32:38 2006 UTC (17 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58m_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58r_post, checkpoint58n_post
Changes since 1.33: +39 -18 lines
put back side-drag (by calling new S/R MOM_W_SIDEDRAG)

1 jmc 1.34 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.33 2006/07/13 14:26:50 jmc Exp $
2 cnh 1.9 C !DESCRIPTION: \bv
3 jmc 1.7 C $Name: $
4 edhill 1.10
5 adcroft 1.1 #include "CPP_OPTIONS.h"
6    
7 cnh 1.9 CBOP
8     C !ROUTINE: CALC_GW
9     C !INTERFACE:
10 jmc 1.25 SUBROUTINE CALC_GW(
11 jmc 1.28 I bi, bj, KappaRU, KappaRV,
12     I myTime, myIter, myThid )
13 cnh 1.9 C !DESCRIPTION: \bv
14     C *==========================================================*
15 jmc 1.25 C | S/R CALC_GW
16 jmc 1.30 C | o Calculate vertical velocity tendency terms
17     C | ( Non-Hydrostatic only )
18 cnh 1.9 C *==========================================================*
19 jmc 1.30 C | In NH, the vertical momentum tendency must be
20 jmc 1.25 C | calculated explicitly and included as a source term
21 cnh 1.9 C | for a 3d pressure eqn. Calculate that term here.
22     C | This routine is not used in HYD calculations.
23     C *==========================================================*
24 jmc 1.25 C \ev
25 cnh 1.9
26     C !USES:
27 adcroft 1.1 IMPLICIT NONE
28     C == Global variables ==
29     #include "SIZE.h"
30     #include "EEPARAMS.h"
31     #include "PARAMS.h"
32     #include "GRID.h"
33 jmc 1.34 #include "SURFACE.h"
34 jmc 1.22 #include "DYNVARS.h"
35     #include "NH_VARS.h"
36 adcroft 1.1
37 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
38 adcroft 1.1 C == Routine arguments ==
39 jmc 1.28 C bi,bj :: current tile indices
40     C KappaRU :: vertical viscosity at U points
41     C KappaRV :: vertical viscosity at V points
42     C myTime :: Current time in simulation
43     C myIter :: Current iteration number in simulation
44     C myThid :: Thread number for this instance of the routine.
45     INTEGER bi,bj
46 baylor 1.27 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
47     _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
48 jmc 1.20 _RL myTime
49     INTEGER myIter
50 adcroft 1.1 INTEGER myThid
51    
52 adcroft 1.3 #ifdef ALLOW_NONHYDROSTATIC
53    
54 cnh 1.9 C !LOCAL VARIABLES:
55 adcroft 1.1 C == Local variables ==
56 jmc 1.30 C iMin,iMax
57     C jMin,jMax
58     C xA :: W-Cell face area normal to X
59     C yA :: W-Cell face area normal to Y
60     C rThickC_W :: thickness (in r-units) of W-Cell at Western Edge
61     C rThickC_S :: thickness (in r-units) of W-Cell at Southern Edge
62 jmc 1.34 C rThickC_C :: thickness (in r-units) of W-Cell (centered on W pt)
63 jmc 1.30 C recip_rThickC :: reciprol thickness of W-Cell (centered on W-point)
64     C flx_NS :: vertical momentum flux, meridional direction
65     C flx_EW :: vertical momentum flux, zonal direction
66     C flxAdvUp :: vertical mom. advective flux, vertical direction (@ level k-1)
67     C flxDisUp :: vertical mom. dissipation flux, vertical direction (@ level k-1)
68     C flx_Dn :: vertical momentum flux, vertical direction (@ level k)
69     C gwDiss :: vertical momentum dissipation tendency
70     C i,j,k :: Loop counters
71 jmc 1.28 INTEGER iMin,iMax,jMin,jMax
72 jmc 1.30 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74     _RL rThickC_W (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75     _RL rThickC_S (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 jmc 1.34 _RL rThickC_C (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 jmc 1.30 _RL recip_rThickC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 jmc 1.28 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79     _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80     _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 jmc 1.30 _RL flxAdvUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82     _RL flxDisUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
83     _RL gwDiss(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RL gwAdd (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85 jmc 1.28 _RL del2w (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     INTEGER i,j,k, kp1
87 adcroft 1.1 _RL wOverride
88 jmc 1.30 _RL tmp_WbarZ
89     _RL uTrans, vTrans, rTrans
90 jmc 1.28 _RL viscLoc
91 jmc 1.30 _RL halfRL
92     _RS halfRS, zeroRS
93     PARAMETER( halfRL = 0.5D0 )
94     PARAMETER( halfRS = 0.5 , zeroRS = 0. )
95 cnh 1.9 CEOP
96 edhill 1.10
97 jmc 1.25 C Catch barotropic mode
98     IF ( Nr .LT. 2 ) RETURN
99    
100 mlosch 1.14 iMin = 1
101     iMax = sNx
102     jMin = 1
103     jMax = sNy
104    
105 jmc 1.28 C-- Initialise gW to zero
106     DO k=1,Nr
107     DO j=1-OLy,sNy+OLy
108     DO i=1-OLx,sNx+OLx
109 adcroft 1.1 gW(i,j,k,bi,bj) = 0.
110     ENDDO
111     ENDDO
112 jmc 1.28 ENDDO
113 jmc 1.30 C- Initialise gwDiss to zero
114     DO j=1-OLy,sNy+OLy
115     DO i=1-OLx,sNx+OLx
116     gwDiss(i,j) = 0.
117     ENDDO
118     ENDDO
119 adcroft 1.1
120 jmc 1.28 C-- Boundaries condition at top
121 jmc 1.30 DO j=1-OLy,sNy+OLy
122     DO i=1-OLx,sNx+OLx
123     flxAdvUp(i,j) = 0.
124     flxDisUp(i,j) = 0.
125 jmc 1.28 ENDDO
126     ENDDO
127 adcroft 1.1
128 jmc 1.28 C--- Sweep down column
129     DO k=2,Nr
130     kp1=k+1
131     wOverRide=1.
132     IF (k.EQ.Nr) THEN
133     kp1=Nr
134 adcroft 1.1 wOverRide=0.
135 jmc 1.28 ENDIF
136 jmc 1.30 C-- Compute grid factor arround a W-point:
137     DO j=1-Oly,sNy+Oly
138     DO i=1-Olx,sNx+Olx
139     C- note: assume fluid @ smaller k than bottom: does not work in p-coordinate !
140     IF ( maskC(i,j,k,bi,bj).EQ.0. ) THEN
141     recip_rThickC(i,j) = 0.
142     ELSE
143     recip_rThickC(i,j) = 1. _d 0 /
144 jmc 1.33 & ( drF(k-1)*halfRS
145 jmc 1.30 & + drF( k )*MIN( _hFacC(i,j, k ,bi,bj), halfRS )
146     & )
147     ENDIF
148 jmc 1.34 c IF (momViscosity) THEN
149     #ifdef NONLIN_FRSURF
150     rThickC_C(i,j) =
151     & drF(k-1)*MAX( h0FacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
152     & + drF( k )*MIN( h0FacC(i,j,k ,bi,bj), halfRS )
153     #else
154     rThickC_C(i,j) =
155     & drF(k-1)*MAX( _hFacC(i,j,k-1,bi,bj)-halfRS, zeroRS )
156     & + drF( k )*MIN( _hFacC(i,j,k ,bi,bj), halfRS )
157     #endif
158     rThickC_W(i,j) =
159     & drF(k-1)*MAX( _hFacW(i,j,k-1,bi,bj)-halfRS, zeroRS )
160     & + drF( k )*MIN( _hFacW(i,j,k ,bi,bj), halfRS )
161     rThickC_S(i,j) =
162     & drF(k-1)*MAX( _hFacS(i,j,k-1,bi,bj)-halfRS, zeroRS )
163     & + drF( k )*MIN( _hFacS(i,j, k ,bi,bj), halfRS )
164 jmc 1.30 C W-Cell Western face area:
165     xA(i,j) = _dyG(i,j,bi,bj)*rThickC_W(i,j)
166     C W-Cell Southern face area:
167     yA(i,j) = _dxG(i,j,bi,bj)*rThickC_S(i,j)
168 jmc 1.34 c ENDIF
169 jmc 1.30 ENDDO
170     ENDDO
171    
172 jmc 1.28 C-- horizontal bi-harmonic dissipation
173     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
174     C- calculate the horizontal Laplacian of vertical flow
175 mlosch 1.18 C Zonal flux d/dx W
176     DO j=1-Oly,sNy+Oly
177 jmc 1.30 flx_EW(1-Olx,j)=0.
178 mlosch 1.18 DO i=1-Olx+1,sNx+Olx
179 jmc 1.30 flx_EW(i,j) =
180     & (wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
181     & *_recip_dxC(i,j,bi,bj)*xA(i,j)
182 mlosch 1.18 #ifdef COSINEMETH_III
183 jmc 1.30 & *sqcosFacU(j,bi,bj)
184 mlosch 1.18 #endif
185     ENDDO
186 jmc 1.28 ENDDO
187 mlosch 1.18 C Meridional flux d/dy W
188     DO i=1-Olx,sNx+Olx
189 jmc 1.30 flx_NS(i,1-Oly)=0.
190 mlosch 1.18 ENDDO
191     DO j=1-Oly+1,sNy+Oly
192     DO i=1-Olx,sNx+Olx
193 jmc 1.30 flx_NS(i,j) =
194     & (wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
195     & *_recip_dyC(i,j,bi,bj)*yA(i,j)
196 mlosch 1.18 #ifdef ISOTROPIC_COS_SCALING
197     #ifdef COSINEMETH_III
198 jmc 1.30 & *sqCosFacV(j,bi,bj)
199 mlosch 1.18 #endif
200     #endif
201     ENDDO
202     ENDDO
203 jmc 1.25
204 mlosch 1.18 C del^2 W
205     C Difference of zonal fluxes ...
206     DO j=1-Oly,sNy+Oly
207     DO i=1-Olx,sNx+Olx-1
208 jmc 1.30 del2w(i,j)=flx_EW(i+1,j)-flx_EW(i,j)
209 mlosch 1.18 ENDDO
210     del2w(sNx+Olx,j)=0.
211     ENDDO
212    
213     C ... add difference of meridional fluxes and divide by volume
214     DO j=1-Oly,sNy+Oly-1
215     DO i=1-Olx,sNx+Olx
216 jmc 1.30 del2w(i,j) = ( del2w(i,j)
217     & +(flx_NS(i,j+1)-flx_NS(i,j))
218     & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
219 mlosch 1.18 ENDDO
220     ENDDO
221     C-- No-slip BCs impose a drag at walls...
222     CML ************************************************************
223     CML No-slip Boundary conditions for bi-harmonic dissipation
224     CML need to be implemented here!
225     CML ************************************************************
226 jmc 1.28 ELSE
227 jmc 1.21 C- Initialize del2w to zero:
228     DO j=1-Oly,sNy+Oly
229     DO i=1-Olx,sNx+Olx
230     del2w(i,j) = 0. _d 0
231     ENDDO
232     ENDDO
233 jmc 1.28 ENDIF
234 mlosch 1.18
235 jmc 1.30 IF (momViscosity) THEN
236     C Viscous Flux on Western face
237     DO j=jMin,jMax
238     DO i=iMin,iMax+1
239     flx_EW(i,j)=
240     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i-1,j,k,bi,bj))*halfRL
241     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
242 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
243     cOld & *_recip_dxC(i,j,bi,bj)*rThickC_W(i,j)
244 jmc 1.34 & *cosFacU(j,bi,bj)
245 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i-1,j,k,bi,bj))*halfRL
246     & *(del2w(i,j)-del2w(i-1,j))
247 jmc 1.31 & *_recip_dxC(i,j,bi,bj)*xA(i,j)
248     cOld & *_recip_dxC(i,j,bi,bj)*drC(k)
249 mlosch 1.18 #ifdef COSINEMETH_III
250 jmc 1.30 & *sqCosFacU(j,bi,bj)
251 mlosch 1.18 #else
252 jmc 1.34 & *cosFacU(j,bi,bj)
253 mlosch 1.18 #endif
254 jmc 1.30 ENDDO
255     ENDDO
256     C Viscous Flux on Southern face
257     DO j=jMin,jMax+1
258     DO i=iMin,iMax
259     flx_NS(i,j)=
260     & - (viscAh_W(i,j,k,bi,bj)+viscAh_W(i,j-1,k,bi,bj))*halfRL
261     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
262 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
263     cOld & *_recip_dyC(i,j,bi,bj)*rThickC_S(i,j)
264 jmc 1.34 #ifdef ISOTROPIC_COS_SCALING
265     & *cosFacV(j,bi,bj)
266     #endif
267 jmc 1.30 & + (viscA4_W(i,j,k,bi,bj)+viscA4_W(i,j-1,k,bi,bj))*halfRL
268     & *(del2w(i,j)-del2w(i,j-1))
269 jmc 1.31 & *_recip_dyC(i,j,bi,bj)*yA(i,j)
270     cOld & *_recip_dyC(i,j,bi,bj)*drC(k)
271 jmc 1.30 #ifdef ISOTROPIC_COS_SCALING
272 mlosch 1.18 #ifdef COSINEMETH_III
273 jmc 1.30 & *sqCosFacV(j,bi,bj)
274 jmc 1.25 #else
275 jmc 1.34 & *cosFacV(j,bi,bj)
276 jmc 1.30 #endif
277 mlosch 1.18 #endif
278 jmc 1.30 ENDDO
279     ENDDO
280     C Viscous Flux on Lower face of W-Cell (= at tracer-cell center, level k)
281     DO j=jMin,jMax
282     DO i=iMin,iMax
283     C Interpolate vert viscosity to center of tracer-cell (level k):
284     viscLoc = ( KappaRU(i,j,k) +KappaRU(i+1,j,k)
285     & +KappaRU(i,j,kp1)+KappaRU(i+1,j,kp1)
286     & +KappaRV(i,j,k) +KappaRV(i,j+1,k)
287     & +KappaRV(i,j,kp1)+KappaRV(i,j+1,kp1)
288     & )*0.125 _d 0
289     flx_Dn(i,j) =
290     & - viscLoc*( wVel(i,j,kp1,bi,bj)*wOverRide
291     & -wVel(i,j, k ,bi,bj) )*rkSign
292 jmc 1.31 & *recip_drF(k)*rA(i,j,bi,bj)
293     cOld & *recip_drF(k)
294 jmc 1.30 ENDDO
295     ENDDO
296     C Tendency is minus divergence of viscous fluxes:
297     DO j=jMin,jMax
298     DO i=iMin,iMax
299     gwDiss(i,j) =
300 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
301     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
302     & + ( flx_Dn(i,j)-flxDisUp(i,j) )*rkSign
303     & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
304     cOld gwDiss(i,j) =
305     cOld & -(
306     cOld & +_recip_dxF(i,j,bi,bj)*( flx_EW(i+1,j)-flx_EW(i,j) )
307     cOld & +_recip_dyF(i,j,bi,bj)*( flx_NS(i,j+1)-flx_NS(i,j) )
308     cOld & + ( flxDisUp(i,j)-flx_Dn(i,j) )
309     c & )*recip_rThickC(i,j)
310     cOld & )*recip_drC(k)
311 jmc 1.28 C-- prepare for next level (k+1)
312 jmc 1.30 flxDisUp(i,j)=flx_Dn(i,j)
313     ENDDO
314     ENDDO
315     ENDIF
316    
317     IF (no_slip_sides) THEN
318     C- No-slip BCs impose a drag at walls...
319 jmc 1.34 CALL MOM_W_SIDEDRAG(
320     I bi,bj,k,
321     I wVel, del2w,
322     I rThickC_C, recip_rThickC,
323     I viscAh_W, viscA4_W,
324     O gwAdd,
325     I myThid )
326     DO j=jMin,jMax
327     DO i=iMin,iMax
328     gwDiss(i,j) = gwDiss(i,j) + gwAdd(i,j)
329     ENDDO
330     ENDDO
331 jmc 1.30 ENDIF
332    
333     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
334    
335     IF ( momAdvection ) THEN
336     C Advective Flux on Western face
337     DO j=jMin,jMax
338     DO i=iMin,iMax+1
339     C transport through Western face area:
340     uTrans = (
341     & drF(k-1)*_hFacW(i,j,k-1,bi,bj)*uVel(i,j,k-1,bi,bj)
342     & + drF( k )*_hFacW(i,j, k ,bi,bj)*uVel(i,j, k ,bi,bj)
343 jmc 1.31 & )*halfRL*_dyG(i,j,bi,bj)
344     cOld & )*halfRL
345 jmc 1.30 flx_EW(i,j)=
346     & uTrans*(wVel(i,j,k,bi,bj)+wVel(i-1,j,k,bi,bj))*halfRL
347     ENDDO
348     ENDDO
349     C Advective Flux on Southern face
350     DO j=jMin,jMax+1
351     DO i=iMin,iMax
352     C transport through Southern face area:
353     vTrans = (
354     & drF(k-1)*_hFacS(i,j,k-1,bi,bj)*vVel(i,j,k-1,bi,bj)
355     & +drF( k )*_hFacS(i,j, k ,bi,bj)*vVel(i,j, k ,bi,bj)
356 jmc 1.31 & )*halfRL*_dxG(i,j,bi,bj)
357     cOld & )*halfRL
358 jmc 1.30 flx_NS(i,j)=
359     & vTrans*(wVel(i,j,k,bi,bj)+wVel(i,j-1,k,bi,bj))*halfRL
360     ENDDO
361     ENDDO
362     C Advective Flux on Lower face of W-Cell (= at tracer-cell center, level k)
363     DO j=jMin,jMax
364     DO i=iMin,iMax
365     tmp_WbarZ = halfRL*( wVel(i,j, k ,bi,bj)
366     & +wVel(i,j,kp1,bi,bj)*wOverRide )
367     C transport through Lower face area:
368     rTrans = tmp_WbarZ*rA(i,j,bi,bj)
369 jmc 1.31 flx_Dn(i,j) = rTrans*tmp_WbarZ
370     cOld flx_Dn(i,j) = tmp_WbarZ*tmp_WbarZ
371 jmc 1.30 ENDDO
372     ENDDO
373     C Tendency is minus divergence of advective fluxes:
374     DO j=jMin,jMax
375     DO i=iMin,iMax
376     gW(i,j,k,bi,bj) =
377 jmc 1.31 & -( ( flx_EW(i+1,j)-flx_EW(i,j) )
378     & + ( flx_NS(i,j+1)-flx_NS(i,j) )
379     & + ( flx_Dn(i,j)-flxAdvUp(i,j) )*rkSign
380     & )*recip_rA(i,j,bi,bj)*recip_rThickC(i,j)
381     cOld gW(i,j,k,bi,bj) =
382     cOld & -(
383     cOld & +_recip_dxF(i,j,bi,bj)*( flx_EW(i+1,j)-flx_EW(i,j) )
384     cOld & +_recip_dyF(i,j,bi,bj)*( flx_NS(i,j+1)-flx_NS(i,j) )
385     cOld & + ( flxAdvUp(i,j)-flx_Dn(i,j) )
386     c & )*recip_rThickC(i,j)
387     cOld & )*recip_drC(k)
388 jmc 1.30 C-- prepare for next level (k+1)
389     flxAdvUp(i,j)=flx_Dn(i,j)
390     ENDDO
391     ENDDO
392     ENDIF
393    
394     IF ( useNHMTerms ) THEN
395     CALL MOM_W_METRIC_NH(
396     I bi,bj,k,
397     I uVel, vVel,
398     O gwAdd,
399     I myThid )
400     DO j=jMin,jMax
401     DO i=iMin,iMax
402     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
403     ENDDO
404     ENDDO
405     ENDIF
406 jmc 1.32 IF ( use3dCoriolis ) THEN
407 jmc 1.30 CALL MOM_W_CORIOLIS_NH(
408     I bi,bj,k,
409     I uVel, vVel,
410     O gwAdd,
411     I myThid )
412     DO j=jMin,jMax
413     DO i=iMin,iMax
414     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwAdd(i,j)
415     ENDDO
416     ENDDO
417     ENDIF
418 adcroft 1.1
419 jmc 1.25 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
420    
421 jmc 1.30 C-- Dissipation term inside the Adams-Bashforth:
422     IF ( momViscosity .AND. momDissip_In_AB) THEN
423     DO j=jMin,jMax
424     DO i=iMin,iMax
425     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
426     ENDDO
427     ENDDO
428     ENDIF
429    
430 jmc 1.25 C- Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights)
431     C and save gW_[n] into gwNm1 for the next time step.
432     c#ifdef ALLOW_ADAMSBASHFORTH_3
433 jmc 1.28 c CALL ADAMS_BASHFORTH3(
434     c I bi, bj, k,
435     c U gW, gwNm,
436     c I momStartAB, myIter, myThid )
437 jmc 1.25 c#else /* ALLOW_ADAMSBASHFORTH_3 */
438 jmc 1.28 CALL ADAMS_BASHFORTH2(
439     I bi, bj, k,
440     U gW, gwNm1,
441     I myIter, myThid )
442 jmc 1.25 c#endif /* ALLOW_ADAMSBASHFORTH_3 */
443 jmc 1.21
444 jmc 1.30 C-- Dissipation term outside the Adams-Bashforth:
445     IF ( momViscosity .AND. .NOT.momDissip_In_AB ) THEN
446     DO j=jMin,jMax
447     DO i=iMin,iMax
448     gW(i,j,k,bi,bj) = gW(i,j,k,bi,bj)+gwDiss(i,j)
449     ENDDO
450     ENDDO
451     ENDIF
452    
453 jmc 1.25 C- end of the k loop
454 adcroft 1.4 ENDDO
455    
456 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
457    
458     RETURN
459     END

  ViewVC Help
Powered by ViewVC 1.1.22