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

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

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


Revision 1.7 - (hide annotations) (download)
Fri Aug 10 19:45:26 2012 UTC (11 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63r, checkpoint63s
Changes since 1.6: +2 -2 lines
include ECCO_OPTIONS.h instead of COST_CPPOPTIONS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22