/[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.7 - (show annotations) (download)
Tue Mar 6 17:10:29 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.6: +2 -3 lines
remove "include CG2D.h"

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gw.F,v 1.6 2001/02/04 14:38:46 cnh Exp $
2 C $Name: $
3 #include "CPP_OPTIONS.h"
4
5 SUBROUTINE CALC_GW(
6 I myThid)
7 C /==========================================================\
8 C | S/R CALC_GW |
9 C \==========================================================/
10 IMPLICIT NONE
11
12 C == Global variables ==
13 #include "SIZE.h"
14 #include "DYNVARS.h"
15 #include "FFIELDS.h"
16 #include "EEPARAMS.h"
17 #include "PARAMS.h"
18 #include "GRID.h"
19 #include "GW.h"
20 #include "CG3D.h"
21
22 C == Routine arguments ==
23 C myThid - Instance number for this innvocation of CALC_GW
24 INTEGER myThid
25
26 #ifdef ALLOW_NONHYDROSTATIC
27
28 C == Local variables ==
29 INTEGER bi,bj,iMin,iMax,jMin,jMax
30 _RL aF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
31 _RL vF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
32 _RL v4F(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
33 _RL cF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
34 _RL mT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
35 _RL pF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
36 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
37 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
38
39 _RL flx_NS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
40 _RL flx_EW(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
41 _RL flx_Dn(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
42 _RL flx_Up(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
43 C I,J,K - Loop counters
44 INTEGER i,j,k, kP1, kUp
45 _RL wOverride
46 _RS hFacROpen
47 _RS hFacRClosed
48 _RL ab15,ab05
49 _RL tmp_VbarZ, tmp_UbarZ, tmp_WbarZ
50
51 _RL Half
52 PARAMETER(Half=0.5D0)
53
54 #define I0 1
55 #define In sNx
56 #define J0 1
57 #define Jn sNy
58
59 C Adams-Bashforth timestepping weights
60 ab15=1.5+abeps
61 ab05=-0.5-abeps
62
63 DO bj=myByLo(myThid),myByHi(myThid)
64 DO bi=myBxLo(myThid),myBxHi(myThid)
65 DO K=1,Nr
66 DO j=1-OLy,sNy+OLy
67 DO i=1-OLx,sNx+OLx
68 gW(i,j,k,bi,bj) = 0.
69 ENDDO
70 ENDDO
71 ENDDO
72 ENDDO
73 ENDDO
74
75 C Catch barotropic mode
76 IF ( Nr .LT. 2 ) RETURN
77
78 C For each tile
79 DO bj=myByLo(myThid),myByHi(myThid)
80 DO bi=myBxLo(myThid),myBxHi(myThid)
81
82 C Boundaries condition at top
83 DO J=J0,Jn
84 DO I=I0,In
85 Flx_Dn(I,J,bi,bj)=0.
86 ENDDO
87 ENDDO
88
89 C Sweep down column
90 DO K=2,Nr
91 Kp1=K+1
92 wOverRide=1.
93 if (K.EQ.Nr) then
94 Kp1=Nr
95 wOverRide=0.
96 endif
97 C Flux on Southern face
98 DO J=J0,Jn+1
99 DO I=I0,In
100 tmp_VbarZ=Half*(
101 & _hFacS(I,J,K-1,bi,bj)*vVel( I ,J,K-1,bi,bj)
102 & +_hFacS(I,J, K ,bi,bj)*vVel( I ,J, K ,bi,bj))
103 Flx_NS(I,J,bi,bj)=
104 & tmp_VbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I,J-1,K,bi,bj))
105 & -viscAh*_recip_dyC(I,J,bi,bj)*(
106 & wVel(I,J,K,bi,bj)-wVel(I,J-1,K,bi,bj) )
107 ENDDO
108 ENDDO
109 C Flux on Western face
110 DO J=J0,Jn
111 DO I=I0,In+1
112 tmp_UbarZ=Half*(
113 & _hFacW(I,J,K-1,bi,bj)*uVel( I ,J,K-1,bi,bj)
114 & +_hFacW(I,J, K ,bi,bj)*uVel( I ,J, K ,bi,bj))
115 Flx_EW(I,J,bi,bj)=
116 & tmp_UbarZ*Half*(wVel(I,J,K,bi,bj)+wVel(I-1,J,K,bi,bj))
117 & -viscAh*_recip_dxC(I,J,bi,bj)*(
118 & wVel(I,J,K,bi,bj)-wVel(I-1,J,K,bi,bj) )
119 ENDDO
120 ENDDO
121 C Flux on Lower face
122 DO J=J0,Jn
123 DO I=I0,In
124 Flx_Up(I,J,bi,bj)=Flx_Dn(I,J,bi,bj)
125 tmp_WbarZ=Half*(wVel(I,J,K,bi,bj)+wVel(I,J,Kp1,bi,bj))
126 Flx_Dn(I,J,bi,bj)=
127 & tmp_WbarZ*tmp_WbarZ
128 & -viscAr*recip_drF(K)*( wVel(I,J,K,bi,bj)
129 & -wOverRide*wVel(I,J,Kp1,bi,bj) )
130 ENDDO
131 ENDDO
132 C Divergence of fluxes
133 DO J=J0,Jn
134 DO I=I0,In
135 gW(I,J,K,bi,bj) = 0.
136 & -(
137 & +_recip_dxF(I,J,bi,bj)*(
138 & Flx_EW(I+1,J,bi,bj)-Flx_EW(I,J,bi,bj) )
139 & +_recip_dyF(I,J,bi,bj)*(
140 & Flx_NS(I,J+1,bi,bj)-Flx_NS(I,J,bi,bj) )
141 & +recip_drC(K) *(
142 & Flx_Up(I,J,bi,bj) -Flx_Dn(I,J,bi,bj) )
143 & )
144 caja * recip_hFacU(I,J,K,bi,bj)
145 caja NOTE: This should be included
146 caja but we need an hFacUW (above U points)
147 caja and an hFacUS (above V points) too...
148 ENDDO
149 ENDDO
150 ENDDO
151 ENDDO
152 ENDDO
153
154
155 DO bj=myByLo(myThid),myByHi(myThid)
156 DO bi=myBxLo(myThid),myBxHi(myThid)
157 DO K=2,Nr
158 DO j=J0,Jn
159 DO i=I0,In
160 wVel(i,j,k,bi,bj) = wVel(i,j,k,bi,bj)
161 & +deltatMom*( ab15*gW(i,j,k,bi,bj)
162 & +ab05*gWNM1(i,j,k,bi,bj) )
163 IF (hFacC(I,J,K,bi,bj).EQ.0.) wVel(i,j,k,bi,bj)=0.
164 ENDDO
165 ENDDO
166 ENDDO
167 ENDDO
168 ENDDO
169 DO bj=myByLo(myThid),myByHi(myThid)
170 DO bi=myBxLo(myThid),myBxHi(myThid)
171 DO K=1,Nr
172 DO j=J0,Jn
173 DO i=I0,In
174 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
175 ENDDO
176 ENDDO
177 ENDDO
178 ENDDO
179 ENDDO
180
181 #ifdef ALLOW_OBCS
182 IF (useOBCS) THEN
183 C-- This call is aesthetic: it makes the W field
184 C consistent with the OBs but this has no algorithmic
185 C impact. This is purely for diagnostic purposes.
186 DO bj=myByLo(myThid),myByHi(myThid)
187 DO bi=myBxLo(myThid),myBxHi(myThid)
188 DO K=1,Nr
189 CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
190 ENDDO
191 ENDDO
192 ENDDO
193 ENDIF
194 #endif /* ALLOW_OBCS */
195
196 #endif /* ALLOW_NONHYDROSTATIC */
197
198 RETURN
199 END
200

  ViewVC Help
Powered by ViewVC 1.1.22