/[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.6 - (show annotations) (download)
Sun Feb 4 14:38:46 2001 UTC (23 years, 3 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint36, checkpoint35
Changes since 1.5: +2 -0 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

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

  ViewVC Help
Powered by ViewVC 1.1.22