/[MITgcm]/MITgcm/pkg/obcs/obcs_diag_balance.F
ViewVC logotype

Contents of /MITgcm/pkg/obcs/obcs_diag_balance.F

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


Revision 1.6 - (show annotations) (download)
Fri Dec 19 05:41:20 2014 UTC (9 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, HEAD
Changes since 1.5: +7 -9 lines
fix bug (typo: was using OBN_connect intead of OBS_connect)

1 C $Header: /u/gcmpack/MITgcm/pkg/obcs/obcs_diag_balance.F,v 1.5 2014/11/25 01:08:55 jmc Exp $
2 C $Name: $
3
4 #include "OBCS_OPTIONS.h"
5
6 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
7 CBOP
8 C !ROUTINE: OBCS_DIAG_BALANCE
9
10 C !INTERFACE:
11 SUBROUTINE OBCS_DIAG_BALANCE(
12 U div2d,
13 I uTrans, vTrans, k,
14 I myTime, myIter, myThid )
15
16 C !DESCRIPTION:
17 C *==========================================================*
18 C | SUBROUTINE OBCS_DIAG_BALANCE
19 C | o For diagnostics purpose, modify horizontal divergence
20 C | next (but outside) OB to ensure zero net inflow
21 C *==========================================================*
22
23 C !USES:
24 IMPLICIT NONE
25
26 C === Global variables ===
27 #include "SIZE.h"
28 #include "EEPARAMS.h"
29 #include "PARAMS.h"
30 #include "GRID.h"
31 #include "OBCS_PARAMS.h"
32 #include "OBCS_GRID.h"
33 #include "OBCS_FIELDS.h"
34
35 C !INPUT/OUTPUT PARAMETERS:
36 C div2d :: horizontal divergence x grid-cell volume [m^3/s]
37 C uTrans :: horizontal transport to balance [m^3/s]
38 C vTrans :: horizontal transport to balance [m^3/s]
39 C k :: current level index
40 C myTime :: current time of simulation (s)
41 C myIter :: current iteration number
42 C myThid :: my Thread Id number
43 _RL div2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
44 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
45 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
46 INTEGER k
47 _RL myTime
48 INTEGER myIter
49 INTEGER myThid
50 CEOP
51
52 #ifdef ALLOW_OBCS
53 #ifdef ALLOW_DIAGNOSTICS
54
55 C !FUNCTIONS:
56
57 C !LOCAL VARIABLES:
58 C bi, bj :: tile indices
59 C i,j,k :: loop indices
60 C iB, jB :: local index of open boundary
61 C msgBuf :: Informational/error message buffer
62 INTEGER bi, bj, n
63 INTEGER i, j, iB, jB, iBt, jBt
64 CHARACTER*(MAX_LEN_MBUF) msgBuf
65 _RL areaOB(OBCS_maxConnect), tmpA
66 _RL inFlow(OBCS_maxConnect)
67 _RL tileAreaOB(nSx,nSy,OBCS_maxConnect)
68 _RL tileInFlow(nSx,nSy,OBCS_maxConnect)
69
70 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
71
72 #ifdef ALLOW_DEBUG
73 IF (debugMode) CALL DEBUG_ENTER('OBCS_DIAG_BALANCE',myThid)
74 #endif
75
76 C-- Integrate the transport through each OB
77 DO n=1,OB_connectNumber(k)
78 areaOB(n)= 0. _d 0
79 inFlow(n)= 0. _d 0
80 DO bj=myByLo(myThid),myByHi(myThid)
81 DO bi=myBxLo(myThid),myBxHi(myThid)
82 tileAreaOB(bi,bj,n) = 0.
83 tileInFlow(bi,bj,n) = 0.
84 ENDDO
85 ENDDO
86 ENDDO
87
88 #ifdef ALLOW_OBCS_EAST
89 DO bj=myByLo(myThid),myByHi(myThid)
90 DO bi=myBxLo(myThid),myBxHi(myThid)
91 IF ( tileHasOBE(bi,bj) ) THEN
92 DO j=1,sNy
93 iB = OB_Ie(j,bi,bj)
94 IF ( iB.NE.OB_indexNone .AND. iB.GT.1 ) THEN
95 n = OBE_connect(j,k,bi,bj)
96 IF ( n.GE.1 ) THEN
97 tmpA = drF(k)*hFacW(iB,j,k,bi,bj)*dyG(iB,j,bi,bj)
98 & *maskInW(iB,j,bi,bj)
99 tileAreaOB(bi,bj,n) = tileAreaOB(bi,bj,n) + tmpA
100 tileInFlow(bi,bj,n) = tileInFlow(bi,bj,n)
101 & - uTrans(iB,j,bi,bj)
102 ENDIF
103 ENDIF
104 ENDDO
105 ENDIF
106 ENDDO
107 ENDDO
108 #endif /* ALLOW_OBCS_EAST */
109
110 #ifdef ALLOW_OBCS_WEST
111 DO bj=myByLo(myThid),myByHi(myThid)
112 DO bi=myBxLo(myThid),myBxHi(myThid)
113 IF ( tileHasOBW(bi,bj) ) THEN
114 DO j=1,sNy
115 iB = OB_Iw(j,bi,bj)
116 IF ( iB.NE.OB_indexNone .AND. iB.LT.sNx ) THEN
117 n = OBW_connect(j,k,bi,bj)
118 IF ( n.GE.1 ) THEN
119 iB = 1+iB
120 tmpA = drF(k)*hFacW(iB,j,k,bi,bj)*dyG(iB,j,bi,bj)
121 & *maskInW(iB,j,bi,bj)
122 tileAreaOB(bi,bj,n) = tileAreaOB(bi,bj,n) + tmpA
123 tileInFlow(bi,bj,n) = tileInFlow(bi,bj,n)
124 & + uTrans(iB,j,bi,bj)
125 ENDIF
126 ENDIF
127 ENDDO
128 ENDIF
129 ENDDO
130 ENDDO
131 #endif /* ALLOW_OBCS_WEST */
132
133 #ifdef ALLOW_OBCS_NORTH
134 DO bj=myByLo(myThid),myByHi(myThid)
135 DO bi=myBxLo(myThid),myBxHi(myThid)
136 IF ( tileHasOBN(bi,bj) ) THEN
137 DO i=1,sNx
138 jB = OB_Jn(i,bi,bj)
139 IF ( jB.NE.OB_indexNone .AND. jB.GT.1 ) THEN
140 n = OBN_connect(i,k,bi,bj)
141 IF ( n.GE.1 ) THEN
142 tmpA = drF(k)*hFacS(i,jB,k,bi,bj)*dxG(i,jB,bi,bj)
143 & *maskInS(i,jB,bi,bj)
144 tileAreaOB(bi,bj,n) = tileAreaOB(bi,bj,n) + tmpA
145 tileInFlow(bi,bj,n) = tileInFlow(bi,bj,n)
146 & - vTrans(i,jB,bi,bj)
147 ENDIF
148 ENDIF
149 ENDDO
150 ENDIF
151 ENDDO
152 ENDDO
153 #endif /* ALLOW_OBCS_NORTH */
154
155 #ifdef ALLOW_OBCS_SOUTH
156 DO bj=myByLo(myThid),myByHi(myThid)
157 DO bi=myBxLo(myThid),myBxHi(myThid)
158 IF ( tileHasOBS(bi,bj) ) THEN
159 DO i=1,sNx
160 jB = OB_Js(i,bi,bj)
161 IF ( jB.NE.OB_indexNone .AND. jB.LT.sNy ) THEN
162 n = OBS_connect(i,k,bi,bj)
163 IF ( n.GE.1 ) THEN
164 jB = 1+jB
165 tmpA = drF(k)*hFacS(i,jB,k,bi,bj)*dxG(i,jB,bi,bj)
166 & *maskInS(i,jB,bi,bj)
167 tileAreaOB(bi,bj,n) = tileAreaOB(bi,bj,n) + tmpA
168 tileInFlow(bi,bj,n) = tileInFlow(bi,bj,n)
169 & + vTrans(i,jB,bi,bj)
170 ENDIF
171 ENDIF
172 ENDDO
173 ENDIF
174 ENDDO
175 ENDDO
176 #endif /* ALLOW_OBCS_SOUTH */
177
178 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
179 C-- For each set of OB connected points, calculate a unique velocity
180 C correction and correct to the corresponding OB point divergence
181
182 DO n=1,OB_connectNumber(k)
183 CALL GLOBAL_SUM_TILE_RL( tileAreaOB(1,1,n), areaOB(n), myThid )
184 IF ( areaOB(n).GT.0. ) THEN
185 CALL GLOBAL_SUM_TILE_RL( tileInFlow(1,1,n), inFlow(n), myThid )
186 IF ( debugLevel.GE.debLevC ) THEN
187 WRITE(msgBuf,'(A,I9,2I4,A,1P2E16.8)')
188 & 'OBCS_DIAG_BALANCE (it,k,n=', myIter, k, n,
189 & ' ) inFlow:',inFlow(n),inFlow(n)/areaOB(n)
190 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
191 & SQUEEZE_RIGHT, myThid )
192 ENDIF
193 inFlow(n) = inFlow(n) / areaOB(n)
194 ENDIF
195 ENDDO
196
197 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
198 C-- Add correction:
199
200 DO bj=myByLo(myThid),myByHi(myThid)
201 DO bi=myBxLo(myThid),myBxHi(myThid)
202 #ifdef ALLOW_OBCS_EAST
203 IF ( tileHasOBE(bi,bj) ) THEN
204 c DO j=1-OLy,sNy+OLy
205 DO j=1,sNy
206 IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
207 iB = OB_Ie(j,bi,bj)
208 n = OBE_connect(j,k,bi,bj)
209 IF ( n.EQ.0 ) THEN
210 div2d(iB ,j,bi,bj) = div2d(iB ,j,bi,bj)
211 & - uTrans(iB,j,bi,bj)
212 ELSE
213 div2d(iB ,j,bi,bj) = div2d(iB,j,bi,bj)
214 & + inFlow(n)*drF(k)*hFacW(iB,j,k,bi,bj)
215 & *dyG(iB,j,bi,bj)*maskInW(iB,j,bi,bj)
216 ENDIF
217 ENDIF
218 ENDDO
219 ENDIF
220 #endif /* ALLOW_OBCS_EAST */
221
222 #ifdef ALLOW_OBCS_WEST
223 IF ( tileHasOBW(bi,bj) ) THEN
224 c DO j=1-OLy,sNy+OLy
225 DO j=1,sNy
226 IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
227 iBt = OB_Iw(j,bi,bj)
228 iB = 1+iBt
229 n = OBW_connect(j,k,bi,bj)
230 IF ( n.EQ.0 ) THEN
231 div2d(iBt,j,bi,bj) = div2d(iBt,j,bi,bj)
232 & + uTrans(iB,j,bi,bj)
233 ELSE
234 div2d(iBt,j,bi,bj) = div2d(iBt,j,bi,bj)
235 & + inFlow(n)*drF(k)*hFacW(iB,j,k,bi,bj)
236 & *dyG(iB,j,bi,bj)*maskInW(iB,j,bi,bj)
237 ENDIF
238 ENDIF
239 ENDDO
240 ENDIF
241 #endif /* ALLOW_OBCS_WEST */
242
243 #ifdef ALLOW_OBCS_NORTH
244 IF ( tileHasOBN(bi,bj) ) THEN
245 c DO i=1-OLx,sNx+OLx
246 DO i=1,sNx
247 IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
248 jB = OB_Jn(i,bi,bj)
249 n = OBN_connect(i,k,bi,bj)
250 IF ( n.EQ.0 ) THEN
251 div2d(i,jB ,bi,bj) = div2d(i,jB ,bi,bj)
252 & - vTrans(i,jB,bi,bj)
253 ELSE
254 div2d(i,jB ,bi,bj) = div2d(i,jB ,bi,bj)
255 & + inFlow(n)*drF(k)*hFacS(i,jB,k,bi,bj)
256 & *dxG(i,jB,bi,bj)*maskInS(i,jB,bi,bj)
257 ENDIF
258 ENDIF
259 ENDDO
260 ENDIF
261 #endif /* ALLOW_OBCS_NORTH */
262
263 #ifdef ALLOW_OBCS_SOUTH
264 IF ( tileHasOBS(bi,bj) ) THEN
265 c DO i=1-OLx,sNx+OLx
266 DO i=1,sNx
267 IF ( OB_Js(i,bi,bj).NE.OB_indexNone ) THEN
268 jBt = OB_Js(i,bi,bj)
269 jB = 1+jBt
270 n = OBS_connect(i,k,bi,bj)
271 IF ( n.EQ.0 ) THEN
272 div2d(i,jBt,bi,bj) = div2d(i,jBt,bi,bj)
273 & + vTrans(i,jB,bi,bj)
274 ELSE
275 div2d(i,jBt,bi,bj) = div2d(i,jBt,bi,bj)
276 & + inFlow(n)*drF(k)*hFacS(i,jB,k,bi,bj)
277 & *dxG(i,jB,bi,bj)*maskInS(i,jB,bi,bj)
278 ENDIF
279 ENDIF
280 ENDDO
281 ENDIF
282 #endif /* ALLOW_OBCS_SOUTH */
283
284 ENDDO
285 ENDDO
286
287 #ifdef ALLOW_DEBUG
288 IF (debugMode) CALL DEBUG_LEAVE('OBCS_DIAG_BALANCE',myThid)
289 #endif
290
291 #endif /* ALLOW_DIAGNOSTICS */
292 #endif /* ALLOW_OBCS */
293
294 RETURN
295 END

  ViewVC Help
Powered by ViewVC 1.1.22