1 |
C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_averagesfields.F,v 1.1.2.8 2003/06/20 05:19:17 heimbach Exp $ |
2 |
|
3 |
#include "COST_CPPOPTIONS.h" |
4 |
#ifdef ALLOW_OBCS |
5 |
# include "OBCS_OPTIONS.h" |
6 |
#endif |
7 |
|
8 |
|
9 |
subroutine cost_averagesfields( mytime, mythid ) |
10 |
|
11 |
c ================================================================== |
12 |
c SUBROUTINE cost_averagesfields |
13 |
c ================================================================== |
14 |
c |
15 |
c o Compute time averages of etaN, theta, and salt. The counters |
16 |
c are explicitly calculated instead of being incremented. This |
17 |
c reduces dependencies. The latter is useful for the adjoint code |
18 |
c generation. |
19 |
c |
20 |
c started: Christian Eckert eckert@mit.edu 30-Jun-1999 |
21 |
c |
22 |
c changed: Christian Eckert eckert@mit.edu 24-Feb-2000 |
23 |
c |
24 |
c - Restructured the code in order to create a package |
25 |
c for the MITgcmUV. |
26 |
c |
27 |
c heimbach@mit.edu Many changes. |
28 |
c |
29 |
c ================================================================== |
30 |
c SUBROUTINE cost_averagesfields |
31 |
c ================================================================== |
32 |
|
33 |
implicit none |
34 |
|
35 |
c == global variables == |
36 |
|
37 |
#include "EEPARAMS.h" |
38 |
#include "SIZE.h" |
39 |
#include "PARAMS.h" |
40 |
#include "DYNVARS.h" |
41 |
#include "CG2D.h" |
42 |
#ifdef ALLOW_DRIFTW_COST_CONTRIBUTION |
43 |
#include "GW.h" |
44 |
#endif |
45 |
|
46 |
#include "exf_fields.h" |
47 |
#include "optim.h" |
48 |
#include "ecco_cost.h" |
49 |
#include "ctrl_dummy.h" |
50 |
|
51 |
c == routine arguments == |
52 |
|
53 |
_RL mytime |
54 |
cph integer myiter |
55 |
integer mythid |
56 |
|
57 |
c == local variables == |
58 |
|
59 |
integer myiter |
60 |
integer bi,bj |
61 |
integer i,j,k |
62 |
integer itlo,ithi |
63 |
integer jtlo,jthi |
64 |
integer jmin,jmax |
65 |
integer imin,imax |
66 |
|
67 |
logical first |
68 |
logical startofday |
69 |
logical startofmonth |
70 |
logical inday |
71 |
logical inmonth |
72 |
logical last |
73 |
logical endofday |
74 |
logical endofmonth |
75 |
|
76 |
integer ilps, ils,ilt |
77 |
|
78 |
character*(128) fnamepsbar |
79 |
character*(128) fnametbar |
80 |
character*(128) fnamesbar |
81 |
character*(128) fnameubar |
82 |
character*(128) fnamevbar |
83 |
character*(128) fnamewbar |
84 |
character*(128) fnametauxbar |
85 |
character*(128) fnametauybar |
86 |
character*(128) fnamehfluxbar |
87 |
character*(128) fnamesfluxbar |
88 |
|
89 |
c == external functions == |
90 |
|
91 |
integer ilnblnk |
92 |
external ilnblnk |
93 |
|
94 |
c == end of interface == |
95 |
|
96 |
jtlo = mybylo(mythid) |
97 |
jthi = mybyhi(mythid) |
98 |
itlo = mybxlo(mythid) |
99 |
ithi = mybxhi(mythid) |
100 |
jmin = 1 |
101 |
jmax = sny |
102 |
imin = 1 |
103 |
imax = snx |
104 |
|
105 |
myiter = niter0 + INT((mytime-starttime)/deltaTClock+0.5) |
106 |
|
107 |
c-- Get the time flags and record numbers for the time averaging. |
108 |
|
109 |
call cost_averagesflags( |
110 |
I myiter, mytime, mythid, |
111 |
O first, startofday, startofmonth, |
112 |
O inday, inmonth, |
113 |
O last, endofday, endofmonth, |
114 |
O sum1day, dayrec, |
115 |
O sum1mon, monrec |
116 |
& ) |
117 |
|
118 |
#ifdef ALLOW_SSH_COST_CONTRIBUTION |
119 |
c-- First, do the daily averages. |
120 |
if (first .or. startofday) then |
121 |
c-- Assign the first value to the array holding the average. |
122 |
do bj = jtlo,jthi |
123 |
do bi = itlo,ithi |
124 |
do j = jmin,jmax |
125 |
do i = imin,imax |
126 |
psbar(i,j,bi,bj) = etaN(i,j,bi,bj) |
127 |
enddo |
128 |
enddo |
129 |
enddo |
130 |
enddo |
131 |
else if (last .or. endofday) then |
132 |
c-- Add the last value and devide by the number of accumulated |
133 |
c-- records. |
134 |
do bj = jtlo,jthi |
135 |
do bi = itlo,ithi |
136 |
do j = jmin,jmax |
137 |
do i = imin,imax |
138 |
psbar(i,j,bi,bj) = (psbar (i,j,bi,bj) + |
139 |
& etaN(i,j,bi,bj) )/ |
140 |
& float(sum1day) |
141 |
enddo |
142 |
enddo |
143 |
enddo |
144 |
enddo |
145 |
|
146 |
c-- Save psbar on file. |
147 |
if (optimcycle .ge. 0) then |
148 |
ilps=ilnblnk( psbarfile ) |
149 |
write(fnamepsbar,'(2a,i10.10)') |
150 |
& psbarfile(1:ilps), '.', optimcycle |
151 |
endif |
152 |
|
153 |
call active_write_xy( fnamepsbar, psbar, dayrec, optimcycle, |
154 |
& mythid, xx_psbar_mean_dummy ) |
155 |
|
156 |
else if ( ( inday ) .and. |
157 |
& .not. (first .or. startofday) .and. |
158 |
& .not. (last .or. endofday ) ) then |
159 |
c-- Accumulate the array holding the average. |
160 |
do bj = jtlo,jthi |
161 |
do bi = itlo,ithi |
162 |
do j = jmin,jmax |
163 |
do i = imin,imax |
164 |
psbar(i,j,bi,bj) = psbar(i,j,bi,bj) + etaN(i,j,bi,bj) |
165 |
enddo |
166 |
enddo |
167 |
enddo |
168 |
enddo |
169 |
else |
170 |
print*,' cost_averagesfields: Daily flags wrong' |
171 |
stop ' ... stopped in cost_averagesfields; psbar part.' |
172 |
endif |
173 |
#endif |
174 |
|
175 |
#if (defined (ALLOW_THETA_COST_CONTRIBUTION) || \ |
176 |
defined (ALLOW_SST_COST_CONTRIBUTION) || \ |
177 |
defined (ALLOW_CTDT_COST_CONTRIBUTION) || \ |
178 |
defined (ALLOW_CTDTCLIM_COST_CONTRIBUTION) || \ |
179 |
defined (ALLOW_XBT_COST_CONTRIBUTION) || \ |
180 |
defined (ALLOW_ARGO_THETA_COST_CONTRIBUTION) || \ |
181 |
defined (ALLOW_DRIFT_COST_CONTRIBUTION) || \ |
182 |
defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
183 |
c-- Next, do the monthly average for temperature. |
184 |
if (first .or. startofmonth) then |
185 |
c-- Assign the first value to the array holding the average. |
186 |
do bj = jtlo,jthi |
187 |
do bi = itlo,ithi |
188 |
do k = 1,nr |
189 |
do j = jmin,jmax |
190 |
do i = imin,imax |
191 |
tbar(i,j,k,bi,bj) = theta(i,j,k,bi,bj) |
192 |
enddo |
193 |
enddo |
194 |
enddo |
195 |
enddo |
196 |
enddo |
197 |
else if (last .or. endofmonth) then |
198 |
c-- Add the last value and devide by the number of accumulated |
199 |
c-- records. |
200 |
do bj = jtlo,jthi |
201 |
do bi = itlo,ithi |
202 |
do k = 1,nr |
203 |
do j = jmin,jmax |
204 |
do i = imin,imax |
205 |
tbar(i,j,k,bi,bj) = (tbar (i,j,k,bi,bj) + |
206 |
& theta(i,j,k,bi,bj) )/ |
207 |
& float(sum1mon) |
208 |
enddo |
209 |
enddo |
210 |
enddo |
211 |
enddo |
212 |
enddo |
213 |
|
214 |
c-- Save tbar on file. |
215 |
if (optimcycle .ge. 0) then |
216 |
ilt=ilnblnk( tbarfile ) |
217 |
write(fnametbar,'(2a,i10.10)') |
218 |
& tbarfile(1:ilt), '.', optimcycle |
219 |
endif |
220 |
|
221 |
call active_write_xyz( fnametbar, tbar, monrec, optimcycle, |
222 |
& mythid, xx_tbar_mean_dummy ) |
223 |
|
224 |
else if ( ( inmonth ) .and. |
225 |
& .not. (first .or. startofmonth) .and. |
226 |
& .not. (last .or. endofmonth ) ) then |
227 |
c-- Accumulate the array holding the average. |
228 |
do bj = jtlo,jthi |
229 |
do bi = itlo,ithi |
230 |
do k = 1,nr |
231 |
do j = jmin,jmax |
232 |
do i = imin,imax |
233 |
tbar(i,j,k,bi,bj) = tbar (i,j,k,bi,bj) + |
234 |
& theta(i,j,k,bi,bj) |
235 |
enddo |
236 |
enddo |
237 |
enddo |
238 |
enddo |
239 |
enddo |
240 |
else |
241 |
print*,' cost_averagesfields: Daily flags wrong' |
242 |
stop ' ... stopped in cost_averagesfields; tbar part (3d).' |
243 |
endif |
244 |
#endif |
245 |
|
246 |
#if (defined (ALLOW_SALT_COST_CONTRIBUTION) || \ |
247 |
defined (ALLOW_SSS_COST_CONTRIBUTION) || \ |
248 |
defined (ALLOW_CTDS_COST_CONTRIBUTION) || \ |
249 |
defined (ALLOW_CTDSCLIM_COST_CONTRIBUTION) || \ |
250 |
defined (ALLOW_ARGO_SALT_COST_CONTRIBUTION) || \ |
251 |
defined (ALLOW_DRIFT_COST_CONTRIBUTION) || \ |
252 |
defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
253 |
c-- Next, do the monthly averages for salinity. |
254 |
if (first .or. startofmonth) then |
255 |
do bj = jtlo,jthi |
256 |
do bi = itlo,ithi |
257 |
do k = 1,nr |
258 |
do j = jmin,jmax |
259 |
do i = imin,imax |
260 |
sbar(i,j,k,bi,bj) = salt(i,j,k,bi,bj) |
261 |
enddo |
262 |
enddo |
263 |
enddo |
264 |
enddo |
265 |
enddo |
266 |
else if (last .or. endofmonth) then |
267 |
do bj = jtlo,jthi |
268 |
do bi = itlo,ithi |
269 |
do k = 1,nr |
270 |
do j = jmin,jmax |
271 |
do i = imin,imax |
272 |
sbar(i,j,k,bi,bj) = (sbar (i,j,k,bi,bj) + |
273 |
& salt(i,j,k,bi,bj) )/ |
274 |
& float(sum1mon) |
275 |
enddo |
276 |
enddo |
277 |
enddo |
278 |
enddo |
279 |
enddo |
280 |
|
281 |
c-- Save sbar. |
282 |
if (optimcycle .ge. 0) then |
283 |
ils=ilnblnk( sbarfile ) |
284 |
write(fnamesbar,'(2a,i10.10)') |
285 |
& sbarfile(1:ils), '.', optimcycle |
286 |
endif |
287 |
|
288 |
call active_write_xyz( fnamesbar, sbar, monrec, optimcycle, |
289 |
& mythid, xx_sbar_mean_dummy) |
290 |
|
291 |
else if ( ( inmonth ) .and. |
292 |
& .not. (first .or. startofmonth) .and. |
293 |
& .not. (last .or. endofmonth ) ) then |
294 |
c-- Accumulate sbar. |
295 |
do bj = jtlo,jthi |
296 |
do bi = itlo,ithi |
297 |
do k = 1,nr |
298 |
do j = jmin,jmax |
299 |
do i = imin,imax |
300 |
sbar(i,j,k,bi,bj) = sbar (i,j,k,bi,bj) + |
301 |
& salt (i,j,k,bi,bj) |
302 |
enddo |
303 |
enddo |
304 |
enddo |
305 |
enddo |
306 |
enddo |
307 |
else |
308 |
print*,' cost_averagesfields: Daily flags wrong' |
309 |
stop ' ... stopped in cost_averagesfields; sbar part.' |
310 |
endif |
311 |
#endif |
312 |
|
313 |
#ifdef ALLOW_DRIFTW_COST_CONTRIBUTION |
314 |
c-- Next, do the averages for velocitty. |
315 |
if (first .or. startofmonth) then |
316 |
c-- Assign the first value to the array holding the average. |
317 |
do bj = jtlo,jthi |
318 |
do bi = itlo,ithi |
319 |
do k = 1,nr |
320 |
do j = jmin,jmax |
321 |
do i = imin,imax |
322 |
wbar(i,j,k,bi,bj) = wvel(i,j,k,bi,bj) |
323 |
enddo |
324 |
enddo |
325 |
enddo |
326 |
enddo |
327 |
enddo |
328 |
else if (last .or. endofmonth) then |
329 |
c-- Add the last value and devide by the number of accumulated |
330 |
c-- records. |
331 |
do bj = jtlo,jthi |
332 |
do bi = itlo,ithi |
333 |
do k = 1,nr |
334 |
do j = jmin,jmax |
335 |
do i = imin,imax |
336 |
wbar(i,j,k,bi,bj) = (wbar (i,j,k,bi,bj) + |
337 |
& wvel(i,j,k,bi,bj) )/ |
338 |
& float(sum1mon) |
339 |
enddo |
340 |
enddo |
341 |
enddo |
342 |
enddo |
343 |
enddo |
344 |
|
345 |
c-- Save wbar on file. |
346 |
if (optimcycle .ge. 0) then |
347 |
ilt=ilnblnk( wbarfile ) |
348 |
write(fnamewbar,'(2a,i10.10)') |
349 |
& wbarfile(1:ilt), '.', optimcycle |
350 |
endif |
351 |
|
352 |
call active_write_xyz( fnamewbar, wbar, monrec, optimcycle, |
353 |
& mythid, xx_wbar_mean_dummy ) |
354 |
|
355 |
else if ( ( inmonth ) .and. |
356 |
& .not. (first .or. startofmonth) .and. |
357 |
& .not. (last .or. endofmonth ) ) then |
358 |
c-- Accumulate the array holding the average. |
359 |
do bj = jtlo,jthi |
360 |
do bi = itlo,ithi |
361 |
do k = 1,nr |
362 |
do j = jmin,jmax |
363 |
do i = imin,imax |
364 |
wbar(i,j,k,bi,bj) = wbar (i,j,k,bi,bj) + |
365 |
& wvel(i,j,k,bi,bj) |
366 |
enddo |
367 |
enddo |
368 |
enddo |
369 |
enddo |
370 |
enddo |
371 |
else |
372 |
print*,' cost_averagesfields: Daily flags wrong' |
373 |
stop ' ... stopped in cost_averagesfields; tbar part (3d).' |
374 |
endif |
375 |
#endif |
376 |
|
377 |
#if (defined (ALLOW_DRIFTER_COST_CONTRIBUTION) || \ |
378 |
defined (ALLOW_CURMTR_COST_CONTRIBUTION) || \ |
379 |
defined (OBCS_AGEOS_COST_CONTRIBUTION)) |
380 |
cph#if (defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
381 |
cph defined (ALLOW_DRIFTER_COST_CONTRIBUTION) || \ |
382 |
cph There's a mismatch between the cost_drifer and the |
383 |
cph cost_obcs usage of ubar, vbar. |
384 |
cph cost_obcs refers to monthly means, cost_drifer to total mean. |
385 |
cph Needs to be updated for cost_drifer. |
386 |
|
387 |
c-- Next, do the averages for velocitty. |
388 |
if (first .or. startofmonth) then |
389 |
do bj = jtlo,jthi |
390 |
do bi = itlo,ithi |
391 |
do k = 1,nr |
392 |
do j = jmin,jmax |
393 |
do i = imin,imax |
394 |
ubar(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) |
395 |
vbar(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) |
396 |
enddo |
397 |
enddo |
398 |
enddo |
399 |
enddo |
400 |
enddo |
401 |
else if (last .or. endofmonth) then |
402 |
do bj = jtlo,jthi |
403 |
do bi = itlo,ithi |
404 |
do k = 1,nr |
405 |
do j = jmin,jmax |
406 |
do i = imin,imax |
407 |
ubar(i,j,k,bi,bj) = (ubar (i,j,k,bi,bj) + |
408 |
& uVel(i,j,k,bi,bj) )/ |
409 |
& float(sum1mon) |
410 |
vbar(i,j,k,bi,bj) = (vbar (i,j,k,bi,bj) + |
411 |
& vVel(i,j,k,bi,bj) )/ |
412 |
& float(sum1mon) |
413 |
enddo |
414 |
enddo |
415 |
enddo |
416 |
enddo |
417 |
enddo |
418 |
|
419 |
c-- Save ubar and vbar. |
420 |
if (optimcycle .ge. 0) then |
421 |
ils=ilnblnk( ubarfile ) |
422 |
write(fnameubar,'(2a,i10.10)') |
423 |
& ubarfile(1:ils), '.', optimcycle |
424 |
write(fnamevbar,'(2a,i10.10)') |
425 |
& vbarfile(1:ils), '.', optimcycle |
426 |
endif |
427 |
|
428 |
call active_write_xyz( fnameubar, ubar, monrec, optimcycle, |
429 |
& mythid, xx_ubar_mean_dummy) |
430 |
|
431 |
call active_write_xyz( fnamevbar, vbar, monrec, optimcycle, |
432 |
& mythid, xx_vbar_mean_dummy) |
433 |
|
434 |
cph else if ( ( inmonth ) .and. |
435 |
cph & .not. (first .or. startofmonth) .and. |
436 |
cph & .not. (last .or. endofmonth) ) then |
437 |
else if ( .not. (first .or. startofmonth) .and. |
438 |
& .not. (last .or. endofmonth ) ) then |
439 |
c-- Accumulate ubar and vbar. |
440 |
do bj = jtlo,jthi |
441 |
do bi = itlo,ithi |
442 |
do k = 1,nr |
443 |
do j = jmin,jmax |
444 |
do i = imin,imax |
445 |
ubar(i,j,k,bi,bj) = ubar (i,j,k,bi,bj) + |
446 |
& uVel (i,j,k,bi,bj) |
447 |
vbar(i,j,k,bi,bj) = vbar (i,j,k,bi,bj) + |
448 |
& vVel (i,j,k,bi,bj) |
449 |
enddo |
450 |
enddo |
451 |
enddo |
452 |
enddo |
453 |
enddo |
454 |
else |
455 |
print*,' cost_averagesfields: Daily flags wrong' |
456 |
stop ' ... stopped in cost_averagesfields; ubar part.' |
457 |
endif |
458 |
|
459 |
#endif |
460 |
|
461 |
#ifdef ALLOW_SCAT_COST_CONTRIBUTION |
462 |
c-- Next, do the averages for velocitty. |
463 |
if (first.or. startofmonth) then |
464 |
do bj = jtlo,jthi |
465 |
do bi = itlo,ithi |
466 |
do j = jmin,jmax |
467 |
do i = imin,imax |
468 |
tauxbar(i,j,bi,bj) = ustress(i,j,bi,bj) |
469 |
tauybar(i,j,bi,bj) = vstress(i,j,bi,bj) |
470 |
enddo |
471 |
enddo |
472 |
enddo |
473 |
enddo |
474 |
else if (last .or. endofmonth) then |
475 |
do bj = jtlo,jthi |
476 |
do bi = itlo,ithi |
477 |
do j = jmin,jmax |
478 |
do i = imin,imax |
479 |
tauxbar(i,j,bi,bj) = (tauxbar (i,j,bi,bj) + |
480 |
& ustress(i,j,bi,bj) )/ |
481 |
& float(sum1mon) |
482 |
tauybar(i,j,bi,bj) = (tauybar (i,j,bi,bj) + |
483 |
& vstress(i,j,bi,bj) )/ |
484 |
& float(sum1mon) |
485 |
enddo |
486 |
enddo |
487 |
enddo |
488 |
enddo |
489 |
|
490 |
c-- Save ubar and vbar. |
491 |
if (optimcycle .ge. 0) then |
492 |
ils=ilnblnk( tauxbarfile ) |
493 |
write(fnametauxbar,'(2a,i10.10)') |
494 |
& tauxbarfile(1:ils), '.', optimcycle |
495 |
ils=ilnblnk( tauybarfile ) |
496 |
write(fnametauybar,'(2a,i10.10)') |
497 |
& tauybarfile(1:ils), '.', optimcycle |
498 |
endif |
499 |
|
500 |
call active_write_xy( fnametauxbar, tauxbar, monrec, optimcycle, |
501 |
& mythid, xx_taux_mean_dummy) |
502 |
|
503 |
call active_write_xy( fnametauybar, tauybar, monrec, optimcycle, |
504 |
& mythid, xx_tauy_mean_dummy) |
505 |
|
506 |
|
507 |
else if ( .not. (first.or. startofmonth) .and. |
508 |
& .not. (last .or. endofmonth) ) then |
509 |
c-- Accumulate ubar and vbar. |
510 |
do bj = jtlo,jthi |
511 |
do bi = itlo,ithi |
512 |
do j = jmin,jmax |
513 |
do i = imin,imax |
514 |
tauxbar(i,j,bi,bj) = tauxbar (i,j,bi,bj) + |
515 |
& ustress (i,j,bi,bj) |
516 |
tauybar(i,j,bi,bj) = tauybar (i,j,bi,bj) + |
517 |
& vstress (i,j,bi,bj) |
518 |
enddo |
519 |
enddo |
520 |
enddo |
521 |
enddo |
522 |
else |
523 |
print*,' cost_averagesfields: Daily flags wrong' |
524 |
stop ' ... stopped in cost_averagesfields; tauxbar part.' |
525 |
endif |
526 |
|
527 |
|
528 |
#endif |
529 |
|
530 |
#ifdef ALLOW_MEAN_HFLUX_COST_CONTRIBUTION |
531 |
cph Armin added swflux here! |
532 |
c-- Next, do the averages for velocitty. |
533 |
if (first) then |
534 |
do bj = jtlo,jthi |
535 |
do bi = itlo,ithi |
536 |
do j = jmin,jmax |
537 |
do i = imin,imax |
538 |
hfluxbar(i,j,bi,bj) = hflux(i,j,bi,bj) + |
539 |
& swflux(i,j,bi,bj) |
540 |
enddo |
541 |
enddo |
542 |
enddo |
543 |
enddo |
544 |
else if (last) then |
545 |
do bj = jtlo,jthi |
546 |
do bi = itlo,ithi |
547 |
do j = jmin,jmax |
548 |
do i = imin,imax |
549 |
hfluxbar(i,j,bi,bj) = (hfluxbar (i,j,bi,bj) + |
550 |
& hflux(i,j,bi,bj)+swflux(i,j,bi,bj))/ |
551 |
& float(nTimeSteps) |
552 |
enddo |
553 |
enddo |
554 |
enddo |
555 |
enddo |
556 |
|
557 |
c-- Save hfluxbar |
558 |
if (optimcycle .ge. 0) then |
559 |
ils=ilnblnk( hfluxbarfile ) |
560 |
write(fnamehfluxbar,'(2a,i10.10)') |
561 |
& hfluxbarfile(1:ils), '.', optimcycle |
562 |
endif |
563 |
|
564 |
call active_write_xy( fnamehfluxbar, hfluxbar, 1, optimcycle, |
565 |
& mythid, xx_hflux_mean_dummy) |
566 |
|
567 |
else if ( .not. (first) .and. |
568 |
& .not. (last ) ) then |
569 |
c-- Accumulate ubar and vbar. |
570 |
do bj = jtlo,jthi |
571 |
do bi = itlo,ithi |
572 |
do j = jmin,jmax |
573 |
do i = imin,imax |
574 |
hfluxbar(i,j,bi,bj) = hfluxbar (i,j,bi,bj) + |
575 |
& hflux (i,j,bi,bj) + swflux(i,j,bi,bj) |
576 |
enddo |
577 |
enddo |
578 |
enddo |
579 |
enddo |
580 |
else |
581 |
print*,' cost_averagesfields: Daily flags wrong' |
582 |
stop ' ... stopped in cost_averagesfields; hfluxbar part.' |
583 |
endif |
584 |
|
585 |
#endif |
586 |
|
587 |
#ifdef ALLOW_MEAN_SFLUX_COST_CONTRIBUTION |
588 |
c-- Next, do the averages for velocitty. |
589 |
if (first ) then |
590 |
do bj = jtlo,jthi |
591 |
do bi = itlo,ithi |
592 |
do j = jmin,jmax |
593 |
do i = imin,imax |
594 |
sfluxbar(i,j,bi,bj) = sflux(i,j,bi,bj) |
595 |
enddo |
596 |
enddo |
597 |
enddo |
598 |
enddo |
599 |
else if (last .or. endofmonth) then |
600 |
do bj = jtlo,jthi |
601 |
do bi = itlo,ithi |
602 |
do j = jmin,jmax |
603 |
do i = imin,imax |
604 |
sfluxbar(i,j,bi,bj) = (sfluxbar (i,j,bi,bj) + |
605 |
& sflux(i,j,bi,bj) )/ |
606 |
& float(nTimeSteps) |
607 |
enddo |
608 |
enddo |
609 |
enddo |
610 |
enddo |
611 |
|
612 |
|
613 |
c-- Save sfluxbar |
614 |
if (optimcycle .ge. 0) then |
615 |
ils=ilnblnk( sfluxbarfile ) |
616 |
write(fnamesfluxbar,'(2a,i10.10)') |
617 |
& sfluxbarfile(1:ils), '.', optimcycle |
618 |
endif |
619 |
|
620 |
call active_write_xy( fnamesfluxbar, sfluxbar, 1, optimcycle, |
621 |
& mythid, xx_sflux_mean_dummy) |
622 |
|
623 |
else if ( .not. (first) .and. |
624 |
& .not. (last ) ) then |
625 |
c-- Accumulate ubar and vbar. |
626 |
do bj = jtlo,jthi |
627 |
do bi = itlo,ithi |
628 |
do j = jmin,jmax |
629 |
do i = imin,imax |
630 |
sfluxbar(i,j,bi,bj) = sfluxbar (i,j,bi,bj) + |
631 |
& sflux (i,j,bi,bj) |
632 |
enddo |
633 |
enddo |
634 |
enddo |
635 |
enddo |
636 |
else |
637 |
print*,' cost_averagesfields: Daily flags wrong' |
638 |
stop ' ... stopped in cost_averagesfields; sfluxbar part.' |
639 |
endif |
640 |
|
641 |
#endif |
642 |
|
643 |
return |
644 |
end |
645 |
|