New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 10136 for NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC/mpp_lnk_generic.h90 – NEMO

Ignore:
Timestamp:
2018-09-17T15:16:43+02:00 (6 years ago)
Author:
dguibert
Message:

bull: async/datatype

Experimental changes to enable/study/bench various mpi "optimisations":

  • BULL_ASYNC
  • BULL_DATATYPE_VECTOR/SUBARRAY

this has been applied to the nonosc subroutine (only for now).

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC/mpp_lnk_generic.h90

    r9844 r10136  
    11#if defined MULTI 
    2 #   define NAT_IN(k)                cd_nat(k)    
     2#   define NAT_IN(k)                cd_nat(k) 
    33#   define SGN_IN(k)                psgn(k) 
    44#   define F_SIZE(ptab)             kfld 
     
    99#      define K_SIZE(ptab)             1 
    1010#      define L_SIZE(ptab)             1 
     11#      define _INDEX(i,j,k,l)          (i+((j)+((0)+(0)*ipk)*jpj)*jpi) 
    1112#   endif 
    1213#   if defined DIM_3d 
     
    1516#      define K_SIZE(ptab)             SIZE(ptab(1)%pt3d,3) 
    1617#      define L_SIZE(ptab)             1 
     18#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(0)*ipk)*jpj)*jpi) 
    1719#   endif 
    1820#   if defined DIM_4d 
     
    2123#      define K_SIZE(ptab)             SIZE(ptab(1)%pt4d,3) 
    2224#      define L_SIZE(ptab)             SIZE(ptab(1)%pt4d,4) 
     25#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(l)*ipk)*jpj)*jpi) 
    2326#   endif 
    2427#else 
     
    2730#   define SGN_IN(k)                psgn 
    2831#   define F_SIZE(ptab)             1 
    29 #   define OPT_K(k)                  
     32#   define OPT_K(k) 
    3033#   if defined DIM_2d 
    3134#      define ARRAY_IN(i,j,k,l,f)   ptab(i,j) 
    3235#      define K_SIZE(ptab)          1 
    3336#      define L_SIZE(ptab)          1 
     37#      define _INDEX(i,j,k,l)          (i+((j)+((0)+(0)*ipk)*jpj)*jpi) 
    3438#   endif 
    3539#   if defined DIM_3d 
     
    3741#      define K_SIZE(ptab)          SIZE(ptab,3) 
    3842#      define L_SIZE(ptab)          1 
     43#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(0)*ipk)*jpj)*jpi) 
    3944#   endif 
    4045#   if defined DIM_4d 
     
    4247#      define K_SIZE(ptab)          SIZE(ptab,3) 
    4348#      define L_SIZE(ptab)          SIZE(ptab,4) 
     49#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(l)*ipk)*jpj)*jpi) 
    4450#   endif 
    4551#endif 
    4652 
    4753#if defined MULTI 
    48    SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn, kfld, cd_mpp, pval ) 
     54#if defined ASYNC 
     55   SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn, loop_fct, kfld, cd_mpp, pval ) 
    4956      INTEGER                     , INTENT(in   ) ::   kfld        ! number of pt3d arrays 
    5057#else 
    51    SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn      , cd_mpp, pval ) 
     58   SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn          , kfld, cd_mpp, pval ) 
     59      INTEGER                     , INTENT(in   ) ::   kfld        ! number of pt3d arrays 
     60#endif 
     61#else 
     62#if defined ASYNC 
     63   SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn, loop_fct,       cd_mpp, pval ) 
     64#else 
     65   SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn,                 cd_mpp, pval ) 
     66#endif 
     67#endif 
     68#ifdef SCOREP_USER_ENABLE 
     69#include "scorep/SCOREP_User.inc" 
     70#else 
     71#define SCOREP_USER_REGION_BEGIN ! 
     72#define SCOREP_USER_REGION_END   ! 
    5273#endif 
    5374      ARRAY_TYPE(:,:,:,:,:)                                        ! array or pointer of arrays on which the boundary condition is applied 
     
    5677      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp      ! fill the overlap area only 
    5778      REAL(wp)        , OPTIONAL  , INTENT(in   ) ::   pval        ! background value (used at closed boundaries) 
     79#ifdef ASYNC 
     80      interface 
     81        subroutine loop_fct(i0, i1, j0, j1, k0, k1, buf) 
     82          integer, intent(in) :: i0, i1, j0, j1, k0, k1 
     83          ! @BULL_FIXME 
     84          ! lib_mpp.f90(4209): error #6683: A kind type parameter must be a compile-time constant.   [WP] 
     85          !          REAL(wp), dimension(:), optional, intent(out) :: buf 
     86          REAL*8, dimension(:,:,:,:,:,:), optional, intent(out) :: buf 
     87        end subroutine loop_fct 
     88      end interface 
     89#endif 
    5890      CHARACTER(len=*),             INTENT(in   ) ::   rname       ! name of the calling subroutine 
    5991      ! 
     
    6698      REAL(wp), DIMENSION(:,:,:,:,:,:), ALLOCATABLE ::   zt3ns, zt3sn   ! north-south & south-north  halos 
    6799      REAL(wp), DIMENSION(:,:,:,:,:,:), ALLOCATABLE ::   zt3ew, zt3we   ! east -west  & west - east  halos 
     100#ifdef ASYNC 
     101      integer :: iflag, i 
     102      logical :: finished 
     103#if (defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     104      integer :: ml_reqs(8,F_SIZE(ptab)) 
     105#else 
     106      integer :: ml_reqs(8) 
     107#endif 
     108#endif 
     109#ifdef SCOREP_USER_ENABLE 
     110      integer :: ier 
     111      SCOREP_USER_REGION_DEFINE( reg_cb ) 
     112      SCOREP_USER_REGION_DEFINE( reg_cbWhole ) 
     113      SCOREP_USER_REGION_DEFINE( reg_cbWE ) 
     114      SCOREP_USER_REGION_DEFINE( reg_cbNS ) 
     115      SCOREP_USER_REGION_DEFINE( reg_cbCenter ) 
     116      SCOREP_USER_REGION_DEFINE( reg_pack ) 
     117      SCOREP_USER_REGION_DEFINE( reg_unpack ) 
     118#if (defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     119      SCOREP_USER_REGION_DEFINE( reg_datatype ) 
     120#endif 
     121#endif 
     122#ifdef MPI_DATATYPE_VECTOR 
     123      integer :: type_ns, type_ew 
     124#endif 
     125#ifdef MPI_DATATYPE_SUBARRAY 
     126      integer :: ndims 
     127      integer, dimension(4) :: array_of_sizes 
     128      integer, dimension(4) :: array_of_subsizes 
     129      integer, dimension(4) :: array_of_starts 
     130      integer :: type_north_halo, type_north_ghost 
     131      integer :: type_south_halo, type_south_ghost 
     132      integer :: type_west_halo, type_west_ghost 
     133      integer :: type_east_halo, type_east_ghost 
     134#endif 
     135      real*8 :: t0 
     136      real*8, save :: time=0.0 
     137 
     138#ifdef ASYNC 
     139      ml_reqs = MPI_REQUEST_NULL 
     140#endif 
    68141      !!---------------------------------------------------------------------- 
    69142      ! 
     
    72145      ipf = F_SIZE(ptab)   ! 5th    -      use in "multi" case (array of pointers) 
    73146      ! 
     147#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
    74148      ALLOCATE( zt3ns(jpi,nn_hls,ipk,ipl,ipf,2), zt3sn(jpi,nn_hls,ipk,ipl,ipf,2),   & 
    75149         &      zt3ew(jpj,nn_hls,ipk,ipl,ipf,2), zt3we(jpj,nn_hls,ipk,ipl,ipf,2)  ) 
     150#endif 
    76151      ! 
    77152      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value 
     
    79154      ENDIF 
    80155 
     156#ifndef ASYNC 
    81157      ! ------------------------------- ! 
    82158      !   standard boundary treatment   !    ! CAUTION: semi-column notation is often impossible 
     
    133209      ! we play with the neigbours AND the row number because of the periodicity 
    134210      ! 
     211      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
    135212      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
    136213      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case) 
     
    147224         END DO 
    148225      END SELECT 
     226      SCOREP_USER_REGION_END( reg_pack ) 
    149227      ! 
    150228      !                           ! Migrations 
     
    210288      iihom = nlci-nn_hls 
    211289      ! 
     290      SCOREP_USER_REGION_BEGIN( reg_unpack, "unpack", SCOREP_USER_REGION_TYPE_COMMON ) 
    212291      SELECT CASE ( nbondi ) 
    213292      CASE ( -1 ) 
     
    243322         END DO 
    244323      END SELECT 
     324      SCOREP_USER_REGION_END( reg_unpack ) 
    245325 
    246326      ! 3. North and south directions 
     
    248328      ! always closed : we play only with the neigbours 
    249329      ! 
     330      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
    250331      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions 
    251332         ijhom = nlcj-nrecj 
     
    261342         END DO 
    262343      ENDIF 
     344#ifdef SCOREP_USER_ENABLE 
     345      SCOREP_USER_REGION_END( reg_pack ) 
     346#endif 
    263347      ! 
    264348      !                           ! Migrations 
     
    271355         ! 
    272356         CALL tic_tac(.TRUE.) 
    273          !  
     357         ! 
    274358         SELECT CASE ( nbondj ) 
    275359         CASE ( -1 ) 
     
    295379      ! 
    296380      !                           ! Write Dirichlet lateral conditions 
     381      SCOREP_USER_REGION_BEGIN( reg_unpack, "unpack", SCOREP_USER_REGION_TYPE_COMMON ) 
    297382      ijhom = nlcj-nn_hls 
    298383      ! 
     
    330415         END DO 
    331416      END SELECT 
    332  
     417      SCOREP_USER_REGION_END( reg_unpack ) 
     418#else 
     419! ASYNC implementation 
     420 
     421! prepare receptions 
     422   !SUBROUTINE mpprecv( ktyp, pmess, kbytes, ksource ) 
     423   !CALL mpi_recv( pmess, kbytes, mpi_double_precision, use_source, ktyp, mpi_comm_oce, istatus, iflag ) 
     424#ifdef MPI_DATATYPE_VECTOR 
     425!      IF( ln_timing )  t0=MPI_Wtime()  
     426      SCOREP_USER_REGION_BEGIN( reg_datatype, "datatype vector", SCOREP_USER_REGION_TYPE_COMMON ) 
     427      ! int MPI_Type_vector(int count, 
     428      !                 int blocklength, 
     429      !                 int stride, 
     430      !                 MPI_Datatype old_type, 
     431      !                 MPI_Datatype *newtype_p)         
     432#ifdef DIM_2d 
     433      ! NS 
     434      call MPI_Type_contiguous((jpi-2*nn_hls), MPI_DOUBLE_PRECISION, type_ns, iflag) 
     435      call MPI_Type_commit(type_ns, iflag) 
     436      ! EW 
     437      call MPI_Type_vector((jpj-2*nn_hls), nn_hls, jpi, MPI_DOUBLE_PRECISION, type_ew, iflag) 
     438      call MPI_Type_commit(type_ew, iflag) 
     439#endif 
     440#   if (defined DIM_3d || defined DIM_4d) 
     441      ! NS 
     442      call MPI_Type_vector(nn_hls             *ipk*ipl,          (jpi-2*nn_hls), jpi*jpj, MPI_DOUBLE_PRECISION, type_ns, iflag) 
     443      call MPI_Type_commit(type_ns, iflag) 
     444      ! EW 
     445      call MPI_Type_vector(       (jpj-2*nn_hls)*ipk*ipl, nn_hls               , jpi,     MPI_DOUBLE_PRECISION, type_ew, iflag) 
     446      call MPI_Type_commit(type_ew, iflag) 
     447#endif 
     448      SCOREP_USER_REGION_END( reg_datatype ) 
     449!      IF( ln_timing ) time=time+MPI_Wtime()-t0 
     450!      IF( ln_timing ) write(*,*) 'timing datatype vector: ',time 
     451 
     452      DO jf = 1, ipf 
     453         SELECT CASE ( nbondi ) 
     454         CASE ( -1 ) 
     455           call mpi_irecv(ARRAY_IN(1,2,1,1,jf), 1, type_ew, noea, 8*jf+1, mpi_comm_oce, ml_reqs(1,jf), iflag) 
     456         CASE ( 0 ) 
     457           call mpi_irecv(ARRAY_IN(1,2,1,1,jf), 1, type_ew, noea, 8*jf+1, mpi_comm_oce, ml_reqs(1,jf), iflag) 
     458           call mpi_irecv(ARRAY_IN(jpi-nn_hls,2,1,1,jf), 1, type_ew, nowe, 8*jf+2, mpi_comm_oce, ml_reqs(2,jf), iflag) 
     459         CASE ( 1 ) 
     460           call mpi_irecv(ARRAY_IN(jpi-nn_hls,2,1,1,jf), 1, type_ew, nowe, 8*jf+2, mpi_comm_oce, ml_reqs(2,jf), iflag) 
     461         END SELECT 
     462 
     463         SELECT CASE ( nbondj ) 
     464         CASE ( -1 ) 
     465           call mpi_irecv(ARRAY_IN(2,jpj-nn_hls,1,1,jf), 1, type_ns, nono, 8*jf+3, mpi_comm_oce, ml_reqs(3,jf), iflag) 
     466         CASE ( 0 ) 
     467           call mpi_irecv(ARRAY_IN(2,jpj-nn_hls,1,1,jf), 1, type_ns, nono, 8*jf+3, mpi_comm_oce, ml_reqs(3,jf), iflag) 
     468           call mpi_irecv(ARRAY_IN(2,1,1,1,jf), 1, type_ns, noso, 8*jf+4, mpi_comm_oce, ml_reqs(4,jf), iflag) 
     469         CASE ( 1 ) 
     470           call mpi_irecv(ARRAY_IN(2,1,1,1,jf), 1, type_ns, noso, 8*jf+4, mpi_comm_oce, ml_reqs(4,jf), iflag) 
     471         END SELECT 
     472      end do 
     473#endif 
     474 
     475#ifdef MPI_DATATYPE_SUBARRAY 
     476      IF( ln_timing )   CALL timing_start('datatype subarray') 
     477      SCOREP_USER_REGION_BEGIN( reg_datatype, "datatype", SCOREP_USER_REGION_TYPE_COMMON ) 
     478 
     479      array_of_sizes = (/ jpi, jpj, ipk, ipl /) 
     480      array_of_subsizes(3:4) = (/ ipk, ipl /) 
     481      array_of_starts(3:4) = 0 
     482#   if defined DIM_2d 
     483      ndims = 2 
     484#   endif 
     485#   if defined DIM_3d 
     486      ndims = 3 
     487#   endif 
     488#   if defined DIM_4d 
     489      ndims = 4 
     490#   endif 
     491      ! ------------------------------- ! 
     492      !      East and west exchange     ! 
     493      ! ------------------------------- ! 
     494      array_of_subsizes(1:2) = (/ nn_hls, jpj-2*nn_hls /) 
     495 
     496      array_of_starts(1:2) = (/ 1, 1 /) ! zero indexing (as in C) 
     497      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_west_halo, iflag) 
     498      call MPI_Type_commit(type_west_halo, iflag) 
     499      array_of_starts(1:2) = (/ 0, 1 /) ! zero indexing (as in C) 
     500      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_west_ghost, iflag) 
     501      call MPI_Type_commit(type_west_ghost, iflag) 
     502 
     503      array_of_starts(1:2) = (/ jpi-1-nn_hls, 1 /) ! zero indexing (as in C) 
     504      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_east_halo, iflag) 
     505      call MPI_Type_commit(type_east_halo, iflag) 
     506      array_of_starts(1:2) = (/ jpi-nn_hls, 1 /) ! zero indexing (as in C) 
     507      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_east_ghost, iflag) 
     508      call MPI_Type_commit(type_east_ghost, iflag) 
     509 
     510      ! ------------------------------- ! 
     511      !      North and south exchange     ! 
     512      ! ------------------------------- ! 
     513      array_of_subsizes(1:2) = (/ jpi-2*nn_hls, nn_hls /) 
     514 
     515      array_of_starts(1:2) = (/ 1, 1 /) ! zero indexing (as in C) 
     516      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_south_halo, iflag) 
     517      call MPI_Type_commit(type_south_halo, iflag) 
     518      array_of_starts(1:2) = (/ 1, 0 /) ! zero indexing (as in C) 
     519      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_south_ghost, iflag) 
     520      call MPI_Type_commit(type_south_ghost, iflag) 
     521 
     522      array_of_starts(1:2) = (/ 1, jpj-1-nn_hls /) ! zero indexing (as in C) 
     523      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_north_halo, iflag) 
     524      call MPI_Type_commit(type_north_halo, iflag) 
     525      array_of_starts(1:2) = (/ 1, jpj-nn_hls /) ! zero indexing (as in C) 
     526      call MPI_Type_create_subarray( ndims, array_of_sizes, array_of_subsizes, array_of_starts, MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, type_north_ghost, iflag) 
     527      call MPI_Type_commit(type_north_ghost, iflag) 
     528#ifdef SCOREP_USER_ENABLE 
     529      SCOREP_USER_REGION_END( reg_datatype ) 
     530#endif 
     531      IF( ln_timing )   CALL timing_stop('datatype subarray') 
     532 
     533      DO jf = 1, ipf 
     534         SELECT CASE ( nbondi ) 
     535         CASE ( -1 ) 
     536           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_east_ghost, noea, 8*jf+1, mpi_comm_oce, ml_reqs(1,jf), iflag) 
     537         CASE ( 0 ) 
     538           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_east_ghost, noea, 8*jf+1, mpi_comm_oce, ml_reqs(1,jf), iflag) 
     539           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_west_ghost, nowe, 8*jf+2, mpi_comm_oce, ml_reqs(2,jf), iflag) 
     540         CASE ( 1 ) 
     541           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_west_ghost, nowe, 8*jf+2, mpi_comm_oce, ml_reqs(2,jf), iflag) 
     542         END SELECT 
     543 
     544         SELECT CASE ( nbondj ) 
     545         CASE ( -1 ) 
     546           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_north_ghost, nono, 8*jf+3, mpi_comm_oce, ml_reqs(3,jf), iflag) 
     547         CASE ( 0 ) 
     548           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_north_ghost, nono, 8*jf+3, mpi_comm_oce, ml_reqs(3,jf), iflag) 
     549           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_south_ghost, noso, 8*jf+4, mpi_comm_oce, ml_reqs(4,jf), iflag) 
     550         CASE ( 1 ) 
     551           call mpi_irecv(ARRAY_IN(:,:,:,:,jf), 1, type_south_ghost, noso, 8*jf+4, mpi_comm_oce, ml_reqs(4,jf), iflag) 
     552         END SELECT 
     553      end do 
     554#endif 
     555#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     556      !                           ! Migrations 
     557      imigr = nn_hls * jpj * ipk * ipl * ipf 
     558      ! 
     559      SELECT CASE ( nbondi ) 
     560      CASE ( -1 ) 
     561        call mpi_irecv(zt3ew(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, noea, 1, mpi_comm_oce, ml_reqs(1), iflag) 
     562      CASE ( 0 ) 
     563        call mpi_irecv(zt3ew(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, noea, 1, mpi_comm_oce, ml_reqs(1), iflag) 
     564        call mpi_irecv(zt3we(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, nowe, 2, mpi_comm_oce, ml_reqs(2), iflag) 
     565      CASE ( 1 ) 
     566        call mpi_irecv(zt3we(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, nowe, 2, mpi_comm_oce, ml_reqs(2), iflag) 
     567      END SELECT 
     568 
     569      imigr = nn_hls * jpi * ipk * ipl * ipf 
     570      SELECT CASE ( nbondj ) 
     571      CASE ( -1 ) 
     572        call mpi_irecv(zt3ns(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, nono, 3, mpi_comm_oce, ml_reqs(3), iflag) 
     573      CASE ( 0 ) 
     574        call mpi_irecv(zt3ns(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, nono, 3, mpi_comm_oce, ml_reqs(3), iflag) 
     575        call mpi_irecv(zt3sn(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, noso, 4, mpi_comm_oce, ml_reqs(4), iflag) 
     576      CASE ( 1 ) 
     577        call mpi_irecv(zt3sn(1,1,1,1,1,2), imigr, MPI_DOUBLE_PRECISION, noso, 4, mpi_comm_oce, ml_reqs(4), iflag) 
     578      END SELECT 
     579#endif 
     580 
     581! compute West 
     582#define TI 1 
     583#define TJ 1 
     584 
     585#define I0 2 
     586#define I1 jpi-1 
     587#define J0 2 
     588#define J1 jpj-1 
     589 
     590#define FULL_ROWS (I0 == 2 && I1 == jpi-1) 
     591#define FULL_COLUMNS (J0 == 2 && J1 == jpi-1) 
     592#define WHOLE_RANGE (FULL_ROWS && FULL_COLUMNS) 
     593 
     594#if (FULL_ROWS && FULL_COLUMNS) 
     595#warning "BULL: lib_mpp will compute whole cb " 
     596      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     597      SCOREP_USER_REGION_BEGIN( reg_cbWhole, "cb whole", SCOREP_USER_REGION_TYPE_COMMON ) 
     598      call loop_fct( I0, I1 & 
     599                   , J0, J1  & ! stand for 3,jpjm2 
     600                   , 1, jpkm1 & ! TODO check if always jpkm1 
     601                   ) 
     602      SCOREP_USER_REGION_END( reg_cbWhole ) 
     603      SCOREP_USER_REGION_END( reg_cb ) 
     604#endif 
     605 
     606#if !FULL_COLUMNS 
     607#warning "BULL: lib_mpp will compute cb S" 
     608      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     609      SCOREP_USER_REGION_BEGIN( reg_cbNS, "cbns", SCOREP_USER_REGION_TYPE_COMMON ) 
     610! asynchrously send South 
     611      call loop_fct( I0, I1 & 
     612                   , J0-1, J0-1 & 
     613                   , 1, jpkm1 & ! TODO check if always jpkm1 
     614                   ) 
     615      SCOREP_USER_REGION_END( reg_cbNS ) 
     616      SCOREP_USER_REGION_END( reg_cb ) 
     617#endif 
     618      ! 3. South directions 
     619      ! ----------------------------- 
     620      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
     621#ifdef MPI_DATATYPE_SUBARRAY 
     622      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     623      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     624        DO jf = 1, ipf 
     625#ifdef BULL_ISEND 
     626          call mpi_isend(ARRAY_IN(:,:,:,:,jf), 1, type_south_halo, noso, 8*jf+3, mpi_comm_oce, ml_reqs(4+3,jf), iflag) 
     627#else 
     628          call mpi_send(ARRAY_IN(:,:,:,:,jf), 1, type_south_halo, noso, 8*jf+3, mpi_comm_oce, iflag) 
     629#endif 
     630        END DO 
     631      END SELECT 
     632#elif (defined MPI_DATATYPE_VECTOR) 
     633      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     634      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     635        DO jf = 1, ipf 
     636#ifdef BULL_ISEND 
     637          call mpi_isend(ARRAY_IN(2,2,1,1,jf), 1, type_ns, noso, 8*jf+3, mpi_comm_oce, ml_reqs(4+3,jf), iflag) 
     638#else 
     639          call mpi_send(ARRAY_IN(2,2,1,1,jf), 1, type_ns, noso, 8*jf+3, mpi_comm_oce, iflag) 
     640#endif 
     641        END DO 
     642      END SELECT 
     643#else 
     644      ! always closed : we play only with the neigbours 
     645      ! 
     646      imigr = nn_hls * jpi * ipk * ipl * ipf 
     647      SELECT CASE ( nbondj )      ! Read Dirichlet lateral conditions 
     648      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     649         ijhom = nlcj-nrecj 
     650         DO jf = 1, ipf 
     651            DO jl = 1, ipl 
     652               DO jk = 1, ipk 
     653                  DO jh = 1, nn_hls 
     654                     zt3ns(:,jh,jk,jl,jf,1) = ARRAY_IN(:,nn_hls+jh,jk,jl,jf) 
     655                  END DO 
     656               END DO 
     657            END DO 
     658         END DO 
     659         call mpi_isend(zt3ns(1,1,1,1,1,1), imigr, MPI_DOUBLE_PRECISION, noso, 3, mpi_comm_oce, ml_reqs(4+3), iflag) 
     660      END SELECT 
     661#endif 
     662      SCOREP_USER_REGION_END( reg_pack ) 
     663 
     664      ! progress all previous operations 
     665#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     666      call MPI_Testall(8, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     667#else 
     668      call MPI_Testall(8*ipf, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     669#endif 
     670  
     671! compute North 
     672#if !FULL_COLUMNS 
     673#warning "BULL: lib_mpp will compute cb N" 
     674      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     675      SCOREP_USER_REGION_BEGIN( reg_cbNS, "cbns", SCOREP_USER_REGION_TYPE_COMMON ) 
     676      call loop_fct( I0, I1   & 
     677                   , J1+1, J1+1   & 
     678                   , 1, jpkm1 & ! TODO check if always jpkm1 
     679                   ) 
     680      SCOREP_USER_REGION_END( reg_cbNS ) 
     681      SCOREP_USER_REGION_END( reg_cb ) 
     682#endif 
     683      ! 3. North directions 
     684      ! ----------------------------- 
     685      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
     686#ifdef MPI_DATATYPE_SUBARRAY 
     687      ! always closed : we play only with the neigbours 
     688 
     689      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     690      CASE ( -1, 0 )                ! all exept 2 (i.e. close case) 
     691        DO jf = 1, ipf 
     692#ifdef BULL_ISEND 
     693          call mpi_isend(ARRAY_IN(:,:,:,:,jf), 1, type_north_halo, nono, 8*jf+4, mpi_comm_oce, ml_reqs(4+4,jf), iflag) 
     694#else 
     695          call mpi_send(ARRAY_IN(:,:,:,:,jf), 1, type_north_halo, nono, 8*jf+4, mpi_comm_oce, iflag) 
     696#endif 
     697        END DO 
     698      END SELECT 
     699#elif (defined MPI_DATATYPE_VECTOR) 
     700      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     701      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     702        DO jf = 1, ipf 
     703#ifdef BULL_ISEND 
     704          call mpi_isend(ARRAY_IN(2,jpj-nn_hls,1,1,jf), 1, type_ns, nono, 8*jf+4, mpi_comm_oce, ml_reqs(4+3,jf), iflag) 
     705#else 
     706          call mpi_send(ARRAY_IN(2,jpj-nn_hls,1,1,jf), 1, type_ns, nono, 8*jf+4, mpi_comm_oce, iflag) 
     707#endif 
     708        END DO 
     709      END SELECT 
     710#else 
     711      ! 
     712      imigr = nn_hls * jpi * ipk * ipl * ipf 
     713      SELECT CASE ( nbondj )      ! Read Dirichlet lateral conditions 
     714      CASE ( -1, 0 )                ! all exept 2 (i.e. close case) 
     715         ijhom = nlcj-nrecj ! jpj-2*nn_hls 
     716         DO jf = 1, ipf 
     717            DO jl = 1, ipl 
     718               DO jk = 1, ipk 
     719                  DO jh = 1, nn_hls 
     720                     zt3sn(:,jh,jk,jl,jf,1) = ARRAY_IN(:,ijhom +jh,jk,jl,jf) 
     721                  END DO 
     722               END DO 
     723            END DO 
     724         END DO 
     725         call mpi_isend(zt3sn(1,1,1,1,1,1), imigr, MPI_DOUBLE_PRECISION, nono, 4, mpi_comm_oce, ml_reqs(4+4), iflag) 
     726      END SELECT 
     727#endif 
     728      SCOREP_USER_REGION_END( reg_pack ) 
     729 
     730      ! progress all previous operations 
     731#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     732      call MPI_Testall(8, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     733#else 
     734      call MPI_Testall(8*ipf, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     735#endif 
     736 
     737#if !FULL_ROWS 
     738#warning "BULL: lib_mpp will compute cb W" 
     739      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     740      SCOREP_USER_REGION_BEGIN( reg_cbWE, "cbew", SCOREP_USER_REGION_TYPE_COMMON ) 
     741      call loop_fct( I0-1, I0-1 & 
     742                   , J0, J1  & ! stand for 3,jpjm2 
     743                   , 1, jpkm1 & ! TODO check if always jpkm1 
     744                   ) 
     745      SCOREP_USER_REGION_END( reg_cbWE ) 
     746      SCOREP_USER_REGION_END( reg_cb ) 
     747#endif 
     748      ! ------------------------------- ! 
     749      !      West exchange     ! 
     750      ! ------------------------------- ! 
     751      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
     752#ifdef MPI_DATATYPE_SUBARRAY 
     753      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     754      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     755        DO jf = 1, ipf 
     756#ifdef BULL_ISEND 
     757           call mpi_isend(ARRAY_IN(:,:,:,:,jf), 1, type_west_halo, nowe, 8*jf+1, mpi_comm_oce, ml_reqs(4+1,jf), iflag) 
     758#else 
     759           call mpi_send(ARRAY_IN(:,:,:,:,jf), 1, type_west_halo, nowe, 8*jf+1, mpi_comm_oce, iflag) 
     760#endif 
     761        END DO 
     762      END SELECT 
     763#elif (defined MPI_DATATYPE_VECTOR) 
     764      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     765      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     766        DO jf = 1, ipf 
     767#ifdef BULL_ISEND 
     768          call mpi_isend(ARRAY_IN(2,2,1,1,jf), 1, type_ew, nowe, 8*jf+1, mpi_comm_oce, ml_reqs(4+3,jf), iflag) 
     769#else 
     770          call mpi_send(ARRAY_IN(2,2,1,1,jf), 1, type_ew, nowe, 8*jf+1, mpi_comm_oce, iflag) 
     771#endif 
     772        END DO 
     773      END SELECT 
     774#else 
     775      ! we play with the neigbours AND the row number because of the periodicity 
     776      ! 
     777      imigr = nn_hls * jpj * ipk * ipl * ipf 
     778      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     779      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     780        DO jf = 1, ipf 
     781           DO jl = 1, ipl 
     782              DO jk = 1, ipk 
     783                 DO jh = 1, nn_hls 
     784                    zt3ew(:,jh,jk,jl,jf,1) = ARRAY_IN(nn_hls+jh,:,jk,jl,jf) 
     785                 END DO 
     786              END DO 
     787           END DO 
     788        END DO 
     789        call mpi_isend(zt3ew(1,1,1,1,1,1), imigr, MPI_DOUBLE_PRECISION, nowe, 1, mpi_comm_oce, ml_reqs(4+1), iflag) 
     790      END SELECT 
     791#endif 
     792      SCOREP_USER_REGION_END( reg_pack ) 
     793 
     794      ! progress all previous operations 
     795#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     796      call MPI_Testall(8, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     797#else 
     798      call MPI_Testall(8*ipf, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     799#endif 
     800 
     801! compute East 
     802#if !FULL_ROWS 
     803#warning "BULL: lib_mpp will compute cb E" 
     804      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     805      SCOREP_USER_REGION_BEGIN( reg_cbWE, "cbew", SCOREP_USER_REGION_TYPE_COMMON ) 
     806      call loop_fct( I1+1, I1+1 & 
     807                   , J0, J1  & ! stand for 3,jpjm2 
     808                   , 1, jpkm1 & ! TODO check if always jpkm1 
     809                   ) 
     810      SCOREP_USER_REGION_END( reg_cbWE ) 
     811      SCOREP_USER_REGION_END( reg_cb ) 
     812#endif 
     813      ! ------------------------------- ! 
     814      !      East exchange     ! 
     815      ! ------------------------------- ! 
     816      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON ) 
     817#ifdef MPI_DATATYPE_SUBARRAY 
     818      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     819      CASE ( -1, 0 )                ! all exept 2 (i.e. close case) 
     820        DO jf = 1, ipf 
     821#ifdef BULL_ISEND 
     822          call mpi_isend(ARRAY_IN(:,:,:,:,jf), 1, type_east_halo, noea, 8*jf+2, mpi_comm_oce, ml_reqs(4+2,jf), iflag) 
     823#else 
     824          call mpi_send(ARRAY_IN(:,:,:,:,jf), 1, type_east_halo, noea, 8*jf+2, mpi_comm_oce, iflag) 
     825#endif 
     826        END DO 
     827      END SELECT 
     828#elif (defined MPI_DATATYPE_VECTOR) 
     829      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     830      CASE ( 0, 1 )                ! all exept 2 (i.e. close case) 
     831        DO jf = 1, ipf 
     832#ifdef BULL_ISEND 
     833          call mpi_isend(ARRAY_IN(jpi-nn_hls,2,1,1,jf), 1, type_ew, noea, 8*jf+2, mpi_comm_oce, ml_reqs(4+3,jf), iflag) 
     834#else 
     835          call mpi_send(ARRAY_IN(jpi-nn_hls,2,1,1,jf), 1, type_ew, noea, 8*jf+2, mpi_comm_oce, iflag) 
     836#endif 
     837        END DO 
     838      END SELECT 
     839#else 
     840      ! we play with the neigbours AND the row number because of the periodicity 
     841      ! 
     842      imigr = nn_hls * jpj * ipk * ipl * ipf 
     843      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     844      CASE ( -1, 0 )                ! all exept 2 (i.e. close case) 
     845         iihom = nlci-nreci ! jpi-2*nn_hls 
     846         DO jf = 1, ipf 
     847            DO jl = 1, ipl 
     848               DO jk = 1, ipk 
     849                  DO jh = 1, nn_hls 
     850                     zt3we(:,jh,jk,jl,jf,1) = ARRAY_IN(iihom +jh,:,jk,jl,jf) 
     851                  END DO 
     852               END DO 
     853            END DO 
     854         END DO 
     855         call mpi_isend(zt3we(1,1,1,1,1,1), imigr, MPI_DOUBLE_PRECISION, noea, 2, mpi_comm_oce, ml_reqs(4+2), iflag) 
     856      END SELECT 
     857#endif 
     858      SCOREP_USER_REGION_END( reg_pack ) 
     859 
     860      ! progress all previous operations 
     861#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     862      call MPI_Testall(8, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     863#else 
     864      call MPI_Testall(8*ipf, ml_reqs, finished, MPI_STATUSES_IGNORE, iflag) 
     865#endif 
     866 
     867! compute Inner 
     868#if !(FULL_ROWS && FULL_COLUMNS) 
     869#warning "BULL: lib_mpp will compute inner cb" 
     870      SCOREP_USER_REGION_BEGIN( reg_cb, "cb", SCOREP_USER_REGION_TYPE_COMMON ) 
     871      SCOREP_USER_REGION_BEGIN( reg_cbCenter, "cbcenter", SCOREP_USER_REGION_TYPE_COMMON ) 
     872      call loop_fct( I0, I1 & 
     873                   , J0, J1  & ! stand for 3,jpjm2 
     874                   , 1, jpkm1 & ! TODO check if always jpkm1 
     875                   ) 
     876      SCOREP_USER_REGION_END( reg_cbCenter ) 
     877      SCOREP_USER_REGION_END( reg_cb ) 
     878#endif 
     879 
     880      ! ------------------------------- ! 
     881      !   standard boundary treatment   !    ! CAUTION: semi-column notation is often impossible 
     882      ! ------------------------------- ! 
     883      ! 
     884      IF( PRESENT( cd_mpp ) ) THEN     !==  halos filled with inner values  ==! 
     885         ! 
     886         DO jf = 1, ipf                      ! number of arrays to be treated 
     887            ! 
     888            DO jl = 1, ipl                   ! CAUTION: ptab is defined only between nld and nle 
     889               DO jk = 1, ipk 
     890                  DO jj = nlcj+1, jpj                 ! added line(s)   (inner only) 
     891                     ARRAY_IN(nldi  :nlei  ,jj,jk,jl,jf) = ARRAY_IN(nldi:nlei,nlej,jk,jl,jf) 
     892                     ARRAY_IN(1     :nldi-1,jj,jk,jl,jf) = ARRAY_IN(nldi     ,nlej,jk,jl,jf) 
     893                     ARRAY_IN(nlei+1:nlci  ,jj,jk,jl,jf) = ARRAY_IN(     nlei,nlej,jk,jl,jf) 
     894                  END DO 
     895                  DO ji = nlci+1, jpi                 ! added column(s) (full) 
     896                     ARRAY_IN(ji,nldj  :nlej  ,jk,jl,jf) = ARRAY_IN(nlei,nldj:nlej,jk,jl,jf) 
     897                     ARRAY_IN(ji,1     :nldj-1,jk,jl,jf) = ARRAY_IN(nlei,nldj     ,jk,jl,jf) 
     898                     ARRAY_IN(ji,nlej+1:jpj   ,jk,jl,jf) = ARRAY_IN(nlei,     nlej,jk,jl,jf) 
     899                  END DO 
     900               END DO 
     901            END DO 
     902            ! 
     903         END DO 
     904         ! 
     905      ELSE                              !==  standard close or cyclic treatment  ==! 
     906         ! 
     907         DO jf = 1, ipf                      ! number of arrays to be treated 
     908            ! 
     909            !                                ! East-West boundaries 
     910            IF( l_Iperio ) THEN                    !* cyclic 
     911               ARRAY_IN( 1 ,:,:,:,jf) = ARRAY_IN(jpim1,:,:,:,jf) 
     912               ARRAY_IN(jpi,:,:,:,jf) = ARRAY_IN(  2  ,:,:,:,jf) 
     913            ELSE                                   !* closed 
     914               IF( .NOT. NAT_IN(jf) == 'F' )   ARRAY_IN(     1       :nn_hls,:,:,:,jf) = zland    ! east except F-point 
     915                                               ARRAY_IN(nlci-nn_hls+1:jpi   ,:,:,:,jf) = zland    ! west 
     916            ENDIF 
     917            !                                ! North-South boundaries 
     918            IF( l_Jperio ) THEN                    !* cyclic (only with no mpp j-split) 
     919               ARRAY_IN(:, 1 ,:,:,jf) = ARRAY_IN(:, jpjm1,:,:,jf) 
     920               ARRAY_IN(:,jpj,:,:,jf) = ARRAY_IN(:,   2  ,:,:,jf) 
     921            ELSE                                   !* closed 
     922               IF( .NOT. NAT_IN(jf) == 'F' )   ARRAY_IN(:,     1       :nn_hls,:,:,jf) = zland    ! south except F-point 
     923                                               ARRAY_IN(:,nlcj-nn_hls+1:jpj   ,:,:,jf) = zland    ! north 
     924            ENDIF 
     925         END DO 
     926         ! 
     927      ENDIF 
     928 
     929! Wait for any reception (unpack?) 
     930#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     931      call MPI_Waitall(4, ml_reqs, MPI_STATUSES_IGNORE, iflag) 
     932#endif 
     933      !                           ! Write Dirichlet lateral conditions 
     934#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     935      SCOREP_USER_REGION_BEGIN( reg_unpack, "unpack", SCOREP_USER_REGION_TYPE_COMMON ) 
     936      iihom = nlci-nn_hls 
     937      ! 
     938      SELECT CASE ( nbondi ) 
     939      CASE ( -1 ) 
     940         DO jf = 1, ipf 
     941            DO jl = 1, ipl 
     942               DO jk = 1, ipk 
     943                  DO jh = 1, nn_hls 
     944                     ARRAY_IN(iihom+jh,:,jk,jl,jf) = zt3ew(:,jh,jk,jl,jf,2) 
     945                  END DO 
     946               END DO 
     947            END DO 
     948         END DO 
     949      CASE ( 0 ) 
     950         DO jf = 1, ipf 
     951            DO jl = 1, ipl 
     952               DO jk = 1, ipk 
     953                  DO jh = 1, nn_hls 
     954                     ARRAY_IN(jh      ,:,jk,jl,jf) = zt3we(:,jh,jk,jl,jf,2) 
     955                     ARRAY_IN(iihom+jh,:,jk,jl,jf) = zt3ew(:,jh,jk,jl,jf,2) 
     956                  END DO 
     957               END DO 
     958            END DO 
     959         END DO 
     960      CASE ( 1 ) 
     961         DO jf = 1, ipf 
     962            DO jl = 1, ipl 
     963               DO jk = 1, ipk 
     964                  DO jh = 1, nn_hls 
     965                     ARRAY_IN(jh      ,:,jk,jl,jf) = zt3we(:,jh,jk,jl,jf,2) 
     966                  END DO 
     967               END DO 
     968            END DO 
     969         END DO 
     970      END SELECT 
     971      !                           ! Write Dirichlet lateral conditions 
     972      ijhom = nlcj-nn_hls 
     973      ! 
     974      SELECT CASE ( nbondj ) 
     975      CASE ( -1 ) 
     976         DO jf = 1, ipf 
     977            DO jl = 1, ipl 
     978               DO jk = 1, ipk 
     979                  DO jh = 1, nn_hls 
     980                     ARRAY_IN(:,ijhom+jh,jk,jl,jf) = zt3ns(:,jh,jk,jl,jf,2) 
     981                  END DO 
     982               END DO 
     983            END DO 
     984         END DO 
     985      CASE ( 0 ) 
     986         DO jf = 1, ipf 
     987            DO jl = 1, ipl 
     988               DO jk = 1, ipk 
     989                  DO jh = 1, nn_hls 
     990                     ARRAY_IN(:,      jh,jk,jl,jf) = zt3sn(:,jh,jk,jl,jf,2) 
     991                     ARRAY_IN(:,ijhom+jh,jk,jl,jf) = zt3ns(:,jh,jk,jl,jf,2) 
     992                  END DO 
     993               END DO 
     994            END DO 
     995         END DO 
     996      CASE ( 1 ) 
     997         DO jf = 1, ipf 
     998            DO jl = 1, ipl 
     999               DO jk = 1, ipk 
     1000                  DO jh = 1, nn_hls 
     1001                     ARRAY_IN(:,jh,jk,jl,jf) = zt3sn(:,jh,jk,jl,jf,2) 
     1002                  END DO 
     1003               END DO 
     1004            END DO 
     1005         END DO 
     1006      END SELECT 
     1007      SCOREP_USER_REGION_END( reg_unpack ) 
     1008#endif 
     1009#endif 
    3331010      ! 4. north fold treatment 
    3341011      ! ----------------------- 
     
    3431020      ENDIF 
    3441021      ! 
     1022#ifdef ASYNC 
     1023! wait all sending messages 
     1024#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
     1025      call MPI_Waitall(4, ml_reqs(5), MPI_STATUSES_IGNORE, iflag) 
     1026#else 
     1027      call MPI_Waitall(8*ipf, ml_reqs, MPI_STATUSES_IGNORE, iflag) 
     1028#endif 
     1029      !                           ! Write Dirichlet lateral conditions 
     1030#ifdef MPI_DATATYPE_SUBARRAY 
     1031      SCOREP_USER_REGION_BEGIN( reg_datatype, "datatype", SCOREP_USER_REGION_TYPE_COMMON ) 
     1032      call MPI_Type_free(type_north_halo, iflag) 
     1033      call MPI_Type_free(type_south_halo, iflag) 
     1034      call MPI_Type_free(type_east_halo, iflag) 
     1035      call MPI_Type_free(type_west_halo, iflag) 
     1036      call MPI_Type_free(type_north_ghost, iflag) 
     1037      call MPI_Type_free(type_south_ghost, iflag) 
     1038      call MPI_Type_free(type_east_ghost, iflag) 
     1039      call MPI_Type_free(type_west_ghost, iflag) 
     1040      SCOREP_USER_REGION_END( reg_datatype ) 
     1041#endif 
     1042#ifdef MPI_DATATYPE_VECTOR 
     1043      SCOREP_USER_REGION_BEGIN( reg_datatype, "datatype vector", SCOREP_USER_REGION_TYPE_COMMON ) 
     1044      call MPI_Type_free(type_ew, iflag) 
     1045      call MPI_Type_free(type_ns, iflag) 
     1046      SCOREP_USER_REGION_END( reg_datatype ) 
     1047#endif 
     1048#endif 
     1049#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR) 
    3451050      DEALLOCATE( zt3ns, zt3sn, zt3ew, zt3we ) 
     1051#endif 
    3461052      ! 
    3471053   END SUBROUTINE ROUTINE_LNK 
     
    3551061#undef F_SIZE 
    3561062#undef OPT_K 
     1063#undef _INDEX 
Note: See TracChangeset for help on using the changeset viewer.