/[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.10 - (hide annotations) (download)
Mon Oct 20 03:16:12 2014 UTC (9 years, 6 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.9: +4 -1 lines
- CTRL_OPTIONS.h is needed when including ctrl.h, etc

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

  ViewVC Help
Powered by ViewVC 1.1.22