ROMS
Loading...
Searching...
No Matches
tl_set_avg_mod Module Reference

Functions/Subroutines

subroutine, public tl_set_avg (ng, tile)
 
subroutine tl_set_avg_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nout, kout)
 

Function/Subroutine Documentation

◆ tl_set_avg()

subroutine, public tl_set_avg_mod::tl_set_avg ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 26 of file tl_set_avg.F.

27!***********************************************************************
28!
29 USE mod_param
30 USE mod_stepping
31!
32! Imported variable declarations.
33!
34 integer, intent(in) :: ng, tile
35!
36! Local variable declarations.
37!
38 character (len=*), parameter :: MyFile = &
39 & __FILE__
40!
41# include "tile.h"
42!
43# ifdef PROFILE
44 CALL wclock_on (ng, itlm, 5, __line__, myfile)
45# endif
46 CALL tl_set_avg_tile (ng, tile, &
47 & lbi, ubi, lbj, ubj, &
48 & imins, imaxs, jmins, jmaxs, &
49# ifdef SOLVE3D
50 & nout, &
51# endif
52 & kout)
53
54# ifdef WET_DRY
55 CALL set_avg_masks (ng, tile, itlm, &
56 & lbi, ubi, lbj, ubj, &
57 & imins, imaxs, jmins, jmaxs, &
58 & grid(ng) % pmask_avg, &
59 & grid(ng) % rmask_avg, &
60 & grid(ng) % umask_avg, &
61 & grid(ng) % vmask_avg)
62# endif
63
64# ifdef PROFILE
65 CALL wclock_off (ng, itlm, 5, __line__, myfile)
66# endif
67!
68 RETURN
integer, parameter itlm
Definition mod_param.F:663
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_param::itlm, tl_set_avg_tile(), wclock_off(), and wclock_on().

Referenced by rp_main3d(), and tl_main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tl_set_avg_tile()

subroutine tl_set_avg_mod::tl_set_avg_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nout,
integer, intent(in) kout )
private

Definition at line 72 of file tl_set_avg.F.

79!***********************************************************************
80!
81 USE mod_param
82 USE mod_ncparam
83 USE mod_average
84 USE mod_forces
85# ifdef SOLVE3D
86 USE mod_grid
87 USE mod_mixing
88# endif
89 USE mod_ocean
90 USE mod_scalars
91!
92 implicit none
93!
94! Imported variable declarations.
95!
96 integer, intent(in) :: ng, tile
97 integer, intent(in) :: LBi, UBi, LBj, UBj
98 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
99 integer, intent(in) :: Kout
100# ifdef SOLVE3D
101 integer, intent(in) :: Nout
102# endif
103!
104!
105! Local variable declarations.
106!
107 integer :: i, it, j, k
108
109 real(r8) :: fac
110
111 real(r8) :: pfac(IminS:ImaxS,JminS:JmaxS)
112 real(r8) :: rfac(IminS:ImaxS,JminS:JmaxS)
113 real(r8) :: ufac(IminS:ImaxS,JminS:JmaxS)
114 real(r8) :: vfac(IminS:ImaxS,JminS:JmaxS)
115
116# include "set_bounds.h"
117!
118!-----------------------------------------------------------------------
119! Return if time-averaging window is zero.
120!-----------------------------------------------------------------------
121!
122 IF (navg(ng).eq.0) RETURN
123!
124!-----------------------------------------------------------------------
125! Initialize time-averaged arrays when appropriate. Notice that
126! fields are initilized twice during re-start. However, the time-
127! averaged fields are computed correctly.
128!-----------------------------------------------------------------------
129!
130 IF (((iic(ng).gt.ntsavg(ng)).and. &
131 & (mod(iic(ng)-1,navg(ng)).eq.1)).or. &
132 & ((nrrec(ng).gt.0).and.(iic(ng).eq.ntstart(ng)))) THEN
133
134# ifdef WET_DRY
135!
136! If wetting and drying, initialize time dependent counters for wet
137! points. The time averaged field at each point is only accumulated
138! over wet points since its multiplied by the appropriate mask.
139!
140 DO j=jstr,jendr
141 DO i=istr,iendr
142 grid(ng)%pmask_avg(i,j)=max(0.0_r8, &
143 & min(grid(ng)%pmask_full(i,j), &
144 & 1.0_r8))
145 END DO
146 END DO
147 DO j=jstrr,jendr
148 DO i=istrr,iendr
149 grid(ng)%rmask_avg(i,j)=max(0.0_r8, &
150 & min(grid(ng)%rmask_full(i,j), &
151 & 1.0_r8))
152 END DO
153 END DO
154 DO j=jstrr,jendr
155 DO i=istr,iendr
156 grid(ng)%umask_avg(i,j)=max(0.0_r8, &
157 & min(grid(ng)%umask_full(i,j), &
158 & 1.0_r8))
159 END DO
160 END DO
161 DO j=jstr,jendr
162 DO i=istrr,iendr
163 grid(ng)%vmask_avg(i,j)=max(0.0_r8, &
164 & min(grid(ng)%vmask_full(i,j), &
165 & 1.0_r8))
166 END DO
167 END DO
168# endif
169!
170! Initialize adjoint state variables.
171!
172 IF (aout(idfsur,ng)) THEN
173 DO j=jstrr,jendr
174 DO i=istrr,iendr
175 average(ng)%avgzeta(i,j)=ocean(ng)%tl_zeta(i,j,kout)
176# ifdef WET_DRY
177 average(ng)%avgzeta(i,j)=average(ng)%avgzeta(i,j)* &
178 & grid(ng)%rmask_full(i,j)
179# endif
180 END DO
181 END DO
182 END IF
183 IF (aout(idubar,ng)) THEN
184 DO j=jstrr,jendr
185 DO i=istr,iendr
186 average(ng)%avgu2d(i,j)=ocean(ng)%tl_ubar(i,j,kout)
187# ifdef WET_DRY
188 average(ng)%avgu2d(i,j)=average(ng)%avgu2d(i,j)* &
189 & grid(ng)%umask_full(i,j)
190# endif
191 END DO
192 END DO
193 END IF
194 IF (aout(idvbar,ng)) THEN
195 DO j=jstr,jendr
196 DO i=istrr,iendr
197 average(ng)%avgv2d(i,j)=ocean(ng)%tl_vbar(i,j,kout)
198# ifdef WET_DRY
199 average(ng)%avgv2d(i,j)=average(ng)%avgv2d(i,j)* &
200 & grid(ng)%vmask_full(i,j)
201# endif
202 END DO
203 END DO
204 END IF
205
206# ifdef SOLVE3D
207 IF (aout(iduvel,ng)) THEN
208 DO k=1,n(ng)
209 DO j=jstrr,jendr
210 DO i=istr,iendr
211 average(ng)%avgu3d(i,j,k)=ocean(ng)%tl_u(i,j,k,nout)
212# ifdef WET_DRY
213 average(ng)%avgu3d(i,j,k)=average(ng)%avgu3d(i,j,k)* &
214 & grid(ng)%umask_full(i,j)
215# endif
216 END DO
217 END DO
218 END DO
219 END IF
220 IF (aout(idvvel,ng)) THEN
221 DO k=1,n(ng)
222 DO j=jstr,jendr
223 DO i=istrr,iendr
224 average(ng)%avgv3d(i,j,k)=ocean(ng)%tl_v(i,j,k,nout)
225# ifdef WET_DRY
226 average(ng)%avgv3d(i,j,k)=average(ng)%avgv3d(i,j,k)* &
227 & grid(ng)%vmask_full(i,j)
228# endif
229 END DO
230 END DO
231 END DO
232 END IF
233 IF (aout(idovel,ng)) THEN
234 DO k=0,n(ng)
235 DO j=jstrr,jendr
236 DO i=istrr,iendr
237 average(ng)%avgw3d(i,j,k)=ocean(ng)%tl_W(i,j,k)* &
238 & grid(ng)%pm(i,j)* &
239 & grid(ng)%pn(i,j)
240# ifdef WET_DRY
241 average(ng)%avgw3d(i,j,k)=average(ng)%avgw3d(i,j,k)* &
242 & grid(ng)%rmask_full(i,j)
243# endif
244 END DO
245 END DO
246 END DO
247 END IF
248
249 IF (aout(iddano,ng)) THEN
250 DO k=1,n(ng)
251 DO j=jstrr,jendr
252 DO i=istrr,iendr
253 average(ng)%avgrho(i,j,k)=ocean(ng)%tl_rho(i,j,k)
254# ifdef WET_DRY
255 average(ng)%avgrho(i,j,k)=average(ng)%avgrho(i,j,k)* &
256 & grid(ng)%rmask_full(i,j)
257# endif
258 END DO
259 END DO
260 END DO
261 END IF
262 DO it=1,nt(ng)
263 IF (aout(idtvar(it),ng)) THEN
264 DO k=1,n(ng)
265 DO j=jstrr,jendr
266 DO i=istrr,iendr
267 average(ng)%avgt(i,j,k,it)=ocean(ng)%tl_t(i,j,k, &
268 & nout,it)
269# ifdef WET_DRY
270 average(ng)%avgt(i,j,k,it)=average(ng)%avgt(i,j,k,it)*&
271 & grid(ng)%rmask_full(i,j)
272# endif
273 END DO
274 END DO
275 END DO
276 END IF
277 END DO
278
279# if defined LMD_MIXING || defined MY25_MIXING || defined GLS_MIXING
280 IF (aout(idvvis,ng)) THEN
281 DO k=0,n(ng)
282 DO j=jstrr,jendr
283 DO i=istrr,iendr
284 average(ng)%avgAKv(i,j,k)=mixing(ng)%tl_Akv(i,j,k)
285# ifdef WET_DRY
286 average(ng)%avgAKv(i,j,k)=average(ng)%avgAKv(i,j,k)* &
287 & grid(ng)%rmask_full(i,j)
288# endif
289 END DO
290 END DO
291 END DO
292 END IF
293 IF (aout(idtdif,ng)) THEN
294 DO k=0,n(ng)
295 DO j=jstrr,jendr
296 DO i=istrr,iendr
297 average(ng)%avgAKt(i,j,k)=mixing(ng)%tl_Akt(i,j,k,itemp)
298# ifdef WET_DRY
299 average(ng)%avgAKt(i,j,k)=average(ng)%avgAKt(i,j,k)* &
300 & grid(ng)%rmask_full(i,j)
301# endif
302 END DO
303 END DO
304 END DO
305 END IF
306 IF (aout(idsdif,ng)) THEN
307 DO k=0,n(ng)
308 DO j=jstrr,jendr
309 DO i=istrr,iendr
310 average(ng)%avgAKs(i,j,k)=mixing(ng)%tl_Akt(i,j,k,isalt)
311# ifdef WET_DRY
312 average(ng)%avgAKs(i,j,k)=average(ng)%avgAKs(i,j,k)* &
313 & grid(ng)%rmask_full(i,j)
314# endif
315 END DO
316 END DO
317 END DO
318 END IF
319# endif
320# endif
321!
322! Initialize adjoint surface and bottom fluxes.
323!
324 IF (aout(idusms,ng)) THEN
325 DO j=jstrr,jendr
326 DO i=istr,iendr
327 average(ng)%avgsus(i,j)=forces(ng)%tl_sustr(i,j)
328# ifdef WET_DRY
329 average(ng)%avgsus(i,j)=average(ng)%avgsus(i,j)* &
330 & grid(ng)%umask_full(i,j)
331# endif
332 END DO
333 END DO
334 END IF
335 IF (aout(idvsms,ng)) THEN
336 DO j=jstr,jendr
337 DO i=istrr,iendr
338 average(ng)%avgsvs(i,j)=forces(ng)%tl_svstr(i,j)
339# ifdef WET_DRY
340 average(ng)%avgsvs(i,j)=average(ng)%avgsvs(i,j)* &
341 & grid(ng)%vmask_full(i,j)
342# endif
343 END DO
344 END DO
345 END IF
346
347 IF (aout(idubms,ng)) THEN
348 DO j=jstrr,jendr
349 DO i=istr,iendr
350 average(ng)%avgbus(i,j)=forces(ng)%tl_bustr(i,j)
351# ifdef WET_DRY
352 average(ng)%avgbus(i,j)=average(ng)%avgbus(i,j)* &
353 & grid(ng)%umask_full(i,j)
354# endif
355 END DO
356 END DO
357 END IF
358 IF (aout(idvbms,ng)) THEN
359 DO j=jstr,jendr
360 DO i=istrr,iendr
361 average(ng)%avgbvs(i,j)=forces(ng)%tl_bvstr(i,j)
362# ifdef WET_DRY
363 average(ng)%avgbvs(i,j)=average(ng)%avgbvs(i,j)* &
364 & grid(ng)%vmask_full(i,j)
365# endif
366 END DO
367 END DO
368 END IF
369
370# ifdef SOLVE3D
371 IF (aout(idtsur(itemp),ng)) THEN
372 DO j=jstrr,jendr
373 DO i=istrr,iendr
374 average(ng)%avgstf(i,j)=forces(ng)%tl_stflx(i,j,itemp)
375# ifdef WET_DRY
376 average(ng)%avgstf(i,j)=average(ng)%avgstf(i,j)* &
377 & grid(ng)%rmask_full(i,j)
378# endif
379 END DO
380 END DO
381 END IF
382 IF (aout(idtsur(isalt),ng)) THEN
383 DO j=jstrr,jendr
384 DO i=istrr,iendr
385 average(ng)%avgswf(i,j)=forces(ng)%tl_stflx(i,j,isalt)
386# ifdef WET_DRY
387 average(ng)%avgswf(i,j)=average(ng)%avgswf(i,j)* &
388 & grid(ng)%rmask_full(i,j)
389# endif
390 END DO
391 END DO
392 END IF
393# ifdef SHORTWAVE
394 IF (aout(idsrad,ng)) THEN
395 DO j=jstrr,jendr
396 DO i=istrr,iendr
397 average(ng)%avgsrf(i,j)=forces(ng)%tl_srflx(i,j)
398# ifdef WET_DRY
399 average(ng)%avgsrf(i,j)=average(ng)%avgsrf(i,j)* &
400 & grid(ng)%rmask_full(i,j)
401# endif
402 END DO
403 END DO
404 END IF
405# endif
406# ifdef BULK_FLUXES
407 IF (aout(idlhea,ng)) THEN
408 DO j=jstrr,jendr
409 DO i=istrr,iendr
410 average(ng)%avglhf(i,j)=forces(ng)%tl_lhflx(i,j)
411# ifdef WET_DRY
412 average(ng)%avglhf(i,j)=average(ng)%avglhf(i,j)* &
413 & grid(ng)%rmask_full(i,j)
414# endif
415 END DO
416 END DO
417 END IF
418 IF (aout(idlrad,ng)) THEN
419 DO j=jstrr,jendr
420 DO i=istrr,iendr
421 average(ng)%avglrf(i,j)=forces(ng)%tl_lrflx(i,j)
422# ifdef WET_DRY
423 average(ng)%avglrf(i,j)=average(ng)%avglrf(i,j)* &
424 & grid(ng)%rmask_full(i,j)
425# endif
426 END DO
427 END DO
428 END IF
429 IF (aout(idshea,ng)) THEN
430 DO j=jstrr,jendr
431 DO i=istrr,iendr
432 average(ng)%avgshf(i,j)=forces(ng)%tl_shflx(i,j)
433# ifdef WET_DRY
434 average(ng)%avgshf(i,j)=average(ng)%avgshf(i,j)* &
435 & grid(ng)%rmask_full(i,j)
436# endif
437 END DO
438 END DO
439 END IF
440# ifdef EMINUSP
441 IF (aout(idevap,ng)) THEN
442 DO j=jstrr,jendr
443 DO i=istrr,iendr
444 average(ng)%avgevap(i,j)=forces(ng)%tl_evap(i,j)
445# ifdef WET_DRY
446 average(ng)%avgevap(i,j)=average(ng)%avgevap(i,j)* &
447 & grid(ng)%rmask_full(i,j)
448# endif
449 END DO
450 END DO
451 END IF
452# endif
453# endif
454# endif
455!
456! Initialize adjoint of quadratic fields.
457!
458 IF (aout(idzzav,ng)) THEN
459 DO j=jstrr,jendr
460 DO i=istrr,iendr
461 average(ng)%avgZZ(i,j)=ocean(ng)%tl_zeta(i,j,kout)* &
462 & ocean(ng)%tl_zeta(i,j,kout)
463# ifdef WET_DRY
464 average(ng)%avgZZ(i,j)=average(ng)%avgZZ(i,j)* &
465 & grid(ng)%rmask_full(i,j)
466# endif
467 END DO
468 END DO
469 END IF
470 IF (aout(idu2av,ng)) THEN
471 DO j=jstrr,jendr
472 DO i=istr,iendr
473 average(ng)%avgU2(i,j)=ocean(ng)%tl_ubar(i,j,kout)* &
474 & ocean(ng)%tl_ubar(i,j,kout)
475# ifdef WET_DRY
476 average(ng)%avgU2(i,j)=average(ng)%avgU2(i,j)* &
477 & grid(ng)%umask_full(i,j)
478# endif
479 END DO
480 END DO
481 END IF
482 IF (aout(idv2av,ng)) THEN
483 DO j=jstr,jendr
484 DO i=istrr,iendr
485 average(ng)%avgV2(i,j)=ocean(ng)%tl_vbar(i,j,kout)* &
486 & ocean(ng)%tl_vbar(i,j,kout)
487# ifdef WET_DRY
488 average(ng)%avgV2(i,j)=average(ng)%avgV2(i,j)* &
489 & grid(ng)%vmask_full(i,j)
490# endif
491 END DO
492 END DO
493 END IF
494
495# ifdef SOLVE3D
496 IF (aout(iduuav,ng)) THEN
497 DO k=1,n(ng)
498 DO j=jstrr,jendr
499 DO i=istr,iendr
500 average(ng)%avgUU(i,j,k)=ocean(ng)%tl_u(i,j,k,nout)* &
501 & ocean(ng)%tl_u(i,j,k,nout)
502# ifdef WET_DRY
503 average(ng)%avgUU(i,j,k)=average(ng)%avgUU(i,j,k)* &
504 & grid(ng)%umask_full(i,j)
505# endif
506 END DO
507 END DO
508 END DO
509 END IF
510 IF (aout(idvvav,ng)) THEN
511 DO k=1,n(ng)
512 DO j=jstr,jendr
513 DO i=istrr,iendr
514 average(ng)%avgVV(i,j,k)=ocean(ng)%tl_v(i,j,k,nout)* &
515 & ocean(ng)%tl_v(i,j,k,nout)
516# ifdef WET_DRY
517 average(ng)%avgVV(i,j,k)=average(ng)%avgVV(i,j,k)* &
518 & grid(ng)%vmask_full(i,j)
519# endif
520 END DO
521 END DO
522 END DO
523 END IF
524 IF (aout(iduvav,ng)) THEN
525 DO k=1,n(ng)
526 DO j=jstr,jend
527 DO i=istr,iend
528 average(ng)%avgUV(i,j,k)=0.25_r8* &
529 & (ocean(ng)%tl_u(i ,j,k,nout)+ &
530 & ocean(ng)%tl_u(i+1,j,k,nout))*&
531 & (ocean(ng)%tl_v(i,j ,k,nout)+ &
532 & ocean(ng)%tl_v(i,j+1,k,nout))
533# ifdef WET_DRY
534 average(ng)%avgUV(i,j,k)=average(ng)%avgUV(i,j,k)* &
535 & grid(ng)%rmask_full(i,j)
536# endif
537 END DO
538 END DO
539 END DO
540 END IF
541
542 DO it=1,nat
543 IF (aout(idttav(it),ng)) THEN
544 DO k=1,n(ng)
545 DO j=jstrr,jendr
546 DO i=istrr,iendr
547 average(ng)%avgTT(i,j,k,it)=ocean(ng)%tl_t(i,j,k, &
548 & nout,it)* &
549 & ocean(ng)%tl_t(i,j,k, &
550 & nout,it)
551# ifdef WET_DRY
552 average(ng)%avgTT(i,j,k,it)=average(ng)%avgTT(i,j,k, &
553 & it)* &
554 & grid(ng)%rmask_full(i,j)
555# endif
556 END DO
557 END DO
558 END DO
559 END IF
560 IF (aout(idutav(it),ng)) THEN
561 DO k=1,n(ng)
562 DO j=jstrr,jendr
563 DO i=istr,iend
564 average(ng)%avgUT(i,j,k,it)=0.5_r8* &
565 & ocean(ng)%tl_u(i,j,k, &
566 & nout)* &
567 & (ocean(ng)%tl_t(i-1,j,k, &
568 & nout,it)+ &
569 & ocean(ng)%tl_t(i ,j,k, &
570 & nout,it))
571# ifdef WET_DRY
572 average(ng)%avgUT(i,j,k,it)=average(ng)%avgUT(i,j,k, &
573 & it)* &
574 & grid(ng)%umask_full(i,j)
575# endif
576 END DO
577 END DO
578 END DO
579 END IF
580 IF (aout(idvtav(it),ng)) THEN
581 DO k=1,n(ng)
582 DO j=jstr,jend
583 DO i=istrr,iendr
584 average(ng)%avgVT(i,j,k,it)=0.5_r8* &
585 & ocean(ng)%tl_v(i,j,k, &
586 & nout)* &
587 & (ocean(ng)%tl_t(i,j-1,k, &
588 & nout,it)+ &
589 & ocean(ng)%tl_t(i,j ,k, &
590 & nout,it))
591# ifdef WET_DRY
592 average(ng)%avgVT(i,j,k,it)=average(ng)%avgVT(i,j,k, &
593 & it)* &
594 & grid(ng)%vmask_full(i,j)
595# endif
596 END DO
597 END DO
598 END DO
599 END IF
600 END DO
601# endif
602!
603!-----------------------------------------------------------------------
604! Accumulate time-averaged fields.
605!-----------------------------------------------------------------------
606!
607 ELSE IF (iic(ng).gt.ntsavg(ng)) THEN
608
609# ifdef WET_DRY
610!
611! If wetting and drying, accumulate wet points counters.
612! points. The time averaged field at each point is only accumulated
613! over wet points since its multiplied by the appropriate mask.
614!
615 DO j=jstr,jendr
616 DO i=istr,iendr
617 grid(ng)%pmask_avg(i,j)=grid(ng)%pmask_avg(i,j)+ &
618 & max(0.0_r8, &
619 & min(grid(ng)%pmask_full(i,j), &
620 & 1.0_r8))
621 END DO
622 END DO
623 DO j=jstrr,jendr
624 DO i=istrr,iendr
625 grid(ng)%rmask_avg(i,j)=grid(ng)%rmask_avg(i,j)+ &
626 & max(0.0_r8, &
627 & min(grid(ng)%rmask_full(i,j), &
628 & 1.0_r8))
629 END DO
630 END DO
631 DO j=jstrr,jendr
632 DO i=istr,iendr
633 grid(ng)%umask_avg(i,j)=grid(ng)%umask_avg(i,j)+ &
634 & max(0.0_r8, &
635 & min(grid(ng)%umask_full(i,j), &
636 & 1.0_r8))
637 END DO
638 END DO
639 DO j=jstr,jendr
640 DO i=istrr,iendr
641 grid(ng)%vmask_avg(i,j)=grid(ng)%vmask_avg(i,j)+ &
642 & max(0.0_r8, &
643 & min(grid(ng)%vmask_full(i,j), &
644 & 1.0_r8))
645 END DO
646 END DO
647# endif
648!
649! Accumulate adjoint state variables.
650!
651 IF (aout(idfsur,ng)) THEN
652 DO j=jstrr,jendr
653 DO i=istrr,iendr
654 average(ng)%avgzeta(i,j)=average(ng)%avgzeta(i,j)+ &
655# ifdef WET_DRY
656 & grid(ng)%rmask_full(i,j)* &
657# endif
658 & ocean(ng)%tl_zeta(i,j,kout)
659 END DO
660 END DO
661 END IF
662 IF (aout(idubar,ng)) THEN
663 DO j=jstrr,jendr
664 DO i=istr,iendr
665 average(ng)%avgu2d(i,j)=average(ng)%avgu2d(i,j)+ &
666# ifdef WET_DRY
667 & grid(ng)%umask_full(i,j)* &
668# endif
669 & ocean(ng)%tl_ubar(i,j,kout)
670 END DO
671 END DO
672 END IF
673 IF (aout(idvbar,ng)) THEN
674 DO j=jstr,jendr
675 DO i=istrr,iendr
676 average(ng)%avgv2d(i,j)=average(ng)%avgv2d(i,j)+ &
677# ifdef WET_DRY
678 & grid(ng)%vmask_full(i,j)* &
679# endif
680 & ocean(ng)%tl_vbar(i,j,kout)
681 END DO
682 END DO
683 END IF
684
685# ifdef SOLVE3D
686 IF (aout(iduvel,ng)) THEN
687 DO k=1,n(ng)
688 DO j=jstrr,jendr
689 DO i=istr,iendr
690 average(ng)%avgu3d(i,j,k)=average(ng)%avgu3d(i,j,k)+ &
691# ifdef WET_DRY
692 & grid(ng)%umask_full(i,j)* &
693# endif
694 & ocean(ng)%tl_u(i,j,k,nout)
695 END DO
696 END DO
697 END DO
698 END IF
699 IF (aout(idvvel,ng)) THEN
700 DO k=1,n(ng)
701 DO j=jstr,jendr
702 DO i=istrr,iendr
703 average(ng)%avgv3d(i,j,k)=average(ng)%avgv3d(i,j,k)+ &
704# ifdef WET_DRY
705 & grid(ng)%vmask_full(i,j)* &
706# endif
707 & ocean(ng)%tl_v(i,j,k,nout)
708 END DO
709 END DO
710 END DO
711 END IF
712 IF (aout(idovel,ng)) THEN
713 DO k=0,n(ng)
714 DO j=jstrr,jendr
715 DO i=istrr,iendr
716 average(ng)%avgw3d(i,j,k)=average(ng)%avgw3d(i,j,k)+ &
717# ifdef WET_DRY
718 & grid(ng)%rmask_full(i,j)* &
719# endif
720 & ocean(ng)%tl_W(i,j,k)* &
721 & grid(ng)%pm(i,j)* &
722 & grid(ng)%pn(i,j)
723 END DO
724 END DO
725 END DO
726 END IF
727
728 IF (aout(iddano,ng)) THEN
729 DO k=1,n(ng)
730 DO j=jstrr,jendr
731 DO i=istrr,iendr
732 average(ng)%avgrho(i,j,k)=average(ng)%avgrho(i,j,k)+ &
733# ifdef WET_DRY
734 & grid(ng)%rmask_full(i,j)* &
735# endif
736 & ocean(ng)%tl_rho(i,j,k)
737 END DO
738 END DO
739 END DO
740 END IF
741 DO it=1,nt(ng)
742 IF (aout(idtvar(it),ng)) THEN
743 DO k=1,n(ng)
744 DO j=jstrr,jendr
745 DO i=istrr,iendr
746 average(ng)%avgt(i,j,k,it)=average(ng)%avgt(i,j,k,it)+&
747# ifdef WET_DRY
748 & grid(ng)%rmask_full(i,j)* &
749# endif
750 & ocean(ng)%tl_t(i,j,k, &
751 & nout,it)
752 END DO
753 END DO
754 END DO
755 END IF
756 END DO
757
758# if defined LMD_MIXING || defined MY25_MIXING || defined GLS_MIXING
759 IF (aout(idvvis,ng)) THEN
760 DO k=0,n(ng)
761 DO j=jstrr,jendr
762 DO i=istrr,iendr
763 average(ng)%avgAKv(i,j,k)=average(ng)%avgAKv(i,j,k)+ &
764# ifdef WET_DRY
765 & grid(ng)%rmask_full(i,j)* &
766# endif
767 & mixing(ng)%tl_Akv(i,j,k)
768 END DO
769 END DO
770 END DO
771 END IF
772 IF (aout(idtdif,ng)) THEN
773 DO k=0,n(ng)
774 DO j=jstrr,jendr
775 DO i=istrr,iendr
776 average(ng)%avgAKt(i,j,k)=average(ng)%avgAKt(i,j,k)+ &
777# ifdef WET_DRY
778 & grid(ng)%rmask_full(i,j)* &
779# endif
780 & mixing(ng)%tl_Akt(i,j,k,itemp)
781 END DO
782 END DO
783 END DO
784 END IF
785 IF (aout(idsdif,ng)) THEN
786 DO k=0,n(ng)
787 DO j=jstrr,jendr
788 DO i=istrr,iendr
789 average(ng)%avgAKs(i,j,k)=average(ng)%avgAKs(i,j,k)+ &
790# ifdef WET_DRY
791 & grid(ng)%rmask_full(i,j)* &
792# endif
793 & mixing(ng)%tl_Akt(i,j,k,isalt)
794 END DO
795 END DO
796 END DO
797 END IF
798# endif
799# endif
800!
801! Accumulate adjoint surface and bottom fluxes.
802!
803 IF (aout(idusms,ng)) THEN
804 DO j=jstrr,jendr
805 DO i=istr,iendr
806 average(ng)%avgsus(i,j)=average(ng)%avgsus(i,j)+ &
807# ifdef WET_DRY
808 & grid(ng)%umask_full(i,j)* &
809# endif
810 & forces(ng)%tl_sustr(i,j)
811 END DO
812 END DO
813 END IF
814 IF (aout(idvsms,ng)) THEN
815 DO j=jstr,jendr
816 DO i=istrr,iendr
817 average(ng)%avgsvs(i,j)=average(ng)%avgsvs(i,j)+ &
818# ifdef WET_DRY
819 & grid(ng)%vmask_full(i,j)* &
820# endif
821 & forces(ng)%tl_svstr(i,j)
822 END DO
823 END DO
824 END IF
825
826 IF (aout(idubms,ng)) THEN
827 DO j=jstrr,jendr
828 DO i=istr,iendr
829 average(ng)%avgbus(i,j)=average(ng)%avgbus(i,j)+ &
830# ifdef WET_DRY
831 & grid(ng)%umask_full(i,j)* &
832# endif
833 & forces(ng)%tl_bustr(i,j)
834 END DO
835 END DO
836 END IF
837 IF (aout(idvbms,ng)) THEN
838 DO j=jstr,jendr
839 DO i=istrr,iendr
840 average(ng)%avgbvs(i,j)=average(ng)%avgbvs(i,j)+ &
841# ifdef WET_DRY
842 & grid(ng)%vmask_full(i,j)* &
843# endif
844 & forces(ng)%tl_bvstr(i,j)
845 END DO
846 END DO
847 END IF
848
849# ifdef SOLVE3D
850 IF (aout(idtsur(itemp),ng)) THEN
851 DO j=jstrr,jendr
852 DO i=istrr,iendr
853 average(ng)%avgstf(i,j)=average(ng)%avgstf(i,j)+ &
854# ifdef WET_DRY
855 & grid(ng)%rmask_full(i,j)* &
856# endif
857 & forces(ng)%tl_stflx(i,j,itemp)
858 END DO
859 END DO
860 END IF
861 IF (aout(idtsur(isalt),ng)) THEN
862 DO j=jstrr,jendr
863 DO i=istrr,iendr
864 average(ng)%avgswf(i,j)=average(ng)%avgswf(i,j)+ &
865# ifdef WET_DRY
866 & grid(ng)%rmask_full(i,j)* &
867# endif
868 & forces(ng)%tl_stflx(i,j,isalt)
869 END DO
870 END DO
871 END IF
872# ifdef SHORTWAVE
873 IF (aout(idsrad,ng)) THEN
874 DO j=jstrr,jendr
875 DO i=istrr,iendr
876 average(ng)%avgsrf(i,j)=average(ng)%avgsrf(i,j)+ &
877# ifdef WET_DRY
878 & grid(ng)%rmask_full(i,j)* &
879# endif
880 & forces(ng)%tl_srflx(i,j)
881 END DO
882 END DO
883 END IF
884# endif
885# ifdef BULK_FLUXES
886 IF (aout(idlhea,ng)) THEN
887 DO j=jstrr,jendr
888 DO i=istrr,iendr
889 average(ng)%avglhf(i,j)=average(ng)%avglhf(i,j)+ &
890# ifdef WET_DRY
891 & grid(ng)%rmask_full(i,j)* &
892# endif
893 & forces(ng)%tl_lhflx(i,j)
894 END DO
895 END DO
896 END IF
897 IF (aout(idlrad,ng)) THEN
898 DO j=jstrr,jendr
899 DO i=istrr,iendr
900 average(ng)%avglrf(i,j)=average(ng)%avglrf(i,j)+ &
901# ifdef WET_DRY
902 & grid(ng)%rmask_full(i,j)* &
903# endif
904 & forces(ng)%tl_lrflx(i,j)
905 END DO
906 END DO
907 END IF
908 IF (aout(idshea,ng)) THEN
909 DO j=jstrr,jendr
910 DO i=istrr,iendr
911 average(ng)%avgshf(i,j)=average(ng)%avgshf(i,j)+ &
912# ifdef WET_DRY
913 & grid(ng)%rmask_full(i,j)* &
914# endif
915 & forces(ng)%tl_shflx(i,j)
916 END DO
917 END DO
918 END IF
919# ifdef EMINUSP
920 IF (aout(idevap,ng)) THEN
921 DO j=jstrr,jendr
922 DO i=istrr,iendr
923 average(ng)%avgevap(i,j)=average(ng)%avgevap(i,j)+ &
924# ifdef WET_DRY
925 & grid(ng)%rmask_full(i,j)* &
926# endif
927 & forces(ng)%tl_evap(i,j)
928 END DO
929 END DO
930 END IF
931# endif
932# endif
933# endif
934!
935! Accumulate adjoint quadratic fields.
936!
937 IF (aout(idzzav,ng)) THEN
938 DO j=jstrr,jendr
939 DO i=istrr,iendr
940 average(ng)%avgZZ(i,j)=average(ng)%avgZZ(i,j)+ &
941# ifdef WET_DRY
942 & grid(ng)%rmask_full(i,j)* &
943# endif
944 & ocean(ng)%tl_zeta(i,j,kout)* &
945 & ocean(ng)%tl_zeta(i,j,kout)
946 END DO
947 END DO
948 END IF
949 IF (aout(idu2av,ng)) THEN
950 DO j=jstrr,jendr
951 DO i=istr,iendr
952 average(ng)%avgU2(i,j)=average(ng)%avgU2(i,j)+ &
953# ifdef WET_DRY
954 & grid(ng)%umask_full(i,j)* &
955# endif
956 & ocean(ng)%tl_ubar(i,j,kout)* &
957 & ocean(ng)%tl_ubar(i,j,kout)
958 END DO
959 END DO
960 END IF
961 IF (aout(idv2av,ng)) THEN
962 DO j=jstr,jendr
963 DO i=istrr,iendr
964 average(ng)%avgV2(i,j)=average(ng)%avgV2(i,j)+ &
965# ifdef WET_DRY
966 & grid(ng)%vmask_full(i,j)* &
967# endif
968 & ocean(ng)%tl_vbar(i,j,kout)* &
969 & ocean(ng)%tl_vbar(i,j,kout)
970 END DO
971 END DO
972 END IF
973
974# ifdef SOLVE3D
975 IF (aout(iduuav,ng)) THEN
976 DO k=1,n(ng)
977 DO j=jstrr,jendr
978 DO i=istr,iendr
979 average(ng)%avgUU(i,j,k)=average(ng)%avgUU(i,j,k)+ &
980# ifdef WET_DRY
981 & grid(ng)%umask_full(i,j)* &
982# endif
983 & ocean(ng)%tl_u(i,j,k,nout)* &
984 & ocean(ng)%tl_u(i,j,k,nout)
985 END DO
986 END DO
987 END DO
988 END IF
989 IF (aout(idvvav,ng)) THEN
990 DO k=1,n(ng)
991 DO j=jstr,jendr
992 DO i=istrr,iendr
993 average(ng)%avgVV(i,j,k)=average(ng)%avgVV(i,j,k)+ &
994# ifdef WET_DRY
995 & grid(ng)%vmask_full(i,j)* &
996# endif
997 & ocean(ng)%tl_v(i,j,k,nout)* &
998 & ocean(ng)%tl_v(i,j,k,nout)
999 END DO
1000 END DO
1001 END DO
1002 END IF
1003 IF (aout(iduvav,ng)) THEN
1004 DO k=1,n(ng)
1005 DO j=jstr,jend
1006 DO i=istr,iend
1007 average(ng)%avgUV(i,j,k)=average(ng)%avgUV(i,j,k)+ &
1008# ifdef WET_DRY
1009 & grid(ng)%rmask_full(i,j)* &
1010# endif
1011 & 0.25_r8* &
1012 & (ocean(ng)%tl_u(i ,j,k,nout)+ &
1013 & ocean(ng)%tl_u(i+1,j,k,nout))*&
1014 & (ocean(ng)%tl_v(i,j ,k,nout)+ &
1015 & ocean(ng)%tl_v(i,j+1,k,nout))
1016 END DO
1017 END DO
1018 END DO
1019 END IF
1020
1021 DO it=1,nat
1022 IF (aout(idttav(it),ng)) THEN
1023 DO k=1,n(ng)
1024 DO j=jstrr,jendr
1025 DO i=istrr,iendr
1026 average(ng)%avgTT(i,j,k,it)=average(ng)%avgTT(i,j,k, &
1027 & it)+ &
1028# ifdef WET_DRY
1029 & grid(ng)%rmask_full(i,j)* &
1030# endif
1031 & ocean(ng)%tl_t(i,j,k, &
1032 & nout,it)* &
1033 & ocean(ng)%tl_t(i,j,k, &
1034 & nout,it)
1035 END DO
1036 END DO
1037 END DO
1038 END IF
1039 IF (aout(idutav(it),ng)) THEN
1040 DO k=1,n(ng)
1041 DO j=jstrr,jendr
1042 DO i=istr,iend
1043 average(ng)%avgUT(i,j,k,it)=average(ng)%avgUT(i,j,k, &
1044 & it)+ &
1045# ifdef WET_DRY
1046 & grid(ng)%umask_full(i,j)* &
1047# endif
1048 & 0.5_r8* &
1049 & ocean(ng)%tl_u(i,j,k, &
1050 & nout)* &
1051 & (ocean(ng)%tl_t(i-1,j,k, &
1052 & nout,it)+ &
1053 & ocean(ng)%tl_t(i ,j,k, &
1054 & nout,it))
1055 END DO
1056 END DO
1057 END DO
1058 END IF
1059 IF (aout(idvtav(it),ng)) THEN
1060 DO k=1,n(ng)
1061 DO j=jstr,jend
1062 DO i=istrr,iendr
1063 average(ng)%avgVT(i,j,k,it)=average(ng)%avgVT(i,j,k, &
1064 & it)+ &
1065# ifdef WET_DRY
1066 & grid(ng)%vmask_full(i,j)* &
1067# endif
1068 & 0.5_r8* &
1069 & ocean(ng)%tl_v(i,j,k, &
1070 & nout)* &
1071 & (ocean(ng)%tl_t(i,j-1,k, &
1072 & nout,it)+ &
1073 & ocean(ng)%tl_t(i,j ,k, &
1074 & nout,it))
1075 END DO
1076 END DO
1077 END DO
1078 END IF
1079 END DO
1080# endif
1081 END IF
1082!
1083!-----------------------------------------------------------------------
1084! Convert accumulated sums into time-averages, if appropriate.
1085!-----------------------------------------------------------------------
1086!
1087 IF ((iic(ng).gt.ntsavg(ng)).and. &
1088 & (mod(iic(ng)-1,navg(ng)).eq.0).and. &
1089 & ((iic(ng).ne.ntstart(ng)).or.(nrrec(ng).eq.0))) THEN
1090 IF (domain(ng)%SouthWest_Test(tile)) THEN
1091 avgtime(ng)=avgtime(ng)+real(navg(ng),r8)*dt(ng)
1092 END IF
1093!
1094! Set time-averaged factors for each C-grid variable type. Notice that
1095! the I- and J-ranges are all grid types are the same for convinience.
1096# ifdef WET_DRY
1097! In wetting and drying, the sums are devided by the number of times
1098! that each qrid point is wet.
1099# endif
1100!
1101# ifdef WET_DRY
1102 DO j=jstrr,jendr
1103 DO i=istrr,iendr
1104 pfac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%pmask_avg(i,j))
1105 rfac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%rmask_avg(i,j))
1106 ufac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%umask_avg(i,j))
1107 vfac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%vmask_avg(i,j))
1108 END DO
1109 END DO
1110# else
1111 fac=1.0_r8/real(navg(ng),r8)
1112 DO j=jstrr,jendr
1113 DO i=istrr,iendr
1114 pfac(i,j)=fac
1115 rfac(i,j)=fac
1116 ufac(i,j)=fac
1117 vfac(i,j)=fac
1118 END DO
1119 END DO
1120# endif
1121!
1122! Process adjoint state variables.
1123!
1124 IF (aout(idfsur,ng)) THEN
1125 DO j=jstrr,jendr
1126 DO i=istrr,iendr
1127 average(ng)%avgzeta(i,j)=rfac(i,j)* &
1128 & average(ng)%avgzeta(i,j)
1129 END DO
1130 END DO
1131 END IF
1132 IF (aout(idubar,ng)) THEN
1133 DO j=jstrr,jendr
1134 DO i=istr,iendr
1135 average(ng)%avgu2d(i,j)=ufac(i,j)* &
1136 & average(ng)%avgu2d(i,j)
1137 END DO
1138 END DO
1139 END IF
1140 IF (aout(idvbar,ng)) THEN
1141 DO j=jstr,jendr
1142 DO i=istrr,iendr
1143 average(ng)%avgv2d(i,j)=vfac(i,j)* &
1144 & average(ng)%avgv2d(i,j)
1145 END DO
1146 END DO
1147 END IF
1148
1149# ifdef SOLVE3D
1150 IF (aout(iduvel,ng)) THEN
1151 DO k=1,n(ng)
1152 DO j=jstrr,jendr
1153 DO i=istr,iendr
1154 average(ng)%avgu3d(i,j,k)=ufac(i,j)* &
1155 & average(ng)%avgu3d(i,j,k)
1156 END DO
1157 END DO
1158 END DO
1159 END IF
1160 IF (aout(idvvel,ng)) THEN
1161 DO k=1,n(ng)
1162 DO j=jstr,jendr
1163 DO i=istrr,iendr
1164 average(ng)%avgv3d(i,j,k)=vfac(i,j)* &
1165 & average(ng)%avgv3d(i,j,k)
1166 END DO
1167 END DO
1168 END DO
1169 END IF
1170
1171 IF (aout(idovel,ng)) THEN
1172 DO k=0,n(ng)
1173 DO j=jstrr,jendr
1174 DO i=istrr,iendr
1175 average(ng)%avgw3d(i,j,k)=rfac(i,j)* &
1176 & average(ng)%avgw3d(i,j,k)
1177 END DO
1178 END DO
1179 END DO
1180 END IF
1181
1182 IF (aout(iddano,ng)) THEN
1183 DO k=1,n(ng)
1184 DO j=jstrr,jendr
1185 DO i=istrr,iendr
1186 average(ng)%avgrho(i,j,k)=rfac(i,j)* &
1187 & average(ng)%avgrho(i,j,k)
1188 END DO
1189 END DO
1190 END DO
1191 END IF
1192 DO it=1,nt(ng)
1193 IF (aout(idtvar(it),ng)) THEN
1194 DO k=1,n(ng)
1195 DO j=jstrr,jendr
1196 DO i=istrr,iendr
1197 average(ng)%avgt(i,j,k,it)=rfac(i,j)* &
1198 & average(ng)%avgt(i,j,k,it)
1199 END DO
1200 END DO
1201 END DO
1202 END IF
1203 END DO
1204
1205# if defined LMD_MIXING || defined MY25_MIXING || defined GLS_MIXING
1206 IF (aout(idvvis,ng)) THEN
1207 DO k=0,n(ng)
1208 DO j=jstrr,jendr
1209 DO i=istrr,iendr
1210 average(ng)%avgAKv(i,j,k)=rfac(i,j)* &
1211 & average(ng)%avgAKv(i,j,k)
1212 END DO
1213 END DO
1214 END DO
1215 END IF
1216 IF (aout(idtdif,ng)) THEN
1217 DO k=0,n(ng)
1218 DO j=jstrr,jendr
1219 DO i=istrr,iendr
1220 average(ng)%avgAKt(i,j,k)=rfac(i,j)* &
1221 & average(ng)%avgAKt(i,j,k)
1222 END DO
1223 END DO
1224 END DO
1225 END IF
1226 IF (aout(idsdif,ng)) THEN
1227 DO k=0,n(ng)
1228 DO j=jstrr,jendr
1229 DO i=istrr,iendr
1230 average(ng)%avgAKs(i,j,k)=rfac(i,j)* &
1231 & average(ng)%avgAKs(i,j,k)
1232 END DO
1233 END DO
1234 END DO
1235 END IF
1236# endif
1237# endif
1238!
1239! Process adjoint surface and bottom fluxes.
1240!
1241 IF (aout(idusms,ng)) THEN
1242 DO j=jstrr,jendr
1243 DO i=istr,iendr
1244 average(ng)%avgsus(i,j)=ufac(i,j)* &
1245 & average(ng)%avgsus(i,j)
1246 END DO
1247 END DO
1248 END IF
1249 IF (aout(idvsms,ng)) THEN
1250 DO j=jstr,jendr
1251 DO i=istrr,iendr
1252 average(ng)%avgsvs(i,j)=vfac(i,j)* &
1253 & average(ng)%avgsvs(i,j)
1254 END DO
1255 END DO
1256 END IF
1257
1258 IF (aout(idubms,ng)) THEN
1259 DO j=jstrr,jendr
1260 DO i=istr,iendr
1261 average(ng)%avgbus(i,j)=ufac(i,j)* &
1262 & average(ng)%avgbus(i,j)
1263 END DO
1264 END DO
1265 END IF
1266 IF (aout(idvbms,ng)) THEN
1267 DO j=jstr,jendr
1268 DO i=istrr,iendr
1269 average(ng)%avgbvs(i,j)=vfac(i,j)* &
1270 & average(ng)%avgbvs(i,j)
1271 END DO
1272 END DO
1273 END IF
1274
1275# ifdef SOLVE3D
1276 IF (aout(idtsur(itemp),ng)) THEN
1277 DO j=jstrr,jendr
1278 DO i=istrr,iendr
1279 average(ng)%avgstf(i,j)=rfac(i,j)* &
1280 & average(ng)%avgstf(i,j)
1281 END DO
1282 END DO
1283 END IF
1284 IF (aout(idtsur(isalt),ng)) THEN
1285 DO j=jstrr,jendr
1286 DO i=istrr,iendr
1287 average(ng)%avgswf(i,j)=rfac(i,j)* &
1288 & average(ng)%avgswf(i,j)
1289 END DO
1290 END DO
1291 END IF
1292# ifdef SHORTWAVE
1293 IF (aout(idsrad,ng)) THEN
1294 DO j=jstrr,jendr
1295 DO i=istrr,iendr
1296 average(ng)%avgsrf(i,j)=rfac(i,j)* &
1297 & average(ng)%avgsrf(i,j)
1298 END DO
1299 END DO
1300 END IF
1301# endif
1302# ifdef BULK_FLUXES
1303 IF (aout(idlhea,ng)) THEN
1304 DO j=jstrr,jendr
1305 DO i=istrr,iendr
1306 average(ng)%avglhf(i,j)=rfac(i,j)* &
1307 & average(ng)%avglhf(i,j)
1308 END DO
1309 END DO
1310 END IF
1311 IF (aout(idshea,ng)) THEN
1312 DO j=jstrr,jendr
1313 DO i=istrr,iendr
1314 average(ng)%avgshf(i,j)=rfac(i,j)* &
1315 & average(ng)%avgshf(i,j)
1316 END DO
1317 END DO
1318 END IF
1319 IF (aout(idlrad,ng)) THEN
1320 DO j=jstrr,jendr
1321 DO i=istrr,iendr
1322 average(ng)%avglrf(i,j)=rfac(i,j)* &
1323 & average(ng)%avglrf(i,j)
1324 END DO
1325 END DO
1326 END IF
1327# ifdef EMINUSP
1328 IF (aout(idevap,ng)) THEN
1329 DO j=jstrr,jendr
1330 DO i=istrr,iendr
1331 average(ng)%avgevap(i,j)=rfac(i,j)* &
1332 & average(ng)%avgevap(i,j)
1333 END DO
1334 END DO
1335 END IF
1336# endif
1337# endif
1338# endif
1339!
1340! Process adjoint quadratic fields.
1341!
1342 IF (aout(idzzav,ng)) THEN
1343 DO j=jstrr,jendr
1344 DO i=istrr,iendr
1345 average(ng)%avgZZ(i,j)=rfac(i,j)* &
1346 & average(ng)%avgZZ(i,j)
1347 END DO
1348 END DO
1349 END IF
1350 IF (aout(idu2av,ng)) THEN
1351 DO j=jstrr,jendr
1352 DO i=istr,iendr
1353 average(ng)%avgU2(i,j)=ufac(i,j)* &
1354 & average(ng)%avgU2(i,j)
1355 END DO
1356 END DO
1357 END IF
1358 IF (aout(idv2av,ng)) THEN
1359 DO j=jstr,jendr
1360 DO i=istrr,iendr
1361 average(ng)%avgV2(i,j)=vfac(i,j)* &
1362 & average(ng)%avgV2(i,j)
1363 END DO
1364 END DO
1365 END IF
1366
1367# ifdef SOLVE3D
1368 IF (aout(iduuav,ng)) THEN
1369 DO k=1,n(ng)
1370 DO j=jstrr,jendr
1371 DO i=istr,iendr
1372 average(ng)%avgUU(i,j,k)=ufac(i,j)* &
1373 & average(ng)%avgUU(i,j,k)
1374 END DO
1375 END DO
1376 END DO
1377 END IF
1378 IF (aout(idvvav,ng)) THEN
1379 DO k=1,n(ng)
1380 DO j=jstr,jendr
1381 DO i=istrr,iendr
1382 average(ng)%avgVV(i,j,k)=vfac(i,j)* &
1383 & average(ng)%avgVV(i,j,k)
1384 END DO
1385 END DO
1386 END DO
1387 END IF
1388 IF (aout(iduvav,ng)) THEN
1389 DO k=1,n(ng)
1390 DO j=jstr,jend
1391 DO i=istr,iend
1392 average(ng)%avgUV(i,j,k)=rfac(i,j)* &
1393 & average(ng)%avgUV(i,j,k)
1394 END DO
1395 END DO
1396 END DO
1397 END IF
1398
1399 DO it=1,nat
1400 IF (aout(idttav(it),ng)) THEN
1401 DO k=1,n(ng)
1402 DO j=jstrr,jendr
1403 DO i=istrr,iendr
1404 average(ng)%avgTT(i,j,k,it)=rfac(i,j)* &
1405 & average(ng)%avgTT(i,j,k, &
1406 & it)
1407 END DO
1408 END DO
1409 END DO
1410 END IF
1411 IF (aout(idutav(it),ng)) THEN
1412 DO k=1,n(ng)
1413 DO j=jstrr,jendr
1414 DO i=istr,iend
1415 average(ng)%avgUT(i,j,k,it)=ufac(i,j)* &
1416 & average(ng)%avgUT(i,j,k, &
1417 & it)
1418 END DO
1419 END DO
1420 END DO
1421 END IF
1422 IF (aout(idvtav(it),ng)) THEN
1423 DO k=1,n(ng)
1424 DO j=jstr,jend
1425 DO i=istrr,iendr
1426 average(ng)%avgVT(i,j,k,it)=vfac(i,j)* &
1427 & average(ng)%avgVT(i,j,k, &
1428 & it)
1429 END DO
1430 END DO
1431 END DO
1432 END IF
1433 END DO
1434# endif
1435 END IF
1436!
1437 RETURN
type(t_average), dimension(:), allocatable average
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer iddano
integer, dimension(:), allocatable idttav
integer idevap
integer idzzav
integer idu2av
integer idubar
integer idvvel
integer idvsms
integer, dimension(:), allocatable idutav
integer idsdif
integer, dimension(:), allocatable idtsur
integer idtdif
integer idfsur
integer, dimension(:), allocatable idtvar
integer idvbms
integer iduvel
integer idovel
integer idvvav
integer idshea
integer idlrad
integer idv2av
integer idusms
integer idvvis
integer, dimension(:), allocatable idvtav
integer idlhea
integer idubms
integer idsrad
integer iduuav
logical, dimension(:,:), allocatable aout
integer iduvav
integer idvbar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nat
Definition mod_param.F:499
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable nrrec
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
real(dp), dimension(:), allocatable avgtime
integer, dimension(:), allocatable navg
integer isalt
integer itemp
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable ntsavg

References mod_ncparam::aout, mod_average::average, mod_scalars::avgtime, mod_param::domain, mod_scalars::dt, mod_forces::forces, mod_grid::grid, mod_ncparam::iddano, mod_ncparam::idevap, mod_ncparam::idfsur, mod_ncparam::idlhea, mod_ncparam::idlrad, mod_ncparam::idovel, mod_ncparam::idsdif, mod_ncparam::idshea, mod_ncparam::idsrad, mod_ncparam::idtdif, mod_ncparam::idtsur, mod_ncparam::idttav, mod_ncparam::idtvar, mod_ncparam::idu2av, mod_ncparam::idubar, mod_ncparam::idubms, mod_ncparam::idusms, mod_ncparam::idutav, mod_ncparam::iduuav, mod_ncparam::iduvav, mod_ncparam::iduvel, mod_ncparam::idv2av, mod_ncparam::idvbar, mod_ncparam::idvbms, mod_ncparam::idvsms, mod_ncparam::idvtav, mod_ncparam::idvvav, mod_ncparam::idvvel, mod_ncparam::idvvis, mod_ncparam::idzzav, mod_scalars::iic, mod_scalars::isalt, mod_scalars::itemp, mod_mixing::mixing, mod_param::n, mod_param::nat, mod_scalars::navg, mod_scalars::nrrec, mod_param::nt, mod_scalars::ntsavg, mod_scalars::ntstart, and mod_ocean::ocean.

Referenced by tl_set_avg().

Here is the caller graph for this function: