/[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.10 - (show annotations) (download)
Thu Oct 9 04:19:18 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint51o_pre, checkpoint51l_post, checkpoint52l_post, checkpoint52k_post, checkpoint52, checkpoint52f_post, checkpoint51t_post, checkpoint51n_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint52e_pre, checkpoint52e_post, checkpoint51n_pre, checkpoint52b_pre, checkpoint51l_pre, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint51r_post, checkpoint51i_post, checkpoint52d_post, checkpoint52a_pre, checkpoint52i_post, checkpoint52h_pre, checkpoint52j_post, branch-netcdf, checkpoint51o_post, checkpoint52a_post, ecco_c52_e35, checkpoint51m_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, tg2-branch, netcdf-sm0, checkpoint51n_branch
Changes since 1.9: +5 -1 lines
 o first check-in for the "branch-genmake2" merge
 o verification suite as run on shelley (gcc 3.2.2):

Wed Oct  8 23:42:29 EDT 2003
                T           S           U           V
G D M    c        m  s        m  s        m  s        m  s
E p a R  g  m  m  e  .  m  m  e  .  m  m  e  .  m  m  e  .
N n k u  2  i  a  a  d  i  a  a  d  i  a  a  d  i  a  a  d
2 d e n  d  n  x  n  .  n  x  n  .  n  x  n  .  n  x  n  .

OPTFILE=NONE

Y Y Y Y 13 16 16 16  0 16 16 16 16 16 16 16 16 13 12  0  0 pass  adjustment.128x64x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16  0  0 16 16  0  0 pass  adjustment.cs-32x32x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16 22  0 16 16 22  0 pass  adjust_nlfs.cs-32x32x1
Y Y Y Y -- 13 13 16 16 13 13 13 13 16 16 16 16 16 16 16 16 N/O   advect_cs
Y Y Y Y -- 22 16 16 16 16 16 16 13 16 16 16 16 16 16 16 16 N/O   advect_xy
Y Y Y Y -- 13 16 13 16 16 16 16 16 16 16 22 16 16 16 16 16 N/O   advect_xz
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  aim.5l_cs
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 16 16 16 16 13 16 pass  aim.5l_Equatorial_Channel
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 13 16 16 13 13 16 pass  aim.5l_LatLon
Y Y Y Y 13 16 16 16 16 16 16 16 16 16 13 12 13 13 16 13 16 pass  exp0
Y Y Y Y 14 16 16 16 16 16 16 16 22 16 16 16 13 16 16 22 16 pass  exp1
Y Y Y Y 13 13 16 13 16 16 16 16 16 13 13 16 16 13 13 13 13 pass  exp2
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  exp4
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 22 16 16 16 22 16 pass  exp5
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  front_relax
Y Y Y Y 14 16 16 13 13 16 16 13 13 16 13 13 16 12 13 13 16 pass  global_ocean.90x40x15
Y Y Y Y 10 16 16 13 13 16 13 16 16 13 13 13 13 16 16 13 16 FAIL  global_ocean.cs32x15
Y Y Y Y  6 11 12 13 13 12 13 16 13  9  9  9  9 10  9  9 11 FAIL  global_ocean_pressure
Y Y Y Y 14 16 16 13 16 16 16 13 13 13 13 13 16 12 16 13 16 pass  global_with_exf
Y Y Y Y 14 16 16 16 16 16 16 16 16 11 13 22 13 16 16  9 16 pass  hs94.128x64x5
Y Y Y Y 13 16 16 16 16 16 16 16 16 11 16 16 16 13 16 22 13 pass  hs94.1x64x5
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 13 13 16 16 22 13 pass  hs94.cs-32x32x5
Y Y Y Y 10 10 16 13 13 16 16 16 22 16 13 13 13 13 13 22 13 FAIL  ideal_2D_oce
Y Y Y Y  8 16 16 16 16 16 16 16 16 13 13  8 16 16 16 16 16 FAIL  internal_wave
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 13 22 13 13 13 22 16 pass  inverted_barometer
Y Y Y Y 12 16 16 16 16 16 16 16 16 16 13 12 13 13 13 13 13 FAIL  lab_sea
Y Y Y Y 11 16 16 16 16 16 16 16 13 13 13 12 13 16 13 12 13 FAIL  natl_box
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  plume_on_slope
Y Y Y Y 13 16 16 16 16 13 16 16 16 16 16 16 16 13 16 16 16 pass  solid-body.cs-32x32x1

1 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 C !DESCRIPTION: \bv
3 C $Name: $
4
5 #include "PACKAGES_CONFIG.h"
6 #include "CPP_OPTIONS.h"
7
8 CBOP
9 C !ROUTINE: CALC_GW
10 C !INTERFACE:
11 SUBROUTINE CALC_GW(
12 I myThid)
13 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 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 C !INPUT/OUTPUT PARAMETERS:
38 C == Routine arguments ==
39 C myThid - Instance number for this innvocation of CALC_GW
40 INTEGER myThid
41
42 #ifdef ALLOW_NONHYDROSTATIC
43
44 C !LOCAL VARIABLES:
45 C == Local variables ==
46 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 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 INTEGER i,j,k, kP1, kUp
60 _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 CEOP
74
75 ceh3 needs an IF ( useNONHYDROSTATIC ) THEN
76
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 gWNM1(i,j,k,bi,bj) = gW(i,j,k,bi,bj)
87 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 ENDDO
186 ENDDO
187 ENDDO
188
189 #ifdef ALLOW_OBCS
190 IF (useOBCS) THEN
191 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 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 ENDDO
200 ENDDO
201 ENDIF
202 #endif /* ALLOW_OBCS */
203
204 #endif /* ALLOW_NONHYDROSTATIC */
205
206 RETURN
207 END
208

  ViewVC Help
Powered by ViewVC 1.1.22