/[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.12 - (show annotations) (download)
Mon Mar 22 02:16:43 2010 UTC (14 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62t
Changes since 1.11: +9 -9 lines
finish removing unbalanced quote (single or double) in commented line

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

  ViewVC Help
Powered by ViewVC 1.1.22