/[MITgcm]/MITgcm/model/src/calc_gw.F
ViewVC logotype

Contents of /MITgcm/model/src/calc_gw.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.11 - (show annotations) (download)
Mon Mar 29 17:58:40 2004 UTC (20 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.10: +34 -10 lines
modified to test free-slip / no-slip B.C. (Matthew's request)

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.10 2003/10/09 04:19:18 edhill Exp $
2 C !DESCRIPTION: \bv
3 C $Name: $
4
5 #include "PACKAGES_CONFIG.h"
6 #include "CPP_OPTIONS.h"
7
8 CBOP
9 C !ROUTINE: CALC_GW
10 C !INTERFACE:
11 SUBROUTINE CALC_GW(
12 I myThid)
13 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 IMPLICIT NONE
27 C == Global variables ==
28 #include "SIZE.h"
29 #include "DYNVARS.h"
30 #include "FFIELDS.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "GRID.h"
34 #include "GW.h"
35 #include "CG3D.h"
36
37 C !INPUT/OUTPUT PARAMETERS:
38 C == Routine arguments ==
39 C myThid - Instance number for this innvocation of CALC_GW
40 INTEGER myThid
41
42 #ifdef ALLOW_NONHYDROSTATIC
43
44 C !LOCAL VARIABLES:
45 C == Local variables ==
46 C bi, bj, :: Loop counters
47 C iMin, iMax,
48 C jMin, jMax
49 C flx_NS :: Temp. used for fVol meridional terms.
50 C flx_EW :: Temp. used for fVol zonal terms.
51 C flx_Up :: Temp. used for fVol vertical terms.
52 C flx_Dn :: Temp. used for fVol vertical terms.
53 INTEGER bi,bj,iMin,iMax,jMin,jMax
54 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55 _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
56 _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57 _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58 C I,J,K - Loop counters
59 INTEGER i,j,k, kP1, kUp
60 _RL wOverride
61 _RS hFacROpen
62 _RS hFacRClosed
63 _RL ab15,ab05
64 _RL slipSideFac, slipBotFac
65 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
66
67 _RL Half
68 PARAMETER(Half=0.5D0)
69
70 #define I0 1
71 #define In sNx
72 #define J0 1
73 #define Jn sNy
74 CEOP
75
76 ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
77
78 C Adams-Bashforth timestepping weights
79 ab15 = 1.5 _d 0 + abeps
80 ab05 = -0.5 _d 0 - abeps
81
82 C Lateral friction (no-slip, free slip, or half slip):
83 IF ( no_slip_sides ) THEN
84 slipSideFac = -Half
85 ELSE
86 slipSideFac = Half
87 ENDIF
88 C- half slip was used before ; keep it for now.
89 slipSideFac = 0. _d 0
90
91 C Bottom friction (no-slip, free slip, or half slip):
92 IF ( no_slip_bottom ) THEN
93 slipBotFac = -1. _d 0
94 ELSE
95 slipBotFac = 1. _d 0
96 ENDIF
97 C- half slip was used before ; keep it for now.
98 slipBotFac = 0. _d 0
99
100 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 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
106 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 DO J=J0,Jn
122 DO I=I0,In
123 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 C Flux on Southern face
136 DO J=J0,Jn+1
137 DO I=I0,In
138 tmp_VbarZ=Half*(
139 & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
140 & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
141 Flx_NS(I,J,bi,bj)=
142 & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
143 & -viscAh*_recip_dyC(I,J,bi,bj)
144 & *(1. _d 0 + slipSideFac*
145 & (maskS(I,J,K-1,bi,bj)+maskS(I,J,K,bi,bj)-2. _d 0))
146 & *(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj))
147 ENDDO
148 ENDDO
149 C Flux on Western face
150 DO J=J0,Jn
151 DO I=I0,In+1
152 tmp_UbarZ=Half*(
153 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
154 & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
155 Flx_EW(I,J,bi,bj)=
156 & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
157 & -viscAh*_recip_dxC(I,J,bi,bj)
158 & *(1. _d 0 + slipSideFac*
159 & (maskW(I,J,K-1,bi,bj)+maskW(I,J,K,bi,bj)-2. _d 0))
160 & *(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj))
161 ENDDO
162 ENDDO
163 C Flux on Lower face
164 DO J=J0,Jn
165 DO I=I0,In
166 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
167 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj))
168 Flx_Dn(I,J,bi,bj)=
169 & tmp_WbarZ*tmp_WbarZ
170 & -viscAr*recip_drF(K)
171 & *( 1. _d 0 + slipBotFac*
172 & (wOverRide*maskC(I,J,Kp1,bi,bj)-1. _d 0) )
173 & *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) )
174 ENDDO
175 ENDDO
176 C Divergence of fluxes
177 DO J=J0,Jn
178 DO I=I0,In
179 gW(I,J,K,bi,bj) = 0.
180 & -(
181 & +_recip_dxF(I,J,bi,bj)*(
182 & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
183 & +_recip_dyF(I,J,bi,bj)*(
184 & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
185 & +recip_drC(K) *(
186 & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
187 & )
188 caja * recip_hFacU(I,J,K,bi,bj)
189 caja NOTE: This should be included
190 caja but we need an hFacUW (above U points)
191 caja and an hFacUS (above V points) too...
192 ENDDO
193 ENDDO
194 ENDDO
195 ENDDO
196 ENDDO
197
198
199 DO bj=myByLo(myThid),myByHi(myThid)
200 DO bi=myBxLo(myThid),myBxHi(myThid)
201 DO K=2,Nr
202 DO j=J0,Jn
203 DO i=I0,In
204 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
205 & +deltatMom*( ab15*gW(i,j,k,bi,bj)
206 & +ab05*gWNM1(i,j,k,bi,bj) )
207 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
208 ENDDO
209 ENDDO
210 ENDDO
211 ENDDO
212 ENDDO
213
214 #ifdef ALLOW_OBCS
215 IF (useOBCS) THEN
216 C-- This call is aesthetic: it makes the W field
217 C consistent with the OBs but this has no algorithmic
218 C impact. This is purely for diagnostic purposes.
219 DO bj=myByLo(myThid),myByHi(myThid)
220 DO bi=myBxLo(myThid),myBxHi(myThid)
221 DO K=1,Nr
222 CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
223 ENDDO
224 ENDDO
225 ENDDO
226 ENDIF
227 #endif /* ALLOW_OBCS */
228
229 #endif /* ALLOW_NONHYDROSTATIC */
230
231 RETURN
232 END

  ViewVC Help
Powered by ViewVC 1.1.22