/[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.16 - (hide annotations) (download)
Fri Jun 25 13:09:40 2004 UTC (19 years, 11 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint54d_post, checkpoint54e_post, checkpoint55, checkpoint54, checkpoint54f_post, checkpoint55c_post, checkpoint54b_post, checkpoint55d_post, checkpoint54a_pre, checkpoint55d_pre, checkpoint54a_post, checkpoint55b_post, checkpoint53g_post, checkpoint53f_post, checkpoint55a_post, checkpoint55e_post, checkpoint54c_post
Changes since 1.15: +3 -3 lines
 o fix syntax for Intel compiler ("no type" error)

1 edhill 1.16 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.15 2004/05/25 07:48:22 mlosch 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     I 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     C myThid - Instance number for this innvocation of CALC_GW
39     INTEGER myThid
40    
41 adcroft 1.3 #ifdef ALLOW_NONHYDROSTATIC
42    
43 cnh 1.9 C !LOCAL VARIABLES:
44 adcroft 1.1 C == Local variables ==
45 cnh 1.9 C bi, bj, :: Loop counters
46     C iMin, iMax,
47     C jMin, jMax
48     C flx_NS :: Temp. used for fVol meridional terms.
49     C flx_EW :: Temp. used for fVol zonal terms.
50     C flx_Up :: Temp. used for fVol vertical terms.
51     C flx_Dn :: Temp. used for fVol vertical terms.
52 adcroft 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
53     _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
54     _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55     _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
56     _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57     C I,J,K - Loop counters
58 adcroft 1.4 INTEGER i,j,k, kP1, kUp
59 adcroft 1.1 _RL wOverride
60 mlosch 1.15 _RS hFacWtmp
61     _RS hFacStmp
62 adcroft 1.1 _RL ab15,ab05
63 jmc 1.12 _RL slipSideFac
64 adcroft 1.1 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
65    
66     _RL Half
67     PARAMETER(Half=0.5D0)
68 cnh 1.9 CEOP
69 edhill 1.10
70     ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
71 adcroft 1.1
72 mlosch 1.14 iMin = 1
73     iMax = sNx
74     jMin = 1
75     jMax = sNy
76    
77 adcroft 1.1 C Adams-Bashforth timestepping weights
78 jmc 1.11 ab15 = 1.5 _d 0 + abeps
79     ab05 = -0.5 _d 0 - abeps
80    
81     C Lateral friction (no-slip, free slip, or half slip):
82     IF ( no_slip_sides ) THEN
83 mlosch 1.15 slipSideFac = -1. _d 0
84 jmc 1.11 ELSE
85 mlosch 1.15 slipSideFac = 1. _d 0
86 jmc 1.11 ENDIF
87 mlosch 1.15 CML half slip was used before ; keep it for now, but half slip is
88     CML not used anywhere in the code as far as I can see
89 mlosch 1.14 C slipSideFac = 0. _d 0
90 jmc 1.11
91 adcroft 1.1 DO bj=myByLo(myThid),myByHi(myThid)
92     DO bi=myBxLo(myThid),myBxHi(myThid)
93     DO K=1,Nr
94     DO j=1-OLy,sNy+OLy
95     DO i=1-OLx,sNx+OLx
96 jmc 1.8 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
97 adcroft 1.1 gW(i,j,k,bi,bj) = 0.
98     ENDDO
99     ENDDO
100     ENDDO
101     ENDDO
102     ENDDO
103    
104     C Catch barotropic mode
105     IF ( Nr .LT. 2 ) RETURN
106    
107     C For each tile
108     DO bj=myByLo(myThid),myByHi(myThid)
109     DO bi=myBxLo(myThid),myBxHi(myThid)
110    
111     C Boundaries condition at top
112 mlosch 1.14 DO J=jMin,jMax
113     DO I=iMin,iMax
114 adcroft 1.1 Flx_Dn(I,J,bi,bj)=0.
115     ENDDO
116     ENDDO
117    
118     C Sweep down column
119     DO K=2,Nr
120     Kp1=K+1
121     wOverRide=1.
122     if (K.EQ.Nr) then
123     Kp1=Nr
124     wOverRide=0.
125     endif
126     C Flux on Southern face
127 mlosch 1.14 DO J=jMin,jMax+1
128     DO I=iMin,iMax
129 mlosch 1.15 C First compute the fraction of open water for the w-control volume
130     C at the southern face
131 edhill 1.16 hFacStmp=max(hFacS(I,J,K-1,bi,bj)-Half,0. _d 0)
132 mlosch 1.15 & + min(hFacS(I,J,K ,bi,bj),Half)
133 adcroft 1.1 tmp_VbarZ=Half*(
134     & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
135     & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
136     Flx_NS(I,J,bi,bj)=
137     & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
138 jmc 1.11 & -viscAh*_recip_dyC(I,J,bi,bj)
139 mlosch 1.15 & *(hFacStmp*(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
140     & +(1. _d 0 - hFacStmp)*(1. _d 0 - slipSideFac)
141     & *wVel(I,J,K,bi,bj))
142     C The last term is the weighted average of the viscous stress at the open
143     C fraction of the w control volume and at the closed fraction of the
144     C the control volume. A more compact but less intelligible version
145     C of the last three lines is:
146     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacStmp))
147     CML & *wVel(I,J,K,bi,bi) + hFacStmp*wVel(I,J-1,K,bi,bj) )
148 adcroft 1.1 ENDDO
149     ENDDO
150     C Flux on Western face
151 mlosch 1.14 DO J=jMin,jMax
152     DO I=iMin,iMax+1
153 mlosch 1.15 C First compute the fraction of open water for the w-control volume
154     C at the western face
155 edhill 1.16 hFacWtmp=max(hFacW(I,J,K-1,bi,bj)-Half,0. _d 0)
156 mlosch 1.15 & + min(hFacW(I,J,K ,bi,bj),Half)
157     tmp_UbarZ=Half*(
158 adcroft 1.1 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
159     & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
160     Flx_EW(I,J,bi,bj)=
161     & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
162 jmc 1.11 & -viscAh*_recip_dxC(I,J,bi,bj)
163 mlosch 1.15 & *(hFacWtmp*(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
164     & +(1 _d 0 - hFacWtmp)*(1 _d 0 - slipSideFac)
165     & *wVel(I,J,K,bi,bj) )
166     C The last term is the weighted average of the viscous stress at the open
167     C fraction of the w control volume and at the closed fraction of the
168     C the control volume. A more compact but less intelligible version
169     C of the last three lines is:
170     CML & *( (1 _d 0 - slipSideFac*(1 _d 0 - hFacWtmp))
171     CML & *wVel(I,J,K,bi,bi) + hFacWtmp*wVel(I-1,J,K,bi,bj) )
172 adcroft 1.1 ENDDO
173     ENDDO
174     C Flux on Lower face
175 mlosch 1.14 DO J=jMin,jMax
176     DO I=iMin,iMax
177 adcroft 1.1 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
178 mlosch 1.14 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)
179     & +wOverRide*wVel(I,J,Kp1,bi,bj))
180 adcroft 1.1 Flx_Dn(I,J,bi,bj)=
181     & tmp_WbarZ*tmp_WbarZ
182 jmc 1.11 & -viscAr*recip_drF(K)
183     & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
184 adcroft 1.1 ENDDO
185     ENDDO
186     C Divergence of fluxes
187 mlosch 1.14 DO J=jMin,jMax
188     DO I=iMin,iMax
189 adcroft 1.1 gW(I,J,K,bi,bj) = 0.
190     & -(
191     & +_recip_dxF(I,J,bi,bj)*(
192     & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
193     & +_recip_dyF(I,J,bi,bj)*(
194     & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
195     & +recip_drC(K) *(
196     & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
197     & )
198     caja * recip_hFacU(I,J,K,bi,bj)
199     caja NOTE: This should be included
200     caja but we need an hFacUW (above U points)
201     caja and an hFacUS (above V points) too...
202     ENDDO
203     ENDDO
204     ENDDO
205     ENDDO
206     ENDDO
207    
208    
209     DO bj=myByLo(myThid),myByHi(myThid)
210     DO bi=myBxLo(myThid),myBxHi(myThid)
211     DO K=2,Nr
212 mlosch 1.14 DO j=jMin,jMax
213     DO i=iMin,iMax
214 adcroft 1.1 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
215     & +deltatMom*( ab15*gW(i,j,k,bi,bj)
216     & +ab05*gWNM1(i,j,k,bi,bj) )
217     IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
218     ENDDO
219     ENDDO
220 adcroft 1.4 ENDDO
221     ENDDO
222     ENDDO
223    
224 adcroft 1.5 #ifdef ALLOW_OBCS
225     IF (useOBCS) THEN
226 adcroft 1.4 C-- This call is aesthetic: it makes the W field
227     C consistent with the OBs but this has no algorithmic
228     C impact. This is purely for diagnostic purposes.
229 adcroft 1.5 DO bj=myByLo(myThid),myByHi(myThid)
230     DO bi=myBxLo(myThid),myBxHi(myThid)
231     DO K=1,Nr
232     CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
233     ENDDO
234 adcroft 1.1 ENDDO
235     ENDDO
236 adcroft 1.5 ENDIF
237     #endif /* ALLOW_OBCS */
238 adcroft 1.1
239     #endif /* ALLOW_NONHYDROSTATIC */
240    
241     RETURN
242     END

  ViewVC Help
Powered by ViewVC 1.1.22