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

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

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


Revision 1.3 - (show annotations) (download)
Mon Aug 11 22:29:55 2008 UTC (15 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61c, checkpoint61l, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i
Changes since 1.2: +30 -9 lines
- replace calls to "FIND_RHO" by new version "FIND_RHO_2D"
- stop if more than 1 tile / proc (Pb with EXCH & FIND_RHO_2D calls)

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_obcs_ageos.F,v 1.2 2007/10/09 00:02:50 jmc Exp $
2 C $Name: $
3
4 #include "COST_CPPOPTIONS.h"
5 #ifdef ALLOW_OBCS
6 # include "OBCS_OPTIONS.h"
7 #endif
8
9 subroutine cost_obcs_ageos( myiter, mytime, mythid )
10
11 c ==================================================================
12 c SUBROUTINE cost_obcs_ageos
13 c ==================================================================
14 c
15 c o cost function contribution obc -- Ageostrophic boundary flow.
16 c
17 c started: G. Gebbie gebbie@mit.edu 4-Feb-2003
18 c
19 c warning: masks redundantly assume no gradient of topography at
20 c boundary.
21 c
22 c ==================================================================
23 c SUBROUTINE cost_obcs_ageos
24 c ==================================================================
25
26 implicit none
27
28 c == global variables ==
29
30 #include "EEPARAMS.h"
31 #include "SIZE.h"
32 #include "GRID.h"
33 #include "DYNVARS.h"
34 #include "PARAMS.h"
35 #ifdef ALLOW_OBCS
36 # include "OBCS.h"
37 #endif
38
39 #include "cal.h"
40 #include "ecco_cost.h"
41 #include "ctrl.h"
42 #include "ctrl_dummy.h"
43 #include "optim.h"
44
45 c == routine arguments ==
46
47 integer myiter
48 _RL mytime
49 integer mythid
50
51 c == local variables ==
52
53 _RS one_rs
54 parameter( one_rs = 1. )
55
56 integer bi,bj
57 integer i,j,k
58 integer itlo,ithi
59 integer jtlo,jthi
60 integer jmin,jmax
61 integer imin,imax
62 integer irec
63 integer levmon
64 integer levoff
65 integer iltheta
66 integer ilsalt
67 integer iluvel
68 integer ilvvel
69 integer ip1, jp1
70
71 _RL fctile
72 _RL fcthread
73
74 _RL rholoc (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
75 _RL xzgrdrho(1-olx:snx+olx,Nr,nsx,nsy)
76 _RL yzgrdrho(1-oly:sny+oly,Nr,nsx,nsy)
77 _RL xzdvel1 (1-olx:snx+olx,nr,nsx,nsy)
78 _RL xzdvel2 (1-olx:snx+olx,nr,nsx,nsy)
79 _RL yzdvel1 (1-oly:sny+oly,nr,nsx,nsy)
80 _RL yzdvel2 (1-oly:sny+oly,nr,nsx,nsy)
81 _RL maskxzageos (1-olx:snx+olx,nr,nsx,nsy)
82 _RL maskyzageos (1-oly:sny+oly,nr,nsx,nsy)
83 _RL dummy
84
85 character*(80) fnametheta
86 character*(80) fnamesalt
87 character*(80) fnameuvel
88 character*(80) fnamevvel
89
90 logical doglobalread
91 logical ladinit
92
93 character*(MAX_LEN_MBUF) msgbuf
94
95 c == external functions ==
96
97 integer ilnblnk
98 external ilnblnk
99
100 c == end of interface ==
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 #ifdef OBCS_AGEOS_COST_CONTRIBUTION
116
117 #ifdef ECCO_VERBOSE
118 _BEGIN_MASTER( mythid )
119 write(msgbuf,'(a)') ' '
120 call print_message( msgbuf, standardmessageunit,
121 & SQUEEZE_RIGHT , mythid)
122 write(msgbuf,'(a,i8.8)')
123 & ' cost_obcs_ageos: number of records to process = ',nmonsrec
124 call print_message( msgbuf, standardmessageunit,
125 & SQUEEZE_RIGHT , mythid)
126 write(msgbuf,'(a)') ' '
127 call print_message( msgbuf, standardmessageunit,
128 & SQUEEZE_RIGHT , mythid)
129 _END_MASTER( mythid )
130 #endif
131
132 if (optimcycle .ge. 0) then
133 iltheta = ilnblnk( tbarfile )
134 write(fnametheta(1:80),'(2a,i10.10)')
135 & tbarfile(1:iltheta),'.',optimcycle
136 ilsalt = ilnblnk( sbarfile )
137 write(fnamesalt(1:80),'(2a,i10.10)')
138 & sbarfile(1:ilsalt),'.',optimcycle
139 iluvel = ilnblnk( ubarfile )
140 write(fnameuvel(1:80),'(2a,i10.10)')
141 & ubarfile(1:iluvel),'.',optimcycle
142 ilvvel = ilnblnk( vbarfile )
143 write(fnamevvel(1:80),'(2a,i10.10)')
144 & vbarfile(1:ilvvel),'.',optimcycle
145 endif
146
147 fcthread = 0. _d 0
148 fctile = 0. _d 0
149
150 cgg Code safe: always initialize to zero.
151 do bj = jtlo,jthi
152 do bi = itlo,ithi
153 do k = 1,nr
154 do i = 1-olx,snx+olx
155 maskxzageos(i,k,bi,bj) = 0. _d 0
156 xzdvel1(i,k,bi,bj) = 0. _d 0
157 xzdvel2(i,k,bi,bj) = 0. _d 0
158 xzgrdrho(i,k,bi,bj) = 0. _d 0
159 enddo
160 do j = 1-oly,sny+oly
161 maskyzageos(j,k,bi,bj) = 0. _d 0
162 yzdvel1(j,k,bi,bj) = 0. _d 0
163 yzdvel2(j,k,bi,bj) = 0. _d 0
164 yzgrdrho(j,k,bi,bj) = 0. _d 0
165 enddo
166 enddo
167 enddo
168 enddo
169
170 do bj = jtlo,jthi
171 do bi = itlo,ithi
172 do j = 1-oly,sny+oly
173 do i = 1-olx,snx+olx
174 rholoc(i,j,bi,bj) = 0. _d 0
175 enddo
176 enddo
177 enddo
178 enddo
179
180 c-- Loop over records.
181 do irec = 1,nmonsrec
182
183 c-- Read time averages and the monthly mean data.
184 call active_read_xyz( fnametheta, tbar, irec,
185 & doglobalread, ladinit,
186 & optimcycle, mythid,
187 & xx_tbar_mean_dummy )
188
189 c-- Read time averages and the monthly mean data.
190 call active_read_xyz( fnamesalt, sbar, irec,
191 & doglobalread, ladinit,
192 & optimcycle, mythid,
193 & xx_sbar_mean_dummy )
194
195 c-- Read time averages and the monthly mean data.
196 call active_read_xyz( fnameuvel, ubar, irec,
197 & doglobalread, ladinit,
198 & optimcycle, mythid,
199 & xx_ubar_mean_dummy )
200
201 c-- Read time averages and the monthly mean data.
202 call active_read_xyz( fnamevvel, vbar, irec,
203 & doglobalread, ladinit,
204 & optimcycle, mythid,
205 & xx_vbar_mean_dummy )
206
207 cgg Minor problem : grad T,S needs overlap values.
208 _BARRIER
209 _EXCH_XYZ_R8(tbar,myThid)
210 _EXCH_XYZ_R8(sbar,myThid)
211
212
213 do bj = jtlo,jthi
214 do bi = itlo,ithi
215
216 #ifdef ALLOW_OBCSN_CONTROL
217 jp1 = 0
218
219 cgg Make a mask for the velocity shear comparison.
220 do k = 1,nr-1
221 do i = imin, imax
222 j = ob_jn(i,bi,bj)
223 cgg All these points need to be wet.
224 if (j .eq. 0) then
225 maskxzageos(i,k,bi,bj) = 0.
226 else
227 maskxzageos(i,k,bi,bj) =
228 & hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
229 & hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
230 & hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
231 & hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
232 endif
233 enddo
234 enddo
235
236 do k = 1,nr
237
238 C-- jmc: both calls below are wrong if more than 1 tile => stop here
239 IF ( bi.NE.1 .OR. bj.NE.1 )
240 & STOP 'COST_OBCS_AGEOS wrong with more than 1 tile/proc'
241 call find_rho_2d(
242 I iMin, iMax, jMin, jMax, k,
243 I tbar(1-OLx,1-OLy,k,bi,bj),
244 I sbar(1-OLx,1-OLy,k,bi,bj),
245 O rholoc,
246 I k, bi, bj, myThid )
247 _EXCH_XY_R8(rholoc , myThid)
248
249 cgg Compute centered difference horizontal gradient on bdy.
250 do i = imin, imax
251 j = ob_jn(i,bi,bj)
252 xzgrdrho(i,k,bi,bj) =
253 & (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
254 & /(2.*dxc(i,j+jp1,bi,bj))
255 enddo
256 enddo
257
258 cgg Compute vertical shear from geostrophy/thermal wind.
259 cgg Above level 4 needs not be geostrophic.
260 cgg Please get rid of the "4" ASAP. Ridiculous.
261 do k = 4,Nr-1
262 do i = imin,imax
263 j = ob_jn(i,bi,bj)
264 xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
265 & - vbar(i,j+jp1,k+1,bi,bj)
266 xzdvel2(i,k,bi,bj)=((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
267 & (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
268 & * gravity / f0 / rhonil
269
270 fctile = fctile + 100*wcurrent(k,bi,bj) *
271 & maskxzageos(i,k,bi,bj)*
272 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
273 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
274 if (maskxzageos(i,k,bi,bj) .ne. 0) then
275 cgg print*,'XX shear,geos shear',xzdvel1(i,k,bi,bj),xzdvel2(i,k,bi,bj)
276 cgg print*,'i,j,k,fctile N',i,j,k,fctile
277 endif
278 enddo
279 enddo
280 c-- End of loop over layers.
281 #endif /* ALLOW_OBCSN_CONTROL */
282
283 #ifdef ALLOW_OBCSS_CONTROL
284 jp1 = 1
285
286 cgg Make a mask for the velocity shear comparison.
287 do k = 1,nr-1
288 do i = imin, imax
289 j = ob_js(i,bi,bj)
290 if (j .eq. 0) then
291 maskxzageos(i,k,bi,bj) = 0.
292 else
293 cgg All these points need to be wet.
294 maskxzageos(i,k,bi,bj) =
295 & hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
296 & hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
297 & hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
298 & hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
299 endif
300 enddo
301 enddo
302
303 do k = 1,nr
304
305 call find_rho_2d(
306 I iMin, iMax, jMin, jMax, k,
307 I tbar(1-OLx,1-OLy,k,bi,bj),
308 I sbar(1-OLx,1-OLy,k,bi,bj),
309 O rholoc,
310 I k, bi, bj, myThid )
311
312 _EXCH_XY_R8(rholoc , myThid)
313
314 cgg Compute centered difference horizontal gradient on bdy.
315 do i = imin, imax
316 j = ob_js(i,bi,bj)
317 xzgrdrho(i,k,bi,bj) =
318 & (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
319 & /(2.*dxc(i,j+jp1,bi,bj))
320 enddo
321 enddo
322
323 cgg Compute vertical shear from geostrophy/thermal wind.
324 do k = 4,Nr-1
325 do i = imin,imax
326 j = ob_js(i,bi,bj)
327 cgg Retrieve the model vertical shear.
328 xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
329 & - vbar(i,j+jp1,k+1,bi,bj)
330
331 cgg Compute vertical shear from geostrophy/thermal wind.
332 xzdvel2(i,k,bi,bj) =((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
333 & (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
334 & * gravity / f0 /rhonil
335
336 cgg Make a comparison.
337 fctile = fctile + 100*wcurrent(k,bi,bj) *
338 & maskxzageos(i,k,bi,bj)*
339 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
340 & (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
341 cgg print*,'fctile S',fctile
342 enddo
343 enddo
344 c-- End of loop over layers.
345 #endif
346
347 #ifdef ALLOW_OBCSW_CONTROL
348 ip1 = 1
349
350 cgg Make a mask for the velocity shear comparison.
351 do k = 1,nr-1
352 do j = jmin, jmax
353 i = ob_iw(j,bi,bj)
354 cgg All these points need to be wet.
355 if (i.eq. 0) then
356 maskyzageos(j,k,bi,bj) = 0.
357 else
358 maskyzageos(j,k,bi,bj) =
359 & hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
360 & hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
361 & hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
362 & hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
363 endif
364 enddo
365 enddo
366
367 do k = 1,nr
368
369 call find_rho_2d(
370 I iMin, iMax, jMin, jMax, k,
371 I tbar(1-OLx,1-OLy,k,bi,bj),
372 I sbar(1-OLx,1-OLy,k,bi,bj),
373 O rholoc,
374 I k, bi, bj, myThid )
375 _EXCH_XY_R8(rholoc , myThid)
376
377 cgg Compute centered difference horizontal gradient on bdy.
378 do j = jmin, jmax
379 i = ob_iw(j,bi,bj)
380 cgg Negative sign due to geostrophy.
381 yzgrdrho(j,k,bi,bj) =
382 & (rholoc(i+ip1,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
383 & /(2.*dyc(i+ip1,j,bi,bj))
384 enddo
385 enddo
386
387 cgg Compute vertical shear from geostrophy/thermal wind.
388 do k = 4,Nr-1
389 do j = jmin,jmax
390 i = ob_iw(j,bi,bj)
391 cgg Retrieve the model vertical shear.
392 yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
393 & - ubar(i+ip1,j,k+1,bi,bj)
394
395 cgg Compute vertical shear from geostrophy/thermal wind.
396 yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
397 & (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
398 & * gravity / f0 / rhonil
399
400 cgg Make a comparison.
401 fctile = fctile + 100*wcurrent(k,bi,bj) *
402 & maskyzageos(j,k,bi,bj) *
403 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
404 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
405 enddo
406 enddo
407 c-- End of loop over layers.
408 #endif
409
410 #ifdef ALLOW_OBCSE_CONTROL
411 ip1 = 0
412
413 cgg Make a mask for the velocity shear comparison.
414 do k = 1,nr-1
415 do j = jmin, jmax
416 i = ob_ie(j,bi,bj)
417 if (i.eq.0) then
418 maskyzageos(j,k,bi,bj) =0.
419 else
420 cgg All these points need to be wet.
421 maskyzageos(j,k,bi,bj) =
422 & hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
423 & hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
424 & hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
425 & hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
426 endif
427 enddo
428 enddo
429
430 do k = 1,nr
431
432 call find_rho_2d(
433 I iMin, iMax, jMin, jMax, k,
434 I tbar(1-OLx,1-OLy,k,bi,bj),
435 I sbar(1-OLx,1-OLy,k,bi,bj),
436 O rholoc,
437 I k, bi, bj, myThid )
438 _EXCH_XY_R8(rholoc , myThid)
439
440 cgg Compute centered difference horizontal gradient on bdy.
441 do j = jmin, jmax
442 i = ob_ie(j,bi,bj)
443 cgg Negative sign due to geostrophy.
444 yzgrdrho(j,k,bi,bj) =
445 & (rholoc(i+ip1,,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
446 & /(2.*dyc(i+ip1,j,bi,bj))
447 enddo
448 enddo
449
450 cgg Compute vertical shear from geostrophy/thermal wind.
451 do k = 4,Nr-1
452 do j = jmin,jmax
453 i = ob_ie(j,bi,bj)
454 cgg Retrieve the model vertical shear.
455 yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
456 & - ubar(i+ip1,j,k+1,bi,bj)
457
458 cgg Compute vertical shear from geostrophy/thermal wind.
459 yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
460 & (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
461 & * gravity / f0 /rhonil
462
463 cgg Make a comparison.
464 fctile = fctile + 100*wcurrent(k,bi,bj) *
465 & maskyzageos(j,k,bi,bj) *
466 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
467 & (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
468
469 enddo
470 enddo
471 c-- End of loop over layers.
472 #endif
473
474 fcthread = fcthread + fctile
475 objf_ageos(bi,bj) = objf_ageos(bi,bj) + fctile
476
477 #ifdef ECCO_VERBOSE
478 c-- Print cost function for each tile in each thread.
479 write(msgbuf,'(a)') ' '
480 call print_message( msgbuf, standardmessageunit,
481 & SQUEEZE_RIGHT , mythid)
482 write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
483 & ' cost_Theta: irec,bi,bj = ',irec,bi,bj
484 call print_message( msgbuf, standardmessageunit,
485 & SQUEEZE_RIGHT , mythid)
486 write(msgbuf,'(a,d22.15)')
487 & ' cost function (temperature) = ',
488 & fctile
489 call print_message( msgbuf, standardmessageunit,
490 & SQUEEZE_RIGHT , mythid)
491 write(msgbuf,'(a)') ' '
492 call print_message( msgbuf, standardmessageunit,
493 & SQUEEZE_RIGHT , mythid)
494 #endif
495
496 enddo
497 enddo
498
499 #ifdef ECCO_VERBOSE
500 c-- Print cost function for all tiles.
501 _GLOBAL_SUM_R8( fcthread , myThid )
502 write(msgbuf,'(a)') ' '
503 call print_message( msgbuf, standardmessageunit,
504 & SQUEEZE_RIGHT , mythid)
505 write(msgbuf,'(a,i8.8)')
506 & ' cost_Theta: irec = ',irec
507 call print_message( msgbuf, standardmessageunit,
508 & SQUEEZE_RIGHT , mythid)
509 write(msgbuf,'(a,a,d22.15)')
510 & ' global cost function value',
511 & ' (temperature) = ',fcthread
512 call print_message( msgbuf, standardmessageunit,
513 & SQUEEZE_RIGHT , mythid)
514 write(msgbuf,'(a)') ' '
515 call print_message( msgbuf, standardmessageunit,
516 & SQUEEZE_RIGHT , mythid)
517 #endif
518
519 enddo
520 c-- End of loop over records.
521
522 #endif
523
524 return
525 end
526
527
528

  ViewVC Help
Powered by ViewVC 1.1.22