1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_volflux.F,v 1.9 2007/10/09 00:00:01 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_CPPOPTIONS.h" |
5 |
#ifdef ALLOW_OBCS |
6 |
# include "OBCS_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
subroutine ctrl_volflux( |
10 |
I obcsncount, |
11 |
O sumarea, |
12 |
O sumflux, mythid |
13 |
& ) |
14 |
|
15 |
c ================================================================== |
16 |
c SUBROUTINE ctrl_volflux |
17 |
c ================================================================== |
18 |
c |
19 |
c o calculate the o.b. volume flux due to control adjustments. |
20 |
c o Assume the calendar is identical |
21 |
c for all open boundaries. Need to save the barotropic adjustment |
22 |
c velocity so it can be used in all ctrl_getobcs files. |
23 |
c o WARNING: eastern boundary (not defined) filenames have been a |
24 |
c problem in the past. |
25 |
c |
26 |
c - started G. Gebbie, MIT-WHOI, 15-June-2002 |
27 |
c ================================================================== |
28 |
c SUBROUTINE ctrl_obcsvol |
29 |
c ================================================================== |
30 |
|
31 |
implicit none |
32 |
|
33 |
c == global variables == |
34 |
|
35 |
#include "EEPARAMS.h" |
36 |
#include "SIZE.h" |
37 |
#include "PARAMS.h" |
38 |
#include "GRID.h" |
39 |
#include "DYNVARS.h" |
40 |
#ifdef ALLOW_OBCS |
41 |
# include "OBCS.h" |
42 |
#endif |
43 |
|
44 |
#include "ctrl.h" |
45 |
#include "ctrl_dummy.h" |
46 |
#include "optim.h" |
47 |
|
48 |
c == routine arguments == |
49 |
|
50 |
integer obcsncount |
51 |
_RL sumflux |
52 |
_RL sumarea |
53 |
integer mythid |
54 |
|
55 |
|
56 |
#ifdef BALANCE_CONTROL_VOLFLUX_GLOBAL |
57 |
c == local variables == |
58 |
|
59 |
integer bi,bj |
60 |
integer i,j,k |
61 |
integer itlo,ithi |
62 |
integer jtlo,jthi |
63 |
integer jmin,jmax |
64 |
integer imin,imax |
65 |
integer irec |
66 |
integer il |
67 |
integer iobcs |
68 |
integer ip1 |
69 |
integer jp1 |
70 |
integer nrec |
71 |
integer ilfld |
72 |
integer igg |
73 |
|
74 |
_RL tmpflux |
75 |
_RL tmparea |
76 |
_RL dummy |
77 |
_RL gg |
78 |
_RL tmpx |
79 |
_RL tmpy |
80 |
_RL obcsnfac |
81 |
character*(80) fnamefldn |
82 |
character*(80) fnameflds |
83 |
character*(80) fnamefldw |
84 |
character*(80) fnameflde |
85 |
|
86 |
logical doglobalread |
87 |
logical ladinit |
88 |
|
89 |
#ifdef ECCO_VERBOSE |
90 |
character*(MAX_LEN_MBUF) msgbuf |
91 |
#endif |
92 |
|
93 |
c == external functions == |
94 |
|
95 |
integer ilnblnk |
96 |
external ilnblnk |
97 |
|
98 |
c == end of interface == |
99 |
|
100 |
jtlo = mybylo(mythid) |
101 |
jthi = mybyhi(mythid) |
102 |
itlo = mybxlo(mythid) |
103 |
ithi = mybxhi(mythid) |
104 |
jmin = 1 |
105 |
jmax = sny |
106 |
imin = 1 |
107 |
imax = snx |
108 |
|
109 |
c-- Read tiled data. |
110 |
doglobalread = .false. |
111 |
ladinit = .false. |
112 |
|
113 |
cgg Assume the number of records is the same for |
114 |
cgg all boundaries. Needs to be improved someday. |
115 |
|
116 |
#if (defined (ALLOW_OBCS_CONTROL) || \ |
117 |
defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
118 |
|
119 |
tmpflux = 0. d 0 |
120 |
tmparea = 0. d 0 |
121 |
sumarea = 0. d 0 |
122 |
sumflux = 0. d 0 |
123 |
|
124 |
#ifdef ECCO_VERBOSE |
125 |
_BEGIN_MASTER( mythid ) |
126 |
write(msgbuf,'(a)') ' ' |
127 |
call print_message( msgbuf, standardmessageunit, |
128 |
& SQUEEZE_RIGHT , mythid) |
129 |
write(msgbuf,'(a)') ' ' |
130 |
call print_message( msgbuf, standardmessageunit, |
131 |
& SQUEEZE_RIGHT , mythid) |
132 |
write(msgbuf,'(a,i9.8)') |
133 |
& ' ctrl_volflux: number of records to process: ',nrec |
134 |
call print_message( msgbuf, standardmessageunit, |
135 |
& SQUEEZE_RIGHT , mythid) |
136 |
write(msgbuf,'(a)') ' ' |
137 |
call print_message( msgbuf, standardmessageunit, |
138 |
& SQUEEZE_RIGHT , mythid) |
139 |
_END_MASTER( mythid ) |
140 |
#endif |
141 |
|
142 |
if (optimcycle .ge. 0) then |
143 |
c |
144 |
#ifdef ALLOW_OBCSN_CONTROL |
145 |
ilfld=ilnblnk( xx_obcsn_file ) |
146 |
write(fnamefldn(1:80),'(2a,i10.10)') |
147 |
& xx_obcsn_file(1:ilfld),'.', optimcycle |
148 |
#endif |
149 |
#ifdef ALLOW_OBCSS_CONTROL |
150 |
ilfld=ilnblnk( xx_obcss_file ) |
151 |
write(fnameflds(1:80),'(2a,i10.10)') |
152 |
& xx_obcss_file(1:ilfld),'.',optimcycle |
153 |
#endif |
154 |
#ifdef ALLOW_OBCSW_CONTROL |
155 |
ilfld=ilnblnk( xx_obcsw_file ) |
156 |
write(fnamefldw(1:80),'(2a,i10.10)') |
157 |
& xx_obcsw_file(1:ilfld),'.',optimcycle |
158 |
#endif |
159 |
#ifdef ALLOW_OBCSE_CONTROL |
160 |
ilfld=ilnblnk( xx_obcse_file ) |
161 |
write(fnameflde(1:80),'(2a,i10.10)') |
162 |
& xx_obcse_file(1:ilfld),'.',optimcycle |
163 |
#endif |
164 |
c |
165 |
endif |
166 |
|
167 |
#ifdef ALLOW_OBCSN_CONTROL |
168 |
jp1 = 0 |
169 |
|
170 |
call active_read_xz(fnamefldn,tmpfldxz, |
171 |
& (obcsncount-1)*nobcs+3, doglobalread, |
172 |
& ladinit, optimcycle, mythid |
173 |
& , xx_obcsn_dummy ) |
174 |
|
175 |
c-- Loop over this thread tiles. |
176 |
do bj = jtlo,jthi |
177 |
do bi = itlo,ithi |
178 |
|
179 |
tmpflux = 0. d0 |
180 |
tmparea = 0. d0 |
181 |
|
182 |
do k = 1, Nr |
183 |
do i = imin,imax |
184 |
j = Ob_Jn(I,bi,bj) |
185 |
if (j.ne.0) then |
186 |
cgg -- Alternatively I could read the maskobcs file. But this gives the same result. |
187 |
if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then |
188 |
cgg -- Do not let the corners contribute to the volume flux. |
189 |
if(ob_iw(j,bi,bj).ne.i .and.ob_ie(j,bi,bj).ne.i)then |
190 |
CGG -- Barotropic velocity stored in level 1. |
191 |
tmpx = tmpfldxz(i,1,bi,bj) |
192 |
|
193 |
cgg -- Pick the special point where barotropic velocity loses one degree of freedom. |
194 |
cgg -- Add up the cross-sectional area of this column for later calculations. |
195 |
if (ob_iw(j,bi,bj).eq.(i-1) .and. |
196 |
& ob_iw(j,bi,bj).ne. 0) then |
197 |
tmpx = 0. |
198 |
tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj) |
199 |
print*,'tmparea',tmparea |
200 |
endif |
201 |
cgg -- Positive is flux in. |
202 |
tmpflux = tmpflux -tmpx*delR(k)*dxg(i,j+jp1,bi,bj) |
203 |
endif |
204 |
endif |
205 |
endif |
206 |
enddo |
207 |
enddo |
208 |
|
209 |
sumarea = sumarea+ tmparea |
210 |
sumflux = sumflux + tmpflux |
211 |
enddo |
212 |
enddo |
213 |
#endif |
214 |
|
215 |
#ifdef ALLOW_OBCSS_CONTROL |
216 |
jp1 = 1 |
217 |
|
218 |
call active_read_xz(fnameflds,tmpfldxz, |
219 |
& (obcsncount-1)*nobcs+3, doglobalread, |
220 |
& ladinit, optimcycle, mythid |
221 |
& , xx_obcss_dummy ) |
222 |
|
223 |
c-- Loop over this thread tiles. |
224 |
do bj = jtlo,jthi |
225 |
do bi = itlo,ithi |
226 |
|
227 |
tmpflux = 0. d 0 |
228 |
#ifndef ALLOW_OBCSN_CONTROL |
229 |
tmparea = 0. d 0 |
230 |
#endif |
231 |
do k = 1, Nr |
232 |
do i = imin,imax |
233 |
j = Ob_Js(I,bi,bj) |
234 |
if (j .ne. 0) then |
235 |
if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then |
236 |
cgg -- Do not let the corners contribute to the volume flux. |
237 |
if (ob_iw(j,bi,bj).ne.i.and.ob_ie(j,bi,bj).ne.i)then |
238 |
tmpx = tmpfldxz(i,1,bi,bj) |
239 |
#ifndef ALLOW_OBCSN_CONTROL |
240 |
cgg -- Pick the special point where barotropic velocity loses one degree of freedom. |
241 |
cgg -- Add up the cross-sectional area of this column for later calculations. |
242 |
cgg -- This is just the backup case where the northern boundary does not exist. |
243 |
cgg -- warning: never been tested. |
244 |
if (ob_iw(j,bi,bj).eq.(i-1).and. |
245 |
& ob_iw(j,bi,bj).ne. 0) then |
246 |
tmpx = 0. |
247 |
tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj) |
248 |
print*,'tmparea',tmparea |
249 |
endif |
250 |
#endif |
251 |
cgg -- Positive is flux in. |
252 |
tmpflux = tmpflux +tmpx*delR(k)*dxg(i,j+jp1,bi,bj) |
253 |
endif |
254 |
endif |
255 |
endif |
256 |
enddo |
257 |
enddo |
258 |
#ifndef ALLOW_OBCSN_CONTROL |
259 |
sumarea = sumarea+ tmparea |
260 |
#endif |
261 |
sumflux = sumflux + tmpflux |
262 |
enddo |
263 |
enddo |
264 |
#endif |
265 |
|
266 |
#ifdef ALLOW_OBCSW_CONTROL |
267 |
ip1 = 1 |
268 |
|
269 |
call active_read_yz( fnamefldw, tmpfldyz, |
270 |
& (obcsncount-1)*nobcs+3, doglobalread, |
271 |
& ladinit, optimcycle, mythid |
272 |
& , xx_obcsw_dummy ) |
273 |
|
274 |
c-- Loop over this thread tiles. |
275 |
do bj = jtlo,jthi |
276 |
do bi = itlo,ithi |
277 |
|
278 |
c-- Determine the weights to be used. |
279 |
tmpflux = 0. d 0 |
280 |
#ifndef ALLOW_OBCSN_CONTROL |
281 |
#ifndef ALLOW_OBCSS_CONTROL |
282 |
tmparea = 0. d 0 |
283 |
#endif |
284 |
#endif |
285 |
do k = 1, Nr |
286 |
do j = jmin,jmax |
287 |
i = ob_iw(j,bi,bj) |
288 |
if ( i .ne. 0) then |
289 |
if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then |
290 |
cgg -- Do not let the corners contribute to the volume flux. |
291 |
if (ob_jn(i,bi,bj).ne.j.and. ob_js(i,bi,bj).ne.j)then |
292 |
tmpy = tmpfldyz(j,1,bi,bj) |
293 |
|
294 |
#ifndef ALLOW_OBCSN_CONTROL |
295 |
#ifndef ALLOW_OBCSS_CONTROL |
296 |
cgg -- Pick the special point where barotropic velocity loses one degree of freedom. |
297 |
cgg -- Add up the cross-sectional area of this column for later calculations. |
298 |
cgg -- This is an untested backup case. |
299 |
if (ob_jn(i,bi,bj).eq.(j+1) .and. |
300 |
& ob_jn(i,bi,bj).ne. 0) then |
301 |
tmpy = 0. |
302 |
tmparea = tmparea + delR(k) * dyg(i+ip1,j,bi,bj) |
303 |
print*,'tmparea',tmparea |
304 |
endif |
305 |
#endif |
306 |
#endif |
307 |
cgg -- Positive is flux in. |
308 |
tmpflux = tmpflux + tmpy* delR(k)*dyg(i+ip1,j,bi,bj) |
309 |
endif |
310 |
endif |
311 |
endif |
312 |
enddo |
313 |
enddo |
314 |
#ifndef ALLOW_OBCSN_CONTROL |
315 |
#ifndef ALLOW_OBCSS_CONTROL |
316 |
sumarea =sumarea + tmparea |
317 |
#endif |
318 |
#endif |
319 |
sumflux = sumflux + tmpflux |
320 |
enddo |
321 |
enddo |
322 |
#endif |
323 |
|
324 |
#ifdef ALLOW_OBCSE_CONTROL |
325 |
ip1 = 0 |
326 |
|
327 |
call active_read_yz( fnameflde, tmpfldyz, |
328 |
(obcsncount-1)*nobcs+3, doglobalread, |
329 |
ladinit, optimcycle, mythid |
330 |
& , xx_obcse_dummy ) |
331 |
|
332 |
c-- Loop over this thread tiles. |
333 |
do bj = jtlo,jthi |
334 |
do bi = itlo,ithi |
335 |
|
336 |
c-- Determine the weights to be used. |
337 |
tmpflux = 0. d 0 |
338 |
|
339 |
#ifndef ALLOW_OBCSN_CONTROL |
340 |
#ifndef ALLOW_OBCSS_CONTROL |
341 |
#ifndef ALLOW_OBCSW_CONTROL |
342 |
tmparea = 0. d 0 |
343 |
#endif |
344 |
#endif |
345 |
#endif |
346 |
|
347 |
do k = 1, Nr |
348 |
do j = jmin,jmax |
349 |
i = ob_ie(j,bi,bj) |
350 |
if ( i .ne. 0) then |
351 |
if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then |
352 |
cgg -- Do not let the corners contribute to the volume flux. |
353 |
if (ob_jn(i,bi,bj).ne.j .and.ob_js(i,bi,bj).ne.j)then |
354 |
tmpy = tmpfldyz(j,1,bi,bj) |
355 |
|
356 |
#ifndef ALLOW_OBCSN_CONTROL |
357 |
#ifndef ALLOW_OBCSS_CONTROL |
358 |
#ifndef ALLOW_OBCSW_CONTROL |
359 |
cgg -- Pick the special point where barotropic velocity loses one degree of freedom. |
360 |
cgg -- Add up the cross-sectional area of this column for later calculations. |
361 |
cgg -- This is an untested backup case. |
362 |
if (ob_jn(i,bi,bj).eq.(j+1) .and. |
363 |
& ob_jn(i,bi,bj).ne. 0) then |
364 |
tmpy = 0. |
365 |
tmparea = tmparea + delR(k) * dyg(i+ip1,j,bi,bj) |
366 |
print*,'tmparea',tmparea |
367 |
endif |
368 |
#endif |
369 |
#endif |
370 |
#endif |
371 |
|
372 |
cgg -- Positive is flux in. |
373 |
tmpflux = tmpflux - tmpy* delR(k)*dyg(i+ip1,j,bi,bj) |
374 |
#ifndef ALLOW_OBCSN_CONTROL |
375 |
#ifndef ALLOW_OBCSS_CONTROL |
376 |
#ifndef ALLOW_OBCSW_CONTROL |
377 |
tmparea = tmparea + delR(k) *dyg(i+ip1,j,bi,bj) |
378 |
#endif |
379 |
#endif |
380 |
#endif |
381 |
endif |
382 |
endif |
383 |
endif |
384 |
enddo |
385 |
enddo |
386 |
|
387 |
#ifndef ALLOW_OBCSN_CONTROL |
388 |
#ifndef ALLOW_OBCSS_CONTROL |
389 |
#ifndef ALLOW_OBCSW_CONTROL |
390 |
sumarea = sumarea+ tmparea |
391 |
#endif |
392 |
#endif |
393 |
#endif |
394 |
sumflux = sumflux + tmpflux |
395 |
enddo |
396 |
enddo |
397 |
#endif |
398 |
|
399 |
#endif |
400 |
|
401 |
#endif /* BALANCE_CONTROL_VOLFLUX_GLOBAL */ |
402 |
|
403 |
return |
404 |
end |
405 |
|
406 |
|
407 |
|
408 |
|
409 |
|
410 |
|
411 |
|