/[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.9 - (show annotations) (download)
Wed Sep 26 18:09:14 2001 UTC (22 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint47e_post, checkpoint44e_post, checkpoint46l_post, checkpoint46g_pre, checkpoint47c_post, release1_p13_pre, checkpoint50c_post, checkpoint46f_post, checkpoint48e_post, checkpoint50c_pre, checkpoint44f_post, checkpoint46b_post, checkpoint43a-release1mods, ecco_c50_e32, ecco_c50_e33, ecco_c50_e30, ecco_c50_e31, release1_p13, checkpoint48i_post, checkpoint46l_pre, chkpt44d_post, checkpoint51, checkpoint50, release1_p8, release1_p9, checkpoint50d_post, release1_p1, release1_p2, release1_p3, release1_p4, release1_p5, release1_p6, release1_p7, checkpoint50b_pre, checkpoint44e_pre, checkpoint51f_post, release1_b1, ecco_c51_e34d, ecco_c51_e34e, ecco_c51_e34f, ecco_c51_e34g, ecco_c51_e34a, ecco_c51_e34b, ecco_c51_e34c, checkpoint48b_post, checkpoint43, checkpoint51d_post, checkpoint48c_pre, checkpoint47d_pre, release1_chkpt44d_post, checkpoint47a_post, checkpoint48d_pre, checkpoint51j_post, checkpoint47i_post, release1_p11, checkpoint47d_post, icebear5, icebear4, icebear3, icebear2, checkpoint46d_pre, checkpoint48d_post, release1-branch_tutorials, checkpoint48f_post, checkpoint45d_post, checkpoint46j_pre, chkpt44a_post, checkpoint44h_pre, checkpoint48h_post, ecco_c50_e29, checkpoint51b_pre, checkpoint46a_post, checkpoint47g_post, checkpoint46j_post, checkpoint51h_pre, checkpoint46k_post, ecco_c50_e28, chkpt44c_pre, checkpoint48a_post, checkpoint45a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, checkpoint47j_post, ecco_c50_e33a, branch-exfmods-tag, checkpoint44g_post, branchpoint-genmake2, checkpoint46e_pre, checkpoint48c_post, checkpoint45b_post, checkpoint46b_pre, release1-branch-end, release1_final_v1, checkpoint51b_post, checkpoint51c_post, checkpoint46c_pre, checkpoint46, checkpoint47b_post, checkpoint44b_post, ecco_c51_e34, checkpoint46h_pre, checkpoint46m_post, checkpoint46a_pre, checkpoint50g_post, checkpoint45c_post, ecco_ice2, ecco_ice1, checkpoint44h_post, checkpoint46g_post, release1_p12_pre, ecco_c44_e22, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, ecco_c44_e25, checkpoint51i_pre, checkpoint47f_post, checkpoint50e_post, chkpt44a_pre, checkpoint46i_post, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, checkpoint46c_post, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, checkpoint50d_pre, checkpoint46e_post, release1_beta1, checkpoint51e_post, checkpoint44b_pre, checkpoint42, checkpoint41, checkpoint47, checkpoint44, checkpoint45, checkpoint48, checkpoint49, checkpoint46h_post, checkpoint51f_pre, chkpt44c_post, checkpoint48g_post, checkpoint47h_post, checkpoint44f_pre, checkpoint51g_post, checkpoint46d_post, checkpoint50b_post, release1-branch_branchpoint, checkpoint51a_post
Branch point for: c24_e25_ice, branch-exfmods-curt, release1_final, release1-branch, branch-genmake2, release1, ecco-branch, release1_50yr, icebear, release1_coupled
Changes since 1.8: +28 -14 lines
Bringing comments up to data and formatting for document extraction.

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

  ViewVC Help
Powered by ViewVC 1.1.22