/[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.23 - (hide annotations) (download)
Tue Dec 13 23:16:52 2005 UTC (18 years, 5 months ago) by jmc
Branch: MAIN
Changes since 1.22: +8 -3 lines
no AB-2 at the 1rst iteration for Gw (like what is done for Gu,Gv,Gt,Gs)

1 jmc 1.23 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.22 2005/11/08 02:14:10 jmc Exp $
2 cnh 1.9 C !DESCRIPTION: \bv
3 jmc 1.7 C $Name: $
4 edhill 1.10
5     #include "PACKAGES_CONFIG.h"
6 adcroft 1.1 #include "CPP_OPTIONS.h"
7    
8 cnh 1.9 CBOP
9     C !ROUTINE: CALC_GW
10     C !INTERFACE:
11 adcroft 1.1 SUBROUTINE CALC_GW(
12 jmc 1.20 I myTime, myIter, myThid )
13 cnh 1.9 C !DESCRIPTION: \bv
14     C *==========================================================*
15     C | S/R CALC_GW
16     C | o Calculate vert. velocity tendency terms ( NH, QH only )
17     C *==========================================================*
18     C | In NH and QH, the vertical momentum tendency must be
19     C | calculated explicitly and included as a source term
20     C | for a 3d pressure eqn. Calculate that term here.
21     C | This routine is not used in HYD calculations.
22     C *==========================================================*
23     C \ev
24    
25     C !USES:
26 adcroft 1.1 IMPLICIT NONE
27     C == Global variables ==
28     #include "SIZE.h"
29     #include "EEPARAMS.h"
30     #include "PARAMS.h"
31     #include "GRID.h"
32 jmc 1.22 #include "DYNVARS.h"
33     #include "NH_VARS.h"
34 adcroft 1.1
35 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
36 adcroft 1.1 C == Routine arguments ==
37 jmc 1.20 C myTime :: Current time in simulation
38     C myIter :: Current iteration number in simulation
39     C myThid :: Thread number for this instance of the routine.
40     _RL myTime
41     INTEGER myIter
42 adcroft 1.1 INTEGER myThid
43    
44 adcroft 1.3 #ifdef ALLOW_NONHYDROSTATIC
45    
46 cnh 1.9 C !LOCAL VARIABLES:
47 adcroft 1.1 C == Local variables ==
48 cnh 1.9 C bi, bj, :: Loop counters
49     C iMin, iMax,
50     C jMin, jMax
51     C flx_NS :: Temp. used for fVol meridional terms.
52     C flx_EW :: Temp. used for fVol zonal terms.
53     C flx_Up :: Temp. used for fVol vertical terms.
54     C flx_Dn :: Temp. used for fVol vertical terms.
55 adcroft 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
56     _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57     _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58     _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
59     _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
60 mlosch 1.18 _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61     _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62     _RL del2w(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63 jmc 1.21 C i,j,k - Loop counters
64     INTEGER i,j,k, kP1
65 adcroft 1.1 _RL wOverride
66 mlosch 1.15 _RS hFacWtmp
67     _RS hFacStmp
68 mlosch 1.18 _RS hFacCtmp
69     _RS recip_hFacCtmp
70 adcroft 1.1 _RL ab15,ab05
71 jmc 1.12 _RL slipSideFac
72 adcroft 1.1 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
73    
74     _RL Half
75     PARAMETER(Half=0.5D0)
76 cnh 1.9 CEOP
77 edhill 1.10
78 mlosch 1.14 iMin = 1
79     iMax = sNx
80     jMin = 1
81     jMax = sNy
82    
83 adcroft 1.1 C Adams-Bashforth timestepping weights
84 jmc 1.23 IF (myIter .EQ. 0) THEN
85     ab15 = 1.0 _d 0
86     ab05 = 0.0 _d 0
87     ELSE
88     ab15 = 1.5 _d 0 + abeps
89     ab05 = -0.5 _d 0 - abeps
90     ENDIF
91 jmc 1.11
92     C Lateral friction (no-slip, free slip, or half slip):
93     IF ( no_slip_sides ) THEN
94 mlosch 1.15 slipSideFac = -1. _d 0
95 jmc 1.11 ELSE
96 mlosch 1.15 slipSideFac = 1. _d 0
97 jmc 1.11 ENDIF
98 mlosch 1.18 CML half slip was used before ; keep the line for now, but half slip is
99     CML not used anywhere in the code as far as I can see.
100 mlosch 1.14 C slipSideFac = 0. _d 0
101 jmc 1.11
102 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
103     DO bi=myBxLo(myThid),myBxHi(myThid)
104     DO K=1,Nr
105     DO j=1-OLy,sNy+OLy
106     DO i=1-OLx,sNx+OLx
107 jmc 1.21 gwNm1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
108 adcroft 1.1 gW(i,j,k,bi,bj) = 0.
109     ENDDO
110     ENDDO
111     ENDDO
112     ENDDO
113     ENDDO
114    
115     C Catch barotropic mode
116     IF ( Nr .LT. 2 ) RETURN
117    
118     C For each tile
119     DO bj=myByLo(myThid),myByHi(myThid)
120     DO bi=myBxLo(myThid),myBxHi(myThid)
121    
122     C Boundaries condition at top
123 mlosch 1.14 DO J=jMin,jMax
124     DO I=iMin,iMax
125 adcroft 1.1 Flx_Dn(I,J,bi,bj)=0.
126     ENDDO
127     ENDDO
128    
129     C Sweep down column
130     DO K=2,Nr
131     Kp1=K+1
132     wOverRide=1.
133     if (K.EQ.Nr) then
134     Kp1=Nr
135     wOverRide=0.
136     endif
137 mlosch 1.18 C horizontal bi-harmonic dissipation
138     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
139     C calculate the horizontal Laplacian of vertical flow
140     C Zonal flux d/dx W
141     DO j=1-Oly,sNy+Oly
142     fZon(1-Olx,j)=0.
143     DO i=1-Olx+1,sNx+Olx
144     fZon(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
145     & *_dyG(i,j,bi,bj)
146     & *_recip_dxC(i,j,bi,bj)
147     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
148     #ifdef COSINEMETH_III
149     & *sqcosFacU(J,bi,bj)
150     #endif
151     ENDDO
152     ENDDO
153     C Meridional flux d/dy W
154     DO i=1-Olx,sNx+Olx
155     fMer(I,1-Oly)=0.
156     ENDDO
157     DO j=1-Oly+1,sNy+Oly
158     DO i=1-Olx,sNx+Olx
159     fMer(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
160     & *_dxG(i,j,bi,bj)
161     & *_recip_dyC(i,j,bi,bj)
162     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
163     #ifdef ISOTROPIC_COS_SCALING
164     #ifdef COSINEMETH_III
165     & *sqCosFacV(j,bi,bj)
166     #endif
167     #endif
168     ENDDO
169     ENDDO
170    
171     C del^2 W
172     C Difference of zonal fluxes ...
173     DO j=1-Oly,sNy+Oly
174     DO i=1-Olx,sNx+Olx-1
175     del2w(i,j)=fZon(i+1,j)-fZon(i,j)
176     ENDDO
177     del2w(sNx+Olx,j)=0.
178     ENDDO
179    
180     C ... add difference of meridional fluxes and divide by volume
181     DO j=1-Oly,sNy+Oly-1
182     DO i=1-Olx,sNx+Olx
183     C First compute the fraction of open water for the w-control volume
184     C at the southern face
185     hFacCtmp=max(hFacC(I,J,K-1,bi,bj)-Half,0. _d 0)
186     & + min(hFacC(I,J,K ,bi,bj),Half)
187     IF (hFacCtmp .GT. 0.) THEN
188     recip_hFacCtmp = 1./hFacCtmp
189     ELSE
190     recip_hFacCtmp = 0. _d 0
191     ENDIF
192     del2w(i,j)=recip_rA(i,j,bi,bj)
193     & *recip_drC(k)*recip_hFacCtmp
194     & *(
195     & del2w(i,j)
196     & +( fMer(i,j+1)-fMer(i,j) )
197     & )
198     ENDDO
199     ENDDO
200     C-- No-slip BCs impose a drag at walls...
201     CML ************************************************************
202     CML No-slip Boundary conditions for bi-harmonic dissipation
203     CML need to be implemented here!
204     CML ************************************************************
205 jmc 1.21 ELSE
206     C- Initialize del2w to zero:
207     DO j=1-Oly,sNy+Oly
208     DO i=1-Olx,sNx+Olx
209     del2w(i,j) = 0. _d 0
210     ENDDO
211     ENDDO
212 mlosch 1.18 ENDIF
213    
214 adcroft 1.1 C Flux on Southern face
215 mlosch 1.14 DO J=jMin,jMax+1
216     DO I=iMin,iMax
217 mlosch 1.15 C First compute the fraction of open water for the w-control volume
218     C at the southern face
219 edhill 1.16 hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
220 mlosch 1.15 & + min(hFacS(I,J,K ,bi,bj),Half)
221 adcroft 1.1 tmp_VbarZ=Half*(
222     & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
223     & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
224     Flx_NS(I,J,bi,bj)=
225     & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
226 mlosch 1.17 & -viscAhW*_recip_dyC(I,J,bi,bj)
227 mlosch 1.15 & *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
228     & +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
229     & *wVel(I,J,K,bi,bj))
230 mlosch 1.18 & +viscA4W*_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
231     #ifdef ISOTROPIC_COS_SCALING
232     #ifdef COSINEMETH_III
233     & *sqCosFacV(j,bi,bj)
234     #else
235     & *CosFacV(j,bi,bj)
236     #endif
237     #endif
238 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
239     C fraction of the w control volume and at the closed fraction of the
240     C the control volume. A more compact but less intelligible version
241     C of the last three lines is:
242     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
243     CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
244 adcroft 1.1 ENDDO
245     ENDDO
246     C Flux on Western face
247 mlosch 1.14 DO J=jMin,jMax
248     DO I=iMin,iMax+1
249 mlosch 1.15 C First compute the fraction of open water for the w-control volume
250     C at the western face
251 edhill 1.16 hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
252 mlosch 1.15 & + min(hFacW(I,J,K ,bi,bj),Half)
253 jmc 1.21 tmp_UbarZ=Half*(
254 adcroft 1.1 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
255     & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
256     Flx_EW(I,J,bi,bj)=
257     & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
258 mlosch 1.17 & -viscAhW*_recip_dxC(I,J,bi,bj)
259 mlosch 1.15 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
260     & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
261     & *wVel(I,J,K,bi,bj) )
262 mlosch 1.18 & +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
263     #ifdef COSINEMETH_III
264     & *sqCosFacU(j,bi,bj)
265     #else
266     & *CosFacU(j,bi,bj)
267     #endif
268 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
269     C fraction of the w control volume and at the closed fraction of the
270     C the control volume. A more compact but less intelligible version
271     C of the last three lines is:
272     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
273     CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
274 adcroft 1.1 ENDDO
275     ENDDO
276     C Flux on Lower face
277 mlosch 1.14 DO J=jMin,jMax
278     DO I=iMin,iMax
279 adcroft 1.1 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
280 mlosch 1.14 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
281     & +wOverRide*wVel(I,J,Kp1,bi,bj))
282 adcroft 1.1 Flx_Dn(I,J,bi,bj)=
283     & tmp_WbarZ*tmp_WbarZ
284 jmc 1.11 & -viscAr*recip_drF(K)
285     & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
286 adcroft 1.1 ENDDO
287     ENDDO
288     C Divergence of fluxes
289 mlosch 1.14 DO J=jMin,jMax
290     DO I=iMin,iMax
291 adcroft 1.1 gW(I,J,K,bi,bj) = 0.
292     & -(
293     & +_recip_dxF(I,J,bi,bj)*(
294     & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
295     & +_recip_dyF(I,J,bi,bj)*(
296     & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
297     & +recip_drC(K) *(
298     & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
299     & )
300     caja * recip_hFacU(I,J,K,bi,bj)
301     caja NOTE: This should be included
302     caja but we need an hFacUW (above U points)
303     caja and an hFacUS (above V points) too...
304     ENDDO
305     ENDDO
306     ENDDO
307     ENDDO
308     ENDDO
309    
310 jmc 1.21
311 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
312     DO bi=myBxLo(myThid),myBxHi(myThid)
313     DO K=2,Nr
314 mlosch 1.14 DO j=jMin,jMax
315     DO i=iMin,iMax
316 adcroft 1.1 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
317 adcroft 1.19 & +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)
318 jmc 1.21 & +ab05*gwNm1(i,j,k,bi,bj) )
319 adcroft 1.1 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
320     ENDDO
321     ENDDO
322 adcroft 1.4 ENDDO
323     ENDDO
324     ENDDO
325    
326 adcroft 1.5 #ifdef ALLOW_OBCS
327     IF (useOBCS) THEN
328 adcroft 1.4 C-- This call is aesthetic: it makes the W field
329     C consistent with the OBs but this has no algorithmic
330     C impact. This is purely for diagnostic purposes.
331 adcroft 1.5 DO bj=myByLo(myThid),myByHi(myThid)
332     DO bi=myBxLo(myThid),myBxHi(myThid)
333     DO K=1,Nr
334     CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
335     ENDDO
336 adcroft 1.1 ENDDO
337     ENDDO
338 adcroft 1.5 ENDIF
339     #endif /* ALLOW_OBCS */
340 adcroft 1.1
341     #endif /* ALLOW_NONHYDROSTATIC */
342    
343     RETURN
344     END

  ViewVC Help
Powered by ViewVC 1.1.22