/[MITgcm]/MITgcm/pkg/ctrl/ctrl_obcsvol.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_obcsvol.F

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


Revision 1.7 - (show annotations) (download)
Fri Sep 10 12:10:28 2004 UTC (19 years, 8 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint55c_post, checkpoint55e_post, checkpoint55d_pre, checkpoint55j_post, checkpoint56b_post, checkpoint55h_post, checkpoint55b_post, checkpoint55, checkpoint55a_post, checkpoint54f_post, checkpoint55g_post, checkpoint55f_post, checkpoint55i_post, checkpoint56, checkpoint56a_post, checkpoint55d_post
Changes since 1.6: +1 -1 lines
 o remove single quote in comments for online code browser

1
2 #include "CTRL_CPPOPTIONS.h"
3 #ifdef ALLOW_OBCS
4 # include "OBCS_OPTIONS.h"
5 #endif
6
7 subroutine ctrl_obcsvol(
8 I mytime,
9 I myiter,
10 I mythid
11 & )
12
13 c ==================================================================
14 c SUBROUTINE ctrl_obcsvol
15 c ==================================================================
16 c
17 c o volumetrically balance the control vector contribution.
18 c o Assume the calendar is identical
19 c for all open boundaries. Need to save the barotropic adjustment
20 c velocity so it can be used in all ctrl_getobcs files.
21 c ==================================================================
22 c SUBROUTINE ctrl_obcsvol
23 c ==================================================================
24
25 implicit none
26
27 c == global variables ==
28
29 #include "EEPARAMS.h"
30 #include "SIZE.h"
31 #include "PARAMS.h"
32 #include "GRID.h"
33 #include "DYNVARS.h"
34 #ifdef ALLOW_OBCS
35 # include "OBCS.h"
36 #endif
37
38 #include "ctrl.h"
39 #include "ctrl_dummy.h"
40 #include "optim.h"
41
42 c == routine arguments ==
43
44 integer myiter
45 _RL mytime
46 integer mythid
47
48 c == local variables ==
49
50 integer bi,bj
51 integer i,j,k
52 integer itlo,ithi
53 integer jtlo,jthi
54 integer jmin,jmax
55 integer imin,imax
56 integer irec
57 integer il
58 integer iobcs
59 integer ip1
60 integer jp1
61 integer nrec
62 integer ilfld
63 integer igg
64
65 _RL sumvol
66 _RL sumarea
67 _RL tmpflux
68 _RL tmparea
69 _RL dummy
70 _RL gg
71 _RL tmpx
72 _RL tmpy
73 character*(80) fnamefldn
74 character*(80) fnameflds
75 character*(80) fnamefldw
76 character*(80) fnameflde
77
78 logical doglobalread
79 logical ladinit
80 logical obcsnfirst, obcsnchanged
81 integer obcsncount0, obcsncount1
82 _RL obcsnfac
83
84 #ifdef ECCO_VERBOSE
85 character*(MAX_LEN_MBUF) msgbuf
86 #endif
87
88 c == external functions ==
89
90 integer ilnblnk
91 external ilnblnk
92
93 c == end of interface ==
94
95 jtlo = mybylo(mythid)
96 jthi = mybyhi(mythid)
97 itlo = mybxlo(mythid)
98 ithi = mybxhi(mythid)
99 jmin = 1
100 jmax = sny
101 imin = 1
102 imax = snx
103
104 c-- Read tiled data.
105 doglobalread = .false.
106 ladinit = .false.
107
108 cgg Assume the number of records is the same for
109 cgg all boundaries.
110
111
112 #ifdef BALANCE_CONTROL_VOLFLUX_GLOBAL
113 tmpflux= 0. d 0
114 tmparea= 0. d 0
115 sumarea= 0. d 0
116 sumvol = 0. d 0
117
118 #ifdef ECCO_VERBOSE
119 _BEGIN_MASTER( mythid )
120 write(msgbuf,'(a)') ' '
121 call print_message( msgbuf, standardmessageunit,
122 & SQUEEZE_RIGHT , mythid)
123 write(msgbuf,'(a)') ' '
124 call print_message( msgbuf, standardmessageunit,
125 & SQUEEZE_RIGHT , mythid)
126 write(msgbuf,'(a,i9.8)')
127 & ' ctrl_obcsvol: number of records to process: ',nrec
128 call print_message( msgbuf, standardmessageunit,
129 & SQUEEZE_RIGHT , mythid)
130 write(msgbuf,'(a)') ' '
131 call print_message( msgbuf, standardmessageunit,
132 & SQUEEZE_RIGHT , mythid)
133 _END_MASTER( mythid )
134 #endif
135
136 c-- Get the counters, flags, and the interpolation factor.
137 call ctrl_GetRec( 'xx_obcsn',
138 O obcsnfac, obcsnfirst, obcsnchanged,
139 O obcsncount0,obcsncount1,
140 I mytime, myiter, mythid )
141
142 if (optimcycle .ge. 0) then
143 c
144 ilfld=ilnblnk( xx_obcsn_file )
145 write(fnamefldn(1:80),'(2a,i10.10)')
146 & xx_obcsn_file(1:ilfld),'.', optimcycle
147 ilfld=ilnblnk( xx_obcss_file )
148 write(fnameflds(1:80),'(2a,i10.10)')
149 & xx_obcss_file(1:ilfld),'.',optimcycle
150 ilfld=ilnblnk( xx_obcsw_file )
151 write(fnamefldw(1:80),'(2a,i10.10)')
152 & xx_obcsw_file(1:ilfld),'.',optimcycle
153 ilfld=ilnblnk( xx_obcse_file )
154 write(fnameflde(1:80),'(2a,i10.10)')
155 & xx_obcse_file(1:ilfld),'.',optimcycle
156 c
157 endif
158
159 c-- Loop over records. For north boundary, we only need V velocity.
160
161 if ( obcsnfirst ) then
162
163 shiftvel(1) = 0. d0
164 shiftvel(2) = 0. d0
165
166 #ifdef ALLOW_OBCSN_CONTROL
167 jp1 = 0
168
169 call active_read_xz_loc(fnamefldn,tmpfldxz,
170 & (obcsncount0-1)*nobcs+4, doglobalread,
171 & ladinit, optimcycle, mythid
172 & , xx_obcsn_dummy )
173
174 c-- Loop over this thread's tiles.
175 do bj = jtlo,jthi
176 do bi = itlo,ithi
177
178 tmpflux = 0. d0
179 tmparea = 0. d0
180
181 do k = 1, Nr
182 do i = imin,imax
183 j = Ob_Jn(I,bi,bj)
184 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
185 cgg -- Do not let the corners contribute to the volume flux.
186 if(ob_iw(j,bi,bj).ne.i .and. ob_ie(j,bi,bj).ne.i) then
187 tmpx = tmpfldxz(i,k,bi,bj)
188 cgg -- Positive is flux in.
189 tmpflux = tmpflux - tmpx* delR(k)*dxg(i,j+jp1,bi,bj)
190 tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj)
191 endif
192 endif
193 enddo
194 enddo
195
196 sumarea = sumarea+ tmparea
197 sumvol = sumvol + tmpflux
198 enddo
199 enddo
200 #endif
201
202 #ifdef ALLOW_OBCSS_CONTROL
203 jp1 = 1
204
205 call active_read_xz_loc(fnameflds,tmpfldxz,
206 & (obcsncount0-1)*nobcs+4, doglobalread,
207 & ladinit, optimcycle, mythid
208 & , xx_obcss_dummy )
209
210 c-- Loop over this thread's tiles.
211 do bj = jtlo,jthi
212 do bi = itlo,ithi
213
214 tmpflux = 0. d 0
215 tmparea = 0. d 0
216
217 do k = 1, Nr
218 do i = imin,imax
219 j = Ob_Js(I,bi,bj)
220 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
221 cgg -- Do not let the corners contribute to the volume flux.
222 if (ob_iw(j,bi,bj).ne.i .and.ob_ie(j,bi,bj).ne.i) then
223 tmpx = tmpfldxz(i,k,bi,bj)
224 cgg -- Positive is flux in.
225 tmpflux = tmpflux + tmpx* delR(k)*dxg(i,j+jp1,bi,bj)
226 tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj)
227 endif
228 endif
229 enddo
230 enddo
231 sumarea = sumarea+ tmparea
232 sumvol = sumvol + tmpflux
233 enddo
234 enddo
235 #endif
236
237 #ifdef ALLOW_OBCSW_CONTROL
238 ip1 = 1
239
240 call active_read_yz_loc( fnamefldw, tmpfldyz,
241 & (obcsncount0-1)*nobcs+3, doglobalread,
242 & ladinit, optimcycle, mythid
243 & , xx_obcsw_dummy )
244
245 c-- Loop over this thread's tiles.
246 do bj = jtlo,jthi
247 do bi = itlo,ithi
248
249 c-- Determine the weights to be used.
250 tmpflux = 0. d 0
251 tmparea = 0. d 0
252
253 do k = 1, Nr
254 do j = jmin,jmax
255 i = ob_iw(j,bi,bj)
256 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
257 cgg -- Do not let the corners contribute to the volume flux.
258 if (ob_jn(i,bi,bj).ne.j .and. ob_js(i,bi,bj).ne.j)then
259 tmpy = tmpfldyz(j,k,bi,bj)
260 cgg -- Positive is flux in.
261 tmpflux = tmpflux + tmpy* delR(k)*dyg(i+ip1,j,bi,bj)
262 tmparea = tmparea + delR(k)*dyg(i+ip1,j,bi,bj)
263 endif
264 endif
265 enddo
266 enddo
267 sumarea =sumarea + tmparea
268 sumvol = sumvol + tmpflux
269 enddo
270 enddo
271 #endif
272
273 #ifdef ALLOW_OBCSE_CONTROL
274 ip1 = 0
275
276 call active_read_yz_loc( fnameflde, tmpfldyz,
277 (obcsncount0-1)*nobcs+3, doglobalread,
278 ladinit, optimcycle, mythid
279 & , xx_obcse_dummy )
280
281 c-- Loop over this thread's tiles.
282 do bj = jtlo,jthi
283 do bi = itlo,ithi
284
285 c-- Determine the weights to be used.
286 tmpflux = 0. d 0
287 tmparea = 0. d 0
288
289 do k = 1, Nr
290 do j = jmin,jmax
291 i = ob_ie(j,bi,bj)
292 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
293 cgg -- Do not let the corners contribute to the volume flux.
294 if (ob_jn(i,bi,bj).ne.j .and. ob_js(i,bi,bj).ne.j)then
295 tmpy = tmpfldyz(j,k,bi,bj)
296 cgg -- Positive is flux in.
297 tmpflux = tmpflux - tmpy* delR(k)*dyg(i+ip1,j,bi,bj)
298 tmparea = tmparea + delR(k) *dyg(i+ip1,j,bi,bj)
299 endif
300 endif
301 enddo
302 enddo
303 sumarea = sumarea+ tmparea
304 sumvol = sumvol + tmpflux
305 enddo
306 enddo
307 #endif
308
309 c-- Do the global summation.
310 _GLOBAL_SUM_R8( sumvol, mythid )
311 _GLOBAL_SUM_R8( sumarea,mythid )
312
313 shiftvel(2) = sumvol /sumarea
314 endif
315 cgg End of the obcsnfirst loop.
316
317 if ( ( obcsnfirst) .or. (obcsnchanged)) then
318
319 cgg Swap the value.
320 shiftvel(1) = shiftvel(2)
321
322 sumvol = 0. d0
323 sumarea= 0. d0
324
325 #ifdef ALLOW_OBCSN_CONTROL
326 jp1 = 0
327
328 call active_read_xz_loc(fnamefldn,tmpfldxz,
329 & (obcsncount1-1)*nobcs+4, doglobalread,
330 & ladinit, optimcycle, mythid
331 & , xx_obcsn_dummy )
332 c-- Loop over this thread's tiles.
333 do bj = jtlo,jthi
334 do bi = itlo,ithi
335
336 tmpflux = 0. d0
337 tmparea = 0. d0
338
339 do k = 1, Nr
340 do i = imin,imax
341 j = Ob_Jn(I,bi,bj)
342 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
343 cgg -- Do not let the corners contribute to the volume flux.
344 if (ob_iw(j,bi,bj).ne.i .and. ob_ie(j,bi,bj).ne.i)then
345 tmpx = tmpfldxz(i,k,bi,bj)
346 cgg -- Positive is flux in.
347 tmpflux = tmpflux - tmpx* delR(k)*dxg(i,j+jp1,bi,bj)
348 tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj)
349 endif
350 endif
351 enddo
352 enddo
353
354 sumarea = sumarea+ tmparea
355 sumvol = sumvol + tmpflux
356 enddo
357 enddo
358 #endif
359
360 #ifdef ALLOW_OBCSS_CONTROL
361 jp1 = 1
362
363 call active_read_xz_loc(fnameflds,tmpfldxz,
364 & (obcsncount1-1)*nobcs+4, doglobalread,
365 & ladinit, optimcycle, mythid
366 & , xx_obcss_dummy )
367
368 c-- Loop over this thread's tiles.
369 do bj = jtlo,jthi
370 do bi = itlo,ithi
371
372 tmpflux = 0. d 0
373 tmparea = 0. d 0
374
375 do k = 1, Nr
376 do i = imin,imax
377 j = Ob_Js(I,bi,bj)
378 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
379 cgg -- Do not let the corners contribute to the volume flux.
380 if (ob_iw(j,bi,bj).ne.i .and. ob_ie(j,bi,bj).ne.i)then
381 tmpx = tmpfldxz(i,k,bi,bj)
382 cgg -- Positive is flux in.
383 tmpflux = tmpflux + tmpx* delR(k)*dxg(i,j+jp1,bi,bj)
384 tmparea = tmparea + delR(k) * dxg(i,j+jp1,bi,bj)
385 endif
386 endif
387 enddo
388 enddo
389 sumarea = sumarea+ tmparea
390 sumvol = sumvol + tmpflux
391 enddo
392 enddo
393 #endif
394
395 #ifdef ALLOW_OBCSW_CONTROL
396 ip1 = 1
397
398 call active_read_yz_loc( fnamefldw, tmpfldyz,
399 & (obcsncount1-1)*nobcs+3, doglobalread,
400 & ladinit, optimcycle, mythid
401 & , xx_obcsw_dummy )
402
403 c-- Loop over this thread's tiles.
404 do bj = jtlo,jthi
405 do bi = itlo,ithi
406
407 c-- Determine the weights to be used.
408 tmpflux = 0. d 0
409 tmparea = 0. d 0
410
411 do k = 1, Nr
412 do j = jmin,jmax
413 i = ob_iw(j,bi,bj)
414 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
415 cgg -- Do not let corners contribute.
416 if (ob_jn(i,bi,bj).ne.j .and. ob_js(i,bi,bj).ne.j)then
417 tmpy = tmpfldyz(j,k,bi,bj)
418 cgg -- Positive is flux in.
419 tmpflux = tmpflux + tmpy* delR(k)*dyg(i+ip1,j,bi,bj)
420 tmparea = tmparea + delR(k)*dyg(i+ip1,j,bi,bj)
421 endif
422 endif
423 enddo
424 enddo
425 sumarea =sumarea + tmparea
426 sumvol = sumvol + tmpflux
427 enddo
428 enddo
429 #endif
430
431 #ifdef ALLOW_OBCSE_CONTROL
432 ip1 = 0
433
434 call active_read_yz_loc( fnameflde, tmpfldyz,
435 (obcsncount1-1)*nobcs+3, doglobalread,
436 ladinit, optimcycle, mythid
437 & , xx_obcse_dummy )
438
439 c-- Loop over this thread's tiles.
440 do bj = jtlo,jthi
441 do bi = itlo,ithi
442
443 c-- Determine the weights to be used.
444 tmpflux = 0. d 0
445 tmparea = 0. d 0
446
447 do k = 1, Nr
448 do j = jmin,jmax
449 i = ob_ie(j,bi,bj)
450 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
451 cgg -- Do not let the corners contribute to the volume flux.
452 if (ob_jn(i,bi,bj).ne.j .and. ob_js(i,bi,bj).ne.j)then
453 tmpy = tmpfldyz(j,k,bi,bj)
454 cgg -- Positive is flux in.
455 tmpflux = tmpflux - tmpy* delR(k)*dyg(i+ip1,j,bi,bj)
456 tmparea = tmparea + delR(k) *dyg(i+ip1,j,bi,bj)
457 endif
458 endif
459 enddo
460 enddo
461 sumarea = sumarea+ tmparea
462 sumvol = sumvol + tmpflux
463 enddo
464 enddo
465 #endif
466
467 c-- Do the global summation.
468 _GLOBAL_SUM_R8( sumvol, mythid )
469 _GLOBAL_SUM_R8( sumarea,mythid )
470
471 shiftvel(2) = sumvol /sumarea
472 endif
473 cgg End of the obcsnfirst, obcsnchanged loop.
474
475 #endif
476
477 return
478 end
479
480
481
482
483
484
485

  ViewVC Help
Powered by ViewVC 1.1.22