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.
lbc_lnk_pt2pt_generic.h90 in NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC – NEMO

source: NEMO/branches/2021/dev_r14312_MPI_Interface/src/OCE/LBC/lbc_lnk_pt2pt_generic.h90 @ 14367

Last change on this file since 14367 was 14367, checked in by smasson, 3 years ago

dev_r14312_MPI_Interface: keep send/recv buffers in memory if smaller than jpi*jpj, #2598

  • Property svn:keywords set to Id
  • Property svn:mime-type set to text/x-fortran
File size: 19.0 KB
Line 
1 
2   SUBROUTINE lbc_lnk_pt2pt_/**/PRECISION( cdname, ptab, cd_nat, psgn, kfld, kfillmode, pfillval, khls, lsend, lrecv, ld4only )
3      CHARACTER(len=*)              , INTENT(in   ) ::   cdname      ! name of the calling subroutine
4      TYPE(PTR_4d_/**/PRECISION),  DIMENSION(:), INTENT(inout) ::   ptab        ! pointer of arrays on which apply the b.c.
5      CHARACTER(len=1), DIMENSION(:), INTENT(in   ) ::   cd_nat      ! nature of array grid-points
6      REAL(PRECISION),  DIMENSION(:), INTENT(in   ) ::   psgn        ! sign used across the north fold boundary
7      INTEGER                       , INTENT(in   ) ::   kfld        ! number of pt3d arrays
8      INTEGER ,             OPTIONAL, INTENT(in   ) ::   kfillmode   ! filling method for halo over land (default = constant)
9      REAL(PRECISION),      OPTIONAL, INTENT(in   ) ::   pfillval    ! background value (used at closed boundaries)
10      INTEGER ,             OPTIONAL, INTENT(in   ) ::   khls        ! halo size, default = nn_hls
11      LOGICAL, DIMENSION(4),OPTIONAL, INTENT(in   ) ::   lsend, lrecv  ! communication with other 4 proc
12      LOGICAL,              OPTIONAL, INTENT(in   ) ::   ld4only     ! if .T., do only 4-neighbour comm (ignore corners)
13      !
14      INTEGER  ::     ji,   jj,  jk,  jl,  jf, jn     ! dummy loop indices
15      INTEGER  ::    ipi,  ipj, ipk, ipl, ipf         ! dimension of the input array
16      INTEGER  ::   ip0i, ip1i, im0i, im1i
17      INTEGER  ::   ip0j, ip1j, im0j, im1j
18      INTEGER  ::   ishti, ishtj, ishti2, ishtj2
19      INTEGER  ::   ifill_nfd, icomm, ierr
20      INTEGER  ::   ihls, idxs, idxr, iszS, iszR
21      INTEGER, DIMENSION(4)  ::   iwewe, issnn
22      INTEGER, DIMENSION(8)  ::   isizei, ishtSi, ishtRi, ishtPi
23      INTEGER, DIMENSION(8)  ::   isizej, ishtSj, ishtRj, ishtPj
24      INTEGER, DIMENSION(8)  ::   ifill, iszall, ishtS, ishtR
25      INTEGER, DIMENSION(8)  ::   ireq             ! mpi_request id
26      INTEGER, DIMENSION(8)  ::   iStag, iRtag     ! Send and Recv mpi_tag id
27      REAL(PRECISION) ::   zland
28      LOGICAL, DIMENSION(8)  ::   llsend, llrecv
29      LOGICAL  ::   ll4only                                        ! default: 8 neighbourgs
30      LOGICAL  ::   ll_IdoNFold
31      !!----------------------------------------------------------------------
32      !
33      ! ----------------------------------------- !
34      !     1. local variables initialization     !
35      ! ----------------------------------------- !
36      !
37      ipi = SIZE(ptab(1)%pt4d,1)
38      ipj = SIZE(ptab(1)%pt4d,2)
39      ipk = SIZE(ptab(1)%pt4d,3)
40      ipl = SIZE(ptab(1)%pt4d,4)
41      ipf = kfld
42      !
43      IF( narea == 1 .AND. numcom == -1 ) CALL mpp_report( cdname, ipk, ipl, ipf, ld_lbc = .TRUE. )
44      !
45      idxs = 1   ! initalize index for send buffer
46      idxr = 1   ! initalize index for recv buffer
47      icomm = mpi_comm_oce        ! shorter name
48      !
49      ! take care of optional parameters
50      !
51      ihls = nn_hls   ! default definition
52      IF( PRESENT( khls ) )   ihls = khls
53      IF( ihls > n_hlsmax ) THEN
54         WRITE(ctmp1,*) TRIM(cdname), ' is calling lbc_lnk with khls > n_hlsmax : ', khls, '>', n_hlsmax
55         CALL ctl_stop( 'STOP', ctmp1 )
56      ENDIF
57      IF( ipi /= Ni_0+2*ihls ) THEN
58         WRITE(ctmp1,*) TRIM(cdname), ' is calling lbc_lnk with an input array which does not match ihls along i: ', ipi, ihls, Ni_0
59         CALL ctl_stop( 'STOP', ctmp1 )
60      ENDIF
61      IF( ipj /= Nj_0+2*ihls ) THEN
62         WRITE(ctmp1,*) TRIM(cdname), ' is calling lbc_lnk with an input array which does not match ihls along j:', ipj, ihls , Nj_0
63         CALL ctl_stop( 'STOP', ctmp1 )
64      ENDIF
65      !
66      ll4only = .FALSE.    ! default definition
67      IF( PRESENT(ld4only) )   ll4only = ld4only
68      !
69      zland = 0._wp                                     ! land filling value: zero by default
70      IF( PRESENT( pfillval ) )   zland = pfillval      ! set land value
71      !
72      ! define llsend and llrecv: logicals which say if mpi-neibourgs for send or receive exist or not.
73      IF     ( PRESENT(lsend) .AND. PRESENT(lrecv) ) THEN   ! localy defined neighbourgs
74         llsend(1:4) = lsend(1:4)   ;   llrecv(1:4) = lrecv(1:4)
75      ELSE IF( PRESENT(lsend) .OR.  PRESENT(lrecv) ) THEN
76         WRITE(ctmp1,*) TRIM(cdname), ' is calling lbc_lnk with only one of the two arguments lsend or lrecv'
77         CALL ctl_stop( 'STOP', ctmp1 )
78      ELSE                                              ! default neighbours
79         llsend(:) = mpiSnei(ihls,:) >= 0
80         IF( ll4only )   llsend(5:8) = .FALSE.          ! exclude corners
81         llrecv(:) = mpiRnei(ihls,:) >= 0
82         IF( ll4only )   llrecv(5:8) = .FALSE.          ! exclude corners
83      ENDIF
84      !
85      ! define ifill: which method should be used to fill each parts (sides+corners) of the halos
86      ! default definition
87      DO jn = 1, 4
88         IF(             llrecv(jn) ) THEN   ;   ifill(jn) = jpfillmpi    ! with an mpi communication
89         ELSEIF(    l_SelfPerio(jn) ) THEN   ;   ifill(jn) = jpfillperio  ! with self-periodicity
90         ELSEIF( PRESENT(kfillmode) ) THEN   ;   ifill(jn) = kfillmode    ! localy defined
91         ELSE                                ;   ifill(jn) = jpfillcst    ! constant value (zland)
92         ENDIF
93      END DO
94      DO jn = 5, 8
95         IF(             llrecv(jn) ) THEN   ;   ifill(jn) = jpfillmpi    ! with an mpi communication
96         ELSE                                ;   ifill(jn) = jpfillnothing! do nothing
97         ENDIF
98      END DO
99         !
100      ! north fold treatment
101      ll_IdoNFold = l_IdoNFold .AND. ifill(jpno) /= jpfillnothing
102      IF( ll_IdoNFold ) THEN
103         ifill_nfd = ifill(jpno)             ! if we are here, this means llrecv(jpno) = .false. and l_SelfPerio(jpno) = .false.
104         ifill( (/jpno/) ) = jpfillnothing   ! we do north fold -> do nothing for northern halo
105      ENDIF
106     
107      ! We first define the localization and size of the parts of the array that will be sent (s), received (r)
108      ! or used for periodocity (p). The localization is defined as "the bottom left corner - 1" in i and j directions.
109      ! This is a shift that will be applied later in the do loops to pick-up the appropriate part of the array
110      !
111      ! all definitions bellow do not refer to N[ij][se]0 so we can use it with any local value of ihls
112      !                   !                       ________________________
113      ip0i =          0   !          im0j = inner |__|__|__________|__|__|
114      ip1i =       ihls   !   im1j = inner - halo |__|__|__________|__|__|
115      im1i = ipi-2*ihls   !                       |  |  |          |  |  |
116      im0i = ipi - ihls   !                       |  |  |          |  |  |
117      ip0j =          0   !                       |  |  |          |  |  |
118      ip1j =       ihls   !                       |__|__|__________|__|__|
119      im1j = ipj-2*ihls   !           ip1j = halo |__|__|__________|__|__|
120      im0j = ipj - ihls   !              ip0j = 0 |__|__|__________|__|__|
121      !                   !                    ip0i ip1i        im1i im0i
122      !
123      iwewe(:) = (/ jpwe,jpea,jpwe,jpea /)   ;   issnn(:) = (/ jpso,jpso,jpno,jpno /)
124      !cd     sides:     west  east south north      ;   corners: so-we, so-ea, no-we, no-ea
125      isizei(1:4) = (/ ihls, ihls,  ipi,  ipi /)   ;   isizei(5:8) = ihls              ! i- count
126      isizej(1:4) = (/  ipj,  ipj, ihls, ihls /)   ;   isizej(5:8) = ihls              ! j- count
127      ishtSi(1:4) = (/ ip1i, im1i, ip0i, ip0i /)   ;   ishtSi(5:8) = ishtSi( iwewe )   ! i- shift send data
128      ishtSj(1:4) = (/ ip0j, ip0j, ip1j, im1j /)   ;   ishtSj(5:8) = ishtSj( issnn )   ! j- shift send data
129      ishtRi(1:4) = (/ ip0i, im0i, ip0i, ip0i /)   ;   ishtRi(5:8) = ishtRi( iwewe )   ! i- shift received data location
130      ishtRj(1:4) = (/ ip0j, ip0j, ip0j, im0j /)   ;   ishtRj(5:8) = ishtRj( issnn )   ! j- shift received data location
131      ishtPi(1:4) = (/ im1i, ip1i, ip0i, ip0i /)   ;   ishtPi(5:8) = ishtPi( iwewe )   ! i- shift data used for periodicity
132      ishtPj(1:4) = (/ ip0j, ip0j, im1j, ip1j /)   ;   ishtPj(5:8) = ishtPj( issnn )   ! j- shift data used for periodicity
133      !
134      ! -------------------------------- !
135      !     2. Prepare MPI exchanges     !
136      ! -------------------------------- !
137      !
138      iStag = (/ 1, 2, 3, 4, 5, 6, 7, 8 /)   ! any value but each one must be different
139      ! define iRtag with the corresponding iStag, e.g. data received at west where sent at east.
140      iRtag(jpwe) = iStag(jpea)   ;   iRtag(jpea) = iStag(jpwe)   ;   iRtag(jpso) = iStag(jpno)   ;   iRtag(jpno) = iStag(jpso)
141      iRtag(jpsw) = iStag(jpne)   ;   iRtag(jpse) = iStag(jpnw)   ;   iRtag(jpnw) = iStag(jpse)   ;   iRtag(jpne) = iStag(jpsw)
142      !
143      iszall(:) = isizei(:) * isizej(:) * ipk * ipl * ipf
144      ishtS(1) = 0
145      DO jn = 2, 8
146         ishtS(jn) = ishtS(jn-1) + iszall(jn-1) * COUNT( (/llsend(jn-1)/) )
147      END DO
148      ishtR(1) = 0
149      DO jn = 2, 8
150         ishtR(jn) = ishtR(jn-1) + iszall(jn-1) * COUNT( (/llrecv(jn-1)/) )
151      END DO
152
153      ! Allocate buffer arrays to be sent/received if needed
154      iszS = SUM(iszall, mask = llsend)                             ! send buffer size
155      IF( ALLOCATED(BUFFSND) ) THEN
156         CALL mpi_waitall(8, nreq_p2p, MPI_STATUSES_IGNORE, ierr)   ! wait for Isend from the PREVIOUS call
157         IF( SIZE(BUFFSND) < iszS )    DEALLOCATE(BUFFSND)          ! send buffer is too small
158      ENDIF
159      IF( .NOT. ALLOCATED(BUFFSND) )   ALLOCATE( BUFFSND(iszS) )
160      iszR = SUM(iszall, mask = llrecv)                             ! recv buffer size
161      IF( ALLOCATED(BUFFRCV) ) THEN
162         IF( SIZE(BUFFRCV) < iszR )    DEALLOCATE(BUFFRCV)          ! recv buffer is too small
163      ENDIF
164      IF( .NOT. ALLOCATED(BUFFRCV) )   ALLOCATE( BUFFRCV(iszR) )
165      !
166      ! default definition when no communication is done. understood by mpi_waitall
167      nreq_p2p(:) = MPI_REQUEST_NULL   ! WARNING: Must be done after the call to mpi_waitall just above
168      !
169      ! ----------------------------------------------- !
170      !     3. Do east and west MPI_Isend if needed     !
171      ! ----------------------------------------------- !
172      !
173      DO jn = 1, 2
174         IF( llsend(jn) ) THEN
175            ishti = ishtSi(jn)
176            ishtj = ishtSj(jn)
177            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
178               BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)
179               idxs = idxs + 1
180            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
181#if ! defined key_mpi_off
182            IF( ln_timing ) CALL tic_tac(.TRUE.)
183            ! non-blocking send of the west/east side using local buffer
184            CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr )
185            IF( ln_timing ) CALL tic_tac(.FALSE.)
186#endif
187         ENDIF
188      END DO
189      !
190      ! ----------------------------------- !
191      !     4. Fill east and west halos     !
192      ! ----------------------------------- !
193      !
194      DO jn = 1, 2
195         ishti = ishtRi(jn)
196         ishtj = ishtRj(jn)
197         SELECT CASE ( ifill(jn) )
198         CASE ( jpfillnothing )               ! no filling
199         CASE ( jpfillmpi   )                 ! fill with data received by MPI
200#if ! defined key_mpi_off
201            IF( ln_timing ) CALL tic_tac(.TRUE.)
202            !                                 ! blocking receive of the west/east halo in local temporary arrays
203            CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )
204            IF( ln_timing ) CALL tic_tac(.FALSE.)
205#endif
206            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
207               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr)
208               idxr = idxr + 1
209            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
210         CASE ( jpfillperio )                 ! use periodicity
211            ishti2 = ishtPi(jn)
212            ishtj2 = ishtPj(jn)
213            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
214               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = ptab(jf)%pt4d(ishti2+ji,ishtj2+jj,jk,jl)
215            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
216         CASE ( jpfillcopy  )                 ! filling with inner domain values
217            ishti2 = ishtSi(jn)
218            ishtj2 = ishtSj(jn)
219            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
220               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = ptab(jf)%pt4d(ishti2+ji,ishtj2+jj,jk,jl)
221            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
222         CASE ( jpfillcst   )                 ! filling with constant value
223            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
224               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zland
225            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
226         END SELECT
227      END DO
228      !
229      ! ------------------------------- !
230      !     5. north fold treatment     !
231      ! ------------------------------- !
232      !
233      ! do it before south directions so concerned processes can do it without waiting for the comm with the sourthern neighbor
234      !
235      IF( ll_IdoNFold ) THEN
236         IF( jpni == 1 )  THEN   ;   CALL lbc_nfd( ptab, cd_nat, psgn                  , ihls, ipf )   ! self NFold
237         ELSE                    ;   CALL mpp_nfd( ptab, cd_nat, psgn, ifill_nfd, zland, ihls, ipf )   ! mpi  NFold
238         ENDIF
239      ENDIF
240      !
241      ! ------------------------------------------------- !
242      !     6. Do north and south MPI_Isend if needed     !
243      ! ------------------------------------------------- !
244      !
245      DO jn = 3, 4
246         IF( llsend(jn) ) THEN
247            ishti = ishtSi(jn)
248            ishtj = ishtSj(jn)
249            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
250               BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)
251               idxs = idxs + 1
252            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
253#if ! defined key_mpi_off
254            IF( ln_timing ) CALL tic_tac(.TRUE.)
255            ! non-blocking send of the south/north side using local buffer
256            CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr )
257            IF( ln_timing ) CALL tic_tac(.FALSE.)
258#endif
259         ENDIF
260      END DO
261      !
262      ! ------------------------------------- !
263      !     7. Fill south and north halos     !
264      ! ------------------------------------- !
265      !
266      DO jn = 3, 4
267         ishti = ishtRi(jn)
268         ishtj = ishtRj(jn)
269         SELECT CASE ( ifill(jn) )
270         CASE ( jpfillnothing )               ! no filling
271         CASE ( jpfillmpi   )                 ! fill with data received by MPI
272#if ! defined key_mpi_off
273            IF( ln_timing ) CALL tic_tac(.TRUE.)
274            !                                 ! blocking receive of the south/north halo in local temporary arrays
275            CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )
276            IF( ln_timing ) CALL tic_tac(.FALSE.)
277#endif
278            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
279               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr)
280               idxr = idxr + 1
281            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
282         CASE ( jpfillperio )                 ! use periodicity
283            ishti2 = ishtPi(jn)
284            ishtj2 = ishtPj(jn)
285            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
286               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = ptab(jf)%pt4d(ishti2+ji,ishtj2+jj,jk,jl)
287            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
288         CASE ( jpfillcopy  )                 ! filling with inner domain values
289            ishti2 = ishtSi(jn)
290            ishtj2 = ishtSj(jn)
291            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
292               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = ptab(jf)%pt4d(ishti2+ji,ishtj2+jj,jk,jl)
293            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
294         CASE ( jpfillcst   )                 ! filling with constant value
295            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
296               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = zland
297            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
298         END SELECT
299      END DO
300      !
301      ! ----------------------------------------------- !
302      !     8. Specific problem in corner treatment     !
303      !              ( very rate case... )              !
304      ! ----------------------------------------------- !
305      !
306      DO jn = 5, 8
307         IF( llsend(jn) ) THEN
308            ishti = ishtSi(jn)
309            ishtj = ishtSj(jn)
310            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
311               BUFFSND(idxs) = ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl)
312               idxs = idxs + 1
313            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
314#if ! defined key_mpi_off
315            IF( ln_timing ) CALL tic_tac(.TRUE.)
316            ! non-blocking send of the corners using local buffer
317            CALL MPI_ISEND( BUFFSND(ishtS(jn)+1), iszall(jn), MPI_TYPE, mpiSnei(ihls,jn), iStag(jn), icomm, nreq_p2p(jn), ierr )
318            IF( ln_timing ) CALL tic_tac(.FALSE.)
319#endif
320         ENDIF
321      END DO
322      DO jn = 5, 8
323         IF( llrecv(jn) ) THEN
324#if ! defined key_mpi_off
325            IF( ln_timing ) CALL tic_tac(.TRUE.)
326            ! blocking receive of the corner halo in local temporary arrays
327            CALL MPI_RECV( BUFFRCV(ishtR(jn)+1), iszall(jn), MPI_TYPE, mpiRnei(ihls,jn), iRtag(jn), icomm, MPI_STATUS_IGNORE, ierr )
328            IF( ln_timing ) CALL tic_tac(.FALSE.)
329#endif
330            ishti = ishtRi(jn)
331            ishtj = ishtRj(jn)
332            DO jf = 1, ipf  ;  DO jl = 1, ipl  ;  DO jk = 1, ipk  ;  DO jj = 1,isizej(jn)  ;  DO ji = 1,isizei(jn)
333               ptab(jf)%pt4d(ishti+ji,ishtj+jj,jk,jl) = BUFFRCV(idxr)
334               idxr = idxr + 1
335            END DO   ;   END DO   ;   END DO   ;   END DO   ;   END DO
336         ENDIF
337      END DO
338      !
339      ! -------------------------------------------- !
340      !     9. deallocate local temporary arrays     !
341      !        if they areg larger than jpi*jpj      !  <- arbitrary max size...
342      ! -------------------------------------------- !
343      !
344      IF( iszR > jpi*jpj )   DEALLOCATE(BUFFRCV)                    ! blocking receive -> can directly deallocate
345      IF( iszS > jpi*jpj ) THEN
346         CALL mpi_waitall(8, nreq_p2p, MPI_STATUSES_IGNORE, ierr)   ! must wait before deallocate send buffer
347         DEALLOCATE(BUFFSND)
348      ENDIF
349      !
350   END SUBROUTINE lbc_lnk_pt2pt_/**/PRECISION
351
Note: See TracBrowser for help on using the repository browser.