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.
mpp_lnk_generic.h90 in NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC – NEMO

source: NEMO/branches/2018/dev_r9759_HPC09_ESIWACE/src/OCE/LBC/mpp_lnk_generic.h90 @ 10136

Last change on this file since 10136 was 10136, checked in by dguibert, 6 years ago

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).

  • Property svn:mime-type set to text/x-fortran
File size: 43.9 KB
Line 
1#if defined MULTI
2#   define NAT_IN(k)                cd_nat(k)
3#   define SGN_IN(k)                psgn(k)
4#   define F_SIZE(ptab)             kfld
5#   define OPT_K(k)                 ,ipf
6#   if defined DIM_2d
7#      define ARRAY_TYPE(i,j,k,l,f)    TYPE(PTR_2D)                , INTENT(inout) ::   ptab(f)
8#      define ARRAY_IN(i,j,k,l,f)      ptab(f)%pt2d(i,j)
9#      define K_SIZE(ptab)             1
10#      define L_SIZE(ptab)             1
11#      define _INDEX(i,j,k,l)          (i+((j)+((0)+(0)*ipk)*jpj)*jpi)
12#   endif
13#   if defined DIM_3d
14#      define ARRAY_TYPE(i,j,k,l,f)    TYPE(PTR_3D)                , INTENT(inout) ::   ptab(f)
15#      define ARRAY_IN(i,j,k,l,f)      ptab(f)%pt3d(i,j,k)
16#      define K_SIZE(ptab)             SIZE(ptab(1)%pt3d,3)
17#      define L_SIZE(ptab)             1
18#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(0)*ipk)*jpj)*jpi)
19#   endif
20#   if defined DIM_4d
21#      define ARRAY_TYPE(i,j,k,l,f)    TYPE(PTR_4D)                , INTENT(inout) ::   ptab(f)
22#      define ARRAY_IN(i,j,k,l,f)      ptab(f)%pt4d(i,j,k,l)
23#      define K_SIZE(ptab)             SIZE(ptab(1)%pt4d,3)
24#      define L_SIZE(ptab)             SIZE(ptab(1)%pt4d,4)
25#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(l)*ipk)*jpj)*jpi)
26#   endif
27#else
28#   define ARRAY_TYPE(i,j,k,l,f)    REAL(wp)                    , INTENT(inout) ::   ARRAY_IN(i,j,k,l,f)
29#   define NAT_IN(k)                cd_nat
30#   define SGN_IN(k)                psgn
31#   define F_SIZE(ptab)             1
32#   define OPT_K(k)
33#   if defined DIM_2d
34#      define ARRAY_IN(i,j,k,l,f)   ptab(i,j)
35#      define K_SIZE(ptab)          1
36#      define L_SIZE(ptab)          1
37#      define _INDEX(i,j,k,l)          (i+((j)+((0)+(0)*ipk)*jpj)*jpi)
38#   endif
39#   if defined DIM_3d
40#      define ARRAY_IN(i,j,k,l,f)   ptab(i,j,k)
41#      define K_SIZE(ptab)          SIZE(ptab,3)
42#      define L_SIZE(ptab)          1
43#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(0)*ipk)*jpj)*jpi)
44#   endif
45#   if defined DIM_4d
46#      define ARRAY_IN(i,j,k,l,f)   ptab(i,j,k,l)
47#      define K_SIZE(ptab)          SIZE(ptab,3)
48#      define L_SIZE(ptab)          SIZE(ptab,4)
49#      define _INDEX(i,j,k,l)          (i+((j)+((k)+(l)*ipk)*jpj)*jpi)
50#   endif
51#endif
52
53#if defined MULTI
54#if defined ASYNC
55   SUBROUTINE ROUTINE_LNK( rname, ptab, cd_nat, psgn, loop_fct, kfld, cd_mpp, pval )
56      INTEGER                     , INTENT(in   ) ::   kfld        ! number of pt3d arrays
57#else
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   !
73#endif
74      ARRAY_TYPE(:,:,:,:,:)                                        ! array or pointer of arrays on which the boundary condition is applied
75      CHARACTER(len=1)            , INTENT(in   ) ::   NAT_IN(:)   ! nature of array grid-points
76      REAL(wp)                    , INTENT(in   ) ::   SGN_IN(:)   ! sign used across the north fold boundary
77      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp      ! fill the overlap area only
78      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
90      CHARACTER(len=*),             INTENT(in   ) ::   rname       ! name of the calling subroutine
91      !
92      INTEGER  ::    ji,  jj,  jk,  jl, jh, jf   ! dummy loop indices
93      INTEGER  ::   ipi, ipj, ipk, ipl, ipf      ! dimension of the input array
94      INTEGER  ::   imigr, iihom, ijhom          ! local integers
95      INTEGER  ::   ml_req1, ml_req2, ml_err     ! for key_mpi_isend
96      REAL(wp) ::   zland
97      INTEGER , DIMENSION(MPI_STATUS_SIZE)      ::   ml_stat        ! for key_mpi_isend
98      REAL(wp), DIMENSION(:,:,:,:,:,:), ALLOCATABLE ::   zt3ns, zt3sn   ! north-south & south-north  halos
99      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
141      !!----------------------------------------------------------------------
142      !
143      ipk = K_SIZE(ptab)   ! 3rd dimension
144      ipl = L_SIZE(ptab)   ! 4th    -
145      ipf = F_SIZE(ptab)   ! 5th    -      use in "multi" case (array of pointers)
146      !
147#if !(defined MPI_DATATYPE_SUBARRAY || defined MPI_DATATYPE_VECTOR)
148      ALLOCATE( zt3ns(jpi,nn_hls,ipk,ipl,ipf,2), zt3sn(jpi,nn_hls,ipk,ipl,ipf,2),   &
149         &      zt3ew(jpj,nn_hls,ipk,ipl,ipf,2), zt3we(jpj,nn_hls,ipk,ipl,ipf,2)  )
150#endif
151      !
152      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value
153      ELSE                         ;   zland = 0._wp     ! zero by default
154      ENDIF
155
156#ifndef ASYNC
157      ! ------------------------------- !
158      !   standard boundary treatment   !    ! CAUTION: semi-column notation is often impossible
159      ! ------------------------------- !
160      !
161      IF( PRESENT( cd_mpp ) ) THEN     !==  halos filled with inner values  ==!
162         !
163         DO jf = 1, ipf                      ! number of arrays to be treated
164            !
165            DO jl = 1, ipl                   ! CAUTION: ptab is defined only between nld and nle
166               DO jk = 1, ipk
167                  DO jj = nlcj+1, jpj                 ! added line(s)   (inner only)
168                     ARRAY_IN(nldi  :nlei  ,jj,jk,jl,jf) = ARRAY_IN(nldi:nlei,nlej,jk,jl,jf)
169                     ARRAY_IN(1     :nldi-1,jj,jk,jl,jf) = ARRAY_IN(nldi     ,nlej,jk,jl,jf)
170                     ARRAY_IN(nlei+1:nlci  ,jj,jk,jl,jf) = ARRAY_IN(     nlei,nlej,jk,jl,jf)
171                  END DO
172                  DO ji = nlci+1, jpi                 ! added column(s) (full)
173                     ARRAY_IN(ji,nldj  :nlej  ,jk,jl,jf) = ARRAY_IN(nlei,nldj:nlej,jk,jl,jf)
174                     ARRAY_IN(ji,1     :nldj-1,jk,jl,jf) = ARRAY_IN(nlei,nldj     ,jk,jl,jf)
175                     ARRAY_IN(ji,nlej+1:jpj   ,jk,jl,jf) = ARRAY_IN(nlei,     nlej,jk,jl,jf)
176                  END DO
177               END DO
178            END DO
179            !
180         END DO
181         !
182      ELSE                              !==  standard close or cyclic treatment  ==!
183         !
184         DO jf = 1, ipf                      ! number of arrays to be treated
185            !
186            !                                ! East-West boundaries
187            IF( l_Iperio ) THEN                    !* cyclic
188               ARRAY_IN( 1 ,:,:,:,jf) = ARRAY_IN(jpim1,:,:,:,jf)
189               ARRAY_IN(jpi,:,:,:,jf) = ARRAY_IN(  2  ,:,:,:,jf)
190            ELSE                                   !* closed
191               IF( .NOT. NAT_IN(jf) == 'F' )   ARRAY_IN(     1       :nn_hls,:,:,:,jf) = zland    ! east except F-point
192                                               ARRAY_IN(nlci-nn_hls+1:jpi   ,:,:,:,jf) = zland    ! west
193            ENDIF
194            !                                ! North-South boundaries
195            IF( l_Jperio ) THEN                    !* cyclic (only with no mpp j-split)
196               ARRAY_IN(:, 1 ,:,:,jf) = ARRAY_IN(:, jpjm1,:,:,jf)
197               ARRAY_IN(:,jpj,:,:,jf) = ARRAY_IN(:,   2  ,:,:,jf)
198            ELSE                                   !* closed
199               IF( .NOT. NAT_IN(jf) == 'F' )   ARRAY_IN(:,     1       :nn_hls,:,:,jf) = zland    ! south except F-point
200                                               ARRAY_IN(:,nlcj-nn_hls+1:jpj   ,:,:,jf) = zland    ! north
201            ENDIF
202         END DO
203         !
204      ENDIF
205
206      ! ------------------------------- !
207      !      East and west exchange     !
208      ! ------------------------------- !
209      ! we play with the neigbours AND the row number because of the periodicity
210      !
211      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON )
212      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions
213      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case)
214         iihom = nlci-nreci
215         DO jf = 1, ipf
216            DO jl = 1, ipl
217               DO jk = 1, ipk
218                  DO jh = 1, nn_hls
219                     zt3ew(:,jh,jk,jl,jf,1) = ARRAY_IN(nn_hls+jh,:,jk,jl,jf)
220                     zt3we(:,jh,jk,jl,jf,1) = ARRAY_IN(iihom +jh,:,jk,jl,jf)
221                  END DO
222               END DO
223            END DO
224         END DO
225      END SELECT
226      SCOREP_USER_REGION_END( reg_pack )
227      !
228      !                           ! Migrations
229      imigr = nn_hls * jpj * ipk * ipl * ipf
230      !
231      IF ( ncom_stp == nit000 ) then
232         n_sequence = n_sequence + 1
233         icomm_sequence(n_sequence,1) = ipk
234         icomm_sequence(n_sequence,2) = ipf
235         ! write(6,'(A,6I4)') 'size comm ', nn_hls, jpi, jpj, ipk, ipl, ipf
236      ELSE IF ( mpprank == 0 .AND. ncom_stp == (nit000+1) ) THEN
237         IF ( l_print_comm_report ) THEN
238            write(6,*) 'Communication pattern report : '
239            write(6,*) ' '
240            write(6,'(A,I3)') ' Exchanged halos : ', n_sequence
241            jj = 0; jk = 0; jf = 0; jh = 0
242            DO ji = 1, n_sequence
243              IF ( icomm_sequence(ji,1) .gt. 1 ) jk = jk + 1
244              IF ( icomm_sequence(ji,2) .gt. 1 ) jf = jf + 1
245              IF ( icomm_sequence(ji,1) .gt. 1 .AND. icomm_sequence(ji,2) .gt. 1 ) jj = jj + 1
246              jh = MAX (jh, icomm_sequence(ji,1)*icomm_sequence(ji,2))
247            END DO
248            write(6,'(A,I3)') ' 3D Exchanged halos : ', jk
249            write(6,'(A,I3)') ' Multi arrays exchanged halos : ', jf
250            write(6,'(A,I3)') '   from which 3D : ', jj
251            write(6,'(A,I10)') ' array max size : ', jh*jpi*jpj
252            write(6,*) ' '
253            l_print_comm_report = .FALSE.
254         END IF
255         write(6,'(A19,A)') 'calling subroutine ', TRIM(rname)
256      END IF
257      !
258      IF (ncom_stp <= ( nit000 + 1 ) .or. mod(ncom_stp,nn_comm_mod) == 0 ) THEN
259         IF ( TRIM(rname) == "simulated_lbc_lnk" ) THEN
260            zt3we = zt3we + 1. ; zt3ew = zt3ew + 1.
261         ENDIF
262         !
263         CALL tic_tac(.TRUE.)
264         !
265         SELECT CASE ( nbondi )
266         CASE ( -1 )
267            CALL mppsend( 2, zt3we(1,1,1,1,1,1), imigr, noea, ml_req1 )
268            CALL mpprecv( 1, zt3ew(1,1,1,1,1,2), imigr, noea )
269            IF(l_isend)   CALL mpi_wait(ml_req1, ml_stat, ml_err)
270         CASE ( 0 )
271            CALL mppsend( 1, zt3ew(1,1,1,1,1,1), imigr, nowe, ml_req1 )
272            CALL mppsend( 2, zt3we(1,1,1,1,1,1), imigr, noea, ml_req2 )
273            CALL mpprecv( 1, zt3ew(1,1,1,1,1,2), imigr, noea )
274            CALL mpprecv( 2, zt3we(1,1,1,1,1,2), imigr, nowe )
275            IF(l_isend)   CALL mpi_wait(ml_req1, ml_stat, ml_err)
276            IF(l_isend)   CALL mpi_wait(ml_req2, ml_stat, ml_err)
277         CASE ( 1 )
278            CALL mppsend( 1, zt3ew(1,1,1,1,1,1), imigr, nowe, ml_req1 )
279            CALL mpprecv( 2, zt3we(1,1,1,1,1,2), imigr, nowe )
280            IF(l_isend)   CALL mpi_wait(ml_req1, ml_stat, ml_err )
281         END SELECT
282         !
283         CALL tic_tac(.FALSE.)
284         !
285      END IF
286      !
287      !                           ! Write Dirichlet lateral conditions
288      iihom = nlci-nn_hls
289      !
290      SCOREP_USER_REGION_BEGIN( reg_unpack, "unpack", SCOREP_USER_REGION_TYPE_COMMON )
291      SELECT CASE ( nbondi )
292      CASE ( -1 )
293         DO jf = 1, ipf
294            DO jl = 1, ipl
295               DO jk = 1, ipk
296                  DO jh = 1, nn_hls
297                     ARRAY_IN(iihom+jh,:,jk,jl,jf) = zt3ew(:,jh,jk,jl,jf,2)
298                  END DO
299               END DO
300            END DO
301         END DO
302      CASE ( 0 )
303         DO jf = 1, ipf
304            DO jl = 1, ipl
305               DO jk = 1, ipk
306                  DO jh = 1, nn_hls
307                     ARRAY_IN(jh      ,:,jk,jl,jf) = zt3we(:,jh,jk,jl,jf,2)
308                     ARRAY_IN(iihom+jh,:,jk,jl,jf) = zt3ew(:,jh,jk,jl,jf,2)
309                  END DO
310               END DO
311            END DO
312         END DO
313      CASE ( 1 )
314         DO jf = 1, ipf
315            DO jl = 1, ipl
316               DO jk = 1, ipk
317                  DO jh = 1, nn_hls
318                     ARRAY_IN(jh      ,:,jk,jl,jf) = zt3we(:,jh,jk,jl,jf,2)
319                  END DO
320               END DO
321            END DO
322         END DO
323      END SELECT
324      SCOREP_USER_REGION_END( reg_unpack )
325
326      ! 3. North and south directions
327      ! -----------------------------
328      ! always closed : we play only with the neigbours
329      !
330      SCOREP_USER_REGION_BEGIN( reg_pack, "pack", SCOREP_USER_REGION_TYPE_COMMON )
331      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions
332         ijhom = nlcj-nrecj
333         DO jf = 1, ipf
334            DO jl = 1, ipl
335               DO jk = 1, ipk
336                  DO jh = 1, nn_hls
337                     zt3sn(:,jh,jk,jl,jf,1) = ARRAY_IN(:,ijhom +jh,jk,jl,jf)
338                     zt3ns(:,jh,jk,jl,jf,1) = ARRAY_IN(:,nn_hls+jh,jk,jl,jf)
339                  END DO
340               END DO
341            END DO
342         END DO
343      ENDIF
344#ifdef SCOREP_USER_ENABLE
345      SCOREP_USER_REGION_END( reg_pack )
346#endif
347      !
348      !                           ! Migrations
349      imigr = nn_hls * jpi * ipk * ipl * ipf
350      !
351      IF (ncom_stp <= ( nit000 + 1 ) .or. mod(ncom_stp,nn_comm_mod) == 0 ) THEN
352         IF ( TRIM(rname) == "simulated_lbc_lnk" ) THEN
353            zt3sn = zt3sn + 1. ; zt3ns = zt3ns + 1.
354         ENDIF
355         !
356         CALL tic_tac(.TRUE.)
357         !
358         SELECT CASE ( nbondj )
359         CASE ( -1 )
360            CALL mppsend( 4, zt3sn(1,1,1,1,1,1), imigr, nono, ml_req1 )
361            CALL mpprecv( 3, zt3ns(1,1,1,1,1,2), imigr, nono )
362            IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err )
363         CASE ( 0 )
364            CALL mppsend( 3, zt3ns(1,1,1,1,1,1), imigr, noso, ml_req1 )
365            CALL mppsend( 4, zt3sn(1,1,1,1,1,1), imigr, nono, ml_req2 )
366            CALL mpprecv( 3, zt3ns(1,1,1,1,1,2), imigr, nono )
367            CALL mpprecv( 4, zt3sn(1,1,1,1,1,2), imigr, noso )
368            IF(l_isend)   CALL mpi_wait(ml_req1, ml_stat, ml_err )
369            IF(l_isend)   CALL mpi_wait(ml_req2, ml_stat, ml_err )
370         CASE ( 1 )
371            CALL mppsend( 3, zt3ns(1,1,1,1,1,1), imigr, noso, ml_req1 )
372            CALL mpprecv( 4, zt3sn(1,1,1,1,1,2), imigr, noso )
373            IF(l_isend)   CALL mpi_wait(ml_req1, ml_stat, ml_err )
374         END SELECT
375         ! imbalance measurement
376         CALL tic_tac(.FALSE.)
377         !
378      END IF
379      !
380      !                           ! Write Dirichlet lateral conditions
381      SCOREP_USER_REGION_BEGIN( reg_unpack, "unpack", SCOREP_USER_REGION_TYPE_COMMON )
382      ijhom = nlcj-nn_hls
383      !
384      SELECT CASE ( nbondj )
385      CASE ( -1 )
386         DO jf = 1, ipf
387            DO jl = 1, ipl
388               DO jk = 1, ipk
389                  DO jh = 1, nn_hls
390                     ARRAY_IN(:,ijhom+jh,jk,jl,jf) = zt3ns(:,jh,jk,jl,jf,2)
391                  END DO
392               END DO
393            END DO
394         END DO
395      CASE ( 0 )
396         DO jf = 1, ipf
397            DO jl = 1, ipl
398               DO jk = 1, ipk
399                  DO jh = 1, nn_hls
400                     ARRAY_IN(:,      jh,jk,jl,jf) = zt3sn(:,jh,jk,jl,jf,2)
401                     ARRAY_IN(:,ijhom+jh,jk,jl,jf) = zt3ns(:,jh,jk,jl,jf,2)
402                  END DO
403               END DO
404            END DO
405         END DO
406      CASE ( 1 )
407         DO jf = 1, ipf
408            DO jl = 1, ipl
409               DO jk = 1, ipk
410                  DO jh = 1, nn_hls
411                     ARRAY_IN(:,jh,jk,jl,jf) = zt3sn(:,jh,jk,jl,jf,2)
412                  END DO
413               END DO
414            END DO
415         END DO
416      END SELECT
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
1010      ! 4. north fold treatment
1011      ! -----------------------
1012      !
1013      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN
1014         !
1015         SELECT CASE ( jpni )
1016         CASE ( 1 )     ;   CALL lbc_nfd( ptab, NAT_IN(:), SGN_IN(:) OPT_K(:) )   ! only 1 northern proc, no mpp
1017         CASE DEFAULT   ;   CALL mpp_nfd( ptab, NAT_IN(:), SGN_IN(:) OPT_K(:) )   ! for all northern procs.
1018         END SELECT
1019         !
1020      ENDIF
1021      !
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)
1050      DEALLOCATE( zt3ns, zt3sn, zt3ew, zt3we )
1051#endif
1052      !
1053   END SUBROUTINE ROUTINE_LNK
1054
1055#undef ARRAY_TYPE
1056#undef NAT_IN
1057#undef SGN_IN
1058#undef ARRAY_IN
1059#undef K_SIZE
1060#undef L_SIZE
1061#undef F_SIZE
1062#undef OPT_K
1063#undef _INDEX
Note: See TracBrowser for help on using the repository browser.