source: codes/icosagcm/devel/src/parallel/transfert_mpi.f90 @ 603

Last change on this file since 603 was 603, checked in by dubos, 7 years ago

devel : output theta and momentum fluxes

File size: 61.0 KB
Line 
1MODULE transfert_mpi_mod
2USE genmod
3USE field_mod
4 
5  TYPE array
6    INTEGER,POINTER :: value(:)
7    INTEGER,POINTER :: sign(:)
8    INTEGER         :: domain
9    INTEGER         :: rank
10    INTEGER         :: tag
11    INTEGER         :: size
12    INTEGER         :: offset
13    INTEGER         :: ireq
14    INTEGER,POINTER :: buffer(:)
15    REAL,POINTER    :: buffer_r(:)
16    INTEGER,POINTER :: src_value(:)
17  END TYPE array
18 
19  TYPE t_buffer
20    REAL,POINTER    :: r(:)
21    INTEGER         :: size
22    INTEGER         :: rank
23  END TYPE t_buffer   
24   
25  TYPE t_request
26    INTEGER :: type_field
27    INTEGER :: max_size
28    INTEGER :: size
29    LOGICAL :: vector
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(:)
37    INTEGER,POINTER :: target_sign(:)
38    INTEGER :: nrecv
39    TYPE(ARRAY),POINTER :: recv(:)
40    INTEGER :: nsend
41    INTEGER :: nreq_mpi
42    INTEGER :: nreq_send
43    INTEGER :: nreq_recv
44    TYPE(ARRAY),POINTER :: send(:)
45  END TYPE t_request
46 
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(:)
50 
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 
61 
62  TYPE t_message
63    CHARACTER(LEN=100) :: name ! for debug
64    TYPE(t_request), POINTER :: request(:)
65    INTEGER :: nreq
66    INTEGER :: nreq_send
67    INTEGER :: nreq_recv
68    TYPE(t_reorder), POINTER :: reorder(:)
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
75    LOGICAL :: open      ! for debug
76    INTEGER :: number
77  END TYPE t_message
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
88 
89CONTAINS
90       
91     
92  SUBROUTINE init_transfert
93  USE domain_mod
94  USE dimensions
95  USE field_mod
96  USE metric
97  USE mpipara
98  USE mpi_mod
99  IMPLICIT NONE
100  INTEGER :: ind,i,j
101  LOGICAL ::ok
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
268  END SUBROUTINE init_transfert
269 
270  SUBROUTINE create_request(type_field,request,vector)
271  USE domain_mod
272  USE field_mod
273  IMPLICIT NONE
274    INTEGER :: type_field
275    TYPE(t_request),POINTER :: request(:)
276    LOGICAL,OPTIONAL :: vector
277   
278    TYPE(t_request),POINTER :: req
279    TYPE(t_domain),POINTER :: d
280    INTEGER :: ind
281    INTEGER :: max_size
282       
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
299      req%vector=.FALSE.
300      IF (PRESENT(vector)) req%vector=vector
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))
308      ALLOCATE(req%target_sign(req%max_size))
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(:)
324    INTEGER,POINTER :: target_sign(:)
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
334    target_sign=>req%target_sign
335
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))
343    ALLOCATE(req%target_sign(req%max_size*2))
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(:)
352    req%target_sign(1:req%max_size)=target_sign(:)
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)
363    DEALLOCATE(target_sign)
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
395        req%target_sign(req%size)=1
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)
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
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)
430        req%target_sign(req%size)=1
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(:)
443    TYPE(t_request),POINTER :: req, req_src
444    INTEGER :: nb_domain_recv(0:mpi_size-1)
445    INTEGER :: nb_domain_send(0:mpi_size-1)
446    INTEGER :: tag_rank(0:mpi_size-1)
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
452    INTEGER :: rank,i,j,pos
453    INTEGER :: size_,ind_glo,ind_loc, ind_src
454    INTEGER :: isend, irecv, ireq, nreq, nsend, nrecv
455    INTEGER, ALLOCATABLE :: mpi_req(:)
456    INTEGER, ALLOCATABLE :: status(:,:)
457    INTEGER, ALLOCATABLE :: rank_list(:)
458    INTEGER, ALLOCATABLE :: offset(:)
459    LOGICAL,PARAMETER :: debug = .FALSE.
460
461 
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
469      tag_rank(:)=0
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))
492          ALLOCATE(req%recv(irecv)%sign(req%recv(irecv)%size))
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
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)
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   
528
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)
535        IF (debug) PRINT *,"Isend ",req%recv(irecv)%domain, "from ", mpi_rank, "to ",req%recv(irecv)%rank, "tag ",0
536      ENDDO
537    ENDDO
538
539    IF (debug) PRINT *,"------------"   
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)
546        IF (debug) PRINT *,"IRecv ",list_domain_send(j), "from ", rank, "to ",mpi_rank, "tag ",0
547      ENDDO
548    ENDDO
549    IF (debug) PRINT *,"------------"   
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
564
565    IF (debug) PRINT *,"------------"   
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)
574       IF (debug) PRINT *,"Isend ",mpi_rank, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
575     ENDDO
576    IF (debug) PRINT *,"------------"   
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)
581       IF (debug) PRINT *,"IRecv ",req%send(isend)%rank, "from ", "xxx", "to ",mpi_rank, "tag ",ind_loc
582     ENDDO
583   ENDDO
584
585   IF (debug) PRINT *,"------------"   
586
587   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
588   CALL MPI_BARRIER(comm_icosa,ierr)
589
590   IF (debug) PRINT *,"------------"   
591
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
598       CALL MPI_ISEND(ind_loc,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr)
599       IF (debug) PRINT *,"Isend ",ind_loc, "from ", mpi_rank, "to ",req%recv(irecv)%rank,"tag ",req%recv(irecv)%domain
600     ENDDO
601
602     IF (debug) PRINT *,"------------"   
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)
607       IF (debug) PRINT *,"IRecv ",req%send(isend)%domain, "from ", req%send(isend)%rank, "to ",mpi_rank, "tag ",ind_loc
608     ENDDO
609   ENDDO
610   IF (debug) PRINT *,"------------"   
611   
612   CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
613   CALL MPI_BARRIER(comm_icosa,ierr)
614   IF (debug) PRINT *,"------------"   
615
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
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
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
651     ENDDO
652     IF (debug) PRINT *,"------------"   
653     
654     DO isend=1,req%nsend
655       ireq=ireq+1
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
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
669       CALL MPI_ISEND(req%recv(irecv)%value,req%recv(irecv)%size,MPI_INTEGER,&
670            req%recv(irecv)%rank,req%recv(irecv)%tag,comm_icosa, mpi_req(ireq),ierr)
671     ENDDO
672     
673     DO isend=1,req%nsend
674       ireq=ireq+1
675       ALLOCATE(req%send(isend)%value(req%send(isend)%size))
676       CALL MPI_IRECV(req%send(isend)%value,req%send(isend)%size,MPI_INTEGER,&
677            req%send(isend)%rank,req%send(isend)%tag,comm_icosa, mpi_req(ireq),ierr)
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(:)
688       req%recv(irecv)%sign(:) =req%recv(irecv)%sign(:)
689       DEALLOCATE(req%recv(irecv)%buffer)
690     ENDDO
691   ENDDO 
692   
693
694! domain is on the same mpi process => copie memory to memory
695   
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
704             IF (req_src%send(isend)%rank==mpi_rank .AND. req_src%send(isend)%tag==req%recv(irecv)%tag) THEN
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
707                 PRINT *,ind_loc, irecv, isend, size(req%recv(irecv)%value), size(req_src%send(isend)%value)
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
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
727   DO ind_loc=1,ndomain
728     req=>request(ind_loc)
729
730     DO isend=1,req%nsend
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
742           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
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
754     ENDDO
755   ENDDO
756
757   DEALLOCATE(rank_list)
758   DEALLOCATE(offset)
759     
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
768     DO irecv=1,req%nrecv
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
780           IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
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
792     ENDDO
793   ENDDO 
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     
816       
817  END SUBROUTINE Finalize_request 
818
819
820  SUBROUTINE init_message_seq(field, request, message, name)
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
830    CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: name
831!$OMP MASTER   
832    message%request=>request
833    IF(PRESENT(name)) THEN
834       message%name = TRIM(name)
835    ELSE
836       message%name = 'unknown'
837    END IF
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   
885
886
887   
888  SUBROUTINE init_message_mpi(field,request, message, name)
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
899    CHARACTER(LEN=*), INTENT(IN),OPTIONAL :: name
900
901    TYPE(ARRAY),POINTER :: recv,send 
902    TYPE(t_request),POINTER :: req
903    INTEGER :: irecv,isend
904    INTEGER :: ireq,nreq, nreq_send
905    INTEGER :: ind
906    INTEGER :: dim3,dim4
907    INTEGER :: i,j
908    INTEGER,SAVE :: message_number=0
909!    TYPE(t_reorder),POINTER :: reorder(:)
910!    TYPE(t_reorder) :: reorder_swap
911
912!$OMP BARRIER
913!$OMP MASTER
914    IF(PRESENT(name)) THEN
915       message%name = TRIM(name)
916    ELSE
917       message%name = 'unknown'
918    END IF
919    message%number=message_number
920    message_number=message_number+1
921    IF (message_number==100) message_number=0
922
923 
924    message%request=>request
925    message%nreq=sum(message%request(:)%nreq_mpi)
926    message%nreq_send=sum(message%request(:)%nreq_send)
927    message%nreq_recv=sum(message%request(:)%nreq_recv)
928    nreq=message%nreq
929
930    ALLOCATE(message%mpi_req(nreq))
931    ALLOCATE(message%buffers(nreq))
932    ALLOCATE(message%status(MPI_STATUS_SIZE,nreq))
933    message%buffers(:)%size=0
934    message%pending=.FALSE.
935    message%completed=.FALSE.
936    message%open=.FALSE.
937
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
957    IF (field(1)%data_type==type_real) THEN
958
959      IF (field(1)%ndim==2) THEN
960     
961        DO ireq=1,message%nreq
962          CALL allocate_mpi_buffer(message%buffers(ireq)%r,message%buffers(ireq)%size)
963        ENDDO
964     
965      ELSE  IF (field(1)%ndim==3) THEN
966     
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)
971        ENDDO
972     
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
982     
983         
984   
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   
1014 
1015!$OMP END MASTER
1016!$OMP BARRIER   
1017
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)
1045        ENDDO
1046     
1047      ELSE  IF (field(1)%ndim==3) THEN
1048
1049        DO ireq=1,message%nreq
1050          CALL free_mpi_buffer(message%buffers(ireq)%r)
1051        ENDDO
1052     
1053      ELSE  IF (field(1)%ndim==4) THEN
1054
1055        DO ireq=1,message%nreq
1056          CALL free_mpi_buffer(message%buffers(ireq)%r)
1057        ENDDO
1058
1059      ENDIF     
1060    ENDIF
1061   
1062    DEALLOCATE(message%mpi_req)
1063    DEALLOCATE(message%buffers)
1064    DEALLOCATE(message%status)
1065
1066!$OMP END MASTER
1067!$OMP BARRIER
1068
1069     
1070  END SUBROUTINE Finalize_message_mpi
1071
1072
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
1105    REAL(rstd),POINTER :: rval2d(:), src_rval2d(:) 
1106    REAL(rstd),POINTER :: rval3d(:,:), src_rval3d(:,:) 
1107    REAL(rstd),POINTER :: rval4d(:,:,:), src_rval4d(:,:,:) 
1108    REAL(rstd),POINTER :: buffer_r(:) 
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
1116    INTEGER :: ireq,nreq
1117    INTEGER :: ind,i,n,l,m
1118    INTEGER :: dim3,dim4,d3,d4
1119    INTEGER,POINTER :: src_value(:)
1120    INTEGER,POINTER :: sign(:)
1121    INTEGER :: offset,msize,rank
1122    INTEGER :: lbegin, lend
1123    INTEGER :: max_req
1124
1125!    CALL trace_start("send_message_mpi")
1126
1127!$OMP BARRIER
1128
1129
1130!$OMP MASTER
1131    IF(message%open) THEN
1132       PRINT *, 'send_message_mpi : message ' // TRIM(message%name) // &
1133            ' is still open, no call to wait_message_mpi after last send_message_mpi'
1134       CALL ABORT
1135    END IF
1136    message%open=.TRUE. ! will be set to .FALSE. by wait_message_mpi
1137
1138    message%field=>field
1139
1140    IF (message%nreq>0) THEN
1141      message%completed=.FALSE.
1142      message%pending=.TRUE.
1143    ELSE
1144      message%completed=.TRUE.
1145      message%pending=.FALSE.
1146    ENDIF
1147!$OMP END MASTER
1148!$OMP BARRIER
1149     
1150    IF (field(1)%data_type==type_real) THEN
1151      IF (field(1)%ndim==2) THEN
1152
1153        DO ind=1,ndomain
1154          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1155         
1156          rval2d=>field(ind)%rval2d
1157       
1158          req=>message%request(ind)
1159          DO isend=1,req%nsend
1160            send=>req%send(isend)
1161            value=>send%value
1162
1163           
1164            IF (send%rank/=mpi_rank) THEN
1165              ireq=send%ireq
1166
1167              buffer_r=>message%buffers(ireq)%r
1168              offset=send%offset
1169                                           
1170              DO n=1,send%size
1171                buffer_r(offset+n)=rval2d(value(n))
1172              ENDDO
1173
1174              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1175                !$OMP CRITICAL           
1176                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1177                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1178                !$OMP END CRITICAL
1179              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1180                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1181                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1182              ENDIF
1183             
1184            ENDIF
1185          ENDDO
1186        ENDDO
1187       
1188        DO ind=1,ndomain
1189          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1190          rval2d=>field(ind)%rval2d
1191          req=>message%request(ind)       
1192
1193          DO irecv=1,req%nrecv
1194            recv=>req%recv(irecv)
1195
1196            IF (recv%rank==mpi_rank) THEN
1197
1198              value=>recv%value
1199              src_value => recv%src_value
1200              src_rval2d=>field(recv%domain)%rval2d
1201              sgn=>recv%sign
1202
1203              DO n=1,recv%size
1204                rval2d(value(n))=src_rval2d(src_value(n))*sgn(n)
1205              ENDDO
1206               
1207                   
1208            ELSE
1209           
1210              ireq=recv%ireq
1211              buffer_r=>message%buffers(ireq)%r
1212              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1213               !$OMP CRITICAL           
1214                CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,               &
1215                  recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1216               !$OMP END CRITICAL
1217              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1218                 CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,              &
1219                   recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1220              ENDIF
1221           
1222            ENDIF
1223          ENDDO
1224       
1225        ENDDO
1226       
1227     
1228      ELSE  IF (field(1)%ndim==3) THEN
1229        max_req=0
1230        DO ind=1,ndomain
1231          req=>message%request(ind)
1232          IF (req%nsend>max_req) max_req=req%nsend
1233        ENDDO
1234             
1235        DO ind=1,ndomain
1236          IF (.NOT. assigned_domain(ind) ) CYCLE
1237
1238          dim3=size(field(ind)%rval3d,2)
1239          CALL distrib_level(1,dim3, lbegin,lend)
1240
1241          rval3d=>field(ind)%rval3d
1242          req=>message%request(ind)
1243 
1244          DO isend=1,req%nsend
1245            send=>req%send(isend)
1246            value=>send%value
1247
1248            IF (send%rank/=mpi_rank) THEN
1249              ireq=send%ireq
1250              buffer_r=>message%buffers(ireq)%r
1251
1252              offset=send%offset*dim3 + (lbegin-1)*send%size
1253         
1254              CALL trace_start("copy_to_buffer")
1255
1256              DO d3=lbegin,lend
1257                DO n=1,send%size
1258                  buffer_r(n+offset)=rval3d(value(n),d3)
1259                ENDDO
1260                offset=offset+send%size
1261              ENDDO
1262              CALL trace_end("copy_to_buffer")
1263
1264              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1265                !$OMP BARRIER
1266              ENDIF
1267             
1268              IF (is_omp_level_master) THEN
1269                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1270                  !$OMP CRITICAL   
1271                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1272                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1273                  !$OMP END CRITICAL
1274                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1275                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1276                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1277                ENDIF
1278              ENDIF
1279            ELSE
1280              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1281                !$OMP BARRIER
1282              ENDIF
1283            ENDIF
1284          ENDDO
1285
1286          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1287            DO isend=req%nsend+1,max_req
1288              !$OMP BARRIER
1289            ENDDO
1290          ENDIF
1291
1292        ENDDO
1293         
1294        DO ind=1,ndomain
1295          IF (.NOT. assigned_domain(ind) ) CYCLE
1296          dim3=size(field(ind)%rval3d,2)
1297          CALL distrib_level(1,dim3, lbegin,lend)
1298          rval3d=>field(ind)%rval3d
1299          req=>message%request(ind)
1300
1301          DO irecv=1,req%nrecv
1302            recv=>req%recv(irecv)
1303
1304            IF (recv%rank==mpi_rank) THEN
1305              value=>recv%value
1306              src_value => recv%src_value
1307              src_rval3d=>field(recv%domain)%rval3d
1308              sgn=>recv%sign
1309
1310              CALL trace_start("copy_data")
1311
1312              DO d3=lbegin,lend
1313                DO n=1,recv%size
1314                  rval3d(value(n),d3)=src_rval3d(src_value(n),d3)*sgn(n)
1315                ENDDO
1316              ENDDO
1317              CALL trace_end("copy_data")
1318
1319            ELSE
1320              ireq=recv%ireq
1321              buffer_r=>message%buffers(ireq)%r
1322 
1323              IF (is_omp_level_master) THEN
1324                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1325                  !$OMP CRITICAL
1326                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1327                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1328                  !$OMP END CRITICAL
1329                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1330                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1331                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1332                ENDIF
1333              ENDIF
1334            ENDIF 
1335          ENDDO
1336       
1337        ENDDO
1338
1339      ELSE  IF (field(1)%ndim==4) THEN
1340
1341        max_req=0
1342        DO ind=1,ndomain
1343          req=>message%request(ind)
1344          IF (req%nsend>max_req) max_req=req%nsend
1345        ENDDO
1346   
1347        DO ind=1,ndomain
1348          IF (.NOT. assigned_domain(ind) ) CYCLE
1349
1350          dim3=size(field(ind)%rval4d,2)
1351          CALL distrib_level(1,dim3, lbegin,lend)
1352          dim4=size(field(ind)%rval4d,3)
1353          rval4d=>field(ind)%rval4d
1354          req=>message%request(ind)
1355
1356          DO isend=1,req%nsend
1357            send=>req%send(isend)
1358            value=>send%value
1359
1360            IF (send%rank/=mpi_rank) THEN
1361
1362              ireq=send%ireq
1363              buffer_r=>message%buffers(ireq)%r
1364
1365              CALL trace_start("copy_to_buffer")
1366              DO d4=1,dim4
1367                offset=send%offset*dim3*dim4 + dim3*send%size*(d4-1) +               &
1368                  (lbegin-1)*send%size
1369
1370                DO d3=lbegin,lend
1371                  DO n=1,send%size
1372                    buffer_r(n+offset)=rval4d(value(n),d3,d4)
1373                  ENDDO
1374                  offset=offset+send%size
1375                ENDDO
1376              ENDDO
1377              CALL trace_end("copy_to_buffer")
1378
1379              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1380                !$OMP BARRIER
1381              ENDIF
1382
1383              IF (is_omp_level_master) THEN
1384                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1385                  !$OMP CRITICAL
1386                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1387                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1388                  !$OMP END CRITICAL
1389                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1390                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1391                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1392                ENDIF
1393              ENDIF
1394            ELSE
1395              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1396                !$OMP BARRIER
1397              ENDIF
1398            ENDIF
1399          ENDDO
1400         
1401          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1402            DO isend=req%nsend+1,max_req
1403              !$OMP BARRIER
1404            ENDDO
1405          ENDIF
1406
1407        ENDDO
1408       
1409        DO ind=1,ndomain
1410          IF (.NOT. assigned_domain(ind) ) CYCLE
1411         
1412          dim3=size(field(ind)%rval4d,2)
1413          CALL distrib_level(1,dim3, lbegin,lend)
1414          dim4=size(field(ind)%rval4d,3)
1415          rval4d=>field(ind)%rval4d
1416          req=>message%request(ind)
1417
1418          DO irecv=1,req%nrecv
1419            recv=>req%recv(irecv)
1420            IF (recv%rank==mpi_rank) THEN
1421              value=>recv%value
1422              src_value => recv%src_value
1423              src_rval4d=>field(recv%domain)%rval4d
1424              sgn=>recv%sign
1425
1426              CALL trace_start("copy_data")
1427              DO d4=1,dim4
1428                DO d3=lbegin,lend
1429                  DO n=1,recv%size
1430                    rval4d(value(n),d3,d4)=src_rval4d(src_value(n),d3,d4)*sgn(n)
1431                  ENDDO
1432                ENDDO
1433              ENDDO
1434               
1435              CALL trace_end("copy_data")
1436                   
1437            ELSE
1438
1439              ireq=recv%ireq
1440              buffer_r=>message%buffers(ireq)%r
1441              IF (is_omp_level_master) THEN
1442                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1443                 !$OMP CRITICAL           
1444                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1445                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1446                  !$OMP END CRITICAL
1447                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1448                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1449                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1450                ENDIF
1451              ENDIF
1452            ENDIF
1453          ENDDO
1454        ENDDO
1455
1456      ENDIF     
1457
1458      IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
1459!$OMP BARRIER
1460!$OMP MASTER       
1461
1462        DO ireq=1,message%nreq_send
1463          buffer_r=>message%buffers(ireq)%r
1464          msize=message%buffers(ireq)%size
1465          rank=message%buffers(ireq)%rank
1466          CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1467            comm_icosa, message%mpi_req(ireq),ierr)
1468        ENDDO
1469
1470        DO ireq=message%nreq_send+1,message%nreq_send+message%nreq_recv
1471          buffer_r=>message%buffers(ireq)%r
1472          msize=message%buffers(ireq)%size
1473          rank=message%buffers(ireq)%rank
1474          CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1475            comm_icosa, message%mpi_req(ireq),ierr)
1476        ENDDO
1477
1478!$OMP END MASTER
1479      ENDIF             
1480    ENDIF
1481   
1482!$OMP BARRIER
1483!    CALL trace_end("send_message_mpi")
1484   
1485  END SUBROUTINE send_message_mpi
1486 
1487  SUBROUTINE test_message_mpi(message)
1488  IMPLICIT NONE
1489    TYPE(t_message) :: message
1490   
1491    INTEGER :: ierr
1492
1493!$OMP MASTER
1494    IF (message%pending .AND. .NOT. message%completed) CALL MPI_TESTALL(message%nreq,&
1495      message%mpi_req,message%completed,message%status,ierr)
1496!$OMP END MASTER
1497  END SUBROUTINE  test_message_mpi
1498 
1499   
1500  SUBROUTINE wait_message_mpi(message)
1501  USE field_mod
1502  USE domain_mod
1503  USE mpi_mod
1504  USE mpipara
1505  USE omp_para
1506  USE trace
1507  IMPLICIT NONE
1508    TYPE(t_message) :: message
1509
1510    TYPE(t_field),POINTER :: field(:)
1511    REAL(rstd),POINTER :: rval2d(:) 
1512    REAL(rstd),POINTER :: rval3d(:,:) 
1513    REAL(rstd),POINTER :: rval4d(:,:,:) 
1514    REAL(rstd),POINTER :: buffer_r(:) 
1515    INTEGER,POINTER :: value(:) 
1516    INTEGER,POINTER :: sgn(:) 
1517    TYPE(ARRAY),POINTER :: recv,send 
1518    TYPE(t_request),POINTER :: req
1519    INTEGER, ALLOCATABLE :: mpi_req(:)
1520    INTEGER, ALLOCATABLE :: status(:,:)
1521    INTEGER :: irecv,isend
1522    INTEGER :: ireq,nreq
1523    INTEGER :: ind,n,l,m,i
1524    INTEGER :: dim3,dim4,d3,d4,lbegin,lend
1525    INTEGER :: offset
1526
1527    message%open=.FALSE.
1528    IF (.NOT. message%pending) RETURN
1529
1530!    CALL trace_start("wait_message_mpi")
1531
1532    field=>message%field
1533    nreq=message%nreq
1534   
1535    IF (field(1)%data_type==type_real) THEN
1536      IF (field(1)%ndim==2) THEN
1537
1538!$OMP MASTER
1539        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1540          message%status,ierr)
1541!$OMP END MASTER
1542!$OMP BARRIER
1543       
1544        DO ind=1,ndomain
1545          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1546         
1547          rval2d=>field(ind)%rval2d
1548          req=>message%request(ind)
1549          DO irecv=1,req%nrecv
1550            recv=>req%recv(irecv)
1551            IF (recv%rank/=mpi_rank) THEN
1552              ireq=recv%ireq
1553              buffer_r=>message%buffers(ireq)%r
1554              value=>recv%value
1555              sgn=>recv%sign
1556              offset=recv%offset
1557              DO n=1,recv%size
1558                rval2d(value(n))=buffer_r(n+offset)*sgn(n) 
1559              ENDDO
1560
1561            ENDIF
1562          ENDDO
1563       
1564        ENDDO
1565     
1566     
1567      ELSE  IF (field(1)%ndim==3) THEN
1568
1569!$OMP MASTER
1570        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1571          message%status,ierr)
1572!$OMP END MASTER
1573!$OMP BARRIER
1574
1575       
1576        DO ind=1,ndomain
1577          IF (.NOT. assigned_domain(ind) ) CYCLE
1578
1579          rval3d=>field(ind)%rval3d
1580          req=>message%request(ind)
1581          DO irecv=1,req%nrecv
1582            recv=>req%recv(irecv)
1583            IF (recv%rank/=mpi_rank) THEN
1584              ireq=recv%ireq
1585              buffer_r=>message%buffers(ireq)%r
1586              value=>recv%value
1587              sgn=>recv%sign
1588             
1589              dim3=size(rval3d,2)
1590   
1591              CALL distrib_level(1,dim3, lbegin,lend)
1592              offset=recv%offset*dim3 + (lbegin-1)*recv%size
1593              CALL trace_start("copy_from_buffer")
1594             
1595              IF (req%vector) THEN
1596                DO d3=lbegin,lend
1597                  DO n=1,recv%size
1598                    rval3d(value(n),d3)=buffer_r(n+offset)*sgn(n) 
1599                  ENDDO
1600                  offset=offset+recv%size
1601                ENDDO
1602              ELSE
1603                DO d3=lbegin,lend
1604                  DO n=1,recv%size
1605                    rval3d(value(n),d3)=buffer_r(n+offset) 
1606                  ENDDO
1607                  offset=offset+recv%size
1608                ENDDO
1609              ENDIF
1610               
1611              CALL trace_end("copy_from_buffer")
1612            ENDIF
1613          ENDDO
1614       
1615        ENDDO
1616
1617      ELSE  IF (field(1)%ndim==4) THEN
1618!$OMP MASTER
1619        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1620          message%status,ierr)
1621!$OMP END MASTER
1622!$OMP BARRIER
1623
1624               
1625        DO ind=1,ndomain
1626          IF (.NOT. assigned_domain(ind) ) CYCLE
1627
1628          rval4d=>field(ind)%rval4d
1629          req=>message%request(ind)
1630          DO irecv=1,req%nrecv
1631            recv=>req%recv(irecv)
1632            IF (recv%rank/=mpi_rank) THEN
1633              ireq=recv%ireq
1634              buffer_r=>message%buffers(ireq)%r
1635              value=>recv%value
1636              sgn=>recv%sign
1637
1638              dim3=size(rval4d,2)
1639              CALL distrib_level(1,dim3, lbegin,lend)
1640              dim4=size(rval4d,3)
1641              CALL trace_start("copy_from_buffer")
1642              DO d4=1,dim4
1643                offset=recv%offset*dim3*dim4 + dim3*recv%size*(d4-1) +               &
1644                  (lbegin-1)*recv%size
1645                DO d3=lbegin,lend
1646                  DO n=1,recv%size
1647                    rval4d(value(n),d3,d4)=buffer_r(n+offset)*sgn(n) 
1648                  ENDDO
1649                  offset=offset+recv%size
1650                ENDDO
1651              ENDDO
1652              CALL trace_end("copy_from_buffer")
1653            ENDIF
1654          ENDDO
1655       
1656        ENDDO
1657     
1658      ENDIF     
1659     
1660    ENDIF
1661
1662!$OMP MASTER
1663    message%pending=.FALSE.
1664!$OMP END MASTER
1665
1666!    CALL trace_end("wait_message_mpi")
1667!$OMP BARRIER
1668   
1669  END SUBROUTINE wait_message_mpi
1670
1671  SUBROUTINE transfert_request_mpi(field,request)
1672  USE field_mod
1673  IMPLICIT NONE
1674    TYPE(t_field),POINTER :: field(:)
1675    TYPE(t_request),POINTER :: request(:)
1676
1677    TYPE(t_message),SAVE :: message
1678   
1679   
1680    CALL init_message_mpi(field,request, message)
1681    CALL transfert_message_mpi(field,message)
1682    CALL finalize_message_mpi(field,message)
1683   
1684  END SUBROUTINE transfert_request_mpi
1685 
1686   
1687   
1688  SUBROUTINE transfert_request_seq(field,request)
1689  USE field_mod
1690  USE domain_mod
1691  IMPLICIT NONE
1692    TYPE(t_field),POINTER :: field(:)
1693    TYPE(t_request),POINTER :: request(:)
1694    REAL(rstd),POINTER :: rval2d(:) 
1695    REAL(rstd),POINTER :: rval3d(:,:) 
1696    REAL(rstd),POINTER :: rval4d(:,:,:) 
1697    INTEGER :: ind
1698    TYPE(t_request),POINTER :: req
1699    INTEGER :: n
1700    REAL(rstd) :: var1,var2
1701   
1702    DO ind=1,ndomain
1703      req=>request(ind)
1704      rval2d=>field(ind)%rval2d
1705      rval3d=>field(ind)%rval3d
1706      rval4d=>field(ind)%rval4d
1707     
1708      IF (field(ind)%data_type==type_real) THEN
1709        IF (field(ind)%ndim==2) THEN
1710          DO n=1,req%size
1711            rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*&
1712              req%target_sign(n)
1713          ENDDO
1714        ELSE IF (field(ind)%ndim==3) THEN
1715          DO n=1,req%size
1716            rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*&
1717              req%target_sign(n)
1718          ENDDO
1719        ELSE IF (field(ind)%ndim==4) THEN
1720          DO n=1,req%size
1721            rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*&
1722              req%target_sign(n)
1723          ENDDO
1724        ENDIF
1725      ENDIF       
1726
1727    ENDDO
1728   
1729  END SUBROUTINE transfert_request_seq
1730 
1731 
1732  SUBROUTINE gather_field(field_loc,field_glo)
1733  USE field_mod
1734  USE domain_mod
1735  USE mpi_mod
1736  USE mpipara
1737  IMPLICIT NONE
1738    TYPE(t_field),POINTER :: field_loc(:)
1739    TYPE(t_field),POINTER :: field_glo(:)
1740    INTEGER, ALLOCATABLE :: mpi_req(:)
1741    INTEGER, ALLOCATABLE :: status(:,:)
1742    INTEGER :: ireq,nreq
1743    INTEGER :: ind_glo,ind_loc   
1744 
1745    IF (.NOT. using_mpi) THEN
1746   
1747      DO ind_loc=1,ndomain
1748        IF (field_loc(ind_loc)%ndim==2) field_glo(ind_loc)%rval2d=field_loc(ind_loc)%rval2d
1749        IF (field_loc(ind_loc)%ndim==3) field_glo(ind_loc)%rval3d=field_loc(ind_loc)%rval3d
1750        IF (field_loc(ind_loc)%ndim==4) field_glo(ind_loc)%rval4d=field_loc(ind_loc)%rval4d
1751      ENDDO
1752   
1753    ELSE
1754         
1755      nreq=ndomain
1756      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
1757      ALLOCATE(mpi_req(nreq))
1758      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
1759   
1760   
1761      ireq=0
1762      IF (mpi_rank==0) THEN
1763        DO ind_glo=1,ndomain_glo
1764          ireq=ireq+1
1765
1766          IF (field_glo(ind_glo)%ndim==2) THEN
1767            CALL MPI_IRECV(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
1768                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1769   
1770          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1771            CALL MPI_IRECV(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
1772                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1773
1774          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1775            CALL MPI_IRECV(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
1776                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1777          ENDIF
1778         
1779        ENDDO
1780      ENDIF
1781 
1782      DO ind_loc=1,ndomain
1783        ireq=ireq+1
1784
1785        IF (field_loc(ind_loc)%ndim==2) THEN
1786          CALL MPI_ISEND(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
1787                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1788        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
1789          CALL MPI_ISEND(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
1790                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1791        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
1792          CALL MPI_ISEND(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
1793                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1794        ENDIF
1795     
1796      ENDDO
1797   
1798      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
1799
1800    ENDIF
1801       
1802  END SUBROUTINE gather_field
1803
1804  SUBROUTINE bcast_field(field_glo)
1805  USE field_mod
1806  USE domain_mod
1807  USE mpi_mod
1808  USE mpipara
1809  IMPLICIT NONE
1810    TYPE(t_field),POINTER :: field_glo(:)
1811    INTEGER :: ind_glo   
1812 
1813    IF (.NOT. using_mpi) THEN
1814   
1815! nothing to do
1816   
1817    ELSE
1818         
1819      DO ind_glo=1,ndomain_glo
1820
1821          IF (field_glo(ind_glo)%ndim==2) THEN
1822            CALL MPI_BCAST(field_glo(ind_glo)%rval2d, size(field_glo(ind_glo)%rval2d) , MPI_REAL8, 0, comm_icosa, ierr)
1823          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1824            CALL MPI_BCAST(field_glo(ind_glo)%rval3d, size(field_glo(ind_glo)%rval3d) , MPI_REAL8, 0, comm_icosa, ierr)
1825          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1826            CALL MPI_BCAST(field_glo(ind_glo)%rval4d, size(field_glo(ind_glo)%rval4d) , MPI_REAL8, 0, comm_icosa, ierr)
1827          ENDIF
1828         
1829        ENDDO
1830      ENDIF
1831       
1832  END SUBROUTINE bcast_field
1833
1834  SUBROUTINE scatter_field(field_glo,field_loc)
1835  USE field_mod
1836  USE domain_mod
1837  USE mpi_mod
1838  USE mpipara
1839  IMPLICIT NONE
1840    TYPE(t_field),POINTER :: field_glo(:)
1841    TYPE(t_field),POINTER :: field_loc(:)
1842    INTEGER, ALLOCATABLE :: mpi_req(:)
1843    INTEGER, ALLOCATABLE :: status(:,:)
1844    INTEGER :: ireq,nreq
1845    INTEGER :: ind_glo,ind_loc   
1846 
1847    IF (.NOT. using_mpi) THEN
1848   
1849      DO ind_loc=1,ndomain
1850        IF (field_loc(ind_loc)%ndim==2) field_loc(ind_loc)%rval2d=field_glo(ind_loc)%rval2d
1851        IF (field_loc(ind_loc)%ndim==3) field_loc(ind_loc)%rval3d=field_glo(ind_loc)%rval3d
1852        IF (field_loc(ind_loc)%ndim==4) field_loc(ind_loc)%rval4d=field_glo(ind_loc)%rval4d
1853      ENDDO
1854   
1855    ELSE
1856         
1857      nreq=ndomain
1858      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
1859      ALLOCATE(mpi_req(nreq))
1860      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
1861   
1862   
1863      ireq=0
1864      IF (mpi_rank==0) THEN
1865        DO ind_glo=1,ndomain_glo
1866          ireq=ireq+1
1867
1868          IF (field_glo(ind_glo)%ndim==2) THEN
1869            CALL MPI_ISEND(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
1870                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1871   
1872          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1873            CALL MPI_ISEND(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
1874                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1875
1876          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1877            CALL MPI_ISEND(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
1878                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1879          ENDIF
1880         
1881        ENDDO
1882      ENDIF
1883 
1884      DO ind_loc=1,ndomain
1885        ireq=ireq+1
1886
1887        IF (field_loc(ind_loc)%ndim==2) THEN
1888          CALL MPI_IRECV(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
1889                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1890        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
1891          CALL MPI_IRECV(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
1892                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1893        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
1894          CALL MPI_IRECV(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
1895                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1896        ENDIF
1897     
1898      ENDDO
1899   
1900      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
1901
1902    ENDIF
1903       
1904  END SUBROUTINE scatter_field
1905 
1906  SUBROUTINE trace_in
1907  USE trace
1908  IMPLICIT NONE
1909 
1910    CALL trace_start("transfert_buffer")
1911  END SUBROUTINE trace_in             
1912
1913  SUBROUTINE trace_out
1914  USE trace
1915  IMPLICIT NONE
1916 
1917    CALL trace_end("transfert_buffer")
1918  END SUBROUTINE trace_out             
1919
1920
1921
1922
1923!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1924!! Definition des Broadcast --> 4D   !!
1925!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1926
1927!! -- Les chaine de charactï¿œre -- !!
1928
1929  SUBROUTINE bcast_mpi_c(var1)
1930  IMPLICIT NONE
1931    CHARACTER(LEN=*),INTENT(INOUT) :: Var1
1932   
1933    CALL bcast_mpi_cgen(Var1,len(Var1))
1934
1935  END SUBROUTINE bcast_mpi_c
1936
1937!! -- Les entiers -- !!
1938 
1939  SUBROUTINE bcast_mpi_i(var)
1940  USE mpipara
1941  IMPLICIT NONE
1942    INTEGER,INTENT(INOUT) :: Var
1943   
1944    INTEGER               :: var_tmp(1)
1945   
1946    IF (is_mpi_master) var_tmp(1)=var
1947    CALL bcast_mpi_igen(Var_tmp,1)
1948    var=var_tmp(1)
1949   
1950  END SUBROUTINE bcast_mpi_i
1951
1952  SUBROUTINE bcast_mpi_i1(var)
1953  IMPLICIT NONE
1954    INTEGER,INTENT(INOUT) :: Var(:)
1955
1956    CALL bcast_mpi_igen(Var,size(Var))
1957   
1958  END SUBROUTINE bcast_mpi_i1
1959
1960  SUBROUTINE bcast_mpi_i2(var)
1961  IMPLICIT NONE
1962    INTEGER,INTENT(INOUT) :: Var(:,:)
1963   
1964    CALL bcast_mpi_igen(Var,size(Var))
1965 
1966  END SUBROUTINE bcast_mpi_i2
1967
1968  SUBROUTINE bcast_mpi_i3(var)
1969  IMPLICIT NONE
1970    INTEGER,INTENT(INOUT) :: Var(:,:,:)
1971   
1972    CALL bcast_mpi_igen(Var,size(Var))
1973
1974  END SUBROUTINE bcast_mpi_i3
1975
1976  SUBROUTINE bcast_mpi_i4(var)
1977  IMPLICIT NONE
1978    INTEGER,INTENT(INOUT) :: Var(:,:,:,:)
1979   
1980    CALL bcast_mpi_igen(Var,size(Var))
1981
1982  END SUBROUTINE bcast_mpi_i4
1983
1984
1985!! -- Les reels -- !!
1986
1987  SUBROUTINE bcast_mpi_r(var)
1988  USE mpipara
1989  IMPLICIT NONE
1990    REAL,INTENT(INOUT) :: Var
1991    REAL               :: var_tmp(1)
1992   
1993    IF (is_mpi_master) var_tmp(1)=var
1994    CALL bcast_mpi_rgen(Var_tmp,1)
1995    var=var_tmp(1)   
1996
1997  END SUBROUTINE bcast_mpi_r
1998
1999  SUBROUTINE bcast_mpi_r1(var)
2000  IMPLICIT NONE
2001    REAL,INTENT(INOUT) :: Var(:)
2002   
2003    CALL bcast_mpi_rgen(Var,size(Var))
2004
2005  END SUBROUTINE bcast_mpi_r1
2006
2007  SUBROUTINE bcast_mpi_r2(var)
2008  IMPLICIT NONE
2009    REAL,INTENT(INOUT) :: Var(:,:)
2010   
2011    CALL bcast_mpi_rgen(Var,size(Var))
2012
2013  END SUBROUTINE bcast_mpi_r2
2014
2015  SUBROUTINE bcast_mpi_r3(var)
2016  IMPLICIT NONE
2017    REAL,INTENT(INOUT) :: Var(:,:,:)
2018   
2019    CALL bcast_mpi_rgen(Var,size(Var))
2020
2021  END SUBROUTINE bcast_mpi_r3
2022
2023  SUBROUTINE bcast_mpi_r4(var)
2024  IMPLICIT NONE
2025    REAL,INTENT(INOUT) :: Var(:,:,:,:)
2026   
2027    CALL bcast_mpi_rgen(Var,size(Var))
2028
2029  END SUBROUTINE bcast_mpi_r4
2030 
2031!! -- Les booleans -- !!
2032
2033  SUBROUTINE bcast_mpi_l(var)
2034  USE mpipara
2035  IMPLICIT NONE
2036    LOGICAL,INTENT(INOUT) :: Var
2037    LOGICAL               :: var_tmp(1)
2038   
2039    IF (is_mpi_master) var_tmp(1)=var
2040    CALL bcast_mpi_lgen(Var_tmp,1)
2041    var=var_tmp(1)   
2042
2043  END SUBROUTINE bcast_mpi_l
2044
2045  SUBROUTINE bcast_mpi_l1(var)
2046  IMPLICIT NONE
2047    LOGICAL,INTENT(INOUT) :: Var(:)
2048   
2049    CALL bcast_mpi_lgen(Var,size(Var))
2050
2051  END SUBROUTINE bcast_mpi_l1
2052
2053  SUBROUTINE bcast_mpi_l2(var)
2054  IMPLICIT NONE
2055    LOGICAL,INTENT(INOUT) :: Var(:,:)
2056   
2057    CALL bcast_mpi_lgen(Var,size(Var))
2058
2059  END SUBROUTINE bcast_mpi_l2
2060
2061  SUBROUTINE bcast_mpi_l3(var)
2062  IMPLICIT NONE
2063    LOGICAL,INTENT(INOUT) :: Var(:,:,:)
2064   
2065    CALL bcast_mpi_lgen(Var,size(Var))
2066
2067  END SUBROUTINE bcast_mpi_l3
2068
2069  SUBROUTINE bcast_mpi_l4(var)
2070  IMPLICIT NONE
2071    LOGICAL,INTENT(INOUT) :: Var(:,:,:,:)
2072   
2073    CALL bcast_mpi_lgen(Var,size(Var))
2074
2075  END SUBROUTINE bcast_mpi_l4
2076 
2077!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2078!! DEFINITION DES FONCTIONS DE TRANSFERT GENERIQUES !
2079!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2080
2081  SUBROUTINE bcast_mpi_cgen(var,nb)
2082    USE mpi_mod
2083    USE mpipara
2084    IMPLICIT NONE
2085   
2086    CHARACTER(LEN=*),INTENT(INOUT) :: Var
2087    INTEGER,INTENT(IN) :: nb
2088
2089    IF (.NOT. using_mpi) RETURN
2090   
2091    CALL MPI_BCAST(Var,nb,MPI_CHARACTER,mpi_master,comm_icosa,ierr)
2092       
2093  END SUBROUTINE bcast_mpi_cgen
2094
2095
2096     
2097  SUBROUTINE bcast_mpi_igen(var,nb)
2098    USE mpi_mod
2099    USE mpipara
2100    IMPLICIT NONE
2101   
2102    INTEGER,DIMENSION(nb),INTENT(INOUT) :: Var
2103    INTEGER,INTENT(IN) :: nb
2104
2105    IF (.NOT. using_mpi) RETURN
2106
2107    CALL MPI_BCAST(Var,nb,MPI_INTEGER,mpi_master,comm_icosa,ierr)
2108       
2109  END SUBROUTINE bcast_mpi_igen
2110
2111
2112
2113 
2114  SUBROUTINE bcast_mpi_rgen(var,nb)
2115    USE mpi_mod
2116    USE mpipara
2117    IMPLICIT NONE
2118   
2119    REAL,DIMENSION(nb),INTENT(INOUT) :: Var
2120    INTEGER,INTENT(IN) :: nb
2121
2122    IF (.NOT. using_mpi) RETURN
2123
2124    CALL MPI_BCAST(Var,nb,MPI_REAL8,mpi_master,comm_icosa,ierr)
2125   
2126  END SUBROUTINE bcast_mpi_rgen
2127 
2128
2129
2130
2131  SUBROUTINE bcast_mpi_lgen(var,nb)
2132    USE mpi_mod
2133    USE mpipara
2134    IMPLICIT NONE
2135   
2136    LOGICAL,DIMENSION(nb),INTENT(INOUT) :: Var
2137    INTEGER,INTENT(IN) :: nb
2138
2139    IF (.NOT. using_mpi) RETURN
2140
2141    CALL MPI_BCAST(Var,nb,MPI_LOGICAL,mpi_master,comm_icosa,ierr)
2142
2143  END SUBROUTINE bcast_mpi_lgen
2144 
2145   
2146END MODULE transfert_mpi_mod
2147     
2148       
2149       
2150       
2151     
Note: See TracBrowser for help on using the repository browser.