/[MITgcm]/MITgcm/pkg/ecco/cost_obcsvol.F
ViewVC logotype

Contents of /MITgcm/pkg/ecco/cost_obcsvol.F

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


Revision 1.11 - (show annotations) (download)
Tue Sep 18 20:16:34 2012 UTC (11 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint64, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65d, checkpoint65e
Changes since 1.10: +5 -3 lines
use new parameter OB_indexNone for null index value (instead of hard-coded 0)

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_obcsvol.F,v 1.10 2012/08/10 19:45:26 jmc Exp $
2 C $Name: $
3
4 #include "ECCO_OPTIONS.h"
5
6 subroutine cost_obcsvol(
7 I myiter,
8 I mytime,
9 I startrec,
10 I endrec,
11 I mythid
12 & )
13
14 c ==================================================================
15 c SUBROUTINE cost_obcsvol
16 c ==================================================================
17 c
18 c o cost function contribution obc -- Volume flux imbalance.
19 c
20 c ==================================================================
21 c SUBROUTINE cost_obcsvol
22 c ==================================================================
23
24 implicit none
25
26 c == global variables ==
27
28 #include "EEPARAMS.h"
29 #include "SIZE.h"
30 #include "PARAMS.h"
31 #include "GRID.h"
32 #ifdef ALLOW_OBCS
33 # include "OBCS_GRID.h"
34 #endif
35
36 #include "cal.h"
37 #include "ecco_cost.h"
38 #include "CTRL_SIZE.h"
39 #include "ctrl.h"
40 #include "ctrl_dummy.h"
41 #include "optim.h"
42
43 c == routine arguments ==
44
45 integer myiter
46 _RL mytime
47 integer mythid
48 integer startrec
49 integer endrec
50
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 iobcs
61 integer nrec
62 integer ilfld
63 integer igg
64
65 _RL fctile
66 _RL sumvol
67 _RL gg
68 _RL tmpx
69 _RL tmpy
70 _RL wobcsvol
71 character*(80) fnamefldn
72 character*(80) fnameflds
73 character*(80) fnamefldw
74 character*(80) fnameflde
75
76 #if (defined ALLOW_OBCSN_CONTROL || defined ALLOW_OBCSS_CONTROL)
77 _RL tmpfldxz (1-olx:snx+olx,nr,nsx,nsy)
78 #endif
79 #if (defined ALLOW_OBCSE_CONTROL || defined ALLOW_OBCSW_CONTROL)
80 _RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy)
81 #endif
82
83 logical doglobalread
84 logical ladinit
85
86 #ifdef ECCO_VERBOSE
87 character*(MAX_LEN_MBUF) msgbuf
88 #endif
89
90 c == external functions ==
91
92 integer ilnblnk
93 external ilnblnk
94
95 c == end of interface ==
96
97 #ifdef OBCS_VOLFLUX_COST_CONTRIBUTION
98 #ifdef BAROTROPIC_OBVEL_CONTROL
99
100 stop 's/r cost_obcsvol needs to be fixed'
101
102 jtlo = mybylo(mythid)
103 jthi = mybyhi(mythid)
104 itlo = mybxlo(mythid)
105 ithi = mybxhi(mythid)
106 jmin = 1
107 jmax = sny
108 imin = 1
109 imax = snx
110
111 c-- Read tiled data.
112 doglobalread = .false.
113 ladinit = .false.
114
115 cgg Assume the number of records is the same for
116 cgg all boundaries.
117 c Number of records to be used.
118 nrec = endrec-startrec+1
119
120 sumvol = 0. _d 0
121 wobcsvol = .01
122 cgg Acceptable volume flux is 10^-3. Corresponds to 5 mm change over a year.
123 cgg Added a factor of 1000 because its very important to me.
124 wobcsvol = 1./(wobcsvol * wobcsvol)
125
126 #ifdef ECCO_VERBOSE
127 _BEGIN_MASTER( mythid )
128 write(msgbuf,'(a)') ' '
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 write(msgbuf,'(a,i9.8)')
135 & ' cost_obcsvol: number of records to process: ',nrec
136 call print_message( msgbuf, standardmessageunit,
137 & SQUEEZE_RIGHT , mythid)
138 write(msgbuf,'(a)') ' '
139 call print_message( msgbuf, standardmessageunit,
140 & SQUEEZE_RIGHT , mythid)
141 _END_MASTER( mythid )
142 #endif
143
144 if (optimcycle .ge. 0) then
145 #ifdef ALLOW_OBCSN_CONTROL
146 ilfld=ilnblnk( xx_obcsn_file )
147 write(fnamefldn(1:80),'(2a,i10.10)')
148 & xx_obcsn_file(1:ilfld),'.', optimcycle
149 #endif
150 #ifdef ALLOW_OBCSS_CONTROL
151 ilfld=ilnblnk( xx_obcss_file )
152 write(fnameflds(1:80),'(2a,i10.10)')
153 & xx_obcss_file(1:ilfld),'.',optimcycle
154 #endif
155 #ifdef ALLOW_OBCSW_CONTROL
156 ilfld=ilnblnk( xx_obcsw_file )
157 write(fnamefldw(1:80),'(2a,i10.10)')
158 & xx_obcsw_file(1:ilfld),'.',optimcycle
159 #endif
160 #ifdef ALLOW_OBCSE_CONTROL
161 ilfld=ilnblnk( xx_obcse_file )
162 write(fnameflde(1:80),'(2a,i10.10)')
163 & xx_obcse_file(1:ilfld),'.',optimcycle
164 #endif
165 else
166 print*
167 print*,' obcs_obcsvol: optimcycle has a wrong value.'
168 print*,' optimcycle = ',optimcycle
169 print*
170 stop ' ... stopped in obcs_obcsvol.'
171 endif
172
173 do irec = 1,nrec
174 c-- Loop over records. For north boundary, we only need V velocity.
175
176 cgg Need to solve for iobcs. Then only keep iobcs=3.or.4.
177 gg = (irec-1)/nobcs
178 igg = int(gg)
179 iobcs = irec - igg*nobcs
180
181 #ifdef ALLOW_OBCSN_CONTROL
182 cgg Assume that nobcs=4, and V velocity is the 4th record. I cannot
183 cgg think of a more general way to do this.
184 if (iobcs.eq.4) then
185 call active_read_xz( fnamefldn, tmpfldxz, irec, doglobalread,
186 & ladinit, optimcycle, mythid
187 & , xx_obcsn_dummy )
188
189 cgg At this point, do not be concerned with the overlap halos.
190 cgg From experience, there is no control contribution in the
191 cgg velocity points outside the boundaries. This has something
192 cgg to do with the computational stencil, and the fact that we
193 cgg are diagonally offset. Could check later by employing both
194 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
195 cgg
196 cgg 25-jan-03 --- no idea what i was talking about ^^^^
197 c-- Loop over this thread tiles.
198 do bj = jtlo,jthi
199 do bi = itlo,ithi
200
201 c-- Determine the weights to be used.
202 fctile = 0. _d 0
203
204 do k = 1, Nr
205 do i = imin,imax
206 j = OB_Jn(I,bi,bj)
207 IF ( j.EQ.OB_indexNone ) j = 1
208 cgg Barotropic velocity is stored in level 1.
209 tmpx = tmpfldxz(i,1,bi,bj)
210 if (maskS(i,j,k,bi,bj) .ne. 0.) then
211 cgg -- Positive is flux in.
212 fctile = fctile - tmpx* drF(k) *dxg(i,j,bi,bj)
213 & * _hFacS(i,j,k,bi,bj)
214 endif
215 enddo
216 enddo
217
218 sumvol = sumvol + fctile
219 enddo
220 enddo
221 endif
222 #endif
223
224 #ifdef ALLOW_OBCSS_CONTROL
225 cgg Assume that nobcs=4, and V velocity is the 4th record. I cannot
226 cgg think of a more general way to do this.
227 if (iobcs.eq.4) then
228 call active_read_xz( fnameflds, tmpfldxz, irec, doglobalread,
229 & ladinit, optimcycle, mythid
230 & , xx_obcss_dummy )
231
232 cgg At this point, do not be concerned with the overlap halos.
233 cgg From experience, there is no control contribution in the
234 cgg velocity points outside the boundaries. This has something
235 cgg to do with the computational stencil, and the fact that we
236 cgg are diagonally offset. Could check later by employing both
237 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
238
239 c-- Loop over this thread tiles.
240 do bj = jtlo,jthi
241 do bi = itlo,ithi
242
243 c-- Determine the weights to be used.
244 fctile = 0. _d 0
245
246 do k = 1, Nr
247 do i = imin,imax
248 j = OB_Js(I,bi,bj)
249 IF ( j.EQ.OB_indexNone ) j = 1
250 cgg Barotropic velocity is stored in level 1.
251 tmpx = tmpfldxz(i,1,bi,bj)
252 if (maskS(i,j+1,k,bi,bj) .ne. 0.) then
253 cgg -- Positive is flux in.
254 fctile = fctile + tmpx* drF(k) *dxg(i,j+1,bi,bj)
255 & * _hFacS(i,j+1,k,bi,bj)
256 endif
257 enddo
258 enddo
259
260 sumvol = sumvol + fctile
261 enddo
262 enddo
263 endif
264
265 #endif
266
267 #ifdef ALLOW_OBCSW_CONTROL
268 cgg Assume that nobcs=4, and V velocity is the 4th record. I cannot
269 cgg think of a more general way to do this.
270 if (iobcs.eq.3) then
271 call active_read_yz( fnamefldw, tmpfldyz, irec, doglobalread,
272 & ladinit, optimcycle, mythid
273 & , xx_obcsw_dummy )
274
275 cgg At this point, do not be concerned with the overlap halos.
276 cgg From experience, there is no control contribution in the
277 cgg velocity points outside the boundaries. This has something
278 cgg to do with the computational stencil, and the fact that we
279 cgg are diagonally offset. Could check later by employing both
280 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
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 fctile = 0. _d 0
288
289 do k = 1, Nr
290 do j = jmin,jmax
291 i = OB_Iw(j,bi,bj)
292 IF ( i.EQ.OB_indexNone ) i = 1
293 cgg Barotropic velocity is stored in the level 1.
294 tmpy = tmpfldyz(j,1,bi,bj)
295 if (maskW(i+1,j,k,bi,bj) .ne. 0.) then
296 cgg -- Positive is flux in.
297 fctile = fctile + tmpy* drF(k) *dyg(i+1,j,bi,bj)
298 & * _hFacW(i+1,j,k,bi,bj)
299 endif
300 enddo
301 enddo
302
303 sumvol = sumvol + fctile
304 enddo
305 enddo
306 endif
307
308 #endif
309
310 #ifdef ALLOW_OBCSE_CONTROL
311 cgg Assume that nobcs=4, and V velocity is the 4th record. I cannot
312 cgg think of a more general way to do this.
313 if (iobcs.eq.3) then
314 call active_read_yz( fnameflde, tmpfldyz, irec, doglobalread,
315 & ladinit, optimcycle, mythid
316 & , xx_obcse_dummy )
317
318 cgg At this point, do not be concerned with the overlap halos.
319 cgg From experience, there is no control contribution in the
320 cgg velocity points outside the boundaries. This has something
321 cgg to do with the computational stencil, and the fact that we
322 cgg are diagonally offset. Could check later by employing both
323 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
324
325 c-- Loop over this thread tiles.
326 do bj = jtlo,jthi
327 do bi = itlo,ithi
328
329 c-- Determine the weights to be used.
330 fctile = 0. _d 0
331
332 do k = 1, Nr
333 do j = jmin,jmax
334 i = OB_Ie(j,bi,bj)
335 IF ( i.EQ.OB_indexNone ) i = 1
336 cgg Barotropic velocity stored in level 1.
337 tmpy = tmpfldyz(j,1,bi,bj)
338 if (maskW(i,j,k,bi,bj) .ne. 0.) then
339 cgg -- Positive is flux in.
340 fctile = fctile - tmpy* drF(k) *dyg(i,j,bi,bj)
341 & * _hFacW(i,j,k,bi,bj)
342 endif
343 enddo
344 enddo
345
346 sumvol = sumvol + fctile
347 enddo
348 enddo
349 endif
350
351 #endif
352
353 enddo
354 c-- End of loop over records.
355
356 c-- Do the global summation.
357 _GLOBAL_SUM_RL( sumvol, mythid )
358 objf_obcsvol = wobcsvol * sumvol* sumvol
359
360 #endif
361 #endif
362
363 return
364 end

  ViewVC Help
Powered by ViewVC 1.1.22