/[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.8 - (hide annotations) (download)
Tue Sep 18 20:16:34 2012 UTC (11 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint64, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65d, checkpoint65e
Changes since 1.7: +25 -20 lines
use new parameter OB_indexNone for null index value (instead of hard-coded 0)

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

  ViewVC Help
Powered by ViewVC 1.1.22