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

Last change on this file since 714 was 714, checked in by dubos, 6 years ago

devel : backported from trunk commits r607,r648,r649,r667,r668,r669,r706

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