source: codes/icosagcm/trunk/src/transfert_mpi.f90 @ 378

Last change on this file since 378 was 364, checked in by dubos, 9 years ago

Bugfix : memory leak in transfert_mpi / New : detect send_message not paired with wait_message

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