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