/[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.20 - (hide annotations) (download)
Sat Jul 30 23:38:54 2005 UTC (18 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57p_post
Changes since 1.19: +7 -3 lines
add arguments: myTime & myIter.

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

  ViewVC Help
Powered by ViewVC 1.1.22