1 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_gw.F,v 1.11 2004/03/29 17:58:40 jmc 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 |
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 |
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 |
gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj) |
97 |
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 |
DO J=J0,Jn |
113 |
DO I=I0,In |
114 |
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 |
DO J=J0,Jn+1 |
128 |
DO I=I0,In |
129 |
tmp_VbarZ=Half*( |
130 |
& _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj) |
131 |
& +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj)) |
132 |
Flx_NS(I,J,bi,bj)= |
133 |
& tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj)) |
134 |
& -viscAh*_recip_dyC(I,J,bi,bj) |
135 |
& *(1. _d 0 + slipSideFac* |
136 |
& (maskS(I,J,K-1,bi,bj)+maskS(I,J,K,bi,bj)-2. _d 0)) |
137 |
& *(wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj)) |
138 |
ENDDO |
139 |
ENDDO |
140 |
C Flux on Western face |
141 |
DO J=J0,Jn |
142 |
DO I=I0,In+1 |
143 |
tmp_UbarZ=Half*( |
144 |
& _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj) |
145 |
& +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj)) |
146 |
Flx_EW(I,J,bi,bj)= |
147 |
& tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj)) |
148 |
& -viscAh*_recip_dxC(I,J,bi,bj) |
149 |
& *(1. _d 0 + slipSideFac* |
150 |
& (maskW(I,J,K-1,bi,bj)+maskW(I,J,K,bi,bj)-2. _d 0)) |
151 |
& *(wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj)) |
152 |
ENDDO |
153 |
ENDDO |
154 |
C Flux on Lower face |
155 |
DO J=J0,Jn |
156 |
DO I=I0,In |
157 |
Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj) |
158 |
tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj)) |
159 |
Flx_Dn(I,J,bi,bj)= |
160 |
& tmp_WbarZ*tmp_WbarZ |
161 |
& -viscAr*recip_drF(K) |
162 |
& *( wVel(I,J,K,bi,bj)-wOverRide*wVel(I,J,Kp1,bi,bj) ) |
163 |
ENDDO |
164 |
ENDDO |
165 |
C Divergence of fluxes |
166 |
DO J=J0,Jn |
167 |
DO I=I0,In |
168 |
gW(I,J,K,bi,bj) = 0. |
169 |
& -( |
170 |
& +_recip_dxF(I,J,bi,bj)*( |
171 |
& Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) ) |
172 |
& +_recip_dyF(I,J,bi,bj)*( |
173 |
& Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) ) |
174 |
& +recip_drC(K) *( |
175 |
& Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) ) |
176 |
& ) |
177 |
caja * recip_hFacU(I,J,K,bi,bj) |
178 |
caja NOTE: This should be included |
179 |
caja but we need an hFacUW (above U points) |
180 |
caja and an hFacUS (above V points) too... |
181 |
ENDDO |
182 |
ENDDO |
183 |
ENDDO |
184 |
ENDDO |
185 |
ENDDO |
186 |
|
187 |
|
188 |
DO bj=myByLo(myThid),myByHi(myThid) |
189 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
190 |
DO K=2,Nr |
191 |
DO j=J0,Jn |
192 |
DO i=I0,In |
193 |
wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj) |
194 |
& +deltatMom*( ab15*gW(i,j,k,bi,bj) |
195 |
& +ab05*gWNM1(i,j,k,bi,bj) ) |
196 |
IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0. |
197 |
ENDDO |
198 |
ENDDO |
199 |
ENDDO |
200 |
ENDDO |
201 |
ENDDO |
202 |
|
203 |
#ifdef ALLOW_OBCS |
204 |
IF (useOBCS) THEN |
205 |
C-- This call is aesthetic: it makes the W field |
206 |
C consistent with the OBs but this has no algorithmic |
207 |
C impact. This is purely for diagnostic purposes. |
208 |
DO bj=myByLo(myThid),myByHi(myThid) |
209 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
210 |
DO K=1,Nr |
211 |
CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid ) |
212 |
ENDDO |
213 |
ENDDO |
214 |
ENDDO |
215 |
ENDIF |
216 |
#endif /* ALLOW_OBCS */ |
217 |
|
218 |
#endif /* ALLOW_NONHYDROSTATIC */ |
219 |
|
220 |
RETURN |
221 |
END |