/[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.21 - (hide annotations) (download)
Tue Aug 16 22:42:49 2005 UTC (18 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57s_post, checkpoint57r_post, checkpoint57t_post, checkpoint57v_post, checkpint57u_post, checkpoint57q_post, checkpoint57w_post
Changes since 1.20: +14 -8 lines
fix initialisation Pb (del2w) responsible for getting NANs with g77 on faulks

1 jmc 1.21 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.20 2005/07/30 23:38:54 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 "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 jmc 1.21 C i,j,k - Loop counters
65     INTEGER i,j,k, kP1
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.21 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     IF (hFacCtmp .GT. 0.) THEN
186     recip_hFacCtmp = 1./hFacCtmp
187     ELSE
188     recip_hFacCtmp = 0. _d 0
189     ENDIF
190     del2w(i,j)=recip_rA(i,j,bi,bj)
191     & *recip_drC(k)*recip_hFacCtmp
192     & *(
193     & del2w(i,j)
194     & +( fMer(i,j+1)-fMer(i,j) )
195     & )
196     ENDDO
197     ENDDO
198     C-- No-slip BCs impose a drag at walls...
199     CML ************************************************************
200     CML No-slip Boundary conditions for bi-harmonic dissipation
201     CML need to be implemented here!
202     CML ************************************************************
203 jmc 1.21 ELSE
204     C- Initialize del2w to zero:
205     DO j=1-Oly,sNy+Oly
206     DO i=1-Olx,sNx+Olx
207     del2w(i,j) = 0. _d 0
208     ENDDO
209     ENDDO
210 mlosch 1.18 ENDIF
211    
212 adcroft 1.1 C Flux on Southern face
213 mlosch 1.14 DO J=jMin,jMax+1
214     DO I=iMin,iMax
215 mlosch 1.15 C First compute the fraction of open water for the w-control volume
216     C at the southern face
217 edhill 1.16 hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
218 mlosch 1.15 & + min(hFacS(I,J,K ,bi,bj),Half)
219 adcroft 1.1 tmp_VbarZ=Half*(
220     & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
221     & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
222     Flx_NS(I,J,bi,bj)=
223     & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
224 mlosch 1.17 & -viscAhW*_recip_dyC(I,J,bi,bj)
225 mlosch 1.15 & *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
226     & +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
227     & *wVel(I,J,K,bi,bj))
228 mlosch 1.18 & +viscA4W*_recip_dyC(I,J,bi,bj)*(del2w(I,J)-del2w(I,J-1))
229     #ifdef ISOTROPIC_COS_SCALING
230     #ifdef COSINEMETH_III
231     & *sqCosFacV(j,bi,bj)
232     #else
233     & *CosFacV(j,bi,bj)
234     #endif
235     #endif
236 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
237     C fraction of the w control volume and at the closed fraction of the
238     C the control volume. A more compact but less intelligible version
239     C of the last three lines is:
240     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
241     CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
242 adcroft 1.1 ENDDO
243     ENDDO
244     C Flux on Western face
245 mlosch 1.14 DO J=jMin,jMax
246     DO I=iMin,iMax+1
247 mlosch 1.15 C First compute the fraction of open water for the w-control volume
248     C at the western face
249 edhill 1.16 hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
250 mlosch 1.15 & + min(hFacW(I,J,K ,bi,bj),Half)
251 jmc 1.21 tmp_UbarZ=Half*(
252 adcroft 1.1 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
253     & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
254     Flx_EW(I,J,bi,bj)=
255     & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
256 mlosch 1.17 & -viscAhW*_recip_dxC(I,J,bi,bj)
257 mlosch 1.15 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
258     & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
259     & *wVel(I,J,K,bi,bj) )
260 mlosch 1.18 & +viscA4W*_recip_dxC(I,J,bi,bj)*(del2w(I,J)-del2w(I-1,J))
261     #ifdef COSINEMETH_III
262     & *sqCosFacU(j,bi,bj)
263     #else
264     & *CosFacU(j,bi,bj)
265     #endif
266 mlosch 1.15 C The last term is the weighted average of the viscous stress at the open
267     C fraction of the w control volume and at the closed fraction of the
268     C the control volume. A more compact but less intelligible version
269     C of the last three lines is:
270     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
271     CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
272 adcroft 1.1 ENDDO
273     ENDDO
274     C Flux on Lower face
275 mlosch 1.14 DO J=jMin,jMax
276     DO I=iMin,iMax
277 adcroft 1.1 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
278 mlosch 1.14 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
279     & +wOverRide*wVel(I,J,Kp1,bi,bj))
280 adcroft 1.1 Flx_Dn(I,J,bi,bj)=
281     & tmp_WbarZ*tmp_WbarZ
282 jmc 1.11 & -viscAr*recip_drF(K)
283     & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
284 adcroft 1.1 ENDDO
285     ENDDO
286     C Divergence of fluxes
287 mlosch 1.14 DO J=jMin,jMax
288     DO I=iMin,iMax
289 adcroft 1.1 gW(I,J,K,bi,bj) = 0.
290     & -(
291     & +_recip_dxF(I,J,bi,bj)*(
292     & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
293     & +_recip_dyF(I,J,bi,bj)*(
294     & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
295     & +recip_drC(K) *(
296     & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
297     & )
298     caja * recip_hFacU(I,J,K,bi,bj)
299     caja NOTE: This should be included
300     caja but we need an hFacUW (above U points)
301     caja and an hFacUS (above V points) too...
302     ENDDO
303     ENDDO
304     ENDDO
305     ENDDO
306     ENDDO
307    
308 jmc 1.21
309 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
310     DO bi=myBxLo(myThid),myBxHi(myThid)
311     DO K=2,Nr
312 mlosch 1.14 DO j=jMin,jMax
313     DO i=iMin,iMax
314 adcroft 1.1 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
315 adcroft 1.19 & +deltatMom*nh_Am2*( ab15*gW(i,j,k,bi,bj)
316 jmc 1.21 & +ab05*gwNm1(i,j,k,bi,bj) )
317 adcroft 1.1 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
318     ENDDO
319     ENDDO
320 adcroft 1.4 ENDDO
321     ENDDO
322     ENDDO
323    
324 adcroft 1.5 #ifdef ALLOW_OBCS
325     IF (useOBCS) THEN
326 adcroft 1.4 C-- This call is aesthetic: it makes the W field
327     C consistent with the OBs but this has no algorithmic
328     C impact. This is purely for diagnostic purposes.
329 adcroft 1.5 DO bj=myByLo(myThid),myByHi(myThid)
330     DO bi=myBxLo(myThid),myBxHi(myThid)
331     DO K=1,Nr
332     CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
333     ENDDO
334 adcroft 1.1 ENDDO
335     ENDDO
336 adcroft 1.5 ENDIF
337     #endif /* ALLOW_OBCS */
338 adcroft 1.1
339     #endif /* ALLOW_NONHYDROSTATIC */
340    
341     RETURN
342     END

  ViewVC Help
Powered by ViewVC 1.1.22