ROMS
Loading...
Searching...
No Matches
set_diags.F File Reference
#include "cppdefs.h"
#include "tile.h"
#include "set_bounds.h"
Include dependency graph for set_diags.F:

Go to the source code of this file.

Functions/Subroutines

subroutine set_diags (ng, tile)
 
subroutine set_diags_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kout, nrhs)
 

Function/Subroutine Documentation

◆ set_diags()

subroutine set_diags ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 3 of file set_diags.F.

4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This subroutine accumulates and computes output time-averaged !
13! diagnostic fields. Due to synchronization, the time-averaged !
14! diagnostic fields are computed in delayed mode. All averages !
15! are accumulated at the beginning of the next time-step. !
16! !
17!=======================================================================
18!
19 USE mod_param
20 USE mod_scalars
21 USE mod_stepping
22!
23 implicit none
24!
25! Imported variable declarations.
26!
27 integer, intent(in) :: ng, tile
28!
29! Local variable declarations.
30!
31 character (len=*), parameter :: MyFile = &
32 & __FILE__
33!
34# include "tile.h"
35!
36# ifdef PROFILE
37 CALL wclock_on (ng, inlm, 5, __line__, myfile)
38# endif
39 CALL set_diags_tile (ng, tile, &
40 & lbi, ubi, lbj, ubj, &
41 & imins, imaxs, jmins, jmaxs, &
42# ifdef SOLVE3D
43 & kstp(ng), &
44# else
45 & knew(ng), &
46# endif
47 & nrhs(ng))
48# ifdef PROFILE
49 CALL wclock_off (ng, inlm, 5, __line__, myfile)
50# endif
51!
52 RETURN
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nrhs
subroutine set_diags_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kout, nrhs)
Definition set_diags.F:60
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::inlm, mod_stepping::knew, mod_stepping::kstp, mod_stepping::nrhs, set_diags_tile(), wclock_off(), and wclock_on().

Referenced by main3d().

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

◆ set_diags_tile()

subroutine set_diags_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) kout,
integer, intent(in) nrhs )

Definition at line 56 of file set_diags.F.

60!***********************************************************************
61!
62 USE mod_param
63# ifdef DIAGNOSTICS_BIO
64 USE mod_biology
65# endif
66 USE mod_diags
67 USE mod_grid
68 USE mod_ocean
69 USE mod_scalars
70!
71 USE bc_2d_mod
72# ifdef SOLVE3D
73 USE bc_3d_mod
74# ifdef ECOSIM
75 USE bc_4d_mod
76# endif
77# endif
79# ifdef DISTRIBUTE
82# ifdef SOLVE3D
84# endif
85# endif
86!
87 implicit none
88!
89! Imported variable declarations.
90!
91 integer, intent(in) :: ng, tile
92 integer, intent(in) :: LBi, UBi, LBj, UBj
93 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
94 integer, intent(in) :: kout, nrhs
95!
96! Local variable declarations.
97!
98 integer :: i, it, j, k
99 integer :: iband, idiag
100
101 real(r8) :: fac
102
103 real(r8) :: rfac(IminS:ImaxS,JminS:JmaxS)
104 real(r8) :: ufac(IminS:ImaxS,JminS:JmaxS)
105 real(r8) :: vfac(IminS:ImaxS,JminS:JmaxS)
106
107# include "set_bounds.h"
108!
109!-----------------------------------------------------------------------
110! Return if time-averaging window is zero.
111!-----------------------------------------------------------------------
112!
113 IF (ndia(ng).eq.0) RETURN
114!
115!-----------------------------------------------------------------------
116! Initialize time-averaged diagnostic arrays when appropriate. Notice
117! that fields are initilized twice during re-start. However, the time-
118! averaged fields are computed correctly.
119!-----------------------------------------------------------------------
120!
121 IF (((iic(ng).gt.ntsdia(ng)).and. &
122 & (mod(iic(ng)-1,ndia(ng)).eq.1)).or. &
123 & ((iic(ng).ge.ntsdia(ng)).and.(ndia(ng).eq.1)).or. &
124 & ((nrrec(ng).gt.0).and.(iic(ng).eq.ntstart(ng)))) THEN
125
126# ifdef WET_DRY
127!
128! If wetting and drying, initialize time dependent counters for wet
129! points. The time averaged field at each point is only accumulated
130! over wet points since it is multiplied by the appropriate mask.
131!
132 DO j=jstr,jendr
133 DO i=istr,iendr
134 grid(ng)%pmask_dia(i,j)=max(0.0_r8, &
135 & min(grid(ng)%pmask_full(i,j), &
136 & 1.0_r8))
137 END DO
138 END DO
139 DO j=jstrr,jendr
140 DO i=istrr,iendr
141 grid(ng)%rmask_dia(i,j)=max(0.0_r8, &
142 & min(grid(ng)%rmask_full(i,j), &
143 & 1.0_r8))
144 END DO
145 END DO
146 DO j=jstrr,jendr
147 DO i=istr,iendr
148 grid(ng)%umask_dia(i,j)=max(0.0_r8, &
149 & min(grid(ng)%umask_full(i,j), &
150 & 1.0_r8))
151 END DO
152 END DO
153 DO j=jstr,jendr
154 DO i=istrr,iendr
155 grid(ng)%vmask_dia(i,j)=max(0.0_r8, &
156 & min(grid(ng)%vmask_full(i,j), &
157 & 1.0_r8))
158 END DO
159 END DO
160# endif
161!
162! Initialize.
163!
164 DO j=jstrr,jendr
165 DO i=istrr,iendr
166 diags(ng)%avgzeta(i,j)=ocean(ng)%zeta(i,j,kout)
167# ifdef WET_DRY
168 diags(ng)%avgzeta(i,j)=diags(ng)%avgzeta(i,j)* &
169 & grid(ng)%rmask_full(i,j)
170# endif
171 END DO
172 END DO
173!
174# if defined DIAGNOSTICS_TS || defined DIAGNOSTICS_UV
175# ifdef DIAGNOSTICS_TS
176 DO idiag=1,ndt
177 DO it=1,nt(ng)
178 DO k=1,n(ng)
179 DO j=jstrr,jendr
180 DO i=istrr,iendr
181 diags(ng)%DiaTrc(i,j,k,it,idiag)= &
182# ifdef WET_DRY
183 & grid(ng)%rmask_full(i,j)* &
184# endif
185 & diags(ng)%DiaTwrk(i,j,k,it,idiag)
186 END DO
187 END DO
188 END DO
189 END DO
190 END DO
191# endif
192# ifdef DIAGNOSTICS_UV
193 DO idiag=1,ndm2d
194 DO j=jstrr,jendr
195 DO i=istr,iendr
196 diags(ng)%DiaU2d(i,j,idiag)= &
197# ifdef WET_DRY
198 & grid(ng)%umask_full(i,j)* &
199# endif
200 & diags(ng)%DiaU2wrk(i,j,idiag)
201 END DO
202 END DO
203 DO j=jstr,jendr
204 DO i=istrr,iendr
205 diags(ng)%DiaV2d(i,j,idiag)= &
206# ifdef WET_DRY
207 & grid(ng)%vmask_full(i,j)* &
208# endif
209 & diags(ng)%DiaV2wrk(i,j,idiag)
210 END DO
211 END DO
212 END DO
213# ifdef SOLVE3D
214 DO idiag=1,ndm3d
215 DO k=1,n(ng)
216 DO j=jstrr,jendr
217 DO i=istr,iendr
218 diags(ng)%DiaU3d(i,j,k,idiag)= &
219# ifdef WET_DRY
220 & grid(ng)%umask_full(i,j)* &
221# endif
222 & diags(ng)%DiaU3wrk(i,j,k,idiag)
223 END DO
224 END DO
225 DO j=jstr,jendr
226 DO i=istrr,iendr
227 diags(ng)%DiaV3d(i,j,k,idiag)= &
228# ifdef WET_DRY
229 & grid(ng)%vmask_full(i,j)* &
230# endif
231 & diags(ng)%DiaV3wrk(i,j,k,idiag)
232 END DO
233 END DO
234 END DO
235 END DO
236# endif
237# endif
238# endif
239!
240!-----------------------------------------------------------------------
241! Accumulate time-averaged fields.
242!-----------------------------------------------------------------------
243!
244 ELSE IF (iic(ng).gt.ntsdia(ng)) THEN
245
246# ifdef WET_DRY
247!
248! If wetting and drying, accumulate wet points counters.
249! points. The time averaged field at each point is only accumulated
250! over wet points since its multiplied by the appropriate mask.
251!
252 DO j=jstr,jendr
253 DO i=istr,iendr
254 grid(ng)%pmask_dia(i,j)=grid(ng)%pmask_dia(i,j)+ &
255 & max(0.0_r8, &
256 & min(grid(ng)%pmask_full(i,j), &
257 & 1.0_r8))
258 END DO
259 END DO
260 DO j=jstrr,jendr
261 DO i=istrr,iendr
262 grid(ng)%rmask_dia(i,j)=grid(ng)%rmask_dia(i,j)+ &
263 & max(0.0_r8, &
264 & min(grid(ng)%rmask_full(i,j), &
265 & 1.0_r8))
266 END DO
267 END DO
268 DO j=jstrr,jendr
269 DO i=istr,iendr
270 grid(ng)%umask_dia(i,j)=grid(ng)%umask_dia(i,j)+ &
271 & max(0.0_r8, &
272 & min(grid(ng)%umask_full(i,j), &
273 & 1.0_r8))
274 END DO
275 END DO
276 DO j=jstr,jendr
277 DO i=istrr,iendr
278 grid(ng)%vmask_dia(i,j)=grid(ng)%vmask_dia(i,j)+ &
279 & max(0.0_r8, &
280 & min(grid(ng)%vmask_full(i,j), &
281 & 1.0_r8))
282 END DO
283 END DO
284# endif
285!
286! Accumulate free-surface.
287!
288 DO j=jstrr,jendr
289 DO i=istrr,iendr
290 diags(ng)%avgzeta(i,j)=diags(ng)%avgzeta(i,j)+ &
291# ifdef WET_DRY
292 & grid(ng)%rmask_full(i,j)* &
293# endif
294 & ocean(ng)%zeta(i,j,kout)
295 END DO
296 END DO
297!
298! Accumulate diagnostic terms.
299!
300# if defined DIAGNOSTICS_TS || defined DIAGNOSTICS_UV
301# ifdef DIAGNOSTICS_TS
302 DO idiag=1,ndt
303 DO it=1,nt(ng)
304 DO k=1,n(ng)
305 DO j=jstrr,jendr
306 DO i=istrr,iendr
307 diags(ng)%DiaTrc(i,j,k,it,idiag)= &
308 & diags(ng)%DiaTrc(i,j,k,it,idiag)+ &
309# ifdef WET_DRY
310 & grid(ng)%rmask_full(i,j)* &
311# endif
312 & diags(ng)%DiaTwrk(i,j,k,it,idiag)
313 END DO
314 END DO
315 END DO
316 END DO
317 END DO
318# endif
319# ifdef DIAGNOSTICS_UV
320 DO idiag=1,ndm2d
321 DO j=jstrr,jendr
322 DO i=istr,iendr
323 diags(ng)%DiaU2d(i,j,idiag)=diags(ng)%DiaU2d(i,j,idiag)+ &
324# ifdef WET_DRY
325 & grid(ng)%umask_full(i,j)* &
326# endif
327 & diags(ng)%DiaU2wrk(i,j,idiag)
328 END DO
329 END DO
330 DO j=jstr,jendr
331 DO i=istrr,iendr
332 diags(ng)%DiaV2d(i,j,idiag)=diags(ng)%DiaV2d(i,j,idiag)+ &
333# ifdef WET_DRY
334 & grid(ng)%vmask_full(i,j)* &
335# endif
336 & diags(ng)%DiaV2wrk(i,j,idiag)
337 END DO
338 END DO
339 END DO
340# ifdef SOLVE3D
341 DO idiag=1,ndm3d
342 DO k=1,n(ng)
343 DO j=jstrr,jendr
344 DO i=istr,iendr
345 diags(ng)%DiaU3d(i,j,k,idiag)= &
346 & diags(ng)%DiaU3d(i,j,k,idiag)+ &
347# ifdef WET_DRY
348 & grid(ng)%umask_full(i,j)* &
349# endif
350 & diags(ng)%DiaU3wrk(i,j,k,idiag)
351 END DO
352 END DO
353 DO j=jstr,jendr
354 DO i=istrr,iendr
355 diags(ng)%DiaV3d(i,j,k,idiag)= &
356 & diags(ng)%DiaV3d(i,j,k,idiag)+ &
357# ifdef WET_DRY
358 & grid(ng)%vmask_full(i,j)* &
359# endif
360 & diags(ng)%DiaV3wrk(i,j,k,idiag)
361 END DO
362 END DO
363 END DO
364 END DO
365# endif
366# endif
367# endif
368 END IF
369!
370!-----------------------------------------------------------------------
371! Convert accumulated sums into time-averages, if appropriate
372!-----------------------------------------------------------------------
373!
374 IF (((iic(ng).gt.ntsdia(ng)).and. &
375 & (mod(iic(ng)-1,ndia(ng)).eq.0).and. &
376 & ((iic(ng).ne.ntstart(ng)).or.(nrrec(ng).eq.0))).or. &
377 & ((iic(ng).ge.ntsdia(ng)).and.(ndia(ng).eq.1))) THEN
378 IF (domain(ng)%SouthWest_Test(tile)) THEN
379 IF (ndia(ng).eq.1) THEN
380 diatime(ng)=time(ng)
381 ELSE
382 diatime(ng)=diatime(ng)+real(ndia(ng),r8)*dt(ng)
383 END IF
384 END IF
385!
386! Set time-averaged factors for each C-grid variable type. Notice that
387! the I- and J-ranges are all grid types are the same for convinience.
388# ifdef WET_DRY
389! In wetting and drying, the sums are devided by the number of times
390! that each qrid point is wet.
391# endif
392!
393# ifdef WET_DRY
394 DO j=jstrr,jendr
395 DO i=istrr,iendr
396 rfac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%rmask_dia(i,j))
397 ufac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%umask_dia(i,j))
398 vfac(i,j)=1.0_r8/max(1.0_r8, grid(ng)%vmask_dia(i,j))
399 END DO
400 END DO
401# else
402 fac=1.0_r8/real(ndia(ng),r8)
403 DO j=jstrr,jendr
404 DO i=istrr,iendr
405 rfac(i,j)=fac
406 ufac(i,j)=fac
407 vfac(i,j)=fac
408 END DO
409 END DO
410# endif
411
412# ifdef WET_DRY
413!
414! If last time-step of average window, convert time dependent counters
415! for wet points to time-averaged Land/Sea masks (dry=0, wet=1) for
416! the current average window period. Notice that a grid point is wet
417! if the count is greater than zero for the current time average
418! window.
419!
420 DO j=jstr,jendr
421 DO i=istr,iendr
422 grid(ng)%pmask_dia(i,j)=min(1.0_r8, grid(ng)%pmask_dia(i,j))
423 END DO
424 END DO
425 DO j=jstrr,jendr
426 DO i=istrr,iendr
427 grid(ng)%rmask_dia(i,j)=min(1.0_r8, grid(ng)%rmask_dia(i,j))
428 END DO
429 END DO
430 DO j=jstrr,jendr
431 DO i=istr,iendr
432 grid(ng)%umask_dia(i,j)=min(1.0_r8, grid(ng)%umask_dia(i,j))
433 END DO
434 END DO
435 DO j=jstr,jendr
436 DO i=istrr,iendr
437 grid(ng)%vmask_dia(i,j)=min(1.0_r8, grid(ng)%vmask_dia(i,j))
438 END DO
439 END DO
440!
441! Exchange boundary data.
442!
443 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
444 CALL exchange_p2d_tile (ng, tile, &
445 & lbi, ubi, lbj, ubj, &
446 & grid(ng)%pmask_dia)
447 CALL exchange_r2d_tile (ng, tile, &
448 & lbi, ubi, lbj, ubj, &
449 & grid(ng)%rmask_dia)
450 CALL exchange_u2d_tile (ng, tile, &
451 & lbi, ubi, lbj, ubj, &
452 & grid(ng)%umask_dia)
453 CALL exchange_v2d_tile (ng, tile, &
454 & lbi, ubi, lbj, ubj, &
455 & grid(ng)%vmask_dia)
456 END IF
457
458# ifdef DISTRIBUTE
459 CALL mp_exchange2d (ng, tile, inlm, 4, &
460 & lbi, ubi, lbj, ubj, &
461 & nghostpoints, &
462 & ewperiodic(ng), nsperiodic(ng), &
463 & grid(ng)%pmask_dia, &
464 & grid(ng)%rmask_dia, &
465 & grid(ng)%umask_dia, &
466 & grid(ng)%vmask_dia)
467# endif
468# endif
469!
470! Convert accumulated sums to time averages.
471!
472 DO j=jstrr,jendr
473 DO i=istrr,iendr
474 diags(ng)%avgzeta(i,j)=rfac(i,j)*diags(ng)%avgzeta(i,j)
475 END DO
476 END DO
477
478# ifdef DIAGNOSTICS_TS
479 DO idiag=1,ndt
480 DO it=1,nt(ng)
481 DO k=1,n(ng)
482 DO j=jstrr,jendr
483 DO i=istrr,iendr
484 diags(ng)%DiaTrc(i,j,k,it,idiag)=rfac(i,j)* &
485 & diags(ng)%DiaTrc(i,j,k,it,idiag)
486 END DO
487 END DO
488 END DO
489 END DO
490 END DO
491# endif
492# ifdef DIAGNOSTICS_BIO
493# if defined BIO_FENNEL || defined HYPOXIA_SRM
494 DO idiag=1,ndbio2d
495 IF (idiag.ne.ipco2) THEN
496 DO j=jstrr,jendr
497 DO i=istrr,iendr
498 diags(ng)%DiaBio2d(i,j,idiag)=rfac(i,j)* &
499 & diags(ng)%DiaBio2d(i,j,idiag)
500 END DO
501 END DO
502 END IF
503 END DO
504# endif
505# if defined BIO_FENNEL
506 DO idiag=1,ndbio3d
507 DO k=1,n(ng)
508 DO j=jstrr,jendr
509 DO i=istrr,iendr
510 diags(ng)%DiaBio3d(i,j,k,idiag)=rfac(i,j)* &
511 & diags(ng)%DiaBio3d(i,j,k,idiag)
512 END DO
513 END DO
514 END DO
515 END DO
516# elif defined ECOSIM
517 DO idiag=1,ndbio3d
518 DO k=1,ndbands
519 DO j=jstrr,jendr
520 DO i=istrr,iendr
521 diags(ng)%DiaBio3d(i,j,k,idiag)=rfac(i,j)* &
522 & diags(ng)%DiaBio3d(i,j,k,idiag)
523 END DO
524 END DO
525 END DO
526 END DO
527 DO idiag=1,ndbio4d
528 DO iband=1,ndbands
529 DO k=1,n(ng)
530 DO j=jstrr,jendr
531 DO i=istrr,iendr
532 diags(ng)%DiaBio4d(i,j,k,iband,idiag)=rfac(i,j)* &
533 & diags(ng)%DiaBio4d(i,j,k,iband,idiag)
534 END DO
535 END DO
536 END DO
537 END DO
538 END DO
539# endif
540# endif
541# ifdef DIAGNOSTICS_UV
542 DO idiag=1,ndm2d
543 DO j=jstrr,jendr
544 DO i=istr,iendr
545 diags(ng)%DiaU2d(i,j,idiag)=ufac(i,j)* &
546 & diags(ng)%DiaU2d(i,j,idiag)
547 END DO
548 END DO
549 DO j=jstr,jendr
550 DO i=istrr,iendr
551 diags(ng)%DiaV2d(i,j,idiag)=vfac(i,j)* &
552 & diags(ng)%DiaV2d(i,j,idiag)
553 END DO
554 END DO
555 END DO
556# ifdef SOLVE3D
557 DO idiag=1,ndm3d
558 DO k=1,n(ng)
559 DO j=jstrr,jendr
560 DO i=istr,iendr
561 diags(ng)%DiaU3d(i,j,k,idiag)=ufac(i,j)* &
562 & diags(ng)%DiaU3d(i,j,k,idiag)
563 END DO
564 END DO
565 DO j=jstr,jendr
566 DO i=istrr,iendr
567 diags(ng)%DiaV3d(i,j,k,idiag)=vfac(i,j)* &
568 & diags(ng)%DiaV3d(i,j,k,idiag)
569 END DO
570 END DO
571 END DO
572 END DO
573# endif
574# endif
575
576# if defined DIAGNOSTICS_TS || defined DIAGNOSTICS_UV || \
577 defined diagnostics_bio
578!
579!-----------------------------------------------------------------------
580! Apply periodic or gradient boundary conditions for output purposes.
581!-----------------------------------------------------------------------
582!
583! Time-averaged free surface.
584!
585 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
586 CALL exchange_r2d_tile (ng, tile, &
587 & lbi, ubi, lbj, ubj, &
588 & diags(ng)%avgzeta)
589 END IF
590# ifdef DISTRIBUTE
591 CALL mp_exchange2d (ng, tile, inlm, 1, &
592 & lbi, ubi, lbj, ubj, &
593 & nghostpoints, &
594 & ewperiodic(ng), nsperiodic(ng), &
595 & diags(ng)%avgzeta)
596# endif
597
598# ifdef DIAGNOSTICS_TS
599!
600! 3D tracer diagnostics.
601!
602 DO idiag=1,ndt
603 DO it=1,nt(ng)
604 CALL bc_r3d_tile (ng, tile, &
605 & lbi, ubi, lbj, ubj, 1, n(ng), &
606 & diags(ng)%DiaTrc(:,:,:,it,idiag))
607 END DO
608# ifdef DISTRIBUTE
609 CALL mp_exchange4d (ng, tile, inlm, 1, &
610 & lbi, ubi, lbj, ubj, 1, n(ng), 1, nt(ng), &
611 & nghostpoints, &
612 & ewperiodic(ng), nsperiodic(ng), &
613 & diags(ng)%DiaTrc(:,:,:,:,idiag))
614# endif
615 END DO
616# endif
617# ifdef DIAGNOSTICS_UV
618!
619! 2D momentum diagnostics.
620!
621 DO idiag=1,ndm2d
622 CALL bc_u2d_tile (ng, tile, &
623 & lbi, ubi, lbj, ubj, &
624 & diags(ng)%DiaU2d(:,:,idiag))
625 CALL bc_v2d_tile (ng, tile, &
626 & lbi, ubi, lbj, ubj, &
627 & diags(ng)%DiaV2d(:,:,idiag))
628 END DO
629# ifdef DISTRIBUTE
630 CALL mp_exchange3d (ng, tile, inlm, 2, &
631 & lbi, ubi, lbj, ubj, 1, ndm2d, &
632 & nghostpoints, &
633 & ewperiodic(ng), nsperiodic(ng), &
634 & diags(ng)%DiaU2d, &
635 & diags(ng)%DiaV2d)
636# endif
637# ifdef SOLVE3D
638!
639! 3D momentum diagnostics.
640!
641 DO idiag=1,ndm3d
642 CALL bc_u3d_tile (ng, tile, &
643 & lbi, ubi, lbj, ubj, 1, n(ng), &
644 & diags(ng)%DiaU3d(:,:,:,idiag))
645 CALL bc_v3d_tile (ng, tile, &
646 & lbi, ubi, lbj, ubj, 1, n(ng), &
647 & diags(ng)%DiaV3d(:,:,:,idiag))
648 END DO
649# ifdef DISTRIBUTE
650 CALL mp_exchange4d (ng, tile, inlm, 2, &
651 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ndm3d, &
652 & nghostpoints, &
653 & ewperiodic(ng), nsperiodic(ng), &
654 & diags(ng)%DiaU3d, &
655 & diags(ng)%DiaV3d)
656# endif
657# endif
658# endif
659# ifdef DIAGNOSTICS_BIO
660# if defined BIO_FENNEL || defined HYPOXIA_SRM
661!
662! Biological terms 2D diagnostics.
663!
664 DO idiag=1,ndbio2d
665 CALL bc_r2d_tile (ng, tile, &
666 & lbi, ubi, lbj, ubj, &
667 & diags(ng)%DiaBio2d(:,:,idiag))
668 END DO
669# ifdef DISTRIBUTE
670 CALL mp_exchange3d (ng, tile, inlm, 1, &
671 & lbi, ubi, lbj, ubj, 1, ndbio2d, &
672 & nghostpoints, &
673 & ewperiodic(ng), nsperiodic(ng), &
674 & diags(ng)%DiaBio2d)
675# endif
676# endif
677# if defined BIO_FENNEL
678!
679! Biological terms 3D diagnostics.
680!
681 DO idiag=1,ndbio3d
682 CALL bc_r3d_tile (ng, tile, &
683 & lbi, ubi, lbj, ubj, 1, n(ng), &
684 & diags(ng)%DiaBio3d(:,:,:,idiag))
685 END DO
686# ifdef DISTRIBUTE
687 CALL mp_exchange4d (ng, tile, inlm, 1, &
688 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ndbio3d, &
689 & nghostpoints, &
690 & ewperiodic(ng), nsperiodic(ng), &
691 & diags(ng)%DiaBio3d)
692# endif
693# elif defined ECOSIM
694!
695! Biological terms 3D diagnostics.
696!
697 DO idiag=1,ndbio3d
698 CALL bc_r3d_tile (ng, tile, &
699 & lbi, ubi, lbj, ubj, 1, ndbands, &
700 & diags(ng)%DiaBio3d(:,:,:,idiag))
701 END DO
702# ifdef DISTRIBUTE
703 CALL mp_exchange4d (ng, tile, inlm, 1, &
704 & lbi, ubi, lbj, ubj, 1, ndbands, 1, ndbio3d, &
705 & nghostpoints, &
706 & ewperiodic(ng), nsperiodic(ng), &
707 & diags(ng)%DiaBio3d)
708# endif
709!
710! Biological terms 4D diagnostics.
711!
712 DO idiag=1,ndbio4d
713 CALL bc_r4d_tile (ng, tile, &
714 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ndbands, &
715 & diags(ng)%DiaBio4d(:,:,:,:,idiag))
716# ifdef DISTRIBUTE
717 CALL mp_exchange4d (ng, tile, inlm, 1, &
718 & lbi, ubi, lbj, ubj, 1, n(ng), 1, ndbands, &
719 & nghostpoints, &
720 & ewperiodic(ng), nsperiodic(ng), &
721 & diags(ng)%DiaBio4d(:,:,:,:,idiag))
722# endif
723 END DO
724# endif
725# endif
726# endif
727 END IF
728!
729 RETURN
subroutine bc_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:44
subroutine bc_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:345
subroutine bc_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition bc_2d.F:167
subroutine bc_u3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
Definition bc_3d.F:187
subroutine bc_r3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
Definition bc_3d.F:48
subroutine bc_v3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
Definition bc_3d.F:389
subroutine bc_r4d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, lbl, ubl, a)
Definition bc_4d.F:50
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_p2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition exchange_2d.F:66
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
integer ipco2
Definition fennel_mod.h:109
integer ndbands
Definition ecosim_mod.h:212
type(t_diags), dimension(:), allocatable diags
Definition mod_diags.F:100
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer ndbio2d
Definition mod_param.F:584
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer ndbio4d
Definition mod_param.F:586
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ndm3d
Definition mod_param.F:579
integer ndt
Definition mod_param.F:574
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer ndbio3d
Definition mod_param.F:585
integer ndm2d
Definition mod_param.F:578
integer, dimension(:), allocatable nrrec
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:), allocatable diatime
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable ndia
integer, dimension(:), allocatable ntsdia
subroutine mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, a, b, c)
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)

References bc_2d_mod::bc_r2d_tile(), bc_3d_mod::bc_r3d_tile(), bc_4d_mod::bc_r4d_tile(), bc_2d_mod::bc_u2d_tile(), bc_3d_mod::bc_u3d_tile(), bc_2d_mod::bc_v2d_tile(), bc_3d_mod::bc_v3d_tile(), mod_diags::diags, mod_scalars::diatime, mod_param::domain, mod_scalars::dt, mod_scalars::ewperiodic, exchange_2d_mod::exchange_p2d_tile(), exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_grid::grid, mod_scalars::iic, mod_param::inlm, mod_biology::ipco2, mp_exchange_mod::mp_exchange2d(), mp_exchange_mod::mp_exchange3d(), mp_exchange_mod::mp_exchange4d(), mod_param::n, mod_biology::ndbands, mod_param::ndbio2d, mod_param::ndbio3d, mod_param::ndbio4d, mod_scalars::ndia, mod_param::ndm2d, mod_param::ndm3d, mod_param::ndt, mod_param::nghostpoints, mod_scalars::nrrec, mod_scalars::nsperiodic, mod_param::nt, mod_scalars::ntsdia, mod_scalars::ntstart, mod_ocean::ocean, and mod_scalars::time.

Referenced by set_diags().

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