/[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.22 - (hide annotations) (download)
Tue Nov 8 02:14:10 2005 UTC (18 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57y_post, checkpoint57y_pre, checkpoint57x_post
Changes since 1.21: +3 -6 lines
put all NH variables (formely in DYNVARS.h & GW.h) in NH_VARS.h

1 jmc 1.22 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.21 2005/08/16 22:42:49 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.11 ab15 = 1.5 _d 0 + abeps
85     ab05 = -0.5 _d 0 - abeps
86    
87     C Lateral friction (no-slip, free slip, or half slip):
88     IF ( no_slip_sides ) THEN
89 mlosch 1.15 slipSideFac = -1. _d 0
90 jmc 1.11 ELSE
91 mlosch 1.15 slipSideFac = 1. _d 0
92 jmc 1.11 ENDIF
93 mlosch 1.18 CML half slip was used before ; keep the line for now, but half slip is
94     CML not used anywhere in the code as far as I can see.
95 mlosch 1.14 C slipSideFac = 0. _d 0
96 jmc 1.11
97 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
98     DO bi=myBxLo(myThid),myBxHi(myThid)
99     DO K=1,Nr
100     DO j=1-OLy,sNy+OLy
101     DO i=1-OLx,sNx+OLx
102 jmc 1.21 gwNm1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
103 adcroft 1.1 gW(i,j,k,bi,bj) = 0.
104     ENDDO
105     ENDDO
106     ENDDO
107     ENDDO
108     ENDDO
109    
110     C Catch barotropic mode
111     IF ( Nr .LT. 2 ) RETURN
112    
113     C For each tile
114     DO bj=myByLo(myThid),myByHi(myThid)
115     DO bi=myBxLo(myThid),myBxHi(myThid)
116    
117     C Boundaries condition at top
118 mlosch 1.14 DO J=jMin,jMax
119     DO I=iMin,iMax
120 adcroft 1.1 Flx_Dn(I,J,bi,bj)=0.
121     ENDDO
122     ENDDO
123    
124     C Sweep down column
125     DO K=2,Nr
126     Kp1=K+1
127     wOverRide=1.
128     if (K.EQ.Nr) then
129     Kp1=Nr
130     wOverRide=0.
131     endif
132 mlosch 1.18 C horizontal bi-harmonic dissipation
133     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
134     C calculate the horizontal Laplacian of vertical flow
135     C Zonal flux d/dx W
136     DO j=1-Oly,sNy+Oly
137     fZon(1-Olx,j)=0.
138     DO i=1-Olx+1,sNx+Olx
139     fZon(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
140     & *_dyG(i,j,bi,bj)
141     & *_recip_dxC(i,j,bi,bj)
142     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
143     #ifdef COSINEMETH_III
144     & *sqcosFacU(J,bi,bj)
145     #endif
146     ENDDO
147     ENDDO
148     C Meridional flux d/dy W
149     DO i=1-Olx,sNx+Olx
150     fMer(I,1-Oly)=0.
151     ENDDO
152     DO j=1-Oly+1,sNy+Oly
153     DO i=1-Olx,sNx+Olx
154     fMer(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
155     & *_dxG(i,j,bi,bj)
156     & *_recip_dyC(i,j,bi,bj)
157     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
158     #ifdef ISOTROPIC_COS_SCALING
159     #ifdef COSINEMETH_III
160     & *sqCosFacV(j,bi,bj)
161     #endif
162     #endif
163     ENDDO
164     ENDDO
165    
166     C del^2 W
167     C Difference of zonal fluxes ...
168     DO j=1-Oly,sNy+Oly
169     DO i=1-Olx,sNx+Olx-1
170     del2w(i,j)=fZon(i+1,j)-fZon(i,j)
171     ENDDO
172     del2w(sNx+Olx,j)=0.
173     ENDDO
174    
175     C ... add difference of meridional fluxes and divide by volume
176     DO j=1-Oly,sNy+Oly-1
177     DO i=1-Olx,sNx+Olx
178     C First compute the fraction of open water for the w-control volume
179     C at the southern face
180     hFacCtmp=max(hFacC(I,J,K-1,bi,bj)-Half,0. _d 0)
181     & + min(hFacC(I,J,K ,bi,bj),Half)
182     IF (hFacCtmp .GT. 0.) THEN
183     recip_hFacCtmp = 1./hFacCtmp
184     ELSE
185     recip_hFacCtmp = 0. _d 0
186     ENDIF
187     del2w(i,j)=recip_rA(i,j,bi,bj)
188     & *recip_drC(k)*recip_hFacCtmp
189     & *(
190     & del2w(i,j)
191     & +( fMer(i,j+1)-fMer(i,j) )
192     & )
193     ENDDO
194     ENDDO
195     C-- No-slip BCs impose a drag at walls...
196     CML ************************************************************
197     CML No-slip Boundary conditions for bi-harmonic dissipation
198     CML need to be implemented here!
199     CML ************************************************************
200 jmc 1.21 ELSE
201     C- Initialize del2w to zero:
202     DO j=1-Oly,sNy+Oly
203     DO i=1-Olx,sNx+Olx
204     del2w(i,j) = 0. _d 0
205     ENDDO
206     ENDDO
207 mlosch 1.18 ENDIF
208    
209 adcroft 1.1 C Flux on Southern face
210 mlosch 1.14 DO J=jMin,jMax+1
211     DO I=iMin,iMax
212 mlosch 1.15 C First compute the fraction of open water for the w-control volume
213     C at the southern face
214 edhill 1.16 hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
215 mlosch 1.15 & + min(hFacS(I,J,K ,bi,bj),Half)
216 adcroft 1.1 tmp_VbarZ=Half*(
217     & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
218     & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
219     Flx_NS(I,J,bi,bj)=
220     & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
221 mlosch 1.17 & -viscAhW*_recip_dyC(I,J,bi,bj)
222 mlosch 1.15 & *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
223     & +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
224     & *wVel(I,J,K,bi,bj))
225 mlosch 1.18 & +viscA4W*_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
226     #ifdef ISOTROPIC_COS_SCALING
227     #ifdef COSINEMETH_III
228     & *sqCosFacV(j,bi,bj)
229     #else
230     & *CosFacV(j,bi,bj)
231     #endif
232     #endif
233 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
234     C fraction of the w control volume and at the closed fraction of the
235     C the control volume. A more compact but less intelligible version
236     C of the last three lines is:
237     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
238     CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
239 adcroft 1.1 ENDDO
240     ENDDO
241     C Flux on Western face
242 mlosch 1.14 DO J=jMin,jMax
243     DO I=iMin,iMax+1
244 mlosch 1.15 C First compute the fraction of open water for the w-control volume
245     C at the western face
246 edhill 1.16 hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
247 mlosch 1.15 & + min(hFacW(I,J,K ,bi,bj),Half)
248 jmc 1.21 tmp_UbarZ=Half*(
249 adcroft 1.1 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
250     & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
251     Flx_EW(I,J,bi,bj)=
252     & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
253 mlosch 1.17 & -viscAhW*_recip_dxC(I,J,bi,bj)
254 mlosch 1.15 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
255     & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
256     & *wVel(I,J,K,bi,bj) )
257 mlosch 1.18 & +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
258     #ifdef COSINEMETH_III
259     & *sqCosFacU(j,bi,bj)
260     #else
261     & *CosFacU(j,bi,bj)
262     #endif
263 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
264     C fraction of the w control volume and at the closed fraction of the
265     C the control volume. A more compact but less intelligible version
266     C of the last three lines is:
267     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
268     CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
269 adcroft 1.1 ENDDO
270     ENDDO
271     C Flux on Lower face
272 mlosch 1.14 DO J=jMin,jMax
273     DO I=iMin,iMax
274 adcroft 1.1 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
275 mlosch 1.14 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
276     & +wOverRide*wVel(I,J,Kp1,bi,bj))
277 adcroft 1.1 Flx_Dn(I,J,bi,bj)=
278     & tmp_WbarZ*tmp_WbarZ
279 jmc 1.11 & -viscAr*recip_drF(K)
280     & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
281 adcroft 1.1 ENDDO
282     ENDDO
283     C Divergence of fluxes
284 mlosch 1.14 DO J=jMin,jMax
285     DO I=iMin,iMax
286 adcroft 1.1 gW(I,J,K,bi,bj) = 0.
287     & -(
288     & +_recip_dxF(I,J,bi,bj)*(
289     & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
290     & +_recip_dyF(I,J,bi,bj)*(
291     & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
292     & +recip_drC(K) *(
293     & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
294     & )
295     caja * recip_hFacU(I,J,K,bi,bj)
296     caja NOTE: This should be included
297     caja but we need an hFacUW (above U points)
298     caja and an hFacUS (above V points) too...
299     ENDDO
300     ENDDO
301     ENDDO
302     ENDDO
303     ENDDO
304    
305 jmc 1.21
306 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
307     DO bi=myBxLo(myThid),myBxHi(myThid)
308     DO K=2,Nr
309 mlosch 1.14 DO j=jMin,jMax
310     DO i=iMin,iMax
311 adcroft 1.1 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
312 adcroft 1.19 & +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)
313 jmc 1.21 & +ab05*gwNm1(i,j,k,bi,bj) )
314 adcroft 1.1 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
315     ENDDO
316     ENDDO
317 adcroft 1.4 ENDDO
318     ENDDO
319     ENDDO
320    
321 adcroft 1.5 #ifdef ALLOW_OBCS
322     IF (useOBCS) THEN
323 adcroft 1.4 C-- This call is aesthetic: it makes the W field
324     C consistent with the OBs but this has no algorithmic
325     C impact. This is purely for diagnostic purposes.
326 adcroft 1.5 DO bj=myByLo(myThid),myByHi(myThid)
327     DO bi=myBxLo(myThid),myBxHi(myThid)
328     DO K=1,Nr
329     CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
330     ENDDO
331 adcroft 1.1 ENDDO
332     ENDDO
333 adcroft 1.5 ENDIF
334     #endif /* ALLOW_OBCS */
335 adcroft 1.1
336     #endif /* ALLOW_NONHYDROSTATIC */
337    
338     RETURN
339     END

  ViewVC Help
Powered by ViewVC 1.1.22