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

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

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


Revision 1.8 - (hide annotations) (download)
Mon Mar 12 20:51:03 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint38, checkpoint40pre2, checkpoint40pre4, pre38tag1, c37_adj, pre38-close, checkpoint39, checkpoint37, checkpoint40pre5, checkpoint40
Branch point for: pre38
Changes since 1.7: +2 -12 lines
copy gW -> gwNm1 moved from the end to the beginning of the routine

1 jmc 1.8 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_gw.F,v 1.7 2001/03/06 17:10:29 jmc Exp $
2 jmc 1.7 C $Name: $
3 adcroft 1.1 #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 adcroft 1.3 #ifdef ALLOW_NONHYDROSTATIC
27    
28 adcroft 1.1 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 adcroft 1.4 INTEGER i,j,k, kP1, kUp
45 adcroft 1.1 _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 jmc 1.8 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
69 adcroft 1.1 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 adcroft 1.4 ENDDO
168     ENDDO
169     ENDDO
170    
171 adcroft 1.5 #ifdef ALLOW_OBCS
172     IF (useOBCS) THEN
173 adcroft 1.4 C-- This call is aesthetic: it makes the W field
174     C consistent with the OBs but this has no algorithmic
175     C impact. This is purely for diagnostic purposes.
176 adcroft 1.5 DO bj=myByLo(myThid),myByHi(myThid)
177     DO bi=myBxLo(myThid),myBxHi(myThid)
178     DO K=1,Nr
179     CALL OBCS_APPLY_W( bi, bj, K, wVel, myThid )
180     ENDDO
181 adcroft 1.1 ENDDO
182     ENDDO
183 adcroft 1.5 ENDIF
184     #endif /* ALLOW_OBCS */
185 adcroft 1.1
186     #endif /* ALLOW_NONHYDROSTATIC */
187    
188     RETURN
189     END
190    

  ViewVC Help
Powered by ViewVC 1.1.22