1 |
edhill |
1.10 |
C $Header: /u/u3/gcmpack/MITgcm/model/src/calc_gw.F,v 1.9.20.1 2003/10/02 18:10:45 edhill 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 "FFIELDS.h" |
31 |
|
|
#include "EEPARAMS.h" |
32 |
|
|
#include "PARAMS.h" |
33 |
|
|
#include "GRID.h" |
34 |
|
|
#include "GW.h" |
35 |
|
|
#include "CG3D.h" |
36 |
|
|
|
37 |
cnh |
1.9 |
C !INPUT/OUTPUT PARAMETERS: |
38 |
adcroft |
1.1 |
C == Routine arguments == |
39 |
|
|
C myThid - Instance number for this innvocation of CALC_GW |
40 |
|
|
INTEGER myThid |
41 |
|
|
|
42 |
adcroft |
1.3 |
#ifdef ALLOW_NONHYDROSTATIC |
43 |
|
|
|
44 |
cnh |
1.9 |
C !LOCAL VARIABLES: |
45 |
adcroft |
1.1 |
C == Local variables == |
46 |
cnh |
1.9 |
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 |
adcroft |
1.1 |
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 |
adcroft |
1.4 |
INTEGER i,j,k, kP1, kUp |
60 |
adcroft |
1.1 |
_RL wOverride |
61 |
|
|
_RS hFacROpen |
62 |
|
|
_RS hFacRClosed |
63 |
|
|
_RL ab15,ab05 |
64 |
|
|
_RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ |
65 |
|
|
|
66 |
|
|
_RL Half |
67 |
|
|
PARAMETER(Half=0.5D0) |
68 |
|
|
|
69 |
|
|
#define I0 1 |
70 |
|
|
#define In sNx |
71 |
|
|
#define J0 1 |
72 |
|
|
#define Jn sNy |
73 |
cnh |
1.9 |
CEOP |
74 |
edhill |
1.10 |
|
75 |
|
|
ceh3 needs an IF ( useNONHYDROSTATIC ) THEN |
76 |
adcroft |
1.1 |
|
77 |
|
|
C Adams-Bashforth timestepping weights |
78 |
|
|
ab15=1.5+abeps |
79 |
|
|
ab05=-0.5-abeps |
80 |
|
|
|
81 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
82 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
83 |
|
|
DO K=1,Nr |
84 |
|
|
DO j=1-OLy,sNy+OLy |
85 |
|
|
DO i=1-OLx,sNx+OLx |
86 |
jmc |
1.8 |
gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj) |
87 |
adcroft |
1.1 |
gW(i,j,k,bi,bj) = 0. |
88 |
|
|
ENDDO |
89 |
|
|
ENDDO |
90 |
|
|
ENDDO |
91 |
|
|
ENDDO |
92 |
|
|
ENDDO |
93 |
|
|
|
94 |
|
|
C Catch barotropic mode |
95 |
|
|
IF ( Nr .LT. 2 ) RETURN |
96 |
|
|
|
97 |
|
|
C For each tile |
98 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
99 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
100 |
|
|
|
101 |
|
|
C Boundaries condition at top |
102 |
|
|
DO J=J0,Jn |
103 |
|
|
DO I=I0,In |
104 |
|
|
Flx_Dn(I,J,bi,bj)=0. |
105 |
|
|
ENDDO |
106 |
|
|
ENDDO |
107 |
|
|
|
108 |
|
|
C Sweep down column |
109 |
|
|
DO K=2,Nr |
110 |
|
|
Kp1=K+1 |
111 |
|
|
wOverRide=1. |
112 |
|
|
if (K.EQ.Nr) then |
113 |
|
|
Kp1=Nr |
114 |
|
|
wOverRide=0. |
115 |
|
|
endif |
116 |
|
|
C Flux on Southern face |
117 |
|
|
DO J=J0,Jn+1 |
118 |
|
|
DO I=I0,In |
119 |
|
|
tmp_VbarZ=Half*( |
120 |
|
|
& _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj) |
121 |
|
|
& +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj)) |
122 |
|
|
Flx_NS(I,J,bi,bj)= |
123 |
|
|
& tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj)) |
124 |
|
|
& -viscAh*_recip_dyC(I,J,bi,bj)*( |
125 |
|
|
& wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj) ) |
126 |
|
|
ENDDO |
127 |
|
|
ENDDO |
128 |
|
|
C Flux on Western face |
129 |
|
|
DO J=J0,Jn |
130 |
|
|
DO I=I0,In+1 |
131 |
|
|
tmp_UbarZ=Half*( |
132 |
|
|
& _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj) |
133 |
|
|
& +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj)) |
134 |
|
|
Flx_EW(I,J,bi,bj)= |
135 |
|
|
& tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj)) |
136 |
|
|
& -viscAh*_recip_dxC(I,J,bi,bj)*( |
137 |
|
|
& wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj) ) |
138 |
|
|
ENDDO |
139 |
|
|
ENDDO |
140 |
|
|
C Flux on Lower face |
141 |
|
|
DO J=J0,Jn |
142 |
|
|
DO I=I0,In |
143 |
|
|
Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj) |
144 |
|
|
tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj)) |
145 |
|
|
Flx_Dn(I,J,bi,bj)= |
146 |
|
|
& tmp_WbarZ*tmp_WbarZ |
147 |
|
|
& -viscAr*recip_drF(K)*( wVel(I,J,K,bi,bj) |
148 |
|
|
& -wOverRide*wVel(I,J,Kp1,bi,bj) ) |
149 |
|
|
ENDDO |
150 |
|
|
ENDDO |
151 |
|
|
C Divergence of fluxes |
152 |
|
|
DO J=J0,Jn |
153 |
|
|
DO I=I0,In |
154 |
|
|
gW(I,J,K,bi,bj) = 0. |
155 |
|
|
& -( |
156 |
|
|
& +_recip_dxF(I,J,bi,bj)*( |
157 |
|
|
& Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) ) |
158 |
|
|
& +_recip_dyF(I,J,bi,bj)*( |
159 |
|
|
& Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) ) |
160 |
|
|
& +recip_drC(K) *( |
161 |
|
|
& Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) ) |
162 |
|
|
& ) |
163 |
|
|
caja * recip_hFacU(I,J,K,bi,bj) |
164 |
|
|
caja NOTE: This should be included |
165 |
|
|
caja but we need an hFacUW (above U points) |
166 |
|
|
caja and an hFacUS (above V points) too... |
167 |
|
|
ENDDO |
168 |
|
|
ENDDO |
169 |
|
|
ENDDO |
170 |
|
|
ENDDO |
171 |
|
|
ENDDO |
172 |
|
|
|
173 |
|
|
|
174 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
175 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
176 |
|
|
DO K=2,Nr |
177 |
|
|
DO j=J0,Jn |
178 |
|
|
DO i=I0,In |
179 |
|
|
wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj) |
180 |
|
|
& +deltatMom*( ab15*gW(i,j,k,bi,bj) |
181 |
|
|
& +ab05*gWNM1(i,j,k,bi,bj) ) |
182 |
|
|
IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0. |
183 |
|
|
ENDDO |
184 |
|
|
ENDDO |
185 |
adcroft |
1.4 |
ENDDO |
186 |
|
|
ENDDO |
187 |
|
|
ENDDO |
188 |
|
|
|
189 |
adcroft |
1.5 |
#ifdef ALLOW_OBCS |
190 |
|
|
IF (useOBCS) THEN |
191 |
adcroft |
1.4 |
C-- This call is aesthetic: it makes the W field |
192 |
|
|
C consistent with the OBs but this has no algorithmic |
193 |
|
|
C impact. This is purely for diagnostic purposes. |
194 |
adcroft |
1.5 |
DO bj=myByLo(myThid),myByHi(myThid) |
195 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
196 |
|
|
DO K=1,Nr |
197 |
|
|
CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid ) |
198 |
|
|
ENDDO |
199 |
adcroft |
1.1 |
ENDDO |
200 |
|
|
ENDDO |
201 |
adcroft |
1.5 |
ENDIF |
202 |
|
|
#endif /* ALLOW_OBCS */ |
203 |
adcroft |
1.1 |
|
204 |
|
|
#endif /* ALLOW_NONHYDROSTATIC */ |
205 |
|
|
|
206 |
|
|
RETURN |
207 |
|
|
END |
208 |
|
|
|