/[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.9 - (hide annotations) (download)
Thu Oct 9 00:50:16 2014 UTC (9 years, 7 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65f
Changes since 1.8: +20 -7 lines
- pkg/ecco/ecco_cost.h : rm obcs ctrl variables (now all in CTRL_OBCS.h).
- pkg/ecco/cost_obcs.F, cost_obcs_ageos.F, cost_obcse.F, cost_obcsw.F,
  cost_obcsn.F, cost_obcss.F, cost_obcsvol.F, ecco_cost_init_varia.F,
  ecco_cost_weights.F, ecco_readparms.F, ecco_cost_final.F : add
  CPP brackets and CTRL_OBCS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22