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

Last change on this file since 846 was 846, checked in by dubos, 5 years ago

devel : towards Fortran driver for unstructured/LAM meshes

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