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

Last change on this file since 352 was 327, checked in by ymipsl, 9 years ago

Merge recent developments from saturn branch onto trunk.

  • lmdz generic physics interface
  • performance improvment on mix mpi/openmp
  • asynchrone and overlaping communication
  • best domain distribution between process and threads
  • ....

This version is compatible with the actual saturn version and the both branches are considered merged on dynamico component.

YM

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