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

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

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


Revision 1.4 - (show annotations) (download)
Fri Dec 3 00:48:57 2004 UTC (19 years, 6 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint57t_post, checkpoint58l_post, checkpoint57o_post, checkpoint57m_post, checkpoint58e_post, checkpoint57v_post, checkpoint57f_post, checkpoint57s_post, checkpoint57j_post, checkpoint58b_post, checkpoint58m_post, checkpoint57b_post, checkpoint57f_pre, checkpoint57k_post, checkpoint57d_post, checkpoint57g_post, checkpoint57c_pre, checkpoint58r_post, checkpoint57i_post, checkpoint57y_post, checkpoint58g_post, checkpoint57x_post, checkpoint58n_post, checkpoint58x_post, checkpoint57g_pre, checkpoint58h_post, checkpoint57e_post, checkpoint58w_post, checkpoint56c_post, checkpoint58j_post, checkpoint57h_post, checkpoint57y_pre, checkpoint57a_post, checkpoint58q_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint57, eckpoint57e_pre, checkpoint57h_done, checkpoint58f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58o_post, checkpoint57z_post, checkpoint57c_post, checkpoint58u_post, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint58p_post, checkpoint58t_post, checkpoint57h_pre, checkpoint57l_post
Changes since 1.3: +3 -0 lines
o OBCS as control variables
  - update ad_diff.list
  - remove balance of obcs controls from default
  - fix index bug nobcs in ctrl_init
  - fix dummy fields filen in ctrl_pack
  - add dummy weights for obcs

1
2 #include "CTRL_CPPOPTIONS.h"
3 #ifdef ALLOW_OBCS
4 # include "OBCS_OPTIONS.h"
5 #endif
6
7 subroutine ctrl_obcsbal(
8 I mytime,
9 I myiter,
10 I mythid
11 & )
12
13 c ==================================================================
14 c SUBROUTINE ctrl_obcsbal
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 o WARNING: eastern boundary (not defined) filenames have been a
22 c problem in the past.
23 c
24 c - started G. Gebbie, MIT-WHOI, 15-June-2002
25 c ==================================================================
26 c SUBROUTINE ctrl_obcsvol
27 c ==================================================================
28
29 implicit none
30
31 c == global variables ==
32
33 #include "EEPARAMS.h"
34 #include "SIZE.h"
35 #include "PARAMS.h"
36 #include "GRID.h"
37 #include "DYNVARS.h"
38 #ifdef ALLOW_OBCS
39 # include "OBCS.h"
40 #endif
41
42 #include "ctrl.h"
43 #include "ctrl_dummy.h"
44
45 c == routine arguments ==
46
47 integer myiter
48 _RL mytime
49 integer mythid
50
51 #ifdef BALANCE_CONTROL_VOLFLUX_GLOBAL
52 c == local variables ==
53
54 integer bi,bj
55 integer i,j,k
56 integer itlo,ithi
57 integer jtlo,jthi
58 integer jmin,jmax
59 integer imin,imax
60 integer irec
61 integer il
62 integer iobcs
63 integer ip1
64 integer jp1
65 integer nrec
66 integer ilfld
67 integer igg
68
69 _RL volflux
70 _RL area
71 _RL tmpflux
72 _RL tmparea
73 _RL dummy
74 _RL gg
75 _RL tmpx
76 _RL tmpy
77 _RL obcsnfac
78 character*(80) fnamefldn
79 character*(80) fnameflds
80 character*(80) fnamefldw
81 character*(80) fnameflde
82
83 logical doglobalread
84 logical ladinit
85 logical obcsnfirst, obcsnchanged
86 integer obcsncount0, obcsncount1
87
88 #ifdef ECCO_VERBOSE
89 character*(MAX_LEN_MBUF) msgbuf
90 #endif
91
92 c == external functions ==
93
94 integer ilnblnk
95 external ilnblnk
96
97 c == end of interface ==
98
99 jtlo = mybylo(mythid)
100 jthi = mybyhi(mythid)
101 itlo = mybxlo(mythid)
102 ithi = mybxhi(mythid)
103 jmin = 1
104 jmax = sny
105 imin = 1
106 imax = snx
107
108 c-- Read tiled data.
109 doglobalread = .false.
110 ladinit = .false.
111
112 cgg Assume the number of records is the same for
113 cgg all boundaries. Needs to be improved someday.
114
115 #if (defined (ALLOW_OBCS_CONTROL) || \
116 defined (ALLOW_OBCS_COST_CONTRIBUTION))
117
118 tmpflux= 0. d 0
119 tmparea= 0. d 0
120 area= 0. d 0
121 volflux = 0. d 0
122
123 #ifdef ECCO_VERBOSE
124 _BEGIN_MASTER( mythid )
125 write(msgbuf,'(a)') ' '
126 call print_message( msgbuf, standardmessageunit,
127 & SQUEEZE_RIGHT , mythid)
128 write(msgbuf,'(a)') ' '
129 call print_message( msgbuf, standardmessageunit,
130 & SQUEEZE_RIGHT , mythid)
131 write(msgbuf,'(a,i9.8)')
132 & ' ctrl_obcsvol: number of records to process: ',nrec
133 call print_message( msgbuf, standardmessageunit,
134 & SQUEEZE_RIGHT , mythid)
135 write(msgbuf,'(a)') ' '
136 call print_message( msgbuf, standardmessageunit,
137 & SQUEEZE_RIGHT , mythid)
138 _END_MASTER( mythid )
139 #endif
140
141 c-- Get the counters, flags, and the interpolation factor.
142 call ctrl_GetRec( 'xx_obcsn',
143 O obcsnfac, obcsnfirst, obcsnchanged,
144 O obcsncount0,obcsncount1,
145 I mytime, myiter, mythid )
146
147 c-- Loop over records. For north boundary, we only need V velocity.
148
149 if ( obcsnfirst ) then
150
151 shiftvel(1) = 0. d0
152 shiftvel(2) = 0. d0
153
154 call ctrl_volflux( obcsncount0, area, volflux, mythid)
155
156 c-- Do the global summation.
157 _GLOBAL_SUM_R8( volflux, mythid )
158 _GLOBAL_SUM_R8( area,mythid )
159
160 shiftvel(2) = volflux / area
161 print*,'volflux,area',volflux,area
162 endif
163 cgg End of the obcsnfirst loop.
164
165 if ( ( obcsnfirst) .or. (obcsnchanged)) then
166
167 cgg Swap the value.
168 shiftvel(1) = shiftvel(2)
169
170 volflux = 0. d0
171 area= 0. d0
172
173 call ctrl_volflux( obcsncount1, area, volflux, mythid)
174
175 c-- Do the global summation.
176 _GLOBAL_SUM_R8( volflux, mythid )
177 _GLOBAL_SUM_R8( area,mythid )
178
179 shiftvel(2) = volflux /area
180 print*,'volflux,area',volflux,area
181
182 endif
183 cgg End of the obcsnfirst, obcsnchanged loop.
184
185 #endif
186
187 #endif /* BALANCE_CONTROL_VOLFLUX_GLOBAL */
188
189 return
190 end
191
192
193
194
195
196
197

  ViewVC Help
Powered by ViewVC 1.1.22