60
   61
   63# ifdef DIAGNOSTICS_BIO
   65# endif
   70
   72# ifdef SOLVE3D
   74#  ifdef ECOSIM
   76#  endif
   77# endif
   79# ifdef DISTRIBUTE
   82#  ifdef SOLVE3D
   84#  endif
   85# endif
   86
   87      implicit none
   88
   89
   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
   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
  111
  112
  113      IF (
ndia(ng).eq.0) 
RETURN 
  114
  115
  116
  117
  118
  119
  120
  122     &     (mod(
iic(ng)-1,
ndia(ng)).eq.1)).or.                          &
 
  125 
  126# ifdef WET_DRY
  127
  128
  129
  130
  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
  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
  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
  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
  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
  242
  243
  245 
  246# ifdef WET_DRY
  247
  248
  249
  250
  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
  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
  299
  300# if defined DIAGNOSTICS_TS || defined DIAGNOSTICS_UV
  301#  ifdef DIAGNOSTICS_TS
  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
  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
  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
  372
  373
  375     &     (mod(
iic(ng)-1,
ndia(ng)).eq.0).and.                          &
 
  378        IF (
domain(ng)%SouthWest_Test(tile)) 
THEN 
  379          IF (
ndia(ng).eq.1) 
THEN 
  381          ELSE
  383          END IF
  384        END IF
  385
  386
  387
  388# ifdef WET_DRY
  389
  390
  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
  415
  416
  417
  418
  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
  442
  445     &                            lbi, ubi, lbj, ubj,                   &
  446     &                            
grid(ng)%pmask_dia)
 
  448     &                            lbi, ubi, lbj, ubj,                   &
  449     &                            
grid(ng)%rmask_dia)
 
  451     &                            lbi, ubi, lbj, ubj,                   &
  452     &                            
grid(ng)%umask_dia)
 
  454     &                            lbi, ubi, lbj, ubj,                   &
  455     &                            
grid(ng)%vmask_dia)
 
  456        END IF
  457 
  458#  ifdef DISTRIBUTE
  460     &                      lbi, ubi, lbj, ubj,                         &
  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
  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
  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
  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
  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
  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
  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
  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
  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
  581
  582
  583
  584
  587     &                          lbi, ubi, lbj, ubj,                     &
  589      END IF
  590#  ifdef DISTRIBUTE
  592     &                    lbi, ubi, lbj, ubj,                           &
  596#  endif
  597 
  598#  ifdef DIAGNOSTICS_TS
  599
  600
  601
  605     &                        lbi, ubi, lbj, ubj, 1, 
n(ng),             &
 
  606     &                        
diags(ng)%DiaTrc(:,:,:,it,idiag))
 
  607          END DO
  608#   ifdef DISTRIBUTE
  610     &                        lbi, ubi, lbj, ubj, 1, 
n(ng), 1, 
nt(ng),  &
 
  613     &                        
diags(ng)%DiaTrc(:,:,:,:,idiag))
 
  614#   endif
  615        END DO
  616#  endif
  617#  ifdef DIAGNOSTICS_UV
  618
  619
  620
  623     &                      lbi, ubi, lbj, ubj,                         &
  624     &                      
diags(ng)%DiaU2d(:,:,idiag))
 
  626     &                      lbi, ubi, lbj, ubj,                         &
  627     &                      
diags(ng)%DiaV2d(:,:,idiag))
 
  628        END DO
  629#   ifdef DISTRIBUTE
  631     &                      lbi, ubi, lbj, ubj, 1, 
ndm2d,               &
 
  634     &                      
diags(ng)%DiaU2d,                           &
 
  636#   endif
  637#   ifdef SOLVE3D
  638
  639
  640
  643     &                      lbi, ubi, lbj, ubj, 1, 
n(ng),               &
 
  644     &                      
diags(ng)%DiaU3d(:,:,:,idiag))
 
  646     &                      lbi, ubi, lbj, ubj, 1, 
n(ng),               &
 
  647     &                      
diags(ng)%DiaV3d(:,:,:,idiag))
 
  648        END DO
  649#    ifdef DISTRIBUTE
  651     &                      lbi, ubi, lbj, ubj, 1, 
n(ng), 1, 
ndm3d,     &
 
  654     &                      
diags(ng)%DiaU3d,                           &
 
  656#    endif
  657#   endif
  658#  endif
  659#  ifdef DIAGNOSTICS_BIO
  660#   if defined BIO_FENNEL || defined HYPOXIA_SRM
  661
  662
  663
  666     &                      lbi, ubi, lbj, ubj,                         &
  667     &                      
diags(ng)%DiaBio2d(:,:,idiag))
 
  668        END DO
  669#    ifdef DISTRIBUTE
  671     &                      lbi, ubi, lbj, ubj, 1, 
ndbio2d,             &
 
  674     &                      
diags(ng)%DiaBio2d)
 
  675#    endif
  676#   endif
  677#   if defined BIO_FENNEL
  678
  679
  680
  683     &                      lbi, ubi, lbj, ubj, 1, 
n(ng),               &
 
  684     &                      
diags(ng)%DiaBio3d(:,:,:,idiag))
 
  685        END DO
  686#    ifdef DISTRIBUTE
  688     &                      lbi, ubi, lbj, ubj, 1, 
n(ng), 1, 
ndbio3d,   &
 
  691     &                      
diags(ng)%DiaBio3d)
 
  692#    endif
  693#   elif defined ECOSIM
  694
  695
  696
  699     &                      lbi, ubi, lbj, ubj, 1, 
ndbands,             &
 
  700     &                      
diags(ng)%DiaBio3d(:,:,:,idiag))
 
  701        END DO
  702#    ifdef DISTRIBUTE
  707     &                      
diags(ng)%DiaBio3d)
 
  708#    endif
  709
  710
  711
  714     &                      lbi, ubi, lbj, ubj, 1, 
n(ng), 1, 
ndbands,   &
 
  715     &                      
diags(ng)%DiaBio4d(:,:,:,:,idiag))
 
  716#    ifdef DISTRIBUTE
  718     &                        lbi, ubi, lbj, ubj, 1, 
n(ng), 1, 
ndbands, &
 
  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)
subroutine bc_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine bc_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine bc_u3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
subroutine bc_r3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
subroutine bc_v3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
subroutine bc_r4d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, lbl, ubl, a)
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)
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_diags), dimension(:), allocatable diags
type(t_grid), dimension(:), allocatable grid
type(t_ocean), dimension(:), allocatable ocean
integer, dimension(:), allocatable n
type(t_domain), dimension(:), allocatable domain
integer, dimension(:), allocatable nt
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)