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

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

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


Revision 1.13 - (hide annotations) (download)
Mon Oct 20 03:16:12 2014 UTC (9 years, 7 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65g, HEAD
Changes since 1.12: +4 -1 lines
- CTRL_OPTIONS.h is needed when including ctrl.h, etc

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

  ViewVC Help
Powered by ViewVC 1.1.22