1 |
C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_averagesfields.F,v 1.1.2.4 2002/04/04 10:58:58 heimbach Exp $ |
2 |
|
3 |
#include "COST_CPPOPTIONS.h" |
4 |
#ifdef ALLOW_OBCS |
5 |
# include "OBCS_OPTIONS.h" |
6 |
#endif |
7 |
|
8 |
subroutine cost_AveragesFields( |
9 |
I mytime, |
10 |
I mythid |
11 |
& ) |
12 |
|
13 |
c ================================================================== |
14 |
c SUBROUTINE cost_AveragesFields |
15 |
c ================================================================== |
16 |
c |
17 |
c o Compute time averages of etaN, theta, and salt. The counters |
18 |
c are explicitly calculated instead of being incremented. This |
19 |
c reduces dependencies. The latter is useful for the adjoint code |
20 |
c generation. |
21 |
c |
22 |
c started: Christian Eckert eckert@mit.edu 30-Jun-1999 |
23 |
c |
24 |
c changed: Christian Eckert eckert@mit.edu 24-Feb-2000 |
25 |
c |
26 |
c - Restructured the code in order to create a package |
27 |
c for the MITgcmUV. |
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 |
integer mythid |
55 |
|
56 |
c == local variables == |
57 |
|
58 |
integer myiter |
59 |
integer bi,bj |
60 |
integer i,j,k |
61 |
integer itlo,ithi |
62 |
integer jtlo,jthi |
63 |
integer jmin,jmax |
64 |
integer imin,imax |
65 |
|
66 |
logical first |
67 |
logical startofday |
68 |
logical startofmonth |
69 |
logical inday |
70 |
logical inmonth |
71 |
logical last |
72 |
logical endofday |
73 |
logical endofmonth |
74 |
|
75 |
integer ilps, ils,ilt |
76 |
|
77 |
character*(128) fnamepsbar |
78 |
character*(128) fnametbar |
79 |
character*(128) fnamesbar |
80 |
character*(128) fnameubar |
81 |
character*(128) fnamevbar |
82 |
character*(128) fnamewbar |
83 |
character*(128) fnametauxbar |
84 |
character*(128) fnametauybar |
85 |
character*(128) fnamehfluxbar |
86 |
character*(128) fnamesfluxbar |
87 |
|
88 |
c == external functions == |
89 |
|
90 |
integer ilnblnk |
91 |
external ilnblnk |
92 |
|
93 |
c == end of interface == |
94 |
|
95 |
jtlo = mybylo(mythid) |
96 |
jthi = mybyhi(mythid) |
97 |
itlo = mybxlo(mythid) |
98 |
ithi = mybxhi(mythid) |
99 |
jmin = 1 |
100 |
jmax = sny |
101 |
imin = 1 |
102 |
imax = snx |
103 |
|
104 |
myiter = niter0 + INT((mytime-starttime)/deltaTClock+0.5) |
105 |
|
106 |
c-- Get the time flags and record numbers for the time averaging. |
107 |
|
108 |
call cost_AveragesFlags( |
109 |
I myiter, mytime, mythid, |
110 |
O first, startofday, startofmonth, |
111 |
O inday, inmonth, |
112 |
O last, endofday, endofmonth, |
113 |
O sum1day, dayrec, |
114 |
O sum1mon, monrec |
115 |
& ) |
116 |
|
117 |
#ifdef ALLOW_SSH_COST_CONTRIBUTION |
118 |
c-- First, do the daily averages. |
119 |
if (first .or. startofday) then |
120 |
c-- Assign the first value to the array holding the average. |
121 |
do bj = jtlo,jthi |
122 |
do bi = itlo,ithi |
123 |
do j = jmin,jmax |
124 |
do i = imin,imax |
125 |
psbar(i,j,bi,bj) = etaN(i,j,bi,bj) |
126 |
enddo |
127 |
enddo |
128 |
enddo |
129 |
enddo |
130 |
else if (last .or. endofday) then |
131 |
c-- Add the last value and devide by the number of accumulated |
132 |
c-- records. |
133 |
do bj = jtlo,jthi |
134 |
do bi = itlo,ithi |
135 |
do j = jmin,jmax |
136 |
do i = imin,imax |
137 |
psbar(i,j,bi,bj) = (psbar (i,j,bi,bj) + |
138 |
& etaN(i,j,bi,bj) )/ |
139 |
& float(sum1day) |
140 |
enddo |
141 |
enddo |
142 |
enddo |
143 |
enddo |
144 |
|
145 |
c-- Save psbar on file. |
146 |
if (optimcycle .ge. 0) then |
147 |
ilps=ilnblnk( psbarfile ) |
148 |
write(fnamepsbar,'(2a,i10.10)') |
149 |
& psbarfile(1:ilps), '.', optimcycle |
150 |
endif |
151 |
|
152 |
call active_write_xy( fnamepsbar, psbar, dayrec, optimcycle, |
153 |
& mythid, xx_psbar_mean_dummy ) |
154 |
|
155 |
else if ( ( inday ) .and. |
156 |
& .not. (first .or. startofday) .and. |
157 |
& .not. (last .or. endofday ) ) then |
158 |
c-- Accumulate the array holding the average. |
159 |
do bj = jtlo,jthi |
160 |
do bi = itlo,ithi |
161 |
do j = jmin,jmax |
162 |
do i = imin,imax |
163 |
psbar(i,j,bi,bj) = psbar(i,j,bi,bj) + etaN(i,j,bi,bj) |
164 |
enddo |
165 |
enddo |
166 |
enddo |
167 |
enddo |
168 |
else |
169 |
print* |
170 |
print*,' cost_AveragesFields: Daily flags are set', |
171 |
& ' inappropriately.' |
172 |
print* |
173 |
stop ' ... stopped in cost_AveragesFields; psbar part.' |
174 |
endif |
175 |
#endif |
176 |
|
177 |
#if (defined (ALLOW_THETA_COST_CONTRIBUTION) || \ |
178 |
defined (ALLOW_CTDT_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)') tbarfile(1:ilt),'.',optimcycle |
218 |
endif |
219 |
|
220 |
call active_write_xyz( fnametbar, tbar, monrec, optimcycle, |
221 |
& mythid, xx_tbar_mean_dummy ) |
222 |
|
223 |
else if ( ( inmonth ) .and. |
224 |
& .not. (first .or. startofmonth) .and. |
225 |
& .not. (last .or. endofmonth ) ) then |
226 |
c-- Accumulate the array holding the average. |
227 |
do bj = jtlo,jthi |
228 |
do bi = itlo,ithi |
229 |
do k = 1,nr |
230 |
do j = jmin,jmax |
231 |
do i = imin,imax |
232 |
tbar(i,j,k,bi,bj) = tbar (i,j,k,bi,bj) + |
233 |
& theta(i,j,k,bi,bj) |
234 |
enddo |
235 |
enddo |
236 |
enddo |
237 |
enddo |
238 |
enddo |
239 |
else |
240 |
print* |
241 |
print*,' cost_AveragesFields: Monthly flags a set', |
242 |
& ' inappropriately.' |
243 |
print* |
244 |
stop ' ... stopped in cost_AveragesFields; tbar part (3d).' |
245 |
endif |
246 |
#else |
247 |
#ifdef ALLOW_SST_COST_CONTRIBUTION |
248 |
c-- Next, do the monthly averages for sea surface temperature. |
249 |
if (first .or. startofmonth) then |
250 |
c-- Assign the first value to the array holding the average. |
251 |
do bj = jtlo,jthi |
252 |
do bi = itlo,ithi |
253 |
k = 1 |
254 |
do j = jmin,jmax |
255 |
do i = imin,imax |
256 |
tbar(i,j,bi,bj) = theta(i,j,k,bi,bj) |
257 |
enddo |
258 |
enddo |
259 |
enddo |
260 |
enddo |
261 |
else if (last .or. endofmonth) then |
262 |
c-- Add the last value and devide by the number of accumulated |
263 |
c-- records. |
264 |
do bj = jtlo,jthi |
265 |
do bi = itlo,ithi |
266 |
k = 1 |
267 |
do j = jmin,jmax |
268 |
do i = imin,imax |
269 |
tbar(i,j,bi,bj) = (tbar (i,j,bi,bj) + |
270 |
& theta(i,j,k,bi,bj) )/ |
271 |
& float(sum1mon) |
272 |
enddo |
273 |
enddo |
274 |
enddo |
275 |
enddo |
276 |
|
277 |
c-- Save tbar (2d) on file. |
278 |
if (optimcycle .ge. 0) then |
279 |
ilt=ilnblnk( tbarfile ) |
280 |
write(fnametbar,'(2a,i10.10)') tbarfile(1:ilt),'.',optimcycle |
281 |
endif |
282 |
|
283 |
call active_write_xy( fnametbar, tbar, monrec, optimcycle, |
284 |
& mythid, xx_tbar_mean_dummy) |
285 |
|
286 |
else if ( ( inmonth ) .and. |
287 |
& .not. (first .or. startofmonth) .and. |
288 |
& .not. (last .or. endofmonth ) ) then |
289 |
c-- Accumulate the array holding the average. |
290 |
do bj = jtlo,jthi |
291 |
do bi = itlo,ithi |
292 |
k = 1 |
293 |
do j = jmin,jmax |
294 |
do i = imin,imax |
295 |
tbar(i,j,bi,bj) = tbar (i,j,bi,bj) + |
296 |
& theta(i,j,k,bi,bj) |
297 |
enddo |
298 |
enddo |
299 |
enddo |
300 |
enddo |
301 |
else |
302 |
print* |
303 |
print*,' cost_AveragesFields: Monthly flags a set', |
304 |
& ' inappropriately.' |
305 |
print* |
306 |
stop ' ... stopped in cost_AveragesFields; tbar part (2d).' |
307 |
endif |
308 |
#endif |
309 |
#endif |
310 |
|
311 |
#if (defined (ALLOW_SALT_COST_CONTRIBUTION) || \ |
312 |
defined (ALLOW_CTDS_COST_CONTRIBUTION) || \ |
313 |
defined (ALLOW_ARGO_SALT_COST_CONTRIBUTION) || \ |
314 |
defined (ALLOW_DRIFT_COST_CONTRIBUTION) || \ |
315 |
defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
316 |
c-- Next, do the monthly averages for salinity. |
317 |
if (first .or. startofmonth) then |
318 |
do bj = jtlo,jthi |
319 |
do bi = itlo,ithi |
320 |
do k = 1,nr |
321 |
do j = jmin,jmax |
322 |
do i = imin,imax |
323 |
sbar(i,j,k,bi,bj) = salt(i,j,k,bi,bj) |
324 |
enddo |
325 |
enddo |
326 |
enddo |
327 |
enddo |
328 |
enddo |
329 |
else if (last .or. endofmonth) then |
330 |
do bj = jtlo,jthi |
331 |
do bi = itlo,ithi |
332 |
do k = 1,nr |
333 |
do j = jmin,jmax |
334 |
do i = imin,imax |
335 |
sbar(i,j,k,bi,bj) = (sbar (i,j,k,bi,bj) + |
336 |
& salt(i,j,k,bi,bj) )/ |
337 |
& float(sum1mon) |
338 |
enddo |
339 |
enddo |
340 |
enddo |
341 |
enddo |
342 |
enddo |
343 |
|
344 |
c-- Save sbar. |
345 |
if (optimcycle .ge. 0) then |
346 |
ils=ilnblnk( sbarfile ) |
347 |
write(fnamesbar,'(2a,i10.10)') sbarfile(1:ils),'.', |
348 |
& optimcycle |
349 |
endif |
350 |
|
351 |
call active_write_xyz( fnamesbar, sbar, monrec, optimcycle, |
352 |
& mythid, xx_sbar_mean_dummy) |
353 |
|
354 |
else if ( ( inmonth ) .and. |
355 |
& .not. (first .or. startofmonth) .and. |
356 |
& .not. (last .or. endofmonth ) ) then |
357 |
c-- Accumulate sbar. |
358 |
do bj = jtlo,jthi |
359 |
do bi = itlo,ithi |
360 |
do k = 1,nr |
361 |
do j = jmin,jmax |
362 |
do i = imin,imax |
363 |
sbar(i,j,k,bi,bj) = sbar (i,j,k,bi,bj) + |
364 |
& salt (i,j,k,bi,bj) |
365 |
enddo |
366 |
enddo |
367 |
enddo |
368 |
enddo |
369 |
enddo |
370 |
else |
371 |
print* |
372 |
print*,' cost_AveragesFields: Monthly flags a set', |
373 |
& ' inappropriately.' |
374 |
print* |
375 |
stop ' ... stopped in cost_AveragesFields; sbar part.' |
376 |
endif |
377 |
#else |
378 |
#ifdef ALLOW_SSS_COST_CONTRIBUTION |
379 |
c-- Next, do the monthly averages for sea surface salinity. |
380 |
if (first .or. startofmonth) then |
381 |
c-- Assign the first value to the array holding the average. |
382 |
do bj = jtlo,jthi |
383 |
do bi = itlo,ithi |
384 |
k = 1 |
385 |
do j = jmin,jmax |
386 |
do i = imin,imax |
387 |
sbar(i,j,bi,bj) = salt(i,j,k,bi,bj) |
388 |
enddo |
389 |
enddo |
390 |
enddo |
391 |
enddo |
392 |
else if (last .or. endofmonth) then |
393 |
c-- Add the last value and devide by the number of accumulated |
394 |
c-- records. |
395 |
do bj = jtlo,jthi |
396 |
do bi = itlo,ithi |
397 |
k = 1 |
398 |
do j = jmin,jmax |
399 |
do i = imin,imax |
400 |
sbar(i,j,bi,bj) = (sbar (i,j,bi,bj) + |
401 |
& salt(i,j,k,bi,bj) )/ |
402 |
& float(sum1mon) |
403 |
enddo |
404 |
enddo |
405 |
enddo |
406 |
enddo |
407 |
|
408 |
c-- Save sbar (2d) on file. |
409 |
if (optimcycle .ge. 0) then |
410 |
ilt=ilnblnk( sbarfile ) |
411 |
write(fnamesbar,'(2a,i10.10)') sbarfile(1:ilt),'.',optimcycle |
412 |
endif |
413 |
|
414 |
call active_write_xy( fnamesbar, sbar, monrec, optimcycle, |
415 |
& mythid, xx_tbar_mean_dummy) |
416 |
|
417 |
|
418 |
else if ( ( inmonth ) .and. |
419 |
& .not. (first .or. startofmonth) .and. |
420 |
& .not. (last .or. endofmonth ) ) then |
421 |
c-- Accumulate the array holding the average. |
422 |
do bj = jtlo,jthi |
423 |
do bi = itlo,ithi |
424 |
k = 1 |
425 |
do j = jmin,jmax |
426 |
do i = imin,imax |
427 |
sbar(i,j,bi,bj) = sbar (i,j,bi,bj) + |
428 |
& salt(i,j,k,bi,bj) |
429 |
enddo |
430 |
enddo |
431 |
enddo |
432 |
enddo |
433 |
else |
434 |
print* |
435 |
print*,' cost_AveragesFields: Monthly flags a set', |
436 |
& ' inappropriately.' |
437 |
print* |
438 |
stop ' ... stopped in cost_AveragesFields; sbar part (2d).' |
439 |
endif |
440 |
#endif |
441 |
#endif |
442 |
|
443 |
#ifdef ALLOW_DRIFTW_COST_CONTRIBUTION |
444 |
c-- Next, do the averages for velocitty. |
445 |
if (first .or. startofmonth) then |
446 |
c-- Assign the first value to the array holding the average. |
447 |
do bj = jtlo,jthi |
448 |
do bi = itlo,ithi |
449 |
do k = 1,nr |
450 |
do j = jmin,jmax |
451 |
do i = imin,imax |
452 |
wbar(i,j,k,bi,bj) = wvel(i,j,k,bi,bj) |
453 |
enddo |
454 |
enddo |
455 |
enddo |
456 |
enddo |
457 |
enddo |
458 |
else if (last .or. endofmonth) then |
459 |
c-- Add the last value and devide by the number of accumulated |
460 |
c-- records. |
461 |
do bj = jtlo,jthi |
462 |
do bi = itlo,ithi |
463 |
do k = 1,nr |
464 |
do j = jmin,jmax |
465 |
do i = imin,imax |
466 |
wbar(i,j,k,bi,bj) = (wbar (i,j,k,bi,bj) + |
467 |
& wvel(i,j,k,bi,bj) )/ |
468 |
& float(sum1mon) |
469 |
enddo |
470 |
enddo |
471 |
enddo |
472 |
enddo |
473 |
enddo |
474 |
|
475 |
c-- Save wbar on file. |
476 |
if (optimcycle .ge. 0) then |
477 |
ilt=ilnblnk( wbarfile ) |
478 |
write(fnamewbar,'(2a,i10.10)') wbarfile(1:ilt),'.',optimcycle |
479 |
endif |
480 |
|
481 |
call active_write_xyz( fnamewbar, wbar, monrec, optimcycle, |
482 |
& mythid, xx_wbar_mean_dummy ) |
483 |
|
484 |
else if ( ( inmonth ) .and. |
485 |
& .not. (first .or. startofmonth) .and. |
486 |
& .not. (last .or. endofmonth ) ) then |
487 |
c-- Accumulate the array holding the average. |
488 |
do bj = jtlo,jthi |
489 |
do bi = itlo,ithi |
490 |
do k = 1,nr |
491 |
do j = jmin,jmax |
492 |
do i = imin,imax |
493 |
wbar(i,j,k,bi,bj) = wbar (i,j,k,bi,bj) + |
494 |
& wvel(i,j,k,bi,bj) |
495 |
enddo |
496 |
enddo |
497 |
enddo |
498 |
enddo |
499 |
enddo |
500 |
else |
501 |
print* |
502 |
print*,' cost_AveragesFields: Monthly flags a set', |
503 |
& ' inappropriately.' |
504 |
print* |
505 |
stop ' ... stopped in cost_AveragesFields; tbar part (3d).' |
506 |
endif |
507 |
#endif |
508 |
|
509 |
#if (defined (ALLOW_DRIFTER_COST_CONTRIBUTION) || \ |
510 |
defined (ALLOW_OBCS_COST_CONTRIBUTION)) |
511 |
cph There's a mismatch between the cost_drifer and the |
512 |
cph cost_obcs usage of ubar, vbar. |
513 |
cph cost_obcs refers to monthly means, cost_drifer to total mean. |
514 |
cph Needs to be updated for cost_obcs!!!. |
515 |
c-- Next, do the averages for velocitty. |
516 |
if (first.or.startofmonth) then |
517 |
do bj = jtlo,jthi |
518 |
do bi = itlo,ithi |
519 |
do k = 1,nr |
520 |
do j = jmin,jmax |
521 |
do i = imin,imax |
522 |
ubar(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) |
523 |
vbar(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) |
524 |
enddo |
525 |
enddo |
526 |
enddo |
527 |
enddo |
528 |
enddo |
529 |
else if (last .or. endofmonth) then |
530 |
do bj = jtlo,jthi |
531 |
do bi = itlo,ithi |
532 |
do k = 1,nr |
533 |
do j = jmin,jmax |
534 |
do i = imin,imax |
535 |
ubar(i,j,k,bi,bj) = (ubar (i,j,k,bi,bj) + |
536 |
& uVel(i,j,k,bi,bj) )/ |
537 |
& float(sum1mon) |
538 |
vbar(i,j,k,bi,bj) = (vbar (i,j,k,bi,bj) + |
539 |
& vVel(i,j,k,bi,bj) )/ |
540 |
& float(sum1mon) |
541 |
enddo |
542 |
enddo |
543 |
enddo |
544 |
enddo |
545 |
enddo |
546 |
|
547 |
c-- Save ubar and vbar. |
548 |
if (optimcycle .ge. 0) then |
549 |
ils=ilnblnk( ubarfile ) |
550 |
write(fnameubar,'(2a,i10.10)') ubarfile(1:ils),'.', |
551 |
& optimcycle |
552 |
write(fnamevbar,'(2a,i10.10)') vbarfile(1:ils),'.', |
553 |
& optimcycle |
554 |
endif |
555 |
|
556 |
call active_write_xyz( fnameubar, ubar, monrec, optimcycle, |
557 |
& mythid, xx_ubar_mean_dummy) |
558 |
|
559 |
call active_write_xyz( fnamevbar, vbar, monrec, optimcycle, |
560 |
& mythid, xx_vbar_mean_dummy) |
561 |
|
562 |
ce , myiter, mytime ) |
563 |
|
564 |
else if ( ( inmonth ) .and. |
565 |
& .not. (first .or. startofmonth) .and. |
566 |
& .not. (last .or. endofmonth ) ) then |
567 |
c-- Accumulate ubar and vbar. |
568 |
do bj = jtlo,jthi |
569 |
do bi = itlo,ithi |
570 |
do k = 1,nr |
571 |
do j = jmin,jmax |
572 |
do i = imin,imax |
573 |
ubar(i,j,k,bi,bj) = ubar (i,j,k,bi,bj) + |
574 |
& uVel (i,j,k,bi,bj) |
575 |
vbar(i,j,k,bi,bj) = vbar (i,j,k,bi,bj) + |
576 |
& vVel (i,j,k,bi,bj) |
577 |
enddo |
578 |
enddo |
579 |
enddo |
580 |
enddo |
581 |
enddo |
582 |
else |
583 |
print* |
584 |
print*,' cost_AveragesFields: Monthly flags a set', |
585 |
& ' inappropriately.' |
586 |
print* |
587 |
stop ' ... stopped in cost_AveragesFields; ubar part.' |
588 |
endif |
589 |
|
590 |
#endif |
591 |
|
592 |
#ifdef ALLOW_SCAT_COST_CONTRIBUTION |
593 |
c-- Next, do the averages for velocitty. |
594 |
if (first.or. startofmonth) then |
595 |
do bj = jtlo,jthi |
596 |
do bi = itlo,ithi |
597 |
do j = jmin,jmax |
598 |
do i = imin,imax |
599 |
tauxbar(i,j,bi,bj) = ustress(i,j,bi,bj) |
600 |
tauybar(i,j,bi,bj) = vstress(i,j,bi,bj) |
601 |
enddo |
602 |
enddo |
603 |
enddo |
604 |
enddo |
605 |
else if (last .or. endofmonth) then |
606 |
do bj = jtlo,jthi |
607 |
do bi = itlo,ithi |
608 |
do j = jmin,jmax |
609 |
do i = imin,imax |
610 |
tauxbar(i,j,bi,bj) = (tauxbar (i,j,bi,bj) + |
611 |
& ustress(i,j,bi,bj) )/ |
612 |
& float(sum1mon) |
613 |
tauybar(i,j,bi,bj) = (tauybar (i,j,bi,bj) + |
614 |
& vstress(i,j,bi,bj) )/ |
615 |
& float(sum1mon) |
616 |
enddo |
617 |
enddo |
618 |
enddo |
619 |
enddo |
620 |
|
621 |
c-- Save ubar and vbar. |
622 |
if (optimcycle .ge. 0) then |
623 |
ils=ilnblnk( tauxbarfile ) |
624 |
write(fnametauxbar,'(2a,i10.10)') tauxbarfile(1:ils),'.', |
625 |
& optimcycle |
626 |
ils=ilnblnk( tauybarfile ) |
627 |
write(fnametauybar,'(2a,i10.10)') tauybarfile(1:ils),'.', |
628 |
& optimcycle |
629 |
endif |
630 |
|
631 |
call active_write_xy( fnametauxbar, tauxbar, monrec, optimcycle, |
632 |
& mythid, xx_taux_mean_dummy) |
633 |
|
634 |
call active_write_xy( fnametauybar, tauybar, monrec, optimcycle, |
635 |
& mythid, xx_tauy_mean_dummy) |
636 |
|
637 |
|
638 |
else if ( .not. (first.or. startofmonth) .and. |
639 |
& .not. (last .or. endofmonth) ) then |
640 |
c-- Accumulate ubar and vbar. |
641 |
do bj = jtlo,jthi |
642 |
do bi = itlo,ithi |
643 |
do j = jmin,jmax |
644 |
do i = imin,imax |
645 |
tauxbar(i,j,bi,bj) = tauxbar (i,j,bi,bj) + |
646 |
& ustress (i,j,bi,bj) |
647 |
tauybar(i,j,bi,bj) = tauybar (i,j,bi,bj) + |
648 |
& vstress (i,j,bi,bj) |
649 |
enddo |
650 |
enddo |
651 |
enddo |
652 |
enddo |
653 |
else |
654 |
print* |
655 |
print*,' cost_AveragesFields: Monthly flags a set', |
656 |
& ' inappropriately.' |
657 |
print* |
658 |
stop ' ... stopped in cost_AveragesFields; tauxbar part.' |
659 |
endif |
660 |
|
661 |
|
662 |
#endif |
663 |
|
664 |
#ifdef ALLOW_MEAN_HFLUX_COST_CONTRIBUTION |
665 |
c-- Next, do the averages for velocitty. |
666 |
if (first) then |
667 |
do bj = jtlo,jthi |
668 |
do bi = itlo,ithi |
669 |
do j = jmin,jmax |
670 |
do i = imin,imax |
671 |
hfluxbar(i,j,bi,bj)=hflux(i,j,bi,bj) |
672 |
$ +swflux(i,j,bi,bj) |
673 |
enddo |
674 |
enddo |
675 |
enddo |
676 |
enddo |
677 |
else if (last) then |
678 |
do bj = jtlo,jthi |
679 |
do bi = itlo,ithi |
680 |
do j = jmin,jmax |
681 |
do i = imin,imax |
682 |
hfluxbar(i,j,bi,bj) = (hfluxbar (i,j,bi,bj) + |
683 |
& hflux(i,j,bi,bj)+swflux(i,j,bi,bj))/ |
684 |
& float(nTimeSteps) |
685 |
enddo |
686 |
enddo |
687 |
enddo |
688 |
enddo |
689 |
|
690 |
c-- Save hfluxbar |
691 |
if (optimcycle .ge. 0) then |
692 |
ils=ilnblnk( hfluxbarfile ) |
693 |
write(fnamehfluxbar,'(2a,i10.10)') hfluxbarfile(1:ils),'.', |
694 |
& optimcycle |
695 |
endif |
696 |
|
697 |
call active_write_xy( fnamehfluxbar, hfluxbar, 1, |
698 |
& optimcycle, mythid, xx_hflux_mean_dummy) |
699 |
|
700 |
else if ( .not. (first) .and. |
701 |
& .not. (last ) ) then |
702 |
c-- Accumulate ubar and vbar. |
703 |
do bj = jtlo,jthi |
704 |
do bi = itlo,ithi |
705 |
do j = jmin,jmax |
706 |
do i = imin,imax |
707 |
hfluxbar(i,j,bi,bj) = hfluxbar (i,j,bi,bj) + |
708 |
& swflux(i,j,bi,bj)+hflux(i,j,bi,bj) |
709 |
enddo |
710 |
enddo |
711 |
enddo |
712 |
enddo |
713 |
else |
714 |
print* |
715 |
print*,' cost_AveragesFields: Monthly flags a set', |
716 |
& ' inappropriately.' |
717 |
print* |
718 |
stop ' ... stopped in cost_AveragesFields; hfluxbar part.' |
719 |
endif |
720 |
|
721 |
#endif |
722 |
|
723 |
#ifdef ALLOW_MEAN_SFLUX_COST_CONTRIBUTION |
724 |
c-- Next, do the averages for velocitty. |
725 |
if (first ) then |
726 |
do bj = jtlo,jthi |
727 |
do bi = itlo,ithi |
728 |
do j = jmin,jmax |
729 |
do i = imin,imax |
730 |
sfluxbar(i,j,bi,bj) = sflux(i,j,bi,bj) |
731 |
enddo |
732 |
enddo |
733 |
enddo |
734 |
enddo |
735 |
else if (last) then |
736 |
do bj = jtlo,jthi |
737 |
do bi = itlo,ithi |
738 |
do j = jmin,jmax |
739 |
do i = imin,imax |
740 |
sfluxbar(i,j,bi,bj) = (sfluxbar (i,j,bi,bj) + |
741 |
& sflux(i,j,bi,bj) )/ |
742 |
& float(nTimeSteps) |
743 |
enddo |
744 |
enddo |
745 |
enddo |
746 |
enddo |
747 |
|
748 |
|
749 |
c-- Save sfluxbar |
750 |
if (optimcycle .ge. 0) then |
751 |
ils=ilnblnk( sfluxbarfile ) |
752 |
write(fnamesfluxbar,'(2a,i10.10)') sfluxbarfile(1:ils),'.', |
753 |
& optimcycle |
754 |
endif |
755 |
|
756 |
call active_write_xy( fnamesfluxbar, sfluxbar, 1, |
757 |
& optimcycle, mythid, xx_sflux_mean_dummy) |
758 |
|
759 |
else if ( .not. (first) .and. |
760 |
& .not. (last ) ) then |
761 |
c-- Accumulate ubar and vbar. |
762 |
do bj = jtlo,jthi |
763 |
do bi = itlo,ithi |
764 |
do j = jmin,jmax |
765 |
do i = imin,imax |
766 |
sfluxbar(i,j,bi,bj) = sfluxbar (i,j,bi,bj) + |
767 |
& sflux (i,j,bi,bj) |
768 |
enddo |
769 |
enddo |
770 |
enddo |
771 |
enddo |
772 |
else |
773 |
print* |
774 |
print*,' cost_AveragesFields: Monthly flags a set', |
775 |
& ' inappropriately.' |
776 |
print* |
777 |
stop ' ... stopped in cost_AveragesFields; sfluxbar part.' |
778 |
endif |
779 |
|
780 |
#endif |
781 |
|
782 |
return |
783 |
end |
784 |
|