/[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.1 - (show annotations) (download)
Thu Nov 6 22:10:07 2003 UTC (20 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint58e_post, checkpoint57v_post, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint58u_post, checkpoint58w_post, checkpoint54a_pre, checkpoint57m_post, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint57s_post, checkpoint54a_post, checkpoint53c_post, checkpoint57k_post, checkpoint55d_pre, checkpoint57d_post, checkpoint57g_post, checkpoint57b_post, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint56b_post, checkpoint57i_post, checkpoint57y_post, hrcube_1, checkpoint57e_post, branch-netcdf, checkpoint52d_pre, checkpoint52l_post, checkpoint55h_post, checkpoint58n_post, checkpoint58x_post, checkpoint52k_post, checkpoint52b_pre, checkpoint57g_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint58t_post, checkpoint58h_post, checkpoint54d_post, checkpoint56c_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55, checkpoint53a_post, checkpoint57f_pre, checkpoint57a_post, checkpoint54, checkpoint58q_post, checkpoint54f_post, checkpoint53b_post, checkpoint55g_post, checkpoint58j_post, checkpoint52a_pre, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint55f_post, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint52d_post, eckpoint57e_pre, checkpoint52a_post, checkpoint57h_done, checkpoint58f_post, checkpoint52b_post, checkpoint53g_post, checkpoint52f_post, checkpoint57x_post, checkpoint57n_post, checkpoint52c_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint57f_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, ecco_c52_e35, hrcube5, checkpoint58o_post, checkpoint57z_post, checkpoint57c_post, checkpoint58y_post, checkpoint55e_post, checkpoint58k_post, checkpoint52i_post, checkpoint52j_pre, checkpoint58v_post, checkpoint53f_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint58s_post, checkpoint58p_post, checkpoint57j_post, checkpoint58b_post, checkpoint57h_pre, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, checkpoint51u_post, checkpoint52h_pre, checkpoint52f_pre, checkpoint57h_post, hrcube_2, hrcube_3, checkpoint56a_post, checkpoint55d_post
Branch point for: netcdf-sm0
o merging from ecco-branch
o pkg/ecco now containes ecco-specific part of cost function
o top level routines the_main_loop, forward_step
  supersede those in model/src/
  previous input data.cost now in data.ecco
  (new namelist ecco_cost_nml)

1
2 #include "COST_CPPOPTIONS.h"
3 #ifdef ALLOW_OBCS
4 # include "OBCS_OPTIONS.h"
5 #endif
6
7 subroutine cost_obcsvol(
8 I myiter,
9 I mytime,
10 I startrec,
11 I endrec,
12 I mythid
13 & )
14
15 c ==================================================================
16 c SUBROUTINE cost_obcsvol
17 c ==================================================================
18 c
19 c o cost function contribution obc -- Volume flux imbalance.
20 c
21 c ==================================================================
22 c SUBROUTINE cost_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 "cal.h"
39 #include "ecco_cost.h"
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 integer startrec
50 integer endrec
51
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 fctile
70 _RL sumvol
71 _RL dummy
72 _RL gg
73 _RL tmpx
74 _RL tmpy
75 _RL wobcsvol
76 character*(80) fnamefldn
77 character*(80) fnameflds
78 character*(80) fnamefldw
79 character*(80) fnameflde
80
81 logical doglobalread
82 logical ladinit
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 c Number of records to be used.
111 nrec = endrec-startrec+1
112
113 #ifdef OBCS_VOLFLUX_COST_CONTRIBUTION
114 #ifdef BAROTROPIC_OBVEL_CONTROL
115
116 sumvol = 0. _d 0
117 wobcsvol = .01
118 cgg Acceptable volume flux is 10^-3. Corresponds to 5 mm change over a year.
119 cgg Added a factor of 1000 because its very important to me.
120 wobcsvol = 1./(wobcsvol * wobcsvol)
121
122 #ifdef ECCO_VERBOSE
123 _BEGIN_MASTER( mythid )
124 write(msgbuf,'(a)') ' '
125 call print_message( msgbuf, standardmessageunit,
126 & SQUEEZE_RIGHT , mythid)
127 write(msgbuf,'(a)') ' '
128 call print_message( msgbuf, standardmessageunit,
129 & SQUEEZE_RIGHT , mythid)
130 write(msgbuf,'(a,i9.8)')
131 & ' cost_obcsvol: number of records to process: ',nrec
132 call print_message( msgbuf, standardmessageunit,
133 & SQUEEZE_RIGHT , mythid)
134 write(msgbuf,'(a)') ' '
135 call print_message( msgbuf, standardmessageunit,
136 & SQUEEZE_RIGHT , mythid)
137 _END_MASTER( mythid )
138 #endif
139
140 if (optimcycle .ge. 0) then
141 #ifdef ALLOW_OBCSN_CONTROL
142 ilfld=ilnblnk( xx_obcsn_file )
143 write(fnamefldn(1:80),'(2a,i10.10)')
144 & xx_obcsn_file(1:ilfld),'.', optimcycle
145 #endif
146 #ifdef ALLOW_OBCSS_CONTROL
147 ilfld=ilnblnk( xx_obcss_file )
148 write(fnameflds(1:80),'(2a,i10.10)')
149 & xx_obcss_file(1:ilfld),'.',optimcycle
150 #endif
151 #ifdef ALLOW_OBCSW_CONTROL
152 ilfld=ilnblnk( xx_obcsw_file )
153 write(fnamefldw(1:80),'(2a,i10.10)')
154 & xx_obcsw_file(1:ilfld),'.',optimcycle
155 #endif
156 #ifdef ALLOW_OBCSE_CONTROL
157 ilfld=ilnblnk( xx_obcse_file )
158 write(fnameflde(1:80),'(2a,i10.10)')
159 & xx_obcse_file(1:ilfld),'.',optimcycle
160 #endif
161 else
162 print*
163 print*,' ctrl_obcsvol: optimcycle has a wrong value.'
164 print*,' optimcycle = ',optimcycle
165 print*
166 stop ' ... stopped in ctrl_obcsvol.'
167 endif
168
169 do irec = 1,nrec
170 c-- Loop over records. For north boundary, we only need V velocity.
171
172 cgg Need to solve for iobcs. Then only keep iobcs=3.or.4.
173 gg = (irec-1)/nobcs
174 igg = int(gg)
175 iobcs = irec - igg*nobcs
176
177 #ifdef ALLOW_OBCSN_CONTROL
178 cgg Assume that nobcs=4, and V velocity is the 4th record. I can't
179 cgg think of a more general way to do this.
180 jp1 = 0
181
182 if (iobcs.eq.4) then
183 call active_read_xz( fnamefldn, tmpfldxz, irec, doglobalread,
184 & ladinit, optimcycle, mythid
185 & , xx_obcsn_dummy )
186
187 cgg At this point, do not be concerned with the overlap halos.
188 cgg From experience, there is no control contribution in the
189 cgg velocity points outside the boundaries. This has something
190 cgg to do with the computational stencil, and the fact that we
191 cgg are diagonally offset. Could check later by employing both
192 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
193 cgg
194 cgg 25-jan-03 --- no idea what i was talking about ^^^^
195 c-- Loop over this thread's tiles.
196 do bj = jtlo,jthi
197 do bi = itlo,ithi
198
199 c-- Determine the weights to be used.
200 fctile = 0. _d 0
201
202 do k = 1, Nr
203 do i = imin,imax
204 j = OB_Jn(I,bi,bj)
205 cgg Barotropic velocity is stored in level 1.
206 tmpx = tmpfldxz(i,1,bi,bj)
207 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
208 cgg -- Positive is flux in.
209 fctile = fctile - tmpx* delZ(k) *dxg(i,j+jp1,bi,bj)
210 endif
211 enddo
212 enddo
213
214 sumvol = sumvol + fctile
215 enddo
216 enddo
217 endif
218 #endif
219
220 #ifdef ALLOW_OBCSS_CONTROL
221 cgg Assume that nobcs=4, and V velocity is the 4th record. I can't
222 cgg think of a more general way to do this.
223 jp1 = 1
224
225 if (iobcs.eq.4) then
226 call active_read_xz( fnameflds, tmpfldxz, irec, doglobalread,
227 & ladinit, optimcycle, mythid
228 & , xx_obcss_dummy )
229
230 cgg At this point, do not be concerned with the overlap halos.
231 cgg From experience, there is no control contribution in the
232 cgg velocity points outside the boundaries. This has something
233 cgg to do with the computational stencil, and the fact that we
234 cgg are diagonally offset. Could check later by employing both
235 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
236
237 c-- Loop over this thread's tiles.
238 do bj = jtlo,jthi
239 do bi = itlo,ithi
240
241 c-- Determine the weights to be used.
242 fctile = 0. _d 0
243
244 do k = 1, Nr
245 do i = imin,imax
246 j = OB_Js(I,bi,bj)
247 cgg Barotropic velocity is stored in level 1.
248 tmpx = tmpfldxz(i,1,bi,bj)
249 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
250 cgg -- Positive is flux in.
251 fctile = fctile + tmpx* delZ(k) *dxg(i,j+jp1,bi,bj)
252 endif
253 enddo
254 enddo
255
256 sumvol = sumvol + fctile
257 enddo
258 enddo
259 endif
260
261 #endif
262
263 #ifdef ALLOW_OBCSW_CONTROL
264 cgg Assume that nobcs=4, and V velocity is the 4th record. I can't
265 cgg think of a more general way to do this.
266 ip1 = 1
267
268 if (iobcs.eq.3) then
269 call active_read_yz( fnamefldw, tmpfldyz, irec, doglobalread,
270 & ladinit, optimcycle, mythid
271 & , xx_obcsw_dummy )
272
273 cgg At this point, do not be concerned with the overlap halos.
274 cgg From experience, there is no control contribution in the
275 cgg velocity points outside the boundaries. This has something
276 cgg to do with the computational stencil, and the fact that we
277 cgg are diagonally offset. Could check later by employing both
278 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
279
280 c-- Loop over this thread's tiles.
281 do bj = jtlo,jthi
282 do bi = itlo,ithi
283
284 c-- Determine the weights to be used.
285 fctile = 0. _d 0
286
287 do k = 1, Nr
288 do j = jmin,jmax
289 i = OB_Iw(j,bi,bj)
290 cgg Barotropic velocity is stored in the level 1.
291 tmpy = tmpfldyz(j,1,bi,bj)
292 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
293 cgg -- Positive is flux in.
294 fctile = fctile + tmpy* delZ(k) *dyg(i+ip1,j,bi,bj)
295 endif
296 enddo
297 enddo
298
299 sumvol = sumvol + fctile
300 enddo
301 enddo
302 endif
303
304 #endif
305
306 #ifdef ALLOW_OBCSE_CONTROL
307 cgg Assume that nobcs=4, and V velocity is the 4th record. I can't
308 cgg think of a more general way to do this.
309 ip1 = 0
310
311 if (iobcs.eq.3) then
312 call active_read_yz( fnameflde, tmpfldyz, irec, doglobalread,
313 & ladinit, optimcycle, mythid
314 & , xx_obcse_dummy )
315
316 cgg At this point, do not be concerned with the overlap halos.
317 cgg From experience, there is no control contribution in the
318 cgg velocity points outside the boundaries. This has something
319 cgg to do with the computational stencil, and the fact that we
320 cgg are diagonally offset. Could check later by employing both
321 cgg BALANCE_CONTROL_VOLFLUX and VOLFLUX_COST_CONTRIBUTION.
322
323 c-- Loop over this thread's tiles.
324 do bj = jtlo,jthi
325 do bi = itlo,ithi
326
327 c-- Determine the weights to be used.
328 fctile = 0. _d 0
329
330 do k = 1, Nr
331 do j = jmin,jmax
332 i = OB_Ie(j,bi,bj)
333 cgg Barotropic velocity stored in level 1.
334 tmpy = tmpfldyz(j,1,bi,bj)
335 if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then
336 cgg -- Positive is flux in.
337 fctile = fctile - tmpy* delZ(k) *dyg(i+ip1,j,bi,bj)
338 endif
339 enddo
340 enddo
341
342 sumvol = sumvol + fctile
343 enddo
344 enddo
345 endif
346
347 #endif
348
349 enddo
350 c-- End of loop over records.
351
352 c-- Do the global summation.
353 _GLOBAL_SUM_R8( sumvol, mythid )
354 objf_obcsvol = wobcsvol * sumvol* sumvol
355
356 #endif
357 #endif
358
359 return
360 end
361
362
363
364
365
366
367
368
369

  ViewVC Help
Powered by ViewVC 1.1.22