source: codes/icosagcm/trunk/src/parallel/transfert_mpi.f90 @ 953

Last change on this file since 953 was 953, checked in by adurocher, 5 years ago

trunk : GPU implementation with OpenACC ( merge from glcp.idris.fr )

File size: 67.8 KB
Line 
1MODULE transfert_mpi_mod
2USE genmod
3USE field_mod
4IMPLICIT NONE
5 
6  TYPE array
7    INTEGER,POINTER :: value(:)=>null()
8    INTEGER,POINTER :: sign(:)=>null()
9    INTEGER         :: domain
10    INTEGER         :: rank
11    INTEGER         :: tag
12    INTEGER         :: size
13    INTEGER         :: offset
14    INTEGER         :: ireq
15    INTEGER,POINTER :: buffer(:)=>null()
16    REAL,POINTER    :: buffer_r(:)=>null()
17    INTEGER,POINTER :: src_value(:)=>null()
18  END TYPE array
19 
20  TYPE t_buffer
21    REAL,POINTER    :: r(:)
22    INTEGER         :: size
23    INTEGER         :: rank
24  END TYPE t_buffer   
25   
26  TYPE t_request
27    INTEGER :: type_field
28    INTEGER :: max_size
29    INTEGER :: size
30    LOGICAL :: vector
31    INTEGER,POINTER :: src_domain(:)
32    INTEGER,POINTER :: src_i(:)
33    INTEGER,POINTER :: src_j(:)
34    INTEGER,POINTER :: src_ind(:)
35    INTEGER,POINTER :: target_ind(:)
36    INTEGER,POINTER :: target_i(:)
37    INTEGER,POINTER :: target_j(:)
38    INTEGER,POINTER :: target_sign(:)
39    INTEGER :: nrecv
40    TYPE(ARRAY),POINTER :: recv(:)
41    INTEGER :: nsend
42    INTEGER :: nreq_mpi
43    INTEGER :: nreq_send
44    INTEGER :: nreq_recv
45    TYPE(ARRAY),POINTER :: send(:)
46  END TYPE t_request
47 
48  TYPE(t_request),SAVE,POINTER :: req_i1(:)
49  TYPE(t_request),SAVE,POINTER :: req_e1_scal(:)
50  TYPE(t_request),SAVE,POINTER :: req_e1_vect(:)
51  TYPE(t_request),SAVE,POINTER :: req_z1_scal(:)
52 
53  TYPE(t_request),SAVE,POINTER :: req_i0(:)
54  TYPE(t_request),SAVE,POINTER :: req_e0_scal(:)
55  TYPE(t_request),SAVE,POINTER :: req_e0_vect(:)
56
57  TYPE t_reorder
58    INTEGER :: ind
59    INTEGER :: rank
60    INTEGER :: tag
61    INTEGER :: isend
62  END TYPE t_reorder 
63 
64  TYPE t_message
65    CHARACTER(LEN=100) :: name ! for debug
66    TYPE(t_request), POINTER :: request(:)
67    INTEGER :: nreq
68    INTEGER :: nreq_send
69    INTEGER :: nreq_recv
70    TYPE(t_reorder), POINTER :: reorder(:)
71    INTEGER, POINTER :: mpi_req(:)
72    INTEGER, POINTER :: status(:,:)
73    TYPE(t_buffer),POINTER :: buffers(:) 
74    TYPE(t_field),POINTER :: field(:)
75    LOGICAL :: completed
76    LOGICAL :: pending
77    LOGICAL :: open      ! for debug
78    INTEGER :: number
79    LOGICAL :: ondevice=.false. !<Ready to transfer ondevice field
80  END TYPE t_message
81
82
83  INTERFACE bcast_mpi
84    MODULE PROCEDURE bcast_mpi_c,                                                     &
85                     bcast_mpi_i,bcast_mpi_i1,bcast_mpi_i2,bcast_mpi_i3,bcast_mpi_i4, &
86                     bcast_mpi_r,bcast_mpi_r1,bcast_mpi_r2,bcast_mpi_r3,bcast_mpi_r4, &
87                     bcast_mpi_l,bcast_mpi_l1,bcast_mpi_l2,bcast_mpi_l3,bcast_mpi_l4
88  END INTERFACE
89
90CONTAINS
91       
92     
93  SUBROUTINE init_transfert
94  USE profiling_mod
95  USE domain_mod
96  USE dimensions
97  USE field_mod
98  USE metric
99  USE mpipara
100  USE mpi_mod
101  IMPLICIT NONE
102  INTEGER :: ind,i,j
103
104    CALL register_id('MPI', id_mpi)
105
106    CALL create_request(field_t,req_i1)
107
108    DO ind=1,ndomain
109      CALL swap_dimensions(ind)
110      DO i=ii_begin,ii_end+1
111        CALL request_add_point(ind,i,jj_begin-1,req_i1)
112      ENDDO
113
114      DO j=jj_begin,jj_end
115        CALL request_add_point(ind,ii_end+1,j,req_i1)
116      ENDDO   
117      DO i=ii_begin,ii_end
118        CALL request_add_point(ind,i,jj_end+1,req_i1)
119      ENDDO   
120
121      DO j=jj_begin,jj_end+1
122        CALL request_add_point(ind,ii_begin-1,j,req_i1)
123      ENDDO   
124   
125    ENDDO
126 
127    CALL finalize_request(req_i1)
128
129
130    CALL create_request(field_t,req_i0)
131
132    DO ind=1,ndomain
133      CALL swap_dimensions(ind)
134   
135      DO i=ii_begin,ii_end
136        CALL request_add_point(ind,i,jj_begin,req_i0)
137      ENDDO
138
139      DO j=jj_begin,jj_end
140        CALL request_add_point(ind,ii_end,j,req_i0)
141      ENDDO   
142   
143      DO i=ii_begin,ii_end
144        CALL request_add_point(ind,i,jj_end,req_i0)
145      ENDDO   
146
147      DO j=jj_begin,jj_end
148        CALL request_add_point(ind,ii_begin,j,req_i0)
149      ENDDO   
150   
151    ENDDO
152 
153    CALL finalize_request(req_i0) 
154
155
156    CALL create_request(field_u,req_e1_scal)
157    DO ind=1,ndomain
158      CALL swap_dimensions(ind)
159      DO i=ii_begin,ii_end
160        CALL request_add_point(ind,i,jj_begin-1,req_e1_scal,rup)
161        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_scal,lup)
162      ENDDO
163
164      DO j=jj_begin,jj_end
165        CALL request_add_point(ind,ii_end+1,j,req_e1_scal,left)
166      ENDDO   
167      DO j=jj_begin,jj_end
168        CALL request_add_point(ind,ii_end+1,j-1,req_e1_scal,lup)
169      ENDDO   
170   
171      DO i=ii_begin,ii_end
172        CALL request_add_point(ind,i,jj_end+1,req_e1_scal,ldown)
173        CALL request_add_point(ind,i-1,jj_end+1,req_e1_scal,rdown)
174      ENDDO   
175
176      DO j=jj_begin,jj_end
177        CALL request_add_point(ind,ii_begin-1,j,req_e1_scal,right)
178      ENDDO   
179      DO j=jj_begin,jj_end
180        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_scal,rdown)
181      ENDDO   
182
183    ENDDO
184
185    CALL finalize_request(req_e1_scal)
186
187
188    CALL create_request(field_u,req_e0_scal)
189    DO ind=1,ndomain
190      CALL swap_dimensions(ind)
191
192
193      DO i=ii_begin+1,ii_end-1
194        CALL request_add_point(ind,i,jj_begin,req_e0_scal,right)
195        CALL request_add_point(ind,i,jj_end,req_e0_scal,right)
196      ENDDO
197   
198      DO j=jj_begin+1,jj_end-1
199        CALL request_add_point(ind,ii_begin,j,req_e0_scal,rup)
200        CALL request_add_point(ind,ii_end,j,req_e0_scal,rup)
201      ENDDO   
202
203      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_scal,left)
204      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_scal,ldown)
205      CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_scal,left)
206      CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_scal,ldown)
207
208    ENDDO
209
210    CALL finalize_request(req_e0_scal)
211
212
213   
214    CALL create_request(field_u,req_e1_vect,.TRUE.)
215    DO ind=1,ndomain
216      CALL swap_dimensions(ind)
217      DO i=ii_begin,ii_end
218        CALL request_add_point(ind,i,jj_begin-1,req_e1_vect,rup)
219        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_vect,lup)
220      ENDDO
221
222      DO j=jj_begin,jj_end
223        CALL request_add_point(ind,ii_end+1,j,req_e1_vect,left)
224      ENDDO   
225      DO j=jj_begin,jj_end
226        CALL request_add_point(ind,ii_end+1,j-1,req_e1_vect,lup)
227      ENDDO   
228   
229      DO i=ii_begin,ii_end
230        CALL request_add_point(ind,i,jj_end+1,req_e1_vect,ldown)
231        CALL request_add_point(ind,i-1,jj_end+1,req_e1_vect,rdown)
232      ENDDO   
233
234      DO j=jj_begin,jj_end
235        CALL request_add_point(ind,ii_begin-1,j,req_e1_vect,right)
236      ENDDO   
237      DO j=jj_begin,jj_end
238        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_vect,rdown)
239      ENDDO   
240
241 
242    ENDDO 
243
244    CALL finalize_request(req_e1_vect)
245   
246   
247    CALL create_request(field_u,req_e0_vect,.TRUE.)
248    DO ind=1,ndomain
249      CALL swap_dimensions(ind)
250 
251      DO i=ii_begin+1,ii_end-1
252        CALL request_add_point(ind,i,jj_begin,req_e0_vect,right)
253        CALL request_add_point(ind,i,jj_end,req_e0_vect,right)
254      ENDDO
255   
256      DO j=jj_begin+1,jj_end-1
257        CALL request_add_point(ind,ii_begin,j,req_e0_vect,rup)
258        CALL request_add_point(ind,ii_end,j,req_e0_vect,rup)
259      ENDDO   
260
261      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_vect,left)
262      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_vect,ldown)
263      CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_vect,left)
264      CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_vect,ldown)
265 
266    ENDDO 
267
268    CALL finalize_request(req_e0_vect)
269
270    CALL create_request(field_z,req_z1_scal)
271    DO ind=1,ndomain
272      CALL swap_dimensions(ind)
273      DO i=ii_begin,ii_end
274        CALL request_add_point(ind,i,jj_begin-1,req_z1_scal,vrup)
275        CALL request_add_point(ind,i+1,jj_begin-1,req_z1_scal,vup)
276      ENDDO
277
278      DO j=jj_begin,jj_end
279        CALL request_add_point(ind,ii_end+1,j,req_z1_scal,vlup)
280      ENDDO   
281      DO j=jj_begin,jj_end
282        CALL request_add_point(ind,ii_end+1,j-1,req_z1_scal,vup)
283      ENDDO   
284   
285      DO i=ii_begin,ii_end
286        CALL request_add_point(ind,i,jj_end+1,req_z1_scal,vdown)
287        CALL request_add_point(ind,i-1,jj_end+1,req_z1_scal,vrdown)
288      ENDDO   
289
290      DO j=jj_begin,jj_end
291        CALL request_add_point(ind,ii_begin-1,j,req_z1_scal,vrup)
292      ENDDO   
293      DO j=jj_begin,jj_end
294        CALL request_add_point(ind,ii_begin-1,j,req_z1_scal,vrdown)
295      ENDDO   
296
297    ENDDO
298
299    CALL finalize_request(req_z1_scal)
300
301  END SUBROUTINE init_transfert
302 
303  SUBROUTINE create_request(type_field,request,vector)
304  USE domain_mod
305  USE field_mod
306  IMPLICIT NONE
307    INTEGER :: type_field
308    TYPE(t_request),POINTER :: request(:)
309    LOGICAL,OPTIONAL :: vector
310   
311    TYPE(t_request),POINTER :: req
312    TYPE(t_domain),POINTER :: d
313    INTEGER :: ind
314    INTEGER :: max_size
315       
316    ALLOCATE(request(ndomain))
317
318    DO ind=1,ndomain
319      req=>request(ind)
320      d=>domain(ind)
321      IF (type_field==field_t) THEN
322        Max_size=2*(d%iim+2)+2*(d%jjm+2)
323      ELSE IF (type_field==field_u) THEN
324        Max_size=3*(2*(d%iim+2)+2*(d%jjm+2))
325      ELSE IF (type_field==field_z) THEN
326        Max_size=2*(2*(d%iim+2)+2*(d%jjm+2))
327      ENDIF
328
329      req%type_field=type_field
330      req%max_size=max_size*2
331      req%size=0
332      req%vector=.FALSE.
333      IF (PRESENT(vector)) req%vector=vector
334      ALLOCATE(req%src_domain(req%max_size))
335      ALLOCATE(req%src_ind(req%max_size))
336      ALLOCATE(req%target_ind(req%max_size))
337      ALLOCATE(req%src_i(req%max_size))
338      ALLOCATE(req%src_j(req%max_size))
339      ALLOCATE(req%target_i(req%max_size))
340      ALLOCATE(req%target_j(req%max_size))
341      ALLOCATE(req%target_sign(req%max_size))
342    ENDDO
343 
344  END SUBROUTINE create_request
345
346  SUBROUTINE reallocate_request(req)
347  IMPLICIT NONE
348    TYPE(t_request),POINTER :: req
349     
350    INTEGER,POINTER :: src_domain(:)
351    INTEGER,POINTER :: src_ind(:)
352    INTEGER,POINTER :: target_ind(:)
353    INTEGER,POINTER :: src_i(:)
354    INTEGER,POINTER :: src_j(:)
355    INTEGER,POINTER :: target_i(:)
356    INTEGER,POINTER :: target_j(:)
357    INTEGER,POINTER :: target_sign(:)
358
359    PRINT *,"REALLOCATE_REQUEST"
360    src_domain=>req%src_domain
361    src_ind=>req%src_ind
362    target_ind=>req%target_ind
363    src_i=>req%src_i
364    src_j=>req%src_j
365    target_i=>req%target_i
366    target_j=>req%target_j
367    target_sign=>req%target_sign
368
369    ALLOCATE(req%src_domain(req%max_size*2))
370    ALLOCATE(req%src_ind(req%max_size*2))
371    ALLOCATE(req%target_ind(req%max_size*2))
372    ALLOCATE(req%src_i(req%max_size*2))
373    ALLOCATE(req%src_j(req%max_size*2))
374    ALLOCATE(req%target_i(req%max_size*2))
375    ALLOCATE(req%target_j(req%max_size*2))
376    ALLOCATE(req%target_sign(req%max_size*2))
377   
378    req%src_domain(1:req%max_size)=src_domain(:)
379    req%src_ind(1:req%max_size)=src_ind(:)
380    req%target_ind(1:req%max_size)=target_ind(:)
381    req%src_i(1:req%max_size)=src_i(:)
382    req%src_j(1:req%max_size)=src_j(:)
383    req%target_i(1:req%max_size)=target_i(:)
384    req%target_j(1:req%max_size)=target_j(:)
385    req%target_sign(1:req%max_size)=target_sign(:)
386   
387    req%max_size=req%max_size*2
388         
389    DEALLOCATE(src_domain)
390    DEALLOCATE(src_ind)
391    DEALLOCATE(target_ind)
392    DEALLOCATE(src_i)
393    DEALLOCATE(src_j)
394    DEALLOCATE(target_i)
395    DEALLOCATE(target_j)
396    DEALLOCATE(target_sign)
397
398  END SUBROUTINE reallocate_request
399
400     
401    SUBROUTINE request_add_point(ind,i,j,request,pos)
402    USE domain_mod
403    USE field_mod
404    IMPLICIT NONE
405      INTEGER,INTENT(IN)            ::  ind
406      INTEGER,INTENT(IN)            :: i
407      INTEGER,INTENT(IN)            :: j
408      TYPE(t_request),POINTER :: request(:)
409      INTEGER,INTENT(IN),OPTIONAL  :: pos
410     
411      INTEGER :: src_domain
412      INTEGER :: src_iim,src_i,src_j,src_n,src_pos,src_delta
413      TYPE(t_request),POINTER :: req
414      TYPE(t_domain),POINTER :: d
415     
416      req=>request(ind)
417      d=>domain(ind)
418     
419      IF (req%max_size==req%size) CALL reallocate_request(req)
420      req%size=req%size+1
421      IF (req%type_field==field_t) THEN
422        src_domain=domain(ind)%assign_domain(i,j)
423        src_iim=domain_glo(src_domain)%iim
424        src_i=domain(ind)%assign_i(i,j)
425        src_j=domain(ind)%assign_j(i,j)
426
427        req%target_ind(req%size)=(j-1)*d%iim+i
428        req%target_sign(req%size)=1
429        req%src_domain(req%size)=src_domain
430        req%src_ind(req%size)=(src_j-1)*src_iim+src_i
431      ELSE IF (req%type_field==field_u) THEN
432        IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme'
433
434        src_domain=domain(ind)%edge_assign_domain(pos-1,i,j)
435        src_iim=domain_glo(src_domain)%iim
436        src_i=domain(ind)%edge_assign_i(pos-1,i,j)
437        src_j=domain(ind)%edge_assign_j(pos-1,i,j)
438        src_n=(src_j-1)*src_iim+src_i
439        src_delta=domain(ind)%delta(i,j)
440        src_pos=domain(ind)%edge_assign_pos(pos-1,i,j)+1
441               
442        req%target_ind(req%size)=(j-1)*d%iim+i+d%u_pos(pos)
443
444        req%target_sign(req%size)= 1
445        IF (req%vector) req%target_sign(req%size)= domain(ind)%edge_assign_sign(pos-1,i,j)
446
447        req%src_domain(req%size)=src_domain
448        req%src_ind(req%size)=src_n+domain_glo(src_domain)%u_pos(src_pos)
449
450      ELSE IF (req%type_field==field_z) THEN
451        IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme'
452
453        src_domain=domain(ind)%vertex_assign_domain(pos-1,i,j)
454        src_iim=domain_glo(src_domain)%iim
455        src_i=domain(ind)%vertex_assign_i(pos-1,i,j)
456        src_j=domain(ind)%vertex_assign_j(pos-1,i,j)
457        src_n=(src_j-1)*src_iim+src_i
458        src_delta=domain(ind)%delta(i,j)
459        src_pos=domain(ind)%vertex_assign_pos(pos-1,i,j)+1
460
461       
462        req%target_ind(req%size)=(j-1)*d%iim+i+d%z_pos(pos)
463        req%target_sign(req%size)=1
464        req%src_domain(req%size)=src_domain
465        req%src_ind(req%size)=src_n+domain_glo(src_domain)%z_pos(src_pos)
466      ENDIF
467  END SUBROUTINE request_add_point
468 
469 
470  SUBROUTINE Finalize_request(request)
471  USE mpipara
472  USE domain_mod
473  USE mpi_mod
474  IMPLICIT NONE
475    TYPE(t_request),POINTER :: request(:)
476    TYPE(t_request),POINTER :: req, req_src
477    INTEGER :: nb_domain_recv(0:mpi_size-1)
478    INTEGER :: nb_domain_send(0:mpi_size-1)
479    INTEGER :: tag_rank(0:mpi_size-1)
480    INTEGER :: nb_data_domain_recv(ndomain_glo)
481    INTEGER :: list_domain_recv(ndomain_glo)
482    INTEGER,ALLOCATABLE :: list_domain_send(:)
483    INTEGER             :: list_domain(ndomain)
484
485    INTEGER :: rank,i,j,pos
486    INTEGER :: size_,ind_glo,ind_loc
487    INTEGER :: isend, irecv, ireq, nreq, nsend, nrecv
488    INTEGER, ALLOCATABLE :: mpi_req(:)
489    INTEGER, ALLOCATABLE :: status(:,:)
490    INTEGER, ALLOCATABLE :: rank_list(:)
491    INTEGER, ALLOCATABLE :: offset(:)
492    LOGICAL,PARAMETER :: debug = .FALSE.
493
494 
495    IF (.NOT. using_mpi) RETURN
496   
497    DO ind_loc=1,ndomain
498      req=>request(ind_loc)
499     
500      nb_data_domain_recv(:) = 0
501      nb_domain_recv(:) = 0
502      tag_rank(:)=0
503     
504      DO i=1,req%size
505        ind_glo=req%src_domain(i)
506        nb_data_domain_recv(ind_glo)=nb_data_domain_recv(ind_glo)+1
507      ENDDO
508 
509      DO ind_glo=1,ndomain_glo
510        IF ( nb_data_domain_recv(ind_glo) > 0 )  nb_domain_recv(domglo_rank(ind_glo))=nb_domain_recv(domglo_rank(ind_glo))+1
511      ENDDO
512
513      req%nrecv=sum(nb_domain_recv(:))
514      ALLOCATE(req%recv(req%nrecv))
515
516      irecv=0
517      DO ind_glo=1,ndomain_glo
518        IF (nb_data_domain_recv(ind_glo)>0) THEN
519          irecv=irecv+1
520          list_domain_recv(ind_glo)=irecv
521          req%recv(irecv)%rank=domglo_rank(ind_glo)
522          req%recv(irecv)%size=nb_data_domain_recv(ind_glo)
523          req%recv(irecv)%domain=domglo_loc_ind(ind_glo)
524          ALLOCATE(req%recv(irecv)%value(req%recv(irecv)%size))
525          ALLOCATE(req%recv(irecv)%sign(req%recv(irecv)%size))
526          ALLOCATE(req%recv(irecv)%buffer(req%recv(irecv)%size))
527        ENDIF
528      ENDDO
529     
530      req%recv(:)%size=0
531      irecv=0
532      DO i=1,req%size
533        irecv=list_domain_recv(req%src_domain(i))
534        req%recv(irecv)%size=req%recv(irecv)%size+1
535        size_=req%recv(irecv)%size
536        req%recv(irecv)%value(size_)=req%src_ind(i)
537        req%recv(irecv)%buffer(size_)=req%target_ind(i)
538        req%recv(irecv)%sign(size_)=req%target_sign(i)
539      ENDDO
540    ENDDO
541
542    nb_domain_recv(:) = 0   
543    DO ind_loc=1,ndomain
544      req=>request(ind_loc)
545     
546      DO irecv=1,req%nrecv
547        rank=req%recv(irecv)%rank
548        nb_domain_recv(rank)=nb_domain_recv(rank)+1
549      ENDDO
550    ENDDO
551   
552    CALL MPI_Alltoall(nb_domain_recv,1,MPI_INTEGER,nb_domain_send,1,MPI_INTEGER,comm_icosa,ierr)     
553   
554
555    ALLOCATE(list_domain_send(sum(nb_domain_send)))
556   
557    nreq=sum(nb_domain_recv(:))+sum(nb_domain_send(:))
558    ALLOCATE(mpi_req(nreq))
559    ALLOCATE(status(MPI_STATUS_SIZE,nreq))
560   
561
562    ireq=0
563    DO ind_loc=1,ndomain
564      req=>request(ind_loc)
565      DO irecv=1,req%nrecv
566        ireq=ireq+1
567        CALL MPI_ISEND(req%recv(irecv)%domain,1,MPI_INTEGER,req%recv(irecv)%rank,0,comm_icosa, mpi_req(ireq),ierr)
568        IF (debug) PRINT *,"Isend ",req%recv(irecv)%domain, "from ", mpi_rank, "to ",req%recv(irecv)%rank, "tag ",0
569      ENDDO
570    ENDDO
571
572    IF (debug) PRINT *,"------------"   
573    j=0
574    DO rank=0,mpi_size-1
575      DO i=1,nb_domain_send(rank)
576        j=j+1
577        ireq=ireq+1
578        CALL MPI_IRECV(list_domain_send(j),1,MPI_INTEGER,rank,0,comm_icosa, mpi_req(ireq),ierr)
579        IF (debug) PRINT *,"IRecv ",list_domain_send(j), "from ", rank, "to ",mpi_rank, "tag ",0
580      ENDDO
581    ENDDO
582    IF (debug) PRINT *,"------------"   
583   
584    CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
585   
586    list_domain(:)=0
587    DO i=1,sum(nb_domain_send)
588      ind_loc=list_domain_send(i)
589      list_domain(ind_loc)=list_domain(ind_loc)+1
590    ENDDO
591   
592    DO ind_loc=1,ndomain
593      req=>request(ind_loc)
594      req%nsend=list_domain(ind_loc)
595      ALLOCATE(req%send(req%nsend))
596    ENDDO
597
598    IF (debug) PRINT *,"------------"   
599   
600   ireq=0 
601   DO ind_loc=1,ndomain
602     req=>request(ind_loc)
603     
604     DO irecv=1,req%nrecv
605       ireq=ireq+1
606       CALL MPI_ISEND(mpi_rank,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
607       IF (debug) PRINT *,"Isend ",mpi_rank, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
608     ENDDO
609    IF (debug) PRINT *,"------------"   
610     
611     DO isend=1,req%nsend
612       ireq=ireq+1
613       CALL MPI_IRECV(req%send(isend)%rank,1,MPI_INTEGER,MPI_ANY_SOURCE,ind_loc,comm_icosa, mpi_req(ireq),ierr)
614       IF (debug) PRINT *,"IRecv ",req%send(isend)%rank, "from ", "xxx", "to ",mpi_rank, "tag ",ind_loc
615     ENDDO
616   ENDDO
617
618   IF (debug) PRINT *,"------------"   
619
620   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
621   CALL MPI_BARRIER(comm_icosa,ierr)
622
623   IF (debug) PRINT *,"------------"   
624
625   ireq=0 
626   DO ind_loc=1,ndomain
627     req=>request(ind_loc)
628     
629     DO irecv=1,req%nrecv
630       ireq=ireq+1
631       CALL MPI_ISEND(ind_loc,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
632       IF (debug) PRINT *,"Isend ",ind_loc, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
633     ENDDO
634
635     IF (debug) PRINT *,"------------"   
636     
637     DO isend=1,req%nsend
638       ireq=ireq+1
639       CALL MPI_IRECV(req%send(isend)%domain,1,MPI_INTEGER,req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr)
640       IF (debug) PRINT *,"IRecv ",req%send(isend)%domain, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
641     ENDDO
642   ENDDO
643   IF (debug) PRINT *,"------------"   
644   
645   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
646   CALL MPI_BARRIER(comm_icosa,ierr)
647   IF (debug) PRINT *,"------------"   
648
649   ireq=0
650   DO ind_loc=1,ndomain
651     req=>request(ind_loc)
652     
653     DO irecv=1,req%nrecv
654       ireq=ireq+1
655       req%recv(irecv)%tag=tag_rank(req%recv(irecv)%rank)
656       tag_rank(req%recv(irecv)%rank)=tag_rank(req%recv(irecv)%rank)+1
657       CALL MPI_ISEND(req%recv(irecv)%tag,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
658       IF (debug) PRINT *,"Isend ",req%recv(irecv)%tag, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
659     ENDDO
660   IF (debug) PRINT *,"------------"   
661     
662     DO isend=1,req%nsend
663       ireq=ireq+1
664       CALL MPI_IRECV(req%send(isend)%tag,1,MPI_INTEGER,req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr)
665       IF (debug) PRINT *,"IRecv ",req%send(isend)%tag, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
666     ENDDO
667   ENDDO
668   IF (debug) PRINT *,"------------"   
669   
670   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
671   CALL MPI_BARRIER(comm_icosa,ierr)
672
673
674   IF (debug) PRINT *,"------------"   
675
676   ireq=0 
677   DO ind_loc=1,ndomain
678     req=>request(ind_loc)
679     
680     DO irecv=1,req%nrecv
681       ireq=ireq+1
682       CALL MPI_ISEND(req%recv(irecv)%size,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
683       IF (debug) PRINT *,"Isend ",req%recv(irecv)%size, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
684     ENDDO
685     IF (debug) PRINT *,"------------"   
686     
687     DO isend=1,req%nsend
688       ireq=ireq+1
689       CALL MPI_IRECV(req%send(isend)%size,1,MPI_INTEGER,req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
690       IF (debug) PRINT *,"IRecv ",req%send(isend)%size, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
691     ENDDO
692   ENDDO
693
694   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
695
696   ireq=0 
697   DO ind_loc=1,ndomain
698     req=>request(ind_loc)
699     
700     DO irecv=1,req%nrecv
701       ireq=ireq+1
702       CALL MPI_ISEND(req%recv(irecv)%value,req%recv(irecv)%size,MPI_INTEGER,&
703            req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
704     ENDDO
705     
706     DO isend=1,req%nsend
707       ireq=ireq+1
708       ALLOCATE(req%send(isend)%value(req%send(isend)%size))
709       CALL MPI_IRECV(req%send(isend)%value,req%send(isend)%size,MPI_INTEGER,&
710            req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
711     ENDDO
712   ENDDO
713
714   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
715
716   DO ind_loc=1,ndomain
717     req=>request(ind_loc)
718     
719     DO irecv=1,req%nrecv
720       req%recv(irecv)%value(:)=req%recv(irecv)%buffer(:)
721       req%recv(irecv)%sign(:) =req%recv(irecv)%sign(:)
722       DEALLOCATE(req%recv(irecv)%buffer)
723     ENDDO
724   ENDDO 
725   
726
727! domain is on the same mpi process => copie memory to memory
728   
729   DO ind_loc=1,ndomain
730     req=>request(ind_loc)
731     
732     DO irecv=1,req%nrecv
733   
734       IF (req%recv(irecv)%rank==mpi_rank) THEN
735           req_src=>request(req%recv(irecv)%domain)
736           DO isend=1,req_src%nsend
737             IF (req_src%send(isend)%rank==mpi_rank .AND. req_src%send(isend)%tag==req%recv(irecv)%tag) THEN
738               req%recv(irecv)%src_value => req_src%send(isend)%value
739               IF ( size(req%recv(irecv)%value) /= size(req_src%send(isend)%value)) THEN
740                 PRINT *,ind_loc, irecv, isend, size(req%recv(irecv)%value), size(req_src%send(isend)%value)
741                 STOP "size(req%recv(irecv)%value) /= size(req_src%send(isend)%value"
742               ENDIF
743             ENDIF
744           ENDDO
745       ENDIF
746     
747     ENDDO
748   ENDDO
749   
750! true number of mpi request
751
752   request(:)%nreq_mpi=0
753   request(:)%nreq_send=0
754   request(:)%nreq_recv=0
755   ALLOCATE(rank_list(sum(request(:)%nsend)))
756   ALLOCATE(offset(sum(request(:)%nsend)))
757   offset(:)=0
758   
759   nsend=0
760   DO ind_loc=1,ndomain
761     req=>request(ind_loc)
762
763     DO isend=1,req%nsend
764       IF (req%send(isend)%rank/=mpi_rank) THEN
765         pos=0
766         DO i=1,nsend
767           IF (req%send(isend)%rank==rank_list(i)) EXIT
768           pos=pos+1
769         ENDDO
770       
771         IF (pos==nsend) THEN
772           nsend=nsend+1
773           req%nreq_mpi=req%nreq_mpi+1
774           req%nreq_send=req%nreq_send+1
775           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
776             rank_list(nsend)=req%send(isend)%rank
777           ELSE
778             rank_list(nsend)=-1
779           ENDIF
780         ENDIF
781         
782         pos=pos+1
783         req%send(isend)%ireq=pos
784         req%send(isend)%offset=offset(pos)
785         offset(pos)=offset(pos)+req%send(isend)%size
786       ENDIF
787     ENDDO
788   ENDDO
789
790   DEALLOCATE(rank_list)
791   DEALLOCATE(offset)
792     
793   ALLOCATE(rank_list(sum(request(:)%nrecv)))
794   ALLOCATE(offset(sum(request(:)%nrecv)))
795   offset(:)=0
796   
797   nrecv=0
798   DO ind_loc=1,ndomain
799     req=>request(ind_loc)
800
801     DO irecv=1,req%nrecv
802       IF (req%recv(irecv)%rank/=mpi_rank) THEN
803         pos=0
804         DO i=1,nrecv
805           IF (req%recv(irecv)%rank==rank_list(i)) EXIT
806           pos=pos+1
807         ENDDO
808       
809         IF (pos==nrecv) THEN
810           nrecv=nrecv+1
811           req%nreq_mpi=req%nreq_mpi+1
812           req%nreq_recv=req%nreq_recv+1
813           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
814             rank_list(nrecv)=req%recv(irecv)%rank
815           ELSE
816             rank_list(nrecv)=-1
817           ENDIF
818         ENDIF
819       
820         pos=pos+1
821         req%recv(irecv)%ireq=nsend+pos
822         req%recv(irecv)%offset=offset(pos)
823         offset(pos)=offset(pos)+req%recv(irecv)%size
824       ENDIF
825     ENDDO
826   ENDDO 
827
828! get the offsets   
829
830   ireq=0 
831   DO ind_loc=1,ndomain
832     req=>request(ind_loc)
833     
834     DO irecv=1,req%nrecv
835       ireq=ireq+1
836       CALL MPI_ISEND(req%recv(irecv)%offset,1,MPI_INTEGER,&
837            req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
838     ENDDO
839     
840     DO isend=1,req%nsend
841       ireq=ireq+1
842       CALL MPI_IRECV(req%send(isend)%offset,1,MPI_INTEGER,&
843            req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
844     ENDDO
845   ENDDO
846
847   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
848     
849       
850  END SUBROUTINE Finalize_request 
851
852
853  SUBROUTINE init_message_seq(field, request, message, name)
854  USE field_mod
855  USE domain_mod
856  USE mpi_mod
857  USE mpipara
858  USE mpi_mod
859  IMPLICIT NONE
860    TYPE(t_field),POINTER :: field(:)
861    TYPE(t_request),POINTER :: request(:)
862    TYPE(t_message) :: message
863    CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: name
864!$OMP MASTER   
865    message%request=>request
866    IF(PRESENT(name)) THEN
867       message%name = TRIM(name)
868    ELSE
869       message%name = 'unknown'
870    END IF
871!$OMP END MASTER   
872!$OMP BARRIER   
873
874  END SUBROUTINE init_message_seq
875
876  SUBROUTINE send_message_seq(field,message)
877  USE field_mod
878  USE domain_mod
879  USE mpi_mod
880  USE mpipara
881  USE omp_para
882  USE trace
883  IMPLICIT NONE
884    TYPE(t_field),POINTER :: field(:)
885    TYPE(t_message) :: message
886
887    CALL transfert_request_seq(field,message%request)
888   
889  END SUBROUTINE send_message_seq
890 
891  SUBROUTINE test_message_seq(message)
892  IMPLICIT NONE
893    TYPE(t_message) :: message
894  END SUBROUTINE  test_message_seq
895 
896   
897  SUBROUTINE wait_message_seq(message)
898  IMPLICIT NONE
899    TYPE(t_message) :: message
900   
901  END SUBROUTINE wait_message_seq   
902
903  SUBROUTINE transfert_message_seq(field,message)
904  USE field_mod
905  USE domain_mod
906  USE mpi_mod
907  USE mpipara
908  USE omp_para
909  USE trace
910  IMPLICIT NONE
911    TYPE(t_field),POINTER :: field(:)
912    TYPE(t_message) :: message
913
914   CALL send_message_seq(field,message)
915   
916  END SUBROUTINE transfert_message_seq   
917   
918
919
920   
921  SUBROUTINE init_message_mpi(field,request, message, name)
922  USE field_mod
923  USE domain_mod
924  USE mpi_mod
925  USE mpipara
926  USE mpi_mod
927  IMPLICIT NONE
928 
929    TYPE(t_field),POINTER :: field(:)
930    TYPE(t_request),POINTER :: request(:)
931    TYPE(t_message) :: message
932    CHARACTER(LEN=*), INTENT(IN),OPTIONAL :: name
933
934    TYPE(t_request),POINTER :: req
935    INTEGER :: irecv,isend
936    INTEGER :: ireq,nreq
937    INTEGER :: ind
938    INTEGER :: dim3,dim4
939    INTEGER,SAVE :: message_number=0
940!    TYPE(t_reorder),POINTER :: reorder(:)
941!    TYPE(t_reorder) :: reorder_swap
942
943!$OMP BARRIER
944!$OMP MASTER
945    !init off-device
946    message%ondevice=.false.
947
948    IF(PRESENT(name)) THEN
949       message%name = TRIM(name)
950    ELSE
951       message%name = 'unknown'
952    END IF
953    message%number=message_number
954    message_number=message_number+1
955    IF (message_number==100) message_number=0
956
957 
958    message%request=>request
959    message%nreq=sum(message%request(:)%nreq_mpi)
960    message%nreq_send=sum(message%request(:)%nreq_send)
961    message%nreq_recv=sum(message%request(:)%nreq_recv)
962    nreq=message%nreq
963
964    ALLOCATE(message%mpi_req(nreq))
965    ALLOCATE(message%buffers(nreq))
966    ALLOCATE(message%status(MPI_STATUS_SIZE,nreq))
967    message%buffers(:)%size=0
968    message%pending=.FALSE.
969    message%completed=.FALSE.
970    message%open=.FALSE.
971
972    DO ind=1,ndomain
973      req=>request(ind)
974      DO isend=1,req%nsend
975        IF (req%send(isend)%rank/=mpi_rank) THEN
976          ireq=req%send(isend)%ireq 
977          message%buffers(ireq)%size=message%buffers(ireq)%size+req%send(isend)%size
978          message%buffers(ireq)%rank=req%send(isend)%rank
979        ENDIF
980      ENDDO
981      DO irecv=1,req%nrecv
982        IF (req%recv(irecv)%rank/=mpi_rank) THEN
983          ireq=req%recv(irecv)%ireq 
984          message%buffers(ireq)%size=message%buffers(ireq)%size+req%recv(irecv)%size
985          message%buffers(ireq)%rank=req%recv(irecv)%rank
986        ENDIF
987      ENDDO
988    ENDDO
989
990
991    IF (field(1)%data_type==type_real) THEN
992
993      IF (field(1)%ndim==2) THEN
994     
995        DO ireq=1,message%nreq
996          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
997        ENDDO
998     
999      ELSE  IF (field(1)%ndim==3) THEN
1000     
1001        dim3=size(field(1)%rval3d,2)
1002        DO ireq=1,message%nreq
1003          message%buffers(ireq)%size=message%buffers(ireq)%size*dim3
1004          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
1005        ENDDO
1006     
1007      ELSE  IF (field(1)%ndim==4) THEN
1008        dim3=size(field(1)%rval4d,2)
1009        dim4=size(field(1)%rval4d,3)
1010        DO ireq=1,message%nreq
1011          message%buffers(ireq)%size=message%buffers(ireq)%size*dim3*dim4
1012          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
1013        ENDDO
1014      ENDIF     
1015    ENDIF
1016     
1017         
1018   
1019! ! Reorder the request, so recv request are done in the same order than send request
1020
1021!    nreq_send=sum(request(:)%nsend) 
1022!    message%nreq_send=nreq_send
1023!    ALLOCATE(message%reorder(nreq_send))
1024!    reorder=>message%reorder
1025!    ireq=0
1026!    DO ind=1,ndomain
1027!      req=>request(ind)
1028!      DO isend=1,req%nsend
1029!        ireq=ireq+1
1030!        reorder(ireq)%ind=ind
1031!        reorder(ireq)%isend=isend
1032!        reorder(ireq)%tag=req%send(isend)%tag
1033!      ENDDO
1034!    ENDDO
1035
1036! ! do a very very bad sort
1037!    DO i=1,nreq_send-1
1038!      DO j=i+1,nreq_send
1039!        IF (reorder(i)%tag > reorder(j)%tag) THEN
1040!          reorder_swap=reorder(i)
1041!          reorder(i)=reorder(j)
1042!          reorder(j)=reorder_swap
1043!        ENDIF
1044!      ENDDO
1045!    ENDDO
1046!    PRINT *,"reorder ",reorder(:)%tag
1047   
1048 
1049!$OMP END MASTER
1050!$OMP BARRIER   
1051
1052  END SUBROUTINE init_message_mpi
1053 
1054  SUBROUTINE Finalize_message_mpi(field,message)
1055  USE field_mod
1056  USE domain_mod
1057  USE mpi_mod
1058  USE mpipara
1059  USE mpi_mod
1060  IMPLICIT NONE
1061    TYPE(t_field),POINTER :: field(:)
1062    TYPE(t_message) :: message
1063
1064    INTEGER :: ireq, ibuff
1065
1066!$OMP BARRIER
1067!$OMP MASTER
1068
1069
1070    IF (field(1)%data_type==type_real) THEN
1071
1072      IF (field(1)%ndim==2) THEN
1073     
1074        DO ireq=1,message%nreq
1075          CALL free_mpi_buffer(message%buffers(ireq)%r)
1076        ENDDO
1077     
1078      ELSE  IF (field(1)%ndim==3) THEN
1079
1080        DO ireq=1,message%nreq
1081          CALL free_mpi_buffer(message%buffers(ireq)%r)
1082        ENDDO
1083     
1084      ELSE  IF (field(1)%ndim==4) THEN
1085
1086        DO ireq=1,message%nreq
1087          CALL free_mpi_buffer(message%buffers(ireq)%r)
1088        ENDDO
1089
1090      ENDIF
1091    ENDIF
1092
1093    !deallocate device data if ondevice
1094    if(message%ondevice) then
1095      do ireq=1, ndomain
1096        do ibuff=1,message%request(ireq)%nsend
1097          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer(:))
1098          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer_r(:))
1099          !$acc exit data delete(message%request(ireq)%send(ibuff)%sign(:))
1100          !$acc exit data delete(message%request(ireq)%send(ibuff)%src_value(:))
1101          !$acc exit data delete(message%request(ireq)%send(ibuff)%value(:))
1102        end do
1103        do ibuff=1,message%request(ireq)%nrecv
1104          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer(:))
1105          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer_r(:))
1106          !$acc exit data delete(message%request(ireq)%recv(ibuff)%sign(:))
1107          !$acc exit data delete(message%request(ireq)%recv(ibuff)%src_value(:))
1108          !$acc exit data delete(message%request(ireq)%recv(ibuff)%value(:))
1109        end do
1110      end do
1111      DO ireq=1,message%nreq
1112        !$acc exit data delete(message%buffers(ireq)%r)
1113      ENDDO
1114      message%ondevice=.false.
1115    end if
1116
1117    DEALLOCATE(message%mpi_req)
1118    DEALLOCATE(message%buffers)
1119    DEALLOCATE(message%status)
1120
1121!$OMP END MASTER
1122!$OMP BARRIER
1123
1124     
1125  END SUBROUTINE Finalize_message_mpi
1126
1127
1128 
1129  SUBROUTINE barrier
1130  USE mpi_mod
1131  USE mpipara
1132  IMPLICIT NONE
1133   
1134    CALL MPI_BARRIER(comm_icosa,ierr)
1135   
1136  END SUBROUTINE barrier 
1137   
1138  SUBROUTINE transfert_message_mpi(field,message)
1139  USE field_mod
1140  IMPLICIT NONE
1141    TYPE(t_field),POINTER :: field(:)
1142    TYPE(t_message) :: message
1143   
1144    CALL send_message_mpi(field,message)
1145    CALL wait_message_mpi(message)
1146   
1147  END SUBROUTINE transfert_message_mpi
1148
1149
1150  !!!Update buffers on device for 'message'
1151  !!! does create_device_message when not already ondevice
1152  SUBROUTINE update_device_message_mpi(message)
1153    USE domain_mod
1154    IMPLICIT NONE
1155    TYPE(t_message), intent(inout) :: message
1156    INTEGER :: ireq, ibuff
1157
1158    !if(.not. message%ondevice) call create_device_message_mpi(message)
1159
1160    do ireq=1, ndomain
1161      do ibuff=1,message%request(ireq)%nsend
1162        !device buffers updated even if pointers not attached :
1163        !non allocated buffers in 'message' must be set to NULL()
1164        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer(:)) async
1165        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer_r(:)) async
1166        !$acc enter data copyin(message%request(ireq)%send(ibuff)%sign(:)) async
1167        !$acc enter data copyin(message%request(ireq)%send(ibuff)%src_value(:)) async
1168        !$acc enter data copyin(message%request(ireq)%send(ibuff)%value(:)) async
1169      end do
1170      do ibuff=1,message%request(ireq)%nrecv
1171        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer(:)) async
1172        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer_r(:)) async
1173        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%sign(:)) async
1174        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%src_value(:)) async
1175        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%value(:)) async
1176      end do
1177    end do
1178    DO ireq=1,message%nreq
1179      !$acc enter data copyin(message%buffers(ireq)%r) async
1180    ENDDO
1181    message%ondevice=.true.
1182  END SUBROUTINE
1183
1184  !TODO : add openacc with multiple process
1185  SUBROUTINE send_message_mpi(field,message)
1186  USE abort_mod
1187  USE profiling_mod
1188  USE field_mod
1189  USE domain_mod
1190  USE mpi_mod
1191  USE mpipara
1192  USE omp_para
1193  USE trace
1194  USE abort_mod
1195  IMPLICIT NONE
1196    TYPE(t_field),POINTER :: field(:)
1197    TYPE(t_message) :: message
1198    REAL(rstd),POINTER :: rval2d(:), src_rval2d(:) 
1199    REAL(rstd),POINTER :: rval3d(:,:), src_rval3d(:,:) 
1200    REAL(rstd),POINTER :: rval4d(:,:,:), src_rval4d(:,:,:) 
1201    REAL(rstd),POINTER :: buffer_r(:) 
1202    INTEGER,POINTER :: value(:) 
1203    INTEGER,POINTER :: sgn(:) 
1204    TYPE(ARRAY),POINTER :: recv,send 
1205    TYPE(t_request),POINTER :: req
1206    INTEGER :: irecv,isend
1207    INTEGER :: ireq
1208    INTEGER :: ind,n
1209    INTEGER :: dim3,dim4,d3,d4
1210    INTEGER,POINTER :: src_value(:)
1211    INTEGER :: offset,msize,moffset,rank
1212    INTEGER :: lbegin, lend
1213    INTEGER :: max_req
1214
1215!    CALL trace_start("send_message_mpi")
1216
1217    CALL enter_profile(id_mpi)
1218
1219    !Prepare 'message' for on-device copies if field is on device
1220    if( field(1)%ondevice .AND. .NOT. message%ondevice ) call update_device_message_mpi(message)
1221
1222!$OMP BARRIER
1223
1224
1225!$OMP MASTER
1226    IF(message%open) THEN
1227       PRINT *, 'send_message_mpi : message ' // TRIM(message%name) // &
1228            ' is still open, no call to wait_message_mpi after last send_message_mpi'
1229       CALL dynamico_abort( "send_message_mpi : message still open" )
1230    END IF
1231    message%open=.TRUE. ! will be set to .FALSE. by wait_message_mpi
1232
1233    message%field=>field
1234
1235    IF (message%nreq>0) THEN
1236      message%completed=.FALSE.
1237      message%pending=.TRUE.
1238    ELSE
1239      message%completed=.TRUE.
1240      message%pending=.FALSE.
1241    ENDIF
1242!$OMP END MASTER
1243!$OMP BARRIER
1244     
1245    IF (field(1)%data_type==type_real) THEN
1246      IF (field(1)%ndim==2) THEN
1247
1248        DO ind=1,ndomain
1249          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1250         
1251          rval2d=>field(ind)%rval2d
1252       
1253          req=>message%request(ind)
1254          DO isend=1,req%nsend
1255            send=>req%send(isend)
1256            value=>send%value
1257
1258           
1259            IF (send%rank/=mpi_rank) THEN
1260              ireq=send%ireq
1261
1262              buffer_r=>message%buffers(ireq)%r
1263              offset=send%offset
1264              msize=send%size
1265              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1266              DO n=1,msize
1267                buffer_r(offset+n)=rval2d(value(n))
1268              ENDDO
1269
1270              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1271                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1272                !$OMP CRITICAL           
1273                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1274                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1275                !$OMP END CRITICAL
1276              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1277                CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1278                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1279                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1280              ENDIF
1281             
1282            ENDIF
1283          ENDDO
1284        ENDDO
1285       
1286        DO ind=1,ndomain
1287          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1288          rval2d=>field(ind)%rval2d
1289          req=>message%request(ind)       
1290
1291          DO irecv=1,req%nrecv
1292            recv=>req%recv(irecv)
1293
1294            IF (recv%rank==mpi_rank) THEN
1295
1296              value=>recv%value
1297              src_value => recv%src_value
1298              src_rval2d=>field(recv%domain)%rval2d
1299              sgn=>recv%sign
1300              msize=recv%size
1301              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1302              DO n=1,msize
1303                rval2d(value(n))=src_rval2d(src_value(n))*sgn(n)
1304              ENDDO
1305               
1306                   
1307            ELSE
1308           
1309              ireq=recv%ireq
1310              buffer_r=>message%buffers(ireq)%r
1311              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1312                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1313               !$OMP CRITICAL           
1314                CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,               &
1315                  recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1316               !$OMP END CRITICAL
1317              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1318                 CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1319                 CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,              &
1320                   recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1321              ENDIF
1322           
1323            ENDIF
1324          ENDDO
1325       
1326        ENDDO
1327       
1328     
1329      ELSE  IF (field(1)%ndim==3) THEN
1330        max_req=0
1331        DO ind=1,ndomain
1332          req=>message%request(ind)
1333          IF (req%nsend>max_req) max_req=req%nsend
1334        ENDDO
1335             
1336        DO ind=1,ndomain
1337          IF (.NOT. assigned_domain(ind) ) CYCLE
1338
1339          dim3=size(field(ind)%rval3d,2)
1340          CALL distrib_level(1,dim3, lbegin,lend)
1341
1342          rval3d=>field(ind)%rval3d
1343          req=>message%request(ind)
1344 
1345          DO isend=1,req%nsend
1346            send=>req%send(isend)
1347            value=>send%value
1348
1349            IF (send%rank/=mpi_rank) THEN
1350              ireq=send%ireq
1351              buffer_r=>message%buffers(ireq)%r
1352
1353              msize=send%size
1354              moffset=send%offset
1355              CALL trace_start("copy_to_buffer")
1356
1357              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1358              DO d3=lbegin,lend
1359                offset=moffset*dim3 + (d3-1)*msize
1360                !$acc loop
1361                DO n=1,msize
1362                  buffer_r(n+offset)=rval3d(value(n),d3)
1363                ENDDO
1364              ENDDO
1365              CALL trace_end("copy_to_buffer")
1366
1367              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1368                !$OMP BARRIER
1369              ENDIF
1370             
1371              IF (is_omp_level_master) THEN
1372                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1373                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1374                  !$OMP CRITICAL   
1375                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1376                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1377                  !$OMP END CRITICAL
1378                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1379                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1380                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1381                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1382                ENDIF
1383              ENDIF
1384            ELSE
1385              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1386                !$OMP BARRIER
1387              ENDIF
1388            ENDIF
1389          ENDDO
1390
1391          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1392            DO isend=req%nsend+1,max_req
1393              !$OMP BARRIER
1394            ENDDO
1395          ENDIF
1396
1397        ENDDO
1398         
1399        DO ind=1,ndomain
1400          IF (.NOT. assigned_domain(ind) ) CYCLE
1401          dim3=size(field(ind)%rval3d,2)
1402          CALL distrib_level(1,dim3, lbegin,lend)
1403          rval3d=>field(ind)%rval3d
1404          req=>message%request(ind)
1405
1406          DO irecv=1,req%nrecv
1407            recv=>req%recv(irecv)
1408
1409            IF (recv%rank==mpi_rank) THEN
1410              value=>recv%value
1411              src_value => recv%src_value
1412              src_rval3d=>field(recv%domain)%rval3d
1413              sgn=>recv%sign
1414              msize=recv%size
1415
1416              CALL trace_start("copy_data")
1417              !$acc parallel loop collapse(2) default(present) async if (field(ind)%ondevice)
1418              DO d3=lbegin,lend
1419                DO n=1,msize
1420                  rval3d(value(n),d3)=src_rval3d(src_value(n),d3)*sgn(n)
1421                ENDDO
1422              ENDDO
1423              CALL trace_end("copy_data")
1424
1425            ELSE
1426              ireq=recv%ireq
1427              buffer_r=>message%buffers(ireq)%r
1428 
1429              IF (is_omp_level_master) THEN
1430                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1431                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1432                  !$OMP CRITICAL
1433                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1434                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1435                  !$OMP END CRITICAL
1436                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1437                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1438                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1439                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1440                ENDIF
1441              ENDIF
1442            ENDIF 
1443          ENDDO
1444       
1445        ENDDO
1446
1447      ELSE  IF (field(1)%ndim==4) THEN
1448
1449        max_req=0
1450        DO ind=1,ndomain
1451          req=>message%request(ind)
1452          IF (req%nsend>max_req) max_req=req%nsend
1453        ENDDO
1454   
1455        DO ind=1,ndomain
1456          IF (.NOT. assigned_domain(ind) ) CYCLE
1457
1458          dim3=size(field(ind)%rval4d,2)
1459          CALL distrib_level(1,dim3, lbegin,lend)
1460          dim4=size(field(ind)%rval4d,3)
1461          rval4d=>field(ind)%rval4d
1462          req=>message%request(ind)
1463
1464          DO isend=1,req%nsend
1465            send=>req%send(isend)
1466            value=>send%value
1467
1468            IF (send%rank/=mpi_rank) THEN
1469
1470              ireq=send%ireq
1471              buffer_r=>message%buffers(ireq)%r
1472              msize=send%size
1473              moffset=send%offset
1474
1475              CALL trace_start("copy_to_buffer")
1476              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice)
1477              DO d4=1,dim4
1478                DO d3=lbegin,lend
1479                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               &
1480                    (d3-1)*msize
1481                  !$acc loop
1482                  DO n=1,msize
1483                    buffer_r(n+offset)=rval4d(value(n),d3,d4)
1484                  ENDDO
1485                ENDDO
1486              ENDDO
1487              CALL trace_end("copy_to_buffer")
1488
1489              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1490                !$OMP BARRIER
1491              ENDIF
1492
1493              IF (is_omp_level_master) THEN
1494                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1495                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1496                  !$OMP CRITICAL
1497                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1498                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1499                  !$OMP END CRITICAL
1500                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1501                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1502                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1503                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1504                ENDIF
1505              ENDIF
1506            ELSE
1507              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1508                !$OMP BARRIER
1509              ENDIF
1510            ENDIF
1511          ENDDO
1512         
1513          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1514            DO isend=req%nsend+1,max_req
1515              !$OMP BARRIER
1516            ENDDO
1517          ENDIF
1518
1519        ENDDO
1520       
1521        DO ind=1,ndomain
1522          IF (.NOT. assigned_domain(ind) ) CYCLE
1523         
1524          dim3=size(field(ind)%rval4d,2)
1525          CALL distrib_level(1,dim3, lbegin,lend)
1526          dim4=size(field(ind)%rval4d,3)
1527          rval4d=>field(ind)%rval4d
1528          req=>message%request(ind)
1529
1530          DO irecv=1,req%nrecv
1531            recv=>req%recv(irecv)
1532            IF (recv%rank==mpi_rank) THEN
1533              value=>recv%value
1534              src_value => recv%src_value
1535              src_rval4d=>field(recv%domain)%rval4d
1536              sgn=>recv%sign
1537              msize=recv%size
1538              CALL trace_start("copy_data")
1539
1540              !$acc parallel loop collapse(3) default(present) async if (field(ind)%ondevice)
1541              DO d4=1,dim4
1542                DO d3=lbegin,lend
1543                  DO n=1,msize
1544                    rval4d(value(n),d3,d4)=src_rval4d(src_value(n),d3,d4)*sgn(n)
1545                  ENDDO
1546                ENDDO
1547              ENDDO
1548               
1549              CALL trace_end("copy_data")
1550                   
1551            ELSE
1552
1553              ireq=recv%ireq
1554              buffer_r=>message%buffers(ireq)%r
1555              IF (is_omp_level_master) THEN
1556                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1557                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1558                 !$OMP CRITICAL           
1559                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1560                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1561                  !$OMP END CRITICAL
1562                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1563                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1564                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1565                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1566                ENDIF
1567              ENDIF
1568            ENDIF
1569          ENDDO
1570        ENDDO
1571
1572      ENDIF     
1573
1574      IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
1575!$acc wait
1576!$OMP BARRIER
1577!$OMP MASTER       
1578
1579        DO ireq=1,message%nreq_send
1580          buffer_r=>message%buffers(ireq)%r
1581          msize=message%buffers(ireq)%size
1582          rank=message%buffers(ireq)%rank
1583          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation
1584          IF (message%ondevice) THEN
1585            !$acc host_data use_device(buffer_r) ! if (message%ondevice)
1586            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1587              comm_icosa, message%mpi_req(ireq),ierr)
1588            !$acc end host_data
1589          ELSE
1590            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1591              comm_icosa, message%mpi_req(ireq),ierr)
1592          ENDIF
1593        ENDDO
1594
1595        DO ireq=message%nreq_send+1,message%nreq_send+message%nreq_recv
1596          buffer_r=>message%buffers(ireq)%r
1597          msize=message%buffers(ireq)%size
1598          rank=message%buffers(ireq)%rank
1599          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation
1600          IF (message%ondevice) THEN
1601            !$acc host_data use_device(buffer_r) ! if (message%ondevice)
1602            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1603              comm_icosa, message%mpi_req(ireq),ierr)
1604            !$acc end host_data
1605          ELSE
1606            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1607              comm_icosa, message%mpi_req(ireq),ierr)
1608          ENDIF
1609        ENDDO
1610
1611!$OMP END MASTER
1612      ENDIF             
1613    ENDIF
1614   
1615!$OMP BARRIER
1616!    CALL trace_end("send_message_mpi")
1617
1618    CALL exit_profile(id_mpi)
1619   
1620  END SUBROUTINE send_message_mpi
1621 
1622  SUBROUTINE test_message_mpi(message)
1623  IMPLICIT NONE
1624    TYPE(t_message) :: message
1625   
1626    INTEGER :: ierr
1627
1628!$OMP MASTER
1629    IF (message%pending .AND. .NOT. message%completed) CALL MPI_TESTALL(message%nreq,&
1630      message%mpi_req,message%completed,message%status,ierr)
1631!$OMP END MASTER
1632  END SUBROUTINE  test_message_mpi
1633 
1634   
1635  SUBROUTINE wait_message_mpi(message)
1636  USE profiling_mod
1637  USE field_mod
1638  USE domain_mod
1639  USE mpi_mod
1640  USE mpipara
1641  USE omp_para
1642  USE trace
1643  IMPLICIT NONE
1644    TYPE(t_message) :: message
1645
1646    TYPE(t_field),POINTER :: field(:)
1647    REAL(rstd),POINTER :: rval2d(:) 
1648    REAL(rstd),POINTER :: rval3d(:,:) 
1649    REAL(rstd),POINTER :: rval4d(:,:,:) 
1650    REAL(rstd),POINTER :: buffer_r(:) 
1651    INTEGER,POINTER :: value(:) 
1652    INTEGER,POINTER :: sgn(:) 
1653    TYPE(ARRAY),POINTER :: recv 
1654    TYPE(t_request),POINTER :: req
1655    INTEGER :: irecv
1656    INTEGER :: ireq,nreq
1657    INTEGER :: ind,n
1658    INTEGER :: dim3,dim4,d3,d4,lbegin,lend
1659    INTEGER :: offset, msize, moffset
1660
1661    message%open=.FALSE.
1662    IF (.NOT. message%pending) RETURN
1663
1664    CALL enter_profile(id_mpi)
1665
1666!    CALL trace_start("wait_message_mpi")
1667
1668    field=>message%field
1669    nreq=message%nreq
1670   
1671    IF (field(1)%data_type==type_real) THEN
1672      IF (field(1)%ndim==2) THEN
1673
1674!$OMP MASTER
1675        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1676          message%status,ierr)
1677!$OMP END MASTER
1678!$OMP BARRIER
1679       
1680        DO ind=1,ndomain
1681          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1682         
1683          rval2d=>field(ind)%rval2d
1684          req=>message%request(ind)
1685          DO irecv=1,req%nrecv
1686            recv=>req%recv(irecv)
1687            IF (recv%rank/=mpi_rank) THEN
1688              ireq=recv%ireq
1689              buffer_r=>message%buffers(ireq)%r
1690              value=>recv%value
1691              sgn=>recv%sign
1692              offset=recv%offset
1693              msize=recv%size
1694              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1695              DO n=1,msize
1696                rval2d(value(n))=buffer_r(n+offset)*sgn(n) 
1697              ENDDO
1698
1699            ENDIF
1700          ENDDO
1701       
1702        ENDDO
1703     
1704     
1705      ELSE  IF (field(1)%ndim==3) THEN
1706
1707!$OMP MASTER
1708        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1709          message%status,ierr)
1710!$OMP END MASTER
1711!$OMP BARRIER
1712
1713       
1714        DO ind=1,ndomain
1715          IF (.NOT. assigned_domain(ind) ) CYCLE
1716
1717          rval3d=>field(ind)%rval3d
1718          req=>message%request(ind)
1719          DO irecv=1,req%nrecv
1720            recv=>req%recv(irecv)
1721            IF (recv%rank/=mpi_rank) THEN
1722              ireq=recv%ireq
1723              buffer_r=>message%buffers(ireq)%r
1724              value=>recv%value
1725              sgn=>recv%sign
1726             
1727              dim3=size(rval3d,2)
1728   
1729              CALL distrib_level(1,dim3, lbegin,lend)
1730              msize=recv%size
1731              moffset=recv%offset
1732              CALL trace_start("copy_from_buffer")
1733             
1734              IF (req%vector) THEN
1735                !$acc parallel loop default(present) async if (field(ind)%ondevice)
1736                DO d3=lbegin,lend
1737                  offset=moffset*dim3 + (d3-1)*msize
1738                  !$acc loop
1739                  DO n=1,msize
1740                    rval3d(value(n),d3)=buffer_r(n+offset)*sgn(n) 
1741                  ENDDO
1742                ENDDO
1743              ELSE
1744                !$acc parallel loop default(present) async if (field(ind)%ondevice)
1745                DO d3=lbegin,lend
1746                  offset=moffset*dim3 + (d3-1)*msize
1747                  !$acc loop
1748                  DO n=1,msize
1749                    rval3d(value(n),d3)=buffer_r(n+offset) 
1750                  ENDDO
1751                ENDDO
1752              ENDIF
1753               
1754              CALL trace_end("copy_from_buffer")
1755            ENDIF
1756          ENDDO
1757       
1758        ENDDO
1759
1760      ELSE  IF (field(1)%ndim==4) THEN
1761!$OMP MASTER
1762        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1763          message%status,ierr)
1764!$OMP END MASTER
1765!$OMP BARRIER
1766
1767               
1768        DO ind=1,ndomain
1769          IF (.NOT. assigned_domain(ind) ) CYCLE
1770
1771          rval4d=>field(ind)%rval4d
1772          req=>message%request(ind)
1773          DO irecv=1,req%nrecv
1774            recv=>req%recv(irecv)
1775            IF (recv%rank/=mpi_rank) THEN
1776              ireq=recv%ireq
1777              buffer_r=>message%buffers(ireq)%r
1778              value=>recv%value
1779              sgn=>recv%sign
1780
1781              dim3=size(rval4d,2)
1782              CALL distrib_level(1,dim3, lbegin,lend)
1783              dim4=size(rval4d,3)
1784              msize=recv%size
1785              moffset=recv%offset
1786              CALL trace_start("copy_from_buffer")
1787              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice)
1788              DO d4=1,dim4
1789                DO d3=lbegin,lend
1790                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               &
1791                    (d3-1)*msize
1792                  !$acc loop
1793                  DO n=1,msize
1794                    rval4d(value(n),d3,d4)=buffer_r(n+offset)*sgn(n) 
1795                  ENDDO
1796                ENDDO
1797              ENDDO
1798              CALL trace_end("copy_from_buffer")
1799            ENDIF
1800          ENDDO
1801       
1802        ENDDO
1803     
1804      ENDIF     
1805     
1806    ENDIF
1807
1808!$OMP MASTER
1809    message%pending=.FALSE.
1810!$OMP END MASTER
1811
1812!    CALL trace_end("wait_message_mpi")
1813!$OMP BARRIER
1814 
1815    CALL exit_profile(id_mpi)
1816
1817  END SUBROUTINE wait_message_mpi
1818
1819  SUBROUTINE transfert_request_mpi(field,request)
1820  USE field_mod
1821  IMPLICIT NONE
1822    TYPE(t_field),POINTER :: field(:)
1823    TYPE(t_request),POINTER :: request(:)
1824
1825    TYPE(t_message),SAVE :: message
1826   
1827   
1828    CALL init_message_mpi(field,request, message)
1829    CALL transfert_message_mpi(field,message)
1830    CALL finalize_message_mpi(field,message)
1831   
1832  END SUBROUTINE transfert_request_mpi
1833 
1834   
1835   
1836  SUBROUTINE transfert_request_seq(field,request)
1837  USE field_mod
1838  USE domain_mod
1839  IMPLICIT NONE
1840    TYPE(t_field),POINTER :: field(:)
1841    TYPE(t_request),POINTER :: request(:)
1842    REAL(rstd),POINTER :: rval2d(:) 
1843    REAL(rstd),POINTER :: rval3d(:,:) 
1844    REAL(rstd),POINTER :: rval4d(:,:,:) 
1845    INTEGER :: ind
1846    TYPE(t_request),POINTER :: req
1847    INTEGER :: n
1848   
1849    DO ind=1,ndomain
1850      req=>request(ind)
1851      rval2d=>field(ind)%rval2d
1852      rval3d=>field(ind)%rval3d
1853      rval4d=>field(ind)%rval4d
1854     
1855      IF (field(ind)%data_type==type_real) THEN
1856        IF (field(ind)%ndim==2) THEN
1857          DO n=1,req%size
1858            rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*&
1859              req%target_sign(n)
1860          ENDDO
1861        ELSE IF (field(ind)%ndim==3) THEN
1862          DO n=1,req%size
1863            rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*&
1864              req%target_sign(n)
1865          ENDDO
1866        ELSE IF (field(ind)%ndim==4) THEN
1867          DO n=1,req%size
1868            rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*&
1869              req%target_sign(n)
1870          ENDDO
1871        ENDIF
1872      ENDIF       
1873
1874    ENDDO
1875   
1876  END SUBROUTINE transfert_request_seq
1877 
1878 
1879  SUBROUTINE gather_field(field_loc,field_glo)
1880  USE field_mod
1881  USE domain_mod
1882  USE mpi_mod
1883  USE mpipara
1884  IMPLICIT NONE
1885    TYPE(t_field),POINTER :: field_loc(:)
1886    TYPE(t_field),POINTER :: field_glo(:)
1887    INTEGER, ALLOCATABLE :: mpi_req(:)
1888    INTEGER, ALLOCATABLE :: status(:,:)
1889    INTEGER :: ireq,nreq
1890    INTEGER :: ind_glo,ind_loc   
1891 
1892    IF (.NOT. using_mpi) THEN
1893   
1894      DO ind_loc=1,ndomain
1895        IF (field_loc(ind_loc)%ndim==2) field_glo(ind_loc)%rval2d=field_loc(ind_loc)%rval2d
1896        IF (field_loc(ind_loc)%ndim==3) field_glo(ind_loc)%rval3d=field_loc(ind_loc)%rval3d
1897        IF (field_loc(ind_loc)%ndim==4) field_glo(ind_loc)%rval4d=field_loc(ind_loc)%rval4d
1898      ENDDO
1899   
1900    ELSE
1901         
1902      nreq=ndomain
1903      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
1904      ALLOCATE(mpi_req(nreq))
1905      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
1906   
1907   
1908      ireq=0
1909      IF (mpi_rank==0) THEN
1910        DO ind_glo=1,ndomain_glo
1911          ireq=ireq+1
1912
1913          IF (field_glo(ind_glo)%ndim==2) THEN
1914            CALL MPI_IRECV(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
1915                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1916   
1917          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1918            CALL MPI_IRECV(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
1919                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1920
1921          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1922            CALL MPI_IRECV(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
1923                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1924          ENDIF
1925         
1926        ENDDO
1927      ENDIF
1928 
1929      DO ind_loc=1,ndomain
1930        ireq=ireq+1
1931
1932        IF (field_loc(ind_loc)%ndim==2) THEN
1933          CALL MPI_ISEND(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
1934                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1935        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
1936          CALL MPI_ISEND(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
1937                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1938        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
1939          CALL MPI_ISEND(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
1940                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1941        ENDIF
1942     
1943      ENDDO
1944   
1945      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
1946
1947    ENDIF
1948       
1949  END SUBROUTINE gather_field
1950
1951  SUBROUTINE bcast_field(field_glo)
1952  USE field_mod
1953  USE domain_mod
1954  USE mpi_mod
1955  USE mpipara
1956  IMPLICIT NONE
1957    TYPE(t_field),POINTER :: field_glo(:)
1958    INTEGER :: ind_glo   
1959 
1960    IF (.NOT. using_mpi) THEN
1961   
1962! nothing to do
1963   
1964    ELSE
1965         
1966      DO ind_glo=1,ndomain_glo
1967
1968          IF (field_glo(ind_glo)%ndim==2) THEN
1969            CALL MPI_BCAST(field_glo(ind_glo)%rval2d, size(field_glo(ind_glo)%rval2d) , MPI_REAL8, 0, comm_icosa, ierr)
1970          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1971            CALL MPI_BCAST(field_glo(ind_glo)%rval3d, size(field_glo(ind_glo)%rval3d) , MPI_REAL8, 0, comm_icosa, ierr)
1972          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1973            CALL MPI_BCAST(field_glo(ind_glo)%rval4d, size(field_glo(ind_glo)%rval4d) , MPI_REAL8, 0, comm_icosa, ierr)
1974          ENDIF
1975         
1976        ENDDO
1977      ENDIF
1978       
1979  END SUBROUTINE bcast_field
1980
1981  SUBROUTINE scatter_field(field_glo,field_loc)
1982  USE field_mod
1983  USE domain_mod
1984  USE mpi_mod
1985  USE mpipara
1986  IMPLICIT NONE
1987    TYPE(t_field),POINTER :: field_glo(:)
1988    TYPE(t_field),POINTER :: field_loc(:)
1989    INTEGER, ALLOCATABLE :: mpi_req(:)
1990    INTEGER, ALLOCATABLE :: status(:,:)
1991    INTEGER :: ireq,nreq
1992    INTEGER :: ind_glo,ind_loc   
1993 
1994    IF (.NOT. using_mpi) THEN
1995   
1996      DO ind_loc=1,ndomain
1997        IF (field_loc(ind_loc)%ndim==2) field_loc(ind_loc)%rval2d=field_glo(ind_loc)%rval2d
1998        IF (field_loc(ind_loc)%ndim==3) field_loc(ind_loc)%rval3d=field_glo(ind_loc)%rval3d
1999        IF (field_loc(ind_loc)%ndim==4) field_loc(ind_loc)%rval4d=field_glo(ind_loc)%rval4d
2000      ENDDO
2001   
2002    ELSE
2003         
2004      nreq=ndomain
2005      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
2006      ALLOCATE(mpi_req(nreq))
2007      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
2008   
2009   
2010      ireq=0
2011      IF (mpi_rank==0) THEN
2012        DO ind_glo=1,ndomain_glo
2013          ireq=ireq+1
2014
2015          IF (field_glo(ind_glo)%ndim==2) THEN
2016            CALL MPI_ISEND(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
2017                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
2018   
2019          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
2020            CALL MPI_ISEND(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
2021                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
2022
2023          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
2024            CALL MPI_ISEND(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
2025                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
2026          ENDIF
2027         
2028        ENDDO
2029      ENDIF
2030 
2031      DO ind_loc=1,ndomain
2032        ireq=ireq+1
2033
2034        IF (field_loc(ind_loc)%ndim==2) THEN
2035          CALL MPI_IRECV(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
2036                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
2037        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
2038          CALL MPI_IRECV(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
2039                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
2040        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
2041          CALL MPI_IRECV(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
2042                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
2043        ENDIF
2044     
2045      ENDDO
2046   
2047      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
2048
2049    ENDIF
2050       
2051  END SUBROUTINE scatter_field
2052 
2053  SUBROUTINE trace_in
2054  USE trace
2055  IMPLICIT NONE
2056 
2057    CALL trace_start("transfert_buffer")
2058  END SUBROUTINE trace_in             
2059
2060  SUBROUTINE trace_out
2061  USE trace
2062  IMPLICIT NONE
2063 
2064    CALL trace_end("transfert_buffer")
2065  END SUBROUTINE trace_out             
2066
2067
2068
2069
2070!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2071!! Definition des Broadcast --> 4D   !!
2072!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2073
2074!! -- Les chaine de charactï¿œre -- !!
2075
2076  SUBROUTINE bcast_mpi_c(var1)
2077  IMPLICIT NONE
2078    CHARACTER(LEN=*),INTENT(INOUT) :: Var1
2079   
2080    CALL bcast_mpi_cgen(Var1,len(Var1))
2081
2082  END SUBROUTINE bcast_mpi_c
2083
2084!! -- Les entiers -- !!
2085 
2086  SUBROUTINE bcast_mpi_i(var)
2087  USE mpipara
2088  IMPLICIT NONE
2089    INTEGER,INTENT(INOUT) :: Var
2090   
2091    INTEGER               :: var_tmp(1)
2092   
2093    IF (is_mpi_master) var_tmp(1)=var
2094    CALL bcast_mpi_igen(Var_tmp,1)
2095    var=var_tmp(1)
2096   
2097  END SUBROUTINE bcast_mpi_i
2098
2099  SUBROUTINE bcast_mpi_i1(var)
2100  IMPLICIT NONE
2101    INTEGER,INTENT(INOUT) :: Var(:)
2102
2103    CALL bcast_mpi_igen(Var,size(Var))
2104   
2105  END SUBROUTINE bcast_mpi_i1
2106
2107  SUBROUTINE bcast_mpi_i2(var)
2108  IMPLICIT NONE
2109    INTEGER,INTENT(INOUT) :: Var(:,:)
2110   
2111    CALL bcast_mpi_igen(Var,size(Var))
2112 
2113  END SUBROUTINE bcast_mpi_i2
2114
2115  SUBROUTINE bcast_mpi_i3(var)
2116  IMPLICIT NONE
2117    INTEGER,INTENT(INOUT) :: Var(:,:,:)
2118   
2119    CALL bcast_mpi_igen(Var,size(Var))
2120
2121  END SUBROUTINE bcast_mpi_i3
2122
2123  SUBROUTINE bcast_mpi_i4(var)
2124  IMPLICIT NONE
2125    INTEGER,INTENT(INOUT) :: Var(:,:,:,:)
2126   
2127    CALL bcast_mpi_igen(Var,size(Var))
2128
2129  END SUBROUTINE bcast_mpi_i4
2130
2131
2132!! -- Les reels -- !!
2133
2134  SUBROUTINE bcast_mpi_r(var)
2135  USE mpipara
2136  IMPLICIT NONE
2137    REAL,INTENT(INOUT) :: Var
2138    REAL               :: var_tmp(1)
2139   
2140    IF (is_mpi_master) var_tmp(1)=var
2141    CALL bcast_mpi_rgen(Var_tmp,1)
2142    var=var_tmp(1)   
2143
2144  END SUBROUTINE bcast_mpi_r
2145
2146  SUBROUTINE bcast_mpi_r1(var)
2147  IMPLICIT NONE
2148    REAL,INTENT(INOUT) :: Var(:)
2149   
2150    CALL bcast_mpi_rgen(Var,size(Var))
2151
2152  END SUBROUTINE bcast_mpi_r1
2153
2154  SUBROUTINE bcast_mpi_r2(var)
2155  IMPLICIT NONE
2156    REAL,INTENT(INOUT) :: Var(:,:)
2157   
2158    CALL bcast_mpi_rgen(Var,size(Var))
2159
2160  END SUBROUTINE bcast_mpi_r2
2161
2162  SUBROUTINE bcast_mpi_r3(var)
2163  IMPLICIT NONE
2164    REAL,INTENT(INOUT) :: Var(:,:,:)
2165   
2166    CALL bcast_mpi_rgen(Var,size(Var))
2167
2168  END SUBROUTINE bcast_mpi_r3
2169
2170  SUBROUTINE bcast_mpi_r4(var)
2171  IMPLICIT NONE
2172    REAL,INTENT(INOUT) :: Var(:,:,:,:)
2173   
2174    CALL bcast_mpi_rgen(Var,size(Var))
2175
2176  END SUBROUTINE bcast_mpi_r4
2177 
2178!! -- Les booleans -- !!
2179
2180  SUBROUTINE bcast_mpi_l(var)
2181  USE mpipara
2182  IMPLICIT NONE
2183    LOGICAL,INTENT(INOUT) :: Var
2184    LOGICAL               :: var_tmp(1)
2185   
2186    IF (is_mpi_master) var_tmp(1)=var
2187    CALL bcast_mpi_lgen(Var_tmp,1)
2188    var=var_tmp(1)   
2189
2190  END SUBROUTINE bcast_mpi_l
2191
2192  SUBROUTINE bcast_mpi_l1(var)
2193  IMPLICIT NONE
2194    LOGICAL,INTENT(INOUT) :: Var(:)
2195   
2196    CALL bcast_mpi_lgen(Var,size(Var))
2197
2198  END SUBROUTINE bcast_mpi_l1
2199
2200  SUBROUTINE bcast_mpi_l2(var)
2201  IMPLICIT NONE
2202    LOGICAL,INTENT(INOUT) :: Var(:,:)
2203   
2204    CALL bcast_mpi_lgen(Var,size(Var))
2205
2206  END SUBROUTINE bcast_mpi_l2
2207
2208  SUBROUTINE bcast_mpi_l3(var)
2209  IMPLICIT NONE
2210    LOGICAL,INTENT(INOUT) :: Var(:,:,:)
2211   
2212    CALL bcast_mpi_lgen(Var,size(Var))
2213
2214  END SUBROUTINE bcast_mpi_l3
2215
2216  SUBROUTINE bcast_mpi_l4(var)
2217  IMPLICIT NONE
2218    LOGICAL,INTENT(INOUT) :: Var(:,:,:,:)
2219   
2220    CALL bcast_mpi_lgen(Var,size(Var))
2221
2222  END SUBROUTINE bcast_mpi_l4
2223 
2224!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2225!! DEFINITION DES FONCTIONS DE TRANSFERT GENERIQUES !
2226!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2227
2228  SUBROUTINE bcast_mpi_cgen(var,nb)
2229    USE mpi_mod
2230    USE mpipara
2231    IMPLICIT NONE
2232   
2233    CHARACTER(LEN=*),INTENT(INOUT) :: Var
2234    INTEGER,INTENT(IN) :: nb
2235
2236    IF (.NOT. using_mpi) RETURN
2237   
2238    CALL MPI_BCAST(Var,nb,MPI_CHARACTER,mpi_master,comm_icosa,ierr)
2239       
2240  END SUBROUTINE bcast_mpi_cgen
2241
2242
2243     
2244  SUBROUTINE bcast_mpi_igen(var,nb)
2245    USE mpi_mod
2246    USE mpipara
2247    IMPLICIT NONE
2248    INTEGER,INTENT(IN) :: nb
2249    INTEGER,DIMENSION(nb),INTENT(INOUT) :: Var
2250   
2251    IF (.NOT. using_mpi) RETURN
2252
2253    CALL MPI_BCAST(Var,nb,MPI_INTEGER,mpi_master,comm_icosa,ierr)
2254       
2255  END SUBROUTINE bcast_mpi_igen
2256
2257
2258
2259 
2260  SUBROUTINE bcast_mpi_rgen(var,nb)
2261    USE mpi_mod
2262    USE mpipara
2263    IMPLICIT NONE
2264    INTEGER,INTENT(IN) :: nb
2265    REAL,DIMENSION(nb),INTENT(INOUT) :: Var
2266
2267    IF (.NOT. using_mpi) RETURN
2268
2269    CALL MPI_BCAST(Var,nb,MPI_REAL8,mpi_master,comm_icosa,ierr)
2270   
2271  END SUBROUTINE bcast_mpi_rgen
2272 
2273
2274
2275
2276  SUBROUTINE bcast_mpi_lgen(var,nb)
2277    USE mpi_mod
2278    USE mpipara
2279    IMPLICIT NONE
2280    INTEGER,INTENT(IN) :: nb
2281    LOGICAL,DIMENSION(nb),INTENT(INOUT) :: Var
2282
2283    IF (.NOT. using_mpi) RETURN
2284
2285    CALL MPI_BCAST(Var,nb,MPI_LOGICAL,mpi_master,comm_icosa,ierr)
2286
2287  END SUBROUTINE bcast_mpi_lgen
2288 
2289   
2290END MODULE transfert_mpi_mod
2291     
2292       
2293       
2294       
2295     
Note: See TracBrowser for help on using the repository browser.