/[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.27 - (hide annotations) (download)
Tue Jun 20 20:57:37 2006 UTC (17 years, 11 months ago) by baylor
Branch: MAIN
CVS Tags: checkpoint58k_post
Changes since 1.26: +21 -8 lines
Pass the variable viscosities on to calc_gw.

1 baylor 1.27 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.26 2006/06/07 01:55:12 heimbach Exp $
2 cnh 1.9 C !DESCRIPTION: \bv
3 jmc 1.7 C $Name: $
4 edhill 1.10
5 jmc 1.25 c #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 jmc 1.25 SUBROUTINE CALC_GW(
12 baylor 1.27 I KappaRU, KappaRV,
13     I myTime, myIter, myThid )
14 cnh 1.9 C !DESCRIPTION: \bv
15     C *==========================================================*
16 jmc 1.25 C | S/R CALC_GW
17     C | o Calculate vert. velocity tendency terms ( NH, QH only )
18 cnh 1.9 C *==========================================================*
19     C | In NH and QH, 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.22 #include "DYNVARS.h"
34     #include "NH_VARS.h"
35 adcroft 1.1
36 cnh 1.9 C !INPUT/OUTPUT PARAMETERS:
37 adcroft 1.1 C == Routine arguments ==
38 baylor 1.27 C KappaRU:: vertical viscosity at U points
39     C KappaRV:: vertical viscosity at V points
40 jmc 1.20 C myTime :: Current time in simulation
41     C myIter :: Current iteration number in simulation
42     C myThid :: Thread number for this instance of the routine.
43 baylor 1.27 _RL KappaRU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
44     _RL KappaRV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
45 jmc 1.20 _RL myTime
46     INTEGER myIter
47 adcroft 1.1 INTEGER myThid
48    
49 adcroft 1.3 #ifdef ALLOW_NONHYDROSTATIC
50    
51 cnh 1.9 C !LOCAL VARIABLES:
52 adcroft 1.1 C == Local variables ==
53 cnh 1.9 C bi, bj, :: Loop counters
54     C iMin, iMax,
55     C jMin, jMax
56     C flx_NS :: Temp. used for fVol meridional terms.
57     C flx_EW :: Temp. used for fVol zonal terms.
58     C flx_Up :: Temp. used for fVol vertical terms.
59     C flx_Dn :: Temp. used for fVol vertical terms.
60 adcroft 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
61     _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
62     _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
63     _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
64     _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
65 mlosch 1.18 _RL fZon(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
66     _RL fMer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
67     _RL del2w(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
68 jmc 1.21 C i,j,k - Loop counters
69     INTEGER i,j,k, kP1
70 adcroft 1.1 _RL wOverride
71 mlosch 1.15 _RS hFacWtmp
72     _RS hFacStmp
73 mlosch 1.18 _RS hFacCtmp
74     _RS recip_hFacCtmp
75 jmc 1.12 _RL slipSideFac
76 adcroft 1.1 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
77 baylor 1.27 _RL vischere
78     _RL visc4here
79 adcroft 1.1 _RL Half
80     PARAMETER(Half=0.5D0)
81 cnh 1.9 CEOP
82 edhill 1.10
83 jmc 1.25 C Catch barotropic mode
84     IF ( Nr .LT. 2 ) RETURN
85    
86 mlosch 1.14 iMin = 1
87     iMax = sNx
88     jMin = 1
89     jMax = sNy
90    
91 jmc 1.11 C Lateral friction (no-slip, free slip, or half slip):
92     IF ( no_slip_sides ) THEN
93 mlosch 1.15 slipSideFac = -1. _d 0
94 jmc 1.11 ELSE
95 jmc 1.25 slipSideFac = 1. _d 0
96 jmc 1.11 ENDIF
97 mlosch 1.18 CML half slip was used before ; keep the line for now, but half slip is
98     CML not used anywhere in the code as far as I can see.
99 mlosch 1.14 C slipSideFac = 0. _d 0
100 jmc 1.11
101 jmc 1.25 C For each tile
102 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
103     DO bi=myBxLo(myThid),myBxHi(myThid)
104 jmc 1.25
105     C Initialise gW to zero
106 adcroft 1.1 DO K=1,Nr
107     DO j=1-OLy,sNy+OLy
108     DO i=1-OLx,sNx+OLx
109     gW(i,j,k,bi,bj) = 0.
110     ENDDO
111     ENDDO
112     ENDDO
113    
114     C Boundaries condition at top
115 mlosch 1.14 DO J=jMin,jMax
116     DO I=iMin,iMax
117 adcroft 1.1 Flx_Dn(I,J,bi,bj)=0.
118     ENDDO
119     ENDDO
120    
121     C Sweep down column
122     DO K=2,Nr
123     Kp1=K+1
124     wOverRide=1.
125     if (K.EQ.Nr) then
126     Kp1=Nr
127     wOverRide=0.
128     endif
129 mlosch 1.18 C horizontal bi-harmonic dissipation
130     IF (momViscosity .AND. viscA4W.NE.0. ) THEN
131     C calculate the horizontal Laplacian of vertical flow
132     C Zonal flux d/dx W
133     DO j=1-Oly,sNy+Oly
134     fZon(1-Olx,j)=0.
135     DO i=1-Olx+1,sNx+Olx
136     fZon(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
137     & *_dyG(i,j,bi,bj)
138     & *_recip_dxC(i,j,bi,bj)
139     & *(wVel(i,j,k,bi,bj)-wVel(i-1,j,k,bi,bj))
140     #ifdef COSINEMETH_III
141     & *sqcosFacU(J,bi,bj)
142     #endif
143     ENDDO
144     ENDDO
145     C Meridional flux d/dy W
146     DO i=1-Olx,sNx+Olx
147     fMer(I,1-Oly)=0.
148     ENDDO
149     DO j=1-Oly+1,sNy+Oly
150     DO i=1-Olx,sNx+Olx
151     fMer(i,j) = drF(k)*_hFacC(i,j,k,bi,bj)
152     & *_dxG(i,j,bi,bj)
153     & *_recip_dyC(i,j,bi,bj)
154     & *(wVel(i,j,k,bi,bj)-wVel(i,j-1,k,bi,bj))
155     #ifdef ISOTROPIC_COS_SCALING
156     #ifdef COSINEMETH_III
157     & *sqCosFacV(j,bi,bj)
158     #endif
159     #endif
160     ENDDO
161     ENDDO
162 jmc 1.25
163 mlosch 1.18 C del^2 W
164     C Difference of zonal fluxes ...
165     DO j=1-Oly,sNy+Oly
166     DO i=1-Olx,sNx+Olx-1
167     del2w(i,j)=fZon(i+1,j)-fZon(i,j)
168     ENDDO
169     del2w(sNx+Olx,j)=0.
170     ENDDO
171    
172     C ... add difference of meridional fluxes and divide by volume
173     DO j=1-Oly,sNy+Oly-1
174     DO i=1-Olx,sNx+Olx
175     C First compute the fraction of open water for the w-control volume
176     C at the southern face
177 heimbach 1.26 hFacCtmp=max( _hFacC(I,J,K-1,bi,bj)-Half,0. _d 0 )
178     & + min( _hFacC(I,J,K ,bi,bj) ,Half )
179 mlosch 1.18 IF (hFacCtmp .GT. 0.) THEN
180     recip_hFacCtmp = 1./hFacCtmp
181     ELSE
182     recip_hFacCtmp = 0. _d 0
183     ENDIF
184     del2w(i,j)=recip_rA(i,j,bi,bj)
185     & *recip_drC(k)*recip_hFacCtmp
186     & *(
187     & del2w(i,j)
188     & +( fMer(i,j+1)-fMer(i,j) )
189     & )
190     ENDDO
191     ENDDO
192     C-- No-slip BCs impose a drag at walls...
193     CML ************************************************************
194     CML No-slip Boundary conditions for bi-harmonic dissipation
195     CML need to be implemented here!
196     CML ************************************************************
197 jmc 1.21 ELSE
198     C- Initialize del2w to zero:
199     DO j=1-Oly,sNy+Oly
200     DO i=1-Olx,sNx+Olx
201     del2w(i,j) = 0. _d 0
202     ENDDO
203     ENDDO
204 mlosch 1.18 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 heimbach 1.26 hFacStmp=max(_hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
212     & + 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 baylor 1.27 & -viscAh_W(I,J,K,bi,bj)*_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 baylor 1.27 & +viscA4_W(I,J,K,bi,bj)
223     & *_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
224 mlosch 1.18 #ifdef ISOTROPIC_COS_SCALING
225     #ifdef COSINEMETH_III
226     & *sqCosFacV(j,bi,bj)
227     #else
228     & *CosFacV(j,bi,bj)
229     #endif
230     #endif
231 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
232     C fraction of the w control volume and at the closed fraction of the
233     C the control volume. A more compact but less intelligible version
234     C of the last three lines is:
235     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
236     CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
237 adcroft 1.1 ENDDO
238     ENDDO
239     C Flux on Western face
240 mlosch 1.14 DO J=jMin,jMax
241     DO I=iMin,iMax+1
242 mlosch 1.15 C First compute the fraction of open water for the w-control volume
243     C at the western face
244 heimbach 1.26 hFacWtmp=max(_hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
245     & + min(_hFacW(I,J,K ,bi,bj),Half)
246 jmc 1.21 tmp_UbarZ=Half*(
247 adcroft 1.1 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
248     & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
249     Flx_EW(I,J,bi,bj)=
250     & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
251 baylor 1.27 & -viscAh_W(I,J,K,bi,bj)*_recip_dxC(I,J,bi,bj)
252 mlosch 1.15 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
253     & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
254     & *wVel(I,J,K,bi,bj) )
255 baylor 1.27 & +viscA4_W(I,J,K,bi,bj)
256     & *_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
257 mlosch 1.18 #ifdef COSINEMETH_III
258     & *sqCosFacU(j,bi,bj)
259 jmc 1.25 #else
260 mlosch 1.18 & *CosFacU(j,bi,bj)
261     #endif
262 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
263     C fraction of the w control volume and at the closed fraction of the
264     C the control volume. A more compact but less intelligible version
265     C of the last three lines is:
266     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
267     CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
268 adcroft 1.1 ENDDO
269     ENDDO
270     C Flux on Lower face
271 mlosch 1.14 DO J=jMin,jMax
272     DO I=iMin,iMax
273 baylor 1.27 C Interpolate vert viscosity to W points
274     vischere=0.125*(kappaRU(I,J,K) +kappaRU(I+1,J,K)
275     & +kappaRU(I,J,Kp1)+kappaRU(I+1,J,Kp1)
276     & +kappaRV(I,J,K) +kappaRV(I,J+1,K)
277     & +kappaRV(I,J,Kp1)+kappaRV(I,J+1,Kp1))
278 adcroft 1.1 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
279 jmc 1.25 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
280 mlosch 1.14 & +wOverRide*wVel(I,J,Kp1,bi,bj))
281 adcroft 1.1 Flx_Dn(I,J,bi,bj)=
282     & tmp_WbarZ*tmp_WbarZ
283 baylor 1.27 & -vischere*recip_drF(K)
284 jmc 1.11 & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
285 adcroft 1.1 ENDDO
286     ENDDO
287     C Divergence of fluxes
288 mlosch 1.14 DO J=jMin,jMax
289     DO I=iMin,iMax
290 adcroft 1.1 gW(I,J,K,bi,bj) = 0.
291     & -(
292     & +_recip_dxF(I,J,bi,bj)*(
293     & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
294     & +_recip_dyF(I,J,bi,bj)*(
295     & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
296     & +recip_drC(K) *(
297     & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
298     & )
299     caja * recip_hFacU(I,J,K,bi,bj)
300 jmc 1.25 caja NOTE: This should be included
301 adcroft 1.1 caja but we need an hFacUW (above U points)
302     caja and an hFacUS (above V points) too...
303     ENDDO
304     ENDDO
305    
306 jmc 1.25 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
307    
308     C- Compute effective gW_[n+1/2] terms (including Adams-Bashforth weights)
309     C and save gW_[n] into gwNm1 for the next time step.
310     c#ifdef ALLOW_ADAMSBASHFORTH_3
311     c CALL ADAMS_BASHFORTH3(
312     c I bi, bj, k,
313     c U gW, gwNm,
314     c I momStartAB, myIter, myThid )
315     c#else /* ALLOW_ADAMSBASHFORTH_3 */
316     CALL ADAMS_BASHFORTH2(
317     I bi, bj, k,
318     U gW, gwNm1,
319     I myIter, myThid )
320     c#endif /* ALLOW_ADAMSBASHFORTH_3 */
321 jmc 1.21
322 jmc 1.25 C- end of the k loop
323 adcroft 1.4 ENDDO
324 jmc 1.25
325     C- end of bi,bj loops
326 adcroft 1.4 ENDDO
327     ENDDO
328    
329 adcroft 1.1 #endif /* ALLOW_NONHYDROSTATIC */
330    
331     RETURN
332     END

  ViewVC Help
Powered by ViewVC 1.1.22