source: codes/icosagcm/trunk/src/parallel/transfert_mpi_legacy.f90

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

Merge 'mpi_rewrite' into trunk

File size: 56.2 KB
Line 
1MODULE transfert_mpi_legacy_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  integer :: profile_mpi_copies, profile_mpi_waitall, profile_mpi_omp_barrier
83
84CONTAINS
85
86
87  SUBROUTINE init_transfert
88  USE profiling_mod
89  USE domain_mod
90  USE dimensions
91  USE field_mod
92  USE metric
93  USE mpipara
94  USE mpi_mod
95  IMPLICIT NONE
96  INTEGER :: ind,i,j
97
98    CALL register_id('MPI', id_mpi)
99    CALL register_id('MPI_copies', profile_mpi_copies)
100    CALL register_id('MPI_waitall', profile_mpi_waitall)
101    CALL register_id('MPI_omp_barrier', profile_mpi_omp_barrier)
102
103    CALL create_request(field_t,req_i1)
104
105    DO ind=1,ndomain
106      CALL swap_dimensions(ind)
107      DO i=ii_begin,ii_end+1
108        CALL request_add_point(ind,i,jj_begin-1,req_i1)
109      ENDDO
110
111      DO j=jj_begin,jj_end
112        CALL request_add_point(ind,ii_end+1,j,req_i1)
113      ENDDO
114      DO i=ii_begin,ii_end
115        CALL request_add_point(ind,i,jj_end+1,req_i1)
116      ENDDO
117
118      DO j=jj_begin,jj_end+1
119        CALL request_add_point(ind,ii_begin-1,j,req_i1)
120      ENDDO
121
122    ENDDO
123
124    CALL finalize_request(req_i1)
125
126
127    CALL create_request(field_t,req_i0)
128
129    DO ind=1,ndomain
130      CALL swap_dimensions(ind)
131
132      DO i=ii_begin,ii_end
133        CALL request_add_point(ind,i,jj_begin,req_i0)
134      ENDDO
135
136      DO j=jj_begin,jj_end
137        CALL request_add_point(ind,ii_end,j,req_i0)
138      ENDDO
139
140      DO i=ii_begin,ii_end
141        CALL request_add_point(ind,i,jj_end,req_i0)
142      ENDDO
143
144      DO j=jj_begin,jj_end
145        CALL request_add_point(ind,ii_begin,j,req_i0)
146      ENDDO
147
148    ENDDO
149
150    CALL finalize_request(req_i0)
151
152
153    CALL create_request(field_u,req_e1_scal)
154    DO ind=1,ndomain
155      CALL swap_dimensions(ind)
156      DO i=ii_begin,ii_end
157        CALL request_add_point(ind,i,jj_begin-1,req_e1_scal,rup)
158        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_scal,lup)
159      ENDDO
160
161      DO j=jj_begin,jj_end
162        CALL request_add_point(ind,ii_end+1,j,req_e1_scal,left)
163      ENDDO
164      DO j=jj_begin,jj_end
165        CALL request_add_point(ind,ii_end+1,j-1,req_e1_scal,lup)
166      ENDDO
167
168      DO i=ii_begin,ii_end
169        CALL request_add_point(ind,i,jj_end+1,req_e1_scal,ldown)
170        CALL request_add_point(ind,i-1,jj_end+1,req_e1_scal,rdown)
171      ENDDO
172
173      DO j=jj_begin,jj_end
174        CALL request_add_point(ind,ii_begin-1,j,req_e1_scal,right)
175      ENDDO
176      DO j=jj_begin,jj_end
177        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_scal,rdown)
178      ENDDO
179
180    ENDDO
181
182    CALL finalize_request(req_e1_scal)
183
184
185    CALL create_request(field_u,req_e0_scal)
186    DO ind=1,ndomain
187      CALL swap_dimensions(ind)
188
189
190      DO i=ii_begin+1,ii_end-1
191        CALL request_add_point(ind,i,jj_begin,req_e0_scal,right)
192        CALL request_add_point(ind,i,jj_end,req_e0_scal,right)
193      ENDDO
194
195      DO j=jj_begin+1,jj_end-1
196        CALL request_add_point(ind,ii_begin,j,req_e0_scal,rup)
197        CALL request_add_point(ind,ii_end,j,req_e0_scal,rup)
198      ENDDO
199
200      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_scal,left)
201      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_scal,ldown)
202      CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_scal,left)
203      CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_scal,ldown)
204
205    ENDDO
206
207    CALL finalize_request(req_e0_scal)
208
209
210
211    CALL create_request(field_u,req_e1_vect,.TRUE.)
212    DO ind=1,ndomain
213      CALL swap_dimensions(ind)
214      DO i=ii_begin,ii_end
215        CALL request_add_point(ind,i,jj_begin-1,req_e1_vect,rup)
216        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_vect,lup)
217      ENDDO
218
219      DO j=jj_begin,jj_end
220        CALL request_add_point(ind,ii_end+1,j,req_e1_vect,left)
221      ENDDO
222      DO j=jj_begin,jj_end
223        CALL request_add_point(ind,ii_end+1,j-1,req_e1_vect,lup)
224      ENDDO
225
226      DO i=ii_begin,ii_end
227        CALL request_add_point(ind,i,jj_end+1,req_e1_vect,ldown)
228        CALL request_add_point(ind,i-1,jj_end+1,req_e1_vect,rdown)
229      ENDDO
230
231      DO j=jj_begin,jj_end
232        CALL request_add_point(ind,ii_begin-1,j,req_e1_vect,right)
233      ENDDO
234      DO j=jj_begin,jj_end
235        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_vect,rdown)
236      ENDDO
237
238
239    ENDDO
240
241    CALL finalize_request(req_e1_vect)
242
243
244    CALL create_request(field_u,req_e0_vect,.TRUE.)
245    DO ind=1,ndomain
246      CALL swap_dimensions(ind)
247
248      DO i=ii_begin+1,ii_end-1
249        CALL request_add_point(ind,i,jj_begin,req_e0_vect,right)
250        CALL request_add_point(ind,i,jj_end,req_e0_vect,right)
251      ENDDO
252
253      DO j=jj_begin+1,jj_end-1
254        CALL request_add_point(ind,ii_begin,j,req_e0_vect,rup)
255        CALL request_add_point(ind,ii_end,j,req_e0_vect,rup)
256      ENDDO
257
258      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_vect,left)
259      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_vect,ldown)
260      CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_vect,left)
261      CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_vect,ldown)
262
263    ENDDO
264
265    CALL finalize_request(req_e0_vect)
266
267    CALL create_request(field_z,req_z1_scal)
268    DO ind=1,ndomain
269      CALL swap_dimensions(ind)
270      DO i=ii_begin,ii_end
271        CALL request_add_point(ind,i,jj_begin-1,req_z1_scal,vrup)
272        CALL request_add_point(ind,i+1,jj_begin-1,req_z1_scal,vup)
273      ENDDO
274
275      DO j=jj_begin,jj_end
276        CALL request_add_point(ind,ii_end+1,j,req_z1_scal,vlup)
277      ENDDO
278      DO j=jj_begin,jj_end
279        CALL request_add_point(ind,ii_end+1,j-1,req_z1_scal,vup)
280      ENDDO
281
282      DO i=ii_begin,ii_end
283        CALL request_add_point(ind,i,jj_end+1,req_z1_scal,vdown)
284        CALL request_add_point(ind,i-1,jj_end+1,req_z1_scal,vrdown)
285      ENDDO
286
287      DO j=jj_begin,jj_end
288        CALL request_add_point(ind,ii_begin-1,j,req_z1_scal,vrup)
289      ENDDO
290      DO j=jj_begin,jj_end
291        CALL request_add_point(ind,ii_begin-1,j,req_z1_scal,vrdown)
292      ENDDO
293
294    ENDDO
295
296    CALL finalize_request(req_z1_scal)
297
298  END SUBROUTINE init_transfert
299
300  SUBROUTINE create_request(type_field,request,vector)
301  USE domain_mod
302  USE field_mod
303  IMPLICIT NONE
304    INTEGER :: type_field
305    TYPE(t_request),POINTER :: request(:)
306    LOGICAL,OPTIONAL :: vector
307
308    TYPE(t_request),POINTER :: req
309    TYPE(t_domain),POINTER :: d
310    INTEGER :: ind
311    INTEGER :: max_size
312
313    ALLOCATE(request(ndomain))
314
315    DO ind=1,ndomain
316      req=>request(ind)
317      d=>domain(ind)
318      IF (type_field==field_t) THEN
319        Max_size=2*(d%iim+2)+2*(d%jjm+2)
320      ELSE IF (type_field==field_u) THEN
321        Max_size=3*(2*(d%iim+2)+2*(d%jjm+2))
322      ELSE IF (type_field==field_z) THEN
323        Max_size=2*(2*(d%iim+2)+2*(d%jjm+2))
324      ENDIF
325
326      req%type_field=type_field
327      req%max_size=max_size*2
328      req%size=0
329      req%vector=.FALSE.
330      IF (PRESENT(vector)) req%vector=vector
331      ALLOCATE(req%src_domain(req%max_size))
332      ALLOCATE(req%src_ind(req%max_size))
333      ALLOCATE(req%target_ind(req%max_size))
334      ALLOCATE(req%src_i(req%max_size))
335      ALLOCATE(req%src_j(req%max_size))
336      ALLOCATE(req%target_i(req%max_size))
337      ALLOCATE(req%target_j(req%max_size))
338      ALLOCATE(req%target_sign(req%max_size))
339    ENDDO
340
341  END SUBROUTINE create_request
342
343  SUBROUTINE reallocate_request(req)
344  IMPLICIT NONE
345    TYPE(t_request),POINTER :: req
346
347    INTEGER,POINTER :: src_domain(:)
348    INTEGER,POINTER :: src_ind(:)
349    INTEGER,POINTER :: target_ind(:)
350    INTEGER,POINTER :: src_i(:)
351    INTEGER,POINTER :: src_j(:)
352    INTEGER,POINTER :: target_i(:)
353    INTEGER,POINTER :: target_j(:)
354    INTEGER,POINTER :: target_sign(:)
355
356    PRINT *,"REALLOCATE_REQUEST"
357    src_domain=>req%src_domain
358    src_ind=>req%src_ind
359    target_ind=>req%target_ind
360    src_i=>req%src_i
361    src_j=>req%src_j
362    target_i=>req%target_i
363    target_j=>req%target_j
364    target_sign=>req%target_sign
365
366    ALLOCATE(req%src_domain(req%max_size*2))
367    ALLOCATE(req%src_ind(req%max_size*2))
368    ALLOCATE(req%target_ind(req%max_size*2))
369    ALLOCATE(req%src_i(req%max_size*2))
370    ALLOCATE(req%src_j(req%max_size*2))
371    ALLOCATE(req%target_i(req%max_size*2))
372    ALLOCATE(req%target_j(req%max_size*2))
373    ALLOCATE(req%target_sign(req%max_size*2))
374
375    req%src_domain(1:req%max_size)=src_domain(:)
376    req%src_ind(1:req%max_size)=src_ind(:)
377    req%target_ind(1:req%max_size)=target_ind(:)
378    req%src_i(1:req%max_size)=src_i(:)
379    req%src_j(1:req%max_size)=src_j(:)
380    req%target_i(1:req%max_size)=target_i(:)
381    req%target_j(1:req%max_size)=target_j(:)
382    req%target_sign(1:req%max_size)=target_sign(:)
383
384    req%max_size=req%max_size*2
385
386    DEALLOCATE(src_domain)
387    DEALLOCATE(src_ind)
388    DEALLOCATE(target_ind)
389    DEALLOCATE(src_i)
390    DEALLOCATE(src_j)
391    DEALLOCATE(target_i)
392    DEALLOCATE(target_j)
393    DEALLOCATE(target_sign)
394
395  END SUBROUTINE reallocate_request
396
397
398    SUBROUTINE request_add_point(ind,i,j,request,pos)
399    USE domain_mod
400    USE field_mod
401    IMPLICIT NONE
402      INTEGER,INTENT(IN)            ::  ind
403      INTEGER,INTENT(IN)            :: i
404      INTEGER,INTENT(IN)            :: j
405      TYPE(t_request),POINTER :: request(:)
406      INTEGER,INTENT(IN),OPTIONAL  :: pos
407
408      INTEGER :: src_domain
409      INTEGER :: src_iim,src_i,src_j,src_n,src_pos,src_delta
410      TYPE(t_request),POINTER :: req
411      TYPE(t_domain),POINTER :: d
412
413      req=>request(ind)
414      d=>domain(ind)
415
416      IF (req%max_size==req%size) CALL reallocate_request(req)
417      req%size=req%size+1
418      IF (req%type_field==field_t) THEN
419        src_domain=domain(ind)%assign_domain(i,j)
420        src_iim=domain_glo(src_domain)%iim
421        src_i=domain(ind)%assign_i(i,j)
422        src_j=domain(ind)%assign_j(i,j)
423
424        req%target_ind(req%size)=(j-1)*d%iim+i
425        req%target_sign(req%size)=1
426        req%src_domain(req%size)=src_domain
427        req%src_ind(req%size)=(src_j-1)*src_iim+src_i
428      ELSE IF (req%type_field==field_u) THEN
429        IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme'
430
431        src_domain=domain(ind)%edge_assign_domain(pos-1,i,j)
432        src_iim=domain_glo(src_domain)%iim
433        src_i=domain(ind)%edge_assign_i(pos-1,i,j)
434        src_j=domain(ind)%edge_assign_j(pos-1,i,j)
435        src_n=(src_j-1)*src_iim+src_i
436        src_delta=domain(ind)%delta(i,j)
437        src_pos=domain(ind)%edge_assign_pos(pos-1,i,j)+1
438
439        req%target_ind(req%size)=(j-1)*d%iim+i+d%u_pos(pos)
440
441        req%target_sign(req%size)= 1
442        IF (req%vector) req%target_sign(req%size)= domain(ind)%edge_assign_sign(pos-1,i,j)
443
444        req%src_domain(req%size)=src_domain
445        req%src_ind(req%size)=src_n+domain_glo(src_domain)%u_pos(src_pos)
446
447      ELSE IF (req%type_field==field_z) THEN
448        IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme'
449
450        src_domain=domain(ind)%vertex_assign_domain(pos-1,i,j)
451        src_iim=domain_glo(src_domain)%iim
452        src_i=domain(ind)%vertex_assign_i(pos-1,i,j)
453        src_j=domain(ind)%vertex_assign_j(pos-1,i,j)
454        src_n=(src_j-1)*src_iim+src_i
455        src_delta=domain(ind)%delta(i,j)
456        src_pos=domain(ind)%vertex_assign_pos(pos-1,i,j)+1
457
458
459        req%target_ind(req%size)=(j-1)*d%iim+i+d%z_pos(pos)
460        req%target_sign(req%size)=1
461        req%src_domain(req%size)=src_domain
462        req%src_ind(req%size)=src_n+domain_glo(src_domain)%z_pos(src_pos)
463      ENDIF
464  END SUBROUTINE request_add_point
465
466
467  SUBROUTINE Finalize_request(request)
468  USE mpipara
469  USE domain_mod
470  USE mpi_mod
471  IMPLICIT NONE
472    TYPE(t_request),POINTER :: request(:)
473    TYPE(t_request),POINTER :: req, req_src
474    INTEGER :: nb_domain_recv(0:mpi_size-1)
475    INTEGER :: nb_domain_send(0:mpi_size-1)
476    INTEGER :: tag_rank(0:mpi_size-1)
477    INTEGER :: nb_data_domain_recv(ndomain_glo)
478    INTEGER :: list_domain_recv(ndomain_glo)
479    INTEGER,ALLOCATABLE :: list_domain_send(:)
480    INTEGER             :: list_domain(ndomain)
481
482    INTEGER :: rank,i,j,pos
483    INTEGER :: size_,ind_glo,ind_loc
484    INTEGER :: isend, irecv, ireq, nreq, nsend, nrecv
485    INTEGER, ALLOCATABLE :: mpi_req(:)
486    INTEGER, ALLOCATABLE :: status(:,:)
487    INTEGER, ALLOCATABLE :: rank_list(:)
488    INTEGER, ALLOCATABLE :: offset(:)
489    LOGICAL,PARAMETER :: debug = .FALSE.
490
491
492    IF (.NOT. using_mpi) RETURN
493
494    DO ind_loc=1,ndomain
495      req=>request(ind_loc)
496
497      nb_data_domain_recv(:) = 0
498      nb_domain_recv(:) = 0
499      tag_rank(:)=0
500
501      DO i=1,req%size
502        ind_glo=req%src_domain(i)
503        nb_data_domain_recv(ind_glo)=nb_data_domain_recv(ind_glo)+1
504      ENDDO
505
506      DO ind_glo=1,ndomain_glo
507        IF ( nb_data_domain_recv(ind_glo) > 0 )  nb_domain_recv(domglo_rank(ind_glo))=nb_domain_recv(domglo_rank(ind_glo))+1
508      ENDDO
509
510      req%nrecv=sum(nb_domain_recv(:))
511      ALLOCATE(req%recv(req%nrecv))
512
513      irecv=0
514      DO ind_glo=1,ndomain_glo
515        IF (nb_data_domain_recv(ind_glo)>0) THEN
516          irecv=irecv+1
517          list_domain_recv(ind_glo)=irecv
518          req%recv(irecv)%rank=domglo_rank(ind_glo)
519          req%recv(irecv)%size=nb_data_domain_recv(ind_glo)
520          req%recv(irecv)%domain=domglo_loc_ind(ind_glo)
521          ALLOCATE(req%recv(irecv)%value(req%recv(irecv)%size))
522          ALLOCATE(req%recv(irecv)%sign(req%recv(irecv)%size))
523          ALLOCATE(req%recv(irecv)%buffer(req%recv(irecv)%size))
524        ENDIF
525      ENDDO
526
527      req%recv(:)%size=0
528      irecv=0
529      DO i=1,req%size
530        irecv=list_domain_recv(req%src_domain(i))
531        req%recv(irecv)%size=req%recv(irecv)%size+1
532        size_=req%recv(irecv)%size
533        req%recv(irecv)%value(size_)=req%src_ind(i)
534        req%recv(irecv)%buffer(size_)=req%target_ind(i)
535        req%recv(irecv)%sign(size_)=req%target_sign(i)
536      ENDDO
537    ENDDO
538
539    nb_domain_recv(:) = 0
540    DO ind_loc=1,ndomain
541      req=>request(ind_loc)
542
543      DO irecv=1,req%nrecv
544        rank=req%recv(irecv)%rank
545        nb_domain_recv(rank)=nb_domain_recv(rank)+1
546      ENDDO
547    ENDDO
548
549    CALL MPI_Alltoall(nb_domain_recv,1,MPI_INTEGER,nb_domain_send,1,MPI_INTEGER,comm_icosa,ierr)
550
551
552    ALLOCATE(list_domain_send(sum(nb_domain_send)))
553
554    nreq=sum(nb_domain_recv(:))+sum(nb_domain_send(:))
555    ALLOCATE(mpi_req(nreq))
556    ALLOCATE(status(MPI_STATUS_SIZE,nreq))
557
558
559    ireq=0
560    DO ind_loc=1,ndomain
561      req=>request(ind_loc)
562      DO irecv=1,req%nrecv
563        ireq=ireq+1
564        CALL MPI_ISEND(req%recv(irecv)%domain,1,MPI_INTEGER,req%recv(irecv)%rank,0,comm_icosa, mpi_req(ireq),ierr)
565        IF (debug) PRINT *,"Isend ",req%recv(irecv)%domain, "from ", mpi_rank, "to ",req%recv(irecv)%rank, "tag ",0
566      ENDDO
567    ENDDO
568
569    IF (debug) PRINT *,"------------"
570    j=0
571    DO rank=0,mpi_size-1
572      DO i=1,nb_domain_send(rank)
573        j=j+1
574        ireq=ireq+1
575        CALL MPI_IRECV(list_domain_send(j),1,MPI_INTEGER,rank,0,comm_icosa, mpi_req(ireq),ierr)
576        IF (debug) PRINT *,"IRecv ",list_domain_send(j), "from ", rank, "to ",mpi_rank, "tag ",0
577      ENDDO
578    ENDDO
579    IF (debug) PRINT *,"------------"
580
581    CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
582
583    list_domain(:)=0
584    DO i=1,sum(nb_domain_send)
585      ind_loc=list_domain_send(i)
586      list_domain(ind_loc)=list_domain(ind_loc)+1
587    ENDDO
588
589    DO ind_loc=1,ndomain
590      req=>request(ind_loc)
591      req%nsend=list_domain(ind_loc)
592      ALLOCATE(req%send(req%nsend))
593    ENDDO
594
595    IF (debug) PRINT *,"------------"
596
597   ireq=0
598   DO ind_loc=1,ndomain
599     req=>request(ind_loc)
600
601     DO irecv=1,req%nrecv
602       ireq=ireq+1
603       CALL MPI_ISEND(mpi_rank,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
604       IF (debug) PRINT *,"Isend ",mpi_rank, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
605     ENDDO
606    IF (debug) PRINT *,"------------"
607
608     DO isend=1,req%nsend
609       ireq=ireq+1
610       CALL MPI_IRECV(req%send(isend)%rank,1,MPI_INTEGER,MPI_ANY_SOURCE,ind_loc,comm_icosa, mpi_req(ireq),ierr)
611       IF (debug) PRINT *,"IRecv ",req%send(isend)%rank, "from ", "xxx", "to ",mpi_rank, "tag ",ind_loc
612     ENDDO
613   ENDDO
614
615   IF (debug) PRINT *,"------------"
616
617   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
618   CALL MPI_BARRIER(comm_icosa,ierr)
619
620   IF (debug) PRINT *,"------------"
621
622   ireq=0
623   DO ind_loc=1,ndomain
624     req=>request(ind_loc)
625
626     DO irecv=1,req%nrecv
627       ireq=ireq+1
628       CALL MPI_ISEND(ind_loc,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
629       IF (debug) PRINT *,"Isend ",ind_loc, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
630     ENDDO
631
632     IF (debug) PRINT *,"------------"
633
634     DO isend=1,req%nsend
635       ireq=ireq+1
636       CALL MPI_IRECV(req%send(isend)%domain,1,MPI_INTEGER,req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr)
637       IF (debug) PRINT *,"IRecv ",req%send(isend)%domain, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
638     ENDDO
639   ENDDO
640   IF (debug) PRINT *,"------------"
641
642   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
643   CALL MPI_BARRIER(comm_icosa,ierr)
644   IF (debug) PRINT *,"------------"
645
646   ireq=0
647   DO ind_loc=1,ndomain
648     req=>request(ind_loc)
649
650     DO irecv=1,req%nrecv
651       ireq=ireq+1
652       req%recv(irecv)%tag=tag_rank(req%recv(irecv)%rank)
653       tag_rank(req%recv(irecv)%rank)=tag_rank(req%recv(irecv)%rank)+1
654       CALL MPI_ISEND(req%recv(irecv)%tag,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
655       IF (debug) PRINT *,"Isend ",req%recv(irecv)%tag, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
656     ENDDO
657   IF (debug) PRINT *,"------------"
658
659     DO isend=1,req%nsend
660       ireq=ireq+1
661       CALL MPI_IRECV(req%send(isend)%tag,1,MPI_INTEGER,req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr)
662       IF (debug) PRINT *,"IRecv ",req%send(isend)%tag, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
663     ENDDO
664   ENDDO
665   IF (debug) PRINT *,"------------"
666
667   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
668   CALL MPI_BARRIER(comm_icosa,ierr)
669
670
671   IF (debug) PRINT *,"------------"
672
673   ireq=0
674   DO ind_loc=1,ndomain
675     req=>request(ind_loc)
676
677     DO irecv=1,req%nrecv
678       ireq=ireq+1
679       CALL MPI_ISEND(req%recv(irecv)%size,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
680       IF (debug) PRINT *,"Isend ",req%recv(irecv)%size, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
681     ENDDO
682     IF (debug) PRINT *,"------------"
683
684     DO isend=1,req%nsend
685       ireq=ireq+1
686       CALL MPI_IRECV(req%send(isend)%size,1,MPI_INTEGER,req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
687       IF (debug) PRINT *,"IRecv ",req%send(isend)%size, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
688     ENDDO
689   ENDDO
690
691   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
692
693   ireq=0
694   DO ind_loc=1,ndomain
695     req=>request(ind_loc)
696
697     DO irecv=1,req%nrecv
698       ireq=ireq+1
699       CALL MPI_ISEND(req%recv(irecv)%value,req%recv(irecv)%size,MPI_INTEGER,&
700            req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
701     ENDDO
702
703     DO isend=1,req%nsend
704       ireq=ireq+1
705       ALLOCATE(req%send(isend)%value(req%send(isend)%size))
706       CALL MPI_IRECV(req%send(isend)%value,req%send(isend)%size,MPI_INTEGER,&
707            req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
708     ENDDO
709   ENDDO
710
711   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
712
713   DO ind_loc=1,ndomain
714     req=>request(ind_loc)
715
716     DO irecv=1,req%nrecv
717       req%recv(irecv)%value(:)=req%recv(irecv)%buffer(:)
718       req%recv(irecv)%sign(:) =req%recv(irecv)%sign(:)
719       DEALLOCATE(req%recv(irecv)%buffer)
720     ENDDO
721   ENDDO
722
723
724! domain is on the same mpi process => copie memory to memory
725
726   DO ind_loc=1,ndomain
727     req=>request(ind_loc)
728
729     DO irecv=1,req%nrecv
730
731       IF (req%recv(irecv)%rank==mpi_rank) THEN
732           req_src=>request(req%recv(irecv)%domain)
733           DO isend=1,req_src%nsend
734             IF (req_src%send(isend)%rank==mpi_rank .AND. req_src%send(isend)%tag==req%recv(irecv)%tag) THEN
735               req%recv(irecv)%src_value => req_src%send(isend)%value
736               IF ( size(req%recv(irecv)%value) /= size(req_src%send(isend)%value)) THEN
737                 PRINT *,ind_loc, irecv, isend, size(req%recv(irecv)%value), size(req_src%send(isend)%value)
738                 STOP "size(req%recv(irecv)%value) /= size(req_src%send(isend)%value"
739               ENDIF
740             ENDIF
741           ENDDO
742       ENDIF
743
744     ENDDO
745   ENDDO
746
747! true number of mpi request
748
749   request(:)%nreq_mpi=0
750   request(:)%nreq_send=0
751   request(:)%nreq_recv=0
752   ALLOCATE(rank_list(sum(request(:)%nsend)))
753   ALLOCATE(offset(sum(request(:)%nsend)))
754   offset(:)=0
755
756   nsend=0
757   DO ind_loc=1,ndomain
758     req=>request(ind_loc)
759
760     DO isend=1,req%nsend
761       IF (req%send(isend)%rank/=mpi_rank) THEN
762         pos=0
763         DO i=1,nsend
764           IF (req%send(isend)%rank==rank_list(i)) EXIT
765           pos=pos+1
766         ENDDO
767
768         IF (pos==nsend) THEN
769           nsend=nsend+1
770           req%nreq_mpi=req%nreq_mpi+1
771           req%nreq_send=req%nreq_send+1
772           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
773             rank_list(nsend)=req%send(isend)%rank
774           ELSE
775             rank_list(nsend)=-1
776           ENDIF
777         ENDIF
778
779         pos=pos+1
780         req%send(isend)%ireq=pos
781         req%send(isend)%offset=offset(pos)
782         offset(pos)=offset(pos)+req%send(isend)%size
783       ENDIF
784     ENDDO
785   ENDDO
786
787   DEALLOCATE(rank_list)
788   DEALLOCATE(offset)
789
790   ALLOCATE(rank_list(sum(request(:)%nrecv)))
791   ALLOCATE(offset(sum(request(:)%nrecv)))
792   offset(:)=0
793
794   nrecv=0
795   DO ind_loc=1,ndomain
796     req=>request(ind_loc)
797
798     DO irecv=1,req%nrecv
799       IF (req%recv(irecv)%rank/=mpi_rank) THEN
800         pos=0
801         DO i=1,nrecv
802           IF (req%recv(irecv)%rank==rank_list(i)) EXIT
803           pos=pos+1
804         ENDDO
805
806         IF (pos==nrecv) THEN
807           nrecv=nrecv+1
808           req%nreq_mpi=req%nreq_mpi+1
809           req%nreq_recv=req%nreq_recv+1
810           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
811             rank_list(nrecv)=req%recv(irecv)%rank
812           ELSE
813             rank_list(nrecv)=-1
814           ENDIF
815         ENDIF
816
817         pos=pos+1
818         req%recv(irecv)%ireq=nsend+pos
819         req%recv(irecv)%offset=offset(pos)
820         offset(pos)=offset(pos)+req%recv(irecv)%size
821       ENDIF
822     ENDDO
823   ENDDO
824
825! get the offsets
826
827   ireq=0
828   DO ind_loc=1,ndomain
829     req=>request(ind_loc)
830
831     DO irecv=1,req%nrecv
832       ireq=ireq+1
833       CALL MPI_ISEND(req%recv(irecv)%offset,1,MPI_INTEGER,&
834            req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
835     ENDDO
836
837     DO isend=1,req%nsend
838       ireq=ireq+1
839       CALL MPI_IRECV(req%send(isend)%offset,1,MPI_INTEGER,&
840            req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
841     ENDDO
842   ENDDO
843
844   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
845
846
847  END SUBROUTINE Finalize_request
848
849
850  SUBROUTINE init_message_seq(field, request, message, name)
851  USE field_mod
852  USE domain_mod
853  USE mpi_mod
854  USE mpipara
855  USE mpi_mod
856  IMPLICIT NONE
857    TYPE(t_field),POINTER :: field(:)
858    TYPE(t_request),POINTER :: request(:)
859    TYPE(t_message) :: message
860    CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: name
861!$OMP MASTER
862    message%request=>request
863    IF(PRESENT(name)) THEN
864       message%name = TRIM(name)
865    ELSE
866       message%name = 'unknown'
867    END IF
868!$OMP END MASTER
869!$OMP BARRIER
870
871  END SUBROUTINE init_message_seq
872
873  SUBROUTINE send_message_seq(field,message)
874  USE field_mod
875  USE domain_mod
876  USE mpi_mod
877  USE mpipara
878  USE omp_para
879  USE trace
880  IMPLICIT NONE
881    TYPE(t_field),POINTER :: field(:)
882    TYPE(t_message) :: message
883
884    CALL transfert_request_seq(field,message%request)
885
886  END SUBROUTINE send_message_seq
887
888  SUBROUTINE test_message_seq(message)
889  IMPLICIT NONE
890    TYPE(t_message) :: message
891  END SUBROUTINE  test_message_seq
892
893
894  SUBROUTINE wait_message_seq(message)
895  IMPLICIT NONE
896    TYPE(t_message) :: message
897
898  END SUBROUTINE wait_message_seq
899
900  SUBROUTINE init_message_mpi(field,request, message, name)
901  USE field_mod
902  USE domain_mod
903  USE mpi_mod
904  USE mpipara
905  USE mpi_mod
906  IMPLICIT NONE
907
908    TYPE(t_field),POINTER :: field(:)
909    TYPE(t_request),POINTER :: request(:)
910    TYPE(t_message) :: message
911    CHARACTER(LEN=*), INTENT(IN),OPTIONAL :: name
912
913    TYPE(t_request),POINTER :: req
914    INTEGER :: irecv,isend
915    INTEGER :: ireq,nreq
916    INTEGER :: ind
917    INTEGER :: dim3,dim4
918    INTEGER,SAVE :: message_number=0
919!    TYPE(t_reorder),POINTER :: reorder(:)
920!    TYPE(t_reorder) :: reorder_swap
921
922!$OMP BARRIER
923!$OMP MASTER
924    !init off-device
925    message%ondevice=.false.
926
927    IF(PRESENT(name)) THEN
928       message%name = TRIM(name)
929    ELSE
930       message%name = 'unknown'
931    END IF
932    message%number=message_number
933    message_number=message_number+1
934    IF (message_number==100) message_number=0
935
936
937    message%request=>request
938    message%nreq=sum(message%request(:)%nreq_mpi)
939    message%nreq_send=sum(message%request(:)%nreq_send)
940    message%nreq_recv=sum(message%request(:)%nreq_recv)
941    nreq=message%nreq
942
943    ALLOCATE(message%mpi_req(nreq))
944    ALLOCATE(message%buffers(nreq))
945    ALLOCATE(message%status(MPI_STATUS_SIZE,nreq))
946    message%buffers(:)%size=0
947    message%pending=.FALSE.
948    message%completed=.FALSE.
949    message%open=.FALSE.
950
951    DO ind=1,ndomain
952      req=>request(ind)
953      DO isend=1,req%nsend
954        IF (req%send(isend)%rank/=mpi_rank) THEN
955          ireq=req%send(isend)%ireq
956          message%buffers(ireq)%size=message%buffers(ireq)%size+req%send(isend)%size
957          message%buffers(ireq)%rank=req%send(isend)%rank
958        ENDIF
959      ENDDO
960      DO irecv=1,req%nrecv
961        IF (req%recv(irecv)%rank/=mpi_rank) THEN
962          ireq=req%recv(irecv)%ireq
963          message%buffers(ireq)%size=message%buffers(ireq)%size+req%recv(irecv)%size
964          message%buffers(ireq)%rank=req%recv(irecv)%rank
965        ENDIF
966      ENDDO
967    ENDDO
968
969
970    IF (field(1)%data_type==type_real) THEN
971
972      IF (field(1)%ndim==2) THEN
973
974        DO ireq=1,message%nreq
975          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
976        ENDDO
977
978      ELSE  IF (field(1)%ndim==3) THEN
979
980        dim3=size(field(1)%rval3d,2)
981        DO ireq=1,message%nreq
982          message%buffers(ireq)%size=message%buffers(ireq)%size*dim3
983          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
984        ENDDO
985
986      ELSE  IF (field(1)%ndim==4) THEN
987        dim3=size(field(1)%rval4d,2)
988        dim4=size(field(1)%rval4d,3)
989        DO ireq=1,message%nreq
990          message%buffers(ireq)%size=message%buffers(ireq)%size*dim3*dim4
991          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
992        ENDDO
993      ENDIF
994    ENDIF
995
996
997
998! ! Reorder the request, so recv request are done in the same order than send request
999
1000!    nreq_send=sum(request(:)%nsend)
1001!    message%nreq_send=nreq_send
1002!    ALLOCATE(message%reorder(nreq_send))
1003!    reorder=>message%reorder
1004!    ireq=0
1005!    DO ind=1,ndomain
1006!      req=>request(ind)
1007!      DO isend=1,req%nsend
1008!        ireq=ireq+1
1009!        reorder(ireq)%ind=ind
1010!        reorder(ireq)%isend=isend
1011!        reorder(ireq)%tag=req%send(isend)%tag
1012!      ENDDO
1013!    ENDDO
1014
1015! ! do a very very bad sort
1016!    DO i=1,nreq_send-1
1017!      DO j=i+1,nreq_send
1018!        IF (reorder(i)%tag > reorder(j)%tag) THEN
1019!          reorder_swap=reorder(i)
1020!          reorder(i)=reorder(j)
1021!          reorder(j)=reorder_swap
1022!        ENDIF
1023!      ENDDO
1024!    ENDDO
1025!    PRINT *,"reorder ",reorder(:)%tag
1026
1027
1028!$OMP END MASTER
1029!$OMP BARRIER
1030
1031  END SUBROUTINE init_message_mpi
1032
1033  SUBROUTINE Finalize_message_mpi(message)
1034  USE field_mod
1035  USE domain_mod
1036  USE mpi_mod
1037  USE mpipara
1038  USE mpi_mod
1039  IMPLICIT NONE
1040    TYPE(t_message) :: message
1041
1042    INTEGER :: ireq, ibuff
1043
1044!$OMP BARRIER
1045!$OMP MASTER
1046
1047
1048    IF (message%field(1)%data_type==type_real) THEN
1049      DO ireq=1,message%nreq
1050        CALL free_mpi_buffer(message%buffers(ireq)%r)
1051      ENDDO
1052    ENDIF
1053
1054    !deallocate device data if ondevice
1055    if(message%ondevice) then
1056      do ireq=1, ndomain
1057        do ibuff=1,message%request(ireq)%nsend
1058          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer(:))
1059          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer_r(:))
1060          !$acc exit data delete(message%request(ireq)%send(ibuff)%sign(:))
1061          !$acc exit data delete(message%request(ireq)%send(ibuff)%src_value(:))
1062          !$acc exit data delete(message%request(ireq)%send(ibuff)%value(:))
1063        end do
1064        do ibuff=1,message%request(ireq)%nrecv
1065          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer(:))
1066          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer_r(:))
1067          !$acc exit data delete(message%request(ireq)%recv(ibuff)%sign(:))
1068          !$acc exit data delete(message%request(ireq)%recv(ibuff)%src_value(:))
1069          !$acc exit data delete(message%request(ireq)%recv(ibuff)%value(:))
1070        end do
1071      end do
1072      DO ireq=1,message%nreq
1073        !$acc exit data delete(message%buffers(ireq)%r)
1074      ENDDO
1075      message%ondevice=.false.
1076    end if
1077
1078    DEALLOCATE(message%mpi_req)
1079    DEALLOCATE(message%buffers)
1080    DEALLOCATE(message%status)
1081
1082!$OMP END MASTER
1083!$OMP BARRIER
1084
1085
1086  END SUBROUTINE Finalize_message_mpi
1087
1088
1089  !!!Update buffers on device for 'message'
1090  !!! does create_device_message when not already ondevice
1091  SUBROUTINE update_device_message_mpi(message)
1092    USE domain_mod
1093    IMPLICIT NONE
1094    TYPE(t_message), intent(inout) :: message
1095    INTEGER :: ireq, ibuff
1096
1097    !if(.not. message%ondevice) call create_device_message_mpi(message)
1098
1099    do ireq=1, ndomain
1100      do ibuff=1,message%request(ireq)%nsend
1101        !device buffers updated even if pointers not attached :
1102        !non allocated buffers in 'message' must be set to NULL()
1103        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer(:)) async
1104        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer_r(:)) async
1105        !$acc enter data copyin(message%request(ireq)%send(ibuff)%sign(:)) async
1106        !$acc enter data copyin(message%request(ireq)%send(ibuff)%src_value(:)) async
1107        !$acc enter data copyin(message%request(ireq)%send(ibuff)%value(:)) async
1108      end do
1109      do ibuff=1,message%request(ireq)%nrecv
1110        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer(:)) async
1111        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer_r(:)) async
1112        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%sign(:)) async
1113        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%src_value(:)) async
1114        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%value(:)) async
1115      end do
1116    end do
1117    DO ireq=1,message%nreq
1118      !$acc enter data copyin(message%buffers(ireq)%r) async
1119    ENDDO
1120    message%ondevice=.true.
1121  END SUBROUTINE
1122
1123  !TODO : add openacc with multiple process
1124  SUBROUTINE send_message_mpi(field,message)
1125  USE abort_mod
1126  USE profiling_mod
1127  USE field_mod
1128  USE domain_mod
1129  USE mpi_mod
1130  USE mpipara
1131  USE omp_para
1132  USE trace
1133  USE abort_mod
1134  IMPLICIT NONE
1135    TYPE(t_field),POINTER :: field(:)
1136    TYPE(t_message) :: message
1137    REAL(rstd),POINTER :: rval2d(:), src_rval2d(:)
1138    REAL(rstd),POINTER :: rval3d(:,:), src_rval3d(:,:)
1139    REAL(rstd),POINTER :: rval4d(:,:,:), src_rval4d(:,:,:)
1140    REAL(rstd),POINTER :: buffer_r(:)
1141    INTEGER,POINTER :: value(:)
1142    INTEGER,POINTER :: sgn(:)
1143    TYPE(ARRAY),POINTER :: recv,send
1144    TYPE(t_request),POINTER :: req
1145    INTEGER :: irecv,isend
1146    INTEGER :: ireq
1147    INTEGER :: ind,n
1148    INTEGER :: dim3,dim4,d3,d4
1149    INTEGER,POINTER :: src_value(:)
1150    INTEGER :: offset,msize,moffset,rank
1151    INTEGER :: lbegin, lend
1152    INTEGER :: max_req
1153
1154!    CALL trace_start("send_message_mpi")
1155
1156    CALL enter_profile(id_mpi)
1157
1158    !Prepare 'message' for on-device copies if field is on device
1159    if( field(1)%ondevice .AND. .NOT. message%ondevice ) call update_device_message_mpi(message)
1160
1161CALL enter_profile(profile_mpi_omp_barrier)
1162!$OMP BARRIER
1163CALL exit_profile(profile_mpi_omp_barrier)
1164
1165
1166!$OMP MASTER
1167    IF(message%open) THEN
1168       PRINT *, 'send_message_mpi : message ' // TRIM(message%name) // &
1169            ' is still open, no call to wait_message_mpi after last send_message_mpi'
1170       CALL dynamico_abort( "send_message_mpi : message still open" )
1171    END IF
1172    message%open=.TRUE. ! will be set to .FALSE. by wait_message_mpi
1173
1174    message%field=>field
1175
1176    IF (message%nreq>0) THEN
1177      message%completed=.FALSE.
1178      message%pending=.TRUE.
1179    ELSE
1180      message%completed=.TRUE.
1181      message%pending=.FALSE.
1182    ENDIF
1183!$OMP END MASTER
1184CALL enter_profile(profile_mpi_omp_barrier)
1185!$OMP BARRIER
1186CALL exit_profile(profile_mpi_omp_barrier)
1187     
1188    IF (field(1)%data_type==type_real) THEN
1189      IF (field(1)%ndim==2) THEN
1190
1191        DO ind=1,ndomain
1192          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1193
1194          rval2d=>field(ind)%rval2d
1195
1196          req=>message%request(ind)
1197          DO isend=1,req%nsend
1198            send=>req%send(isend)
1199            value=>send%value
1200
1201
1202            IF (send%rank/=mpi_rank) THEN
1203              ireq=send%ireq
1204
1205              buffer_r=>message%buffers(ireq)%r
1206              offset=send%offset
1207              msize=send%size
1208              call enter_profile(profile_mpi_copies)
1209              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1210              DO n=1,msize
1211                buffer_r(offset+n)=rval2d(value(n))
1212              ENDDO
1213              call exit_profile(profile_mpi_copies)
1214
1215              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1216                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1217                !$OMP CRITICAL
1218                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1219                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1220                !$OMP END CRITICAL
1221              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1222                CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1223                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1224                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1225              ENDIF
1226
1227            ENDIF
1228          ENDDO
1229        ENDDO
1230
1231        DO ind=1,ndomain
1232          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1233          rval2d=>field(ind)%rval2d
1234          req=>message%request(ind)
1235
1236          DO irecv=1,req%nrecv
1237            recv=>req%recv(irecv)
1238
1239            IF (recv%rank==mpi_rank) THEN
1240
1241              value=>recv%value
1242              src_value => recv%src_value
1243              src_rval2d=>field(recv%domain)%rval2d
1244              sgn=>recv%sign
1245              msize=recv%size
1246              call enter_profile(profile_mpi_copies)
1247              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1248              DO n=1,msize
1249                rval2d(value(n))=src_rval2d(src_value(n))*sgn(n)
1250              ENDDO
1251              call exit_profile(profile_mpi_copies)
1252                   
1253            ELSE
1254
1255              ireq=recv%ireq
1256              buffer_r=>message%buffers(ireq)%r
1257              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1258                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1259               !$OMP CRITICAL
1260                CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,               &
1261                  recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1262               !$OMP END CRITICAL
1263              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1264                 CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1265                 CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,              &
1266                   recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1267              ENDIF
1268
1269            ENDIF
1270          ENDDO
1271
1272        ENDDO
1273
1274
1275      ELSE  IF (field(1)%ndim==3) THEN
1276        max_req=0
1277        DO ind=1,ndomain
1278          req=>message%request(ind)
1279          IF (req%nsend>max_req) max_req=req%nsend
1280        ENDDO
1281
1282        DO ind=1,ndomain
1283          IF (.NOT. assigned_domain(ind) ) CYCLE
1284
1285          dim3=size(field(ind)%rval3d,2)
1286          CALL distrib_level(1,dim3, lbegin,lend)
1287
1288          rval3d=>field(ind)%rval3d
1289          req=>message%request(ind)
1290
1291          DO isend=1,req%nsend
1292            send=>req%send(isend)
1293            value=>send%value
1294
1295            IF (send%rank/=mpi_rank) THEN
1296              ireq=send%ireq
1297              buffer_r=>message%buffers(ireq)%r
1298
1299              msize=send%size
1300              moffset=send%offset
1301              call enter_profile(profile_mpi_copies)
1302
1303              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1304              DO d3=lbegin,lend
1305                offset=moffset*dim3 + (d3-1)*msize
1306                !$acc loop
1307                DO n=1,msize
1308                  buffer_r(n+offset)=rval3d(value(n),d3)
1309                ENDDO
1310              ENDDO
1311              call exit_profile(profile_mpi_copies)
1312
1313              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1314CALL enter_profile(profile_mpi_omp_barrier)
1315!$OMP BARRIER
1316CALL exit_profile(profile_mpi_omp_barrier)
1317
1318              ENDIF
1319
1320              IF (is_omp_level_master) THEN
1321                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1322                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1323                  !$OMP CRITICAL
1324                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1325                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1326                  !$OMP END CRITICAL
1327                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1328                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1329                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1330                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1331                ENDIF
1332              ENDIF
1333            ELSE
1334              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1335CALL enter_profile(profile_mpi_omp_barrier)
1336!$OMP BARRIER
1337CALL exit_profile(profile_mpi_omp_barrier)
1338
1339              ENDIF
1340            ENDIF
1341          ENDDO
1342
1343          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1344            DO isend=req%nsend+1,max_req
1345CALL enter_profile(profile_mpi_omp_barrier)
1346!$OMP BARRIER
1347CALL exit_profile(profile_mpi_omp_barrier)
1348
1349            ENDDO
1350          ENDIF
1351
1352        ENDDO
1353
1354        DO ind=1,ndomain
1355          IF (.NOT. assigned_domain(ind) ) CYCLE
1356          dim3=size(field(ind)%rval3d,2)
1357          CALL distrib_level(1,dim3, lbegin,lend)
1358          rval3d=>field(ind)%rval3d
1359          req=>message%request(ind)
1360
1361          DO irecv=1,req%nrecv
1362            recv=>req%recv(irecv)
1363
1364            IF (recv%rank==mpi_rank) THEN
1365              value=>recv%value
1366              src_value => recv%src_value
1367              src_rval3d=>field(recv%domain)%rval3d
1368              sgn=>recv%sign
1369              msize=recv%size
1370
1371              call enter_profile(profile_mpi_copies)
1372              CALL trace_start("copy_data")
1373              !$acc parallel loop collapse(2) default(present) async if (field(ind)%ondevice)
1374              DO d3=lbegin,lend
1375                DO n=1,msize
1376                  rval3d(value(n),d3)=src_rval3d(src_value(n),d3)*sgn(n)
1377                ENDDO
1378              ENDDO
1379              call exit_profile(profile_mpi_copies)
1380              CALL trace_end("copy_data")
1381
1382            ELSE
1383              ireq=recv%ireq
1384              buffer_r=>message%buffers(ireq)%r
1385
1386              IF (is_omp_level_master) THEN
1387                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1388                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1389                  !$OMP CRITICAL
1390                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1391                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1392                  !$OMP END CRITICAL
1393                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1394                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1395                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1396                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1397                ENDIF
1398              ENDIF
1399            ENDIF
1400          ENDDO
1401
1402        ENDDO
1403
1404      ELSE  IF (field(1)%ndim==4) THEN
1405
1406        max_req=0
1407        DO ind=1,ndomain
1408          req=>message%request(ind)
1409          IF (req%nsend>max_req) max_req=req%nsend
1410        ENDDO
1411
1412        DO ind=1,ndomain
1413          IF (.NOT. assigned_domain(ind) ) CYCLE
1414
1415          dim3=size(field(ind)%rval4d,2)
1416          CALL distrib_level(1,dim3, lbegin,lend)
1417          dim4=size(field(ind)%rval4d,3)
1418          rval4d=>field(ind)%rval4d
1419          req=>message%request(ind)
1420
1421          DO isend=1,req%nsend
1422            send=>req%send(isend)
1423            value=>send%value
1424
1425            IF (send%rank/=mpi_rank) THEN
1426
1427              ireq=send%ireq
1428              buffer_r=>message%buffers(ireq)%r
1429              msize=send%size
1430              moffset=send%offset
1431
1432              call enter_profile(profile_mpi_copies)
1433              CALL trace_start("copy_to_buffer")
1434              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice)
1435              DO d4=1,dim4
1436                DO d3=lbegin,lend
1437                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               &
1438                    (d3-1)*msize
1439                  !$acc loop
1440                  DO n=1,msize
1441                    buffer_r(n+offset)=rval4d(value(n),d3,d4)
1442                  ENDDO
1443                ENDDO
1444              ENDDO
1445              CALL trace_end("copy_to_buffer")
1446              call exit_profile(profile_mpi_copies)
1447
1448              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1449CALL enter_profile(profile_mpi_omp_barrier)
1450!$OMP BARRIER
1451CALL exit_profile(profile_mpi_omp_barrier)
1452
1453              ENDIF
1454
1455              IF (is_omp_level_master) THEN
1456                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1457                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1458                  !$OMP CRITICAL
1459                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1460                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1461                  !$OMP END CRITICAL
1462                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1463                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1464                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1465                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1466                ENDIF
1467              ENDIF
1468            ELSE
1469              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1470CALL enter_profile(profile_mpi_omp_barrier)
1471!$OMP BARRIER
1472CALL exit_profile(profile_mpi_omp_barrier)
1473
1474              ENDIF
1475            ENDIF
1476          ENDDO
1477
1478          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1479            DO isend=req%nsend+1,max_req
1480CALL enter_profile(profile_mpi_omp_barrier)
1481!$OMP BARRIER
1482CALL exit_profile(profile_mpi_omp_barrier)
1483
1484            ENDDO
1485          ENDIF
1486
1487        ENDDO
1488
1489        DO ind=1,ndomain
1490          IF (.NOT. assigned_domain(ind) ) CYCLE
1491
1492          dim3=size(field(ind)%rval4d,2)
1493          CALL distrib_level(1,dim3, lbegin,lend)
1494          dim4=size(field(ind)%rval4d,3)
1495          rval4d=>field(ind)%rval4d
1496          req=>message%request(ind)
1497         
1498          DO irecv=1,req%nrecv
1499            recv=>req%recv(irecv)
1500            IF (recv%rank==mpi_rank) THEN
1501              value=>recv%value
1502              src_value => recv%src_value
1503              src_rval4d=>field(recv%domain)%rval4d
1504              sgn=>recv%sign
1505              msize=recv%size
1506              call enter_profile(profile_mpi_copies)
1507              CALL trace_start("copy_data")
1508              !$acc parallel loop collapse(3) default(present) async if (field(ind)%ondevice)
1509              DO d4=1,dim4
1510                DO d3=lbegin,lend
1511                  DO n=1,msize
1512                    rval4d(value(n),d3,d4)=src_rval4d(src_value(n),d3,d4)*sgn(n)
1513                  ENDDO
1514                ENDDO
1515              ENDDO
1516              call exit_profile(profile_mpi_copies)
1517              CALL trace_end("copy_data")
1518
1519            ELSE
1520
1521              ireq=recv%ireq
1522              buffer_r=>message%buffers(ireq)%r
1523              IF (is_omp_level_master) THEN
1524                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1525                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED")
1526                 !$OMP CRITICAL
1527                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1528                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1529                  !$OMP END CRITICAL
1530                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1531                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE")
1532                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1533                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1534                ENDIF
1535              ENDIF
1536            ENDIF
1537          ENDDO
1538        ENDDO
1539
1540      ENDIF
1541
1542      IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
1543CALL enter_profile(profile_mpi_omp_barrier)
1544!$acc wait
1545!$OMP BARRIER
1546CALL exit_profile(profile_mpi_omp_barrier)
1547!$OMP MASTER       
1548
1549        DO ireq=1,message%nreq_send
1550          buffer_r=>message%buffers(ireq)%r
1551          msize=message%buffers(ireq)%size
1552          rank=message%buffers(ireq)%rank
1553          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation
1554          IF (message%ondevice) THEN
1555            !$acc host_data use_device(buffer_r) ! if (message%ondevice)
1556            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1557              comm_icosa, message%mpi_req(ireq),ierr)
1558            !$acc end host_data
1559          ELSE
1560            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1561              comm_icosa, message%mpi_req(ireq),ierr)
1562          ENDIF
1563        ENDDO
1564
1565        DO ireq=message%nreq_send+1,message%nreq_send+message%nreq_recv
1566          buffer_r=>message%buffers(ireq)%r
1567          msize=message%buffers(ireq)%size
1568          rank=message%buffers(ireq)%rank
1569          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation
1570          IF (message%ondevice) THEN
1571            !$acc host_data use_device(buffer_r) ! if (message%ondevice)
1572            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1573              comm_icosa, message%mpi_req(ireq),ierr)
1574            !$acc end host_data
1575          ELSE
1576            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1577              comm_icosa, message%mpi_req(ireq),ierr)
1578          ENDIF
1579        ENDDO
1580
1581!$OMP END MASTER
1582      ENDIF
1583    ENDIF
1584CALL enter_profile(profile_mpi_omp_barrier)
1585!$OMP BARRIER
1586CALL exit_profile(profile_mpi_omp_barrier)
1587 
1588!    CALL trace_end("send_message_mpi")
1589
1590    CALL exit_profile(id_mpi)
1591
1592  END SUBROUTINE send_message_mpi
1593
1594  SUBROUTINE test_message_mpi(message)
1595  IMPLICIT NONE
1596    TYPE(t_message) :: message
1597
1598    INTEGER :: ierr
1599
1600!$OMP MASTER
1601    IF (message%pending .AND. .NOT. message%completed) CALL MPI_TESTALL(message%nreq,&
1602      message%mpi_req,message%completed,message%status,ierr)
1603!$OMP END MASTER
1604  END SUBROUTINE  test_message_mpi
1605
1606
1607  SUBROUTINE wait_message_mpi(message)
1608  USE profiling_mod
1609  USE field_mod
1610  USE domain_mod
1611  USE mpi_mod
1612  USE mpipara
1613  USE omp_para
1614  USE trace
1615  IMPLICIT NONE
1616    TYPE(t_message) :: message
1617
1618    TYPE(t_field),POINTER :: field(:)
1619    REAL(rstd),POINTER :: rval2d(:)
1620    REAL(rstd),POINTER :: rval3d(:,:)
1621    REAL(rstd),POINTER :: rval4d(:,:,:)
1622    REAL(rstd),POINTER :: buffer_r(:)
1623    INTEGER,POINTER :: value(:)
1624    INTEGER,POINTER :: sgn(:)
1625    TYPE(ARRAY),POINTER :: recv
1626    TYPE(t_request),POINTER :: req
1627    INTEGER :: irecv
1628    INTEGER :: ireq,nreq
1629    INTEGER :: ind,n
1630    INTEGER :: dim3,dim4,d3,d4,lbegin,lend
1631    INTEGER :: offset, msize, moffset
1632
1633    message%open=.FALSE.
1634    IF (.NOT. message%pending) RETURN
1635
1636    CALL enter_profile(id_mpi)
1637
1638!    CALL trace_start("wait_message_mpi")
1639
1640    field=>message%field
1641    nreq=message%nreq
1642
1643    IF (field(1)%data_type==type_real) THEN
1644      IF (field(1)%ndim==2) THEN
1645
1646call enter_profile(profile_mpi_waitall)
1647!$OMP MASTER         
1648         IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1649          message%status,ierr)
1650!$OMP END MASTER
1651!$OMP BARRIER
1652call exit_profile(profile_mpi_waitall)
1653        call enter_profile(profile_mpi_copies)
1654        DO ind=1,ndomain
1655          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1656
1657          rval2d=>field(ind)%rval2d
1658          req=>message%request(ind)
1659          DO irecv=1,req%nrecv
1660            recv=>req%recv(irecv)
1661            IF (recv%rank/=mpi_rank) THEN
1662              ireq=recv%ireq
1663              buffer_r=>message%buffers(ireq)%r
1664              value=>recv%value
1665              sgn=>recv%sign
1666              offset=recv%offset
1667              msize=recv%size
1668              !$acc parallel loop default(present) async if (field(ind)%ondevice)
1669              DO n=1,msize
1670                rval2d(value(n))=buffer_r(n+offset)*sgn(n)
1671              ENDDO
1672
1673            ENDIF
1674          ENDDO
1675
1676        ENDDO
1677        call exit_profile(profile_mpi_copies)
1678     
1679     
1680      ELSE  IF (field(1)%ndim==3) THEN
1681         call enter_profile(profile_mpi_waitall)
1682!$OMP MASTER
1683         IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1684              message%status,ierr)
1685!$OMP END MASTER
1686!$OMP BARRIER
1687        call exit_profile(profile_mpi_waitall)
1688
1689
1690        DO ind=1,ndomain
1691          IF (.NOT. assigned_domain(ind) ) CYCLE
1692
1693          rval3d=>field(ind)%rval3d
1694          req=>message%request(ind)
1695          DO irecv=1,req%nrecv
1696            recv=>req%recv(irecv)
1697            IF (recv%rank/=mpi_rank) THEN
1698              ireq=recv%ireq
1699              buffer_r=>message%buffers(ireq)%r
1700              value=>recv%value
1701              sgn=>recv%sign
1702
1703              dim3=size(rval3d,2)
1704
1705              CALL distrib_level(1,dim3, lbegin,lend)
1706              msize=recv%size
1707              moffset=recv%offset
1708              call enter_profile(profile_mpi_copies)
1709              CALL trace_start("copy_from_buffer")
1710
1711              IF (req%vector) THEN
1712                !$acc parallel loop default(present) async if (field(ind)%ondevice)
1713                DO d3=lbegin,lend
1714                  offset=moffset*dim3 + (d3-1)*msize
1715                  !$acc loop
1716                  DO n=1,msize
1717                    rval3d(value(n),d3)=buffer_r(n+offset)*sgn(n)
1718                  ENDDO
1719                ENDDO
1720              ELSE
1721                !$acc parallel loop default(present) async if (field(ind)%ondevice)
1722                DO d3=lbegin,lend
1723                  offset=moffset*dim3 + (d3-1)*msize
1724                  !$acc loop
1725                  DO n=1,msize
1726                    rval3d(value(n),d3)=buffer_r(n+offset)
1727                  ENDDO
1728                ENDDO
1729              ENDIF
1730
1731              CALL trace_end("copy_from_buffer")
1732              call exit_profile(profile_mpi_copies)
1733            ENDIF
1734          ENDDO
1735
1736        ENDDO
1737
1738      ELSE  IF (field(1)%ndim==4) THEN
1739call enter_profile(profile_mpi_waitall)
1740!$OMP MASTER
1741        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1742          message%status,ierr)
1743!$OMP END MASTER
1744!$OMP BARRIER
1745        call exit_profile(profile_mpi_waitall)
1746
1747
1748        DO ind=1,ndomain
1749          IF (.NOT. assigned_domain(ind) ) CYCLE
1750
1751          rval4d=>field(ind)%rval4d
1752          req=>message%request(ind)
1753          DO irecv=1,req%nrecv
1754            recv=>req%recv(irecv)
1755            IF (recv%rank/=mpi_rank) THEN
1756              ireq=recv%ireq
1757              buffer_r=>message%buffers(ireq)%r
1758              value=>recv%value
1759              sgn=>recv%sign
1760
1761              dim3=size(rval4d,2)
1762              CALL distrib_level(1,dim3, lbegin,lend)
1763              dim4=size(rval4d,3)
1764              msize=recv%size
1765              moffset=recv%offset
1766              call enter_profile(profile_mpi_copies)
1767              CALL trace_start("copy_from_buffer")
1768              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice)
1769              DO d4=1,dim4
1770                DO d3=lbegin,lend
1771                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               &
1772                    (d3-1)*msize
1773                  !$acc loop
1774                  DO n=1,msize
1775                    rval4d(value(n),d3,d4)=buffer_r(n+offset)*sgn(n)
1776                  ENDDO
1777                ENDDO
1778              ENDDO
1779              CALL trace_end("copy_from_buffer")
1780              call exit_profile(profile_mpi_copies)
1781            ENDIF
1782          ENDDO
1783
1784        ENDDO
1785
1786      ENDIF
1787
1788    ENDIF
1789
1790!$OMP MASTER
1791    message%pending=.FALSE.
1792!$OMP END MASTER
1793
1794!    CALL trace_end("wait_message_mpi")
1795!$OMP BARRIER
1796
1797    CALL exit_profile(id_mpi)
1798
1799  END SUBROUTINE wait_message_mpi
1800
1801
1802  SUBROUTINE transfert_request_seq(field,request)
1803  USE field_mod
1804  USE domain_mod
1805  IMPLICIT NONE
1806    TYPE(t_field),POINTER :: field(:)
1807    TYPE(t_request),POINTER :: request(:)
1808    REAL(rstd),POINTER :: rval2d(:)
1809    REAL(rstd),POINTER :: rval3d(:,:)
1810    REAL(rstd),POINTER :: rval4d(:,:,:)
1811    INTEGER :: ind
1812    TYPE(t_request),POINTER :: req
1813    INTEGER :: n
1814
1815    DO ind=1,ndomain
1816      req=>request(ind)
1817      rval2d=>field(ind)%rval2d
1818      rval3d=>field(ind)%rval3d
1819      rval4d=>field(ind)%rval4d
1820
1821      IF (field(ind)%data_type==type_real) THEN
1822        IF (field(ind)%ndim==2) THEN
1823          DO n=1,req%size
1824            rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*&
1825              req%target_sign(n)
1826          ENDDO
1827        ELSE IF (field(ind)%ndim==3) THEN
1828          DO n=1,req%size
1829            rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*&
1830              req%target_sign(n)
1831          ENDDO
1832        ELSE IF (field(ind)%ndim==4) THEN
1833          DO n=1,req%size
1834            rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*&
1835              req%target_sign(n)
1836          ENDDO
1837        ENDIF
1838      ENDIF
1839
1840    ENDDO
1841
1842  END SUBROUTINE transfert_request_seq
1843
1844
1845END MODULE transfert_mpi_legacy_mod
1846
1847
1848
1849
1850
Note: See TracBrowser for help on using the repository browser.