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

Last change on this file was 1026, checked in by dubos, 4 years ago

devel : towards conformity to F2008 standard

File size: 62.6 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  USE grid_param
1145  IMPLICIT NONE
1146    TYPE(t_field),POINTER :: field(:)
1147    TYPE(t_message) :: message
1148    REAL(rstd),POINTER :: rval2d(:), src_rval2d(:) 
1149    REAL(rstd),POINTER :: rval3d(:,:), src_rval3d(:,:) 
1150    REAL(rstd),POINTER :: rval4d(:,:,:), src_rval4d(:,:,:) 
1151    REAL(rstd),POINTER :: buffer_r(:) 
1152    INTEGER,POINTER :: value(:) 
1153    INTEGER,POINTER :: sgn(:) 
1154    TYPE(ARRAY),POINTER :: recv,send 
1155    TYPE(t_request),POINTER :: req
1156    INTEGER, ALLOCATABLE :: mpi_req(:)
1157    INTEGER, ALLOCATABLE :: status(:,:)
1158    INTEGER :: irecv,isend
1159    INTEGER :: ireq,nreq
1160    INTEGER :: ind,i,n,l,m
1161    INTEGER :: dim3,dim4,d3,d4
1162    INTEGER,POINTER :: src_value(:)
1163    INTEGER,POINTER :: sign(:)
1164    INTEGER :: offset,msize,rank
1165    INTEGER :: lbegin, lend
1166    INTEGER :: max_req
1167
1168!    CALL trace_start("send_message_mpi")
1169
1170    IF(grid_type == grid_unst) RETURN
1171
1172    CALL enter_profile(id_mpi)
1173
1174!$OMP BARRIER
1175
1176
1177!$OMP MASTER
1178    IF(message%open) THEN
1179       PRINT *, 'send_message_mpi : message ' // TRIM(message%name) // &
1180            ' is still open, no call to wait_message_mpi after last send_message_mpi'
1181       CALL ABORT
1182    END IF
1183    message%open=.TRUE. ! will be set to .FALSE. by wait_message_mpi
1184
1185    message%field=>field
1186
1187    IF (message%nreq>0) THEN
1188      message%completed=.FALSE.
1189      message%pending=.TRUE.
1190    ELSE
1191      message%completed=.TRUE.
1192      message%pending=.FALSE.
1193    ENDIF
1194!$OMP END MASTER
1195!$OMP BARRIER
1196     
1197    IF (field(1)%data_type==type_real) THEN
1198      IF (field(1)%ndim==2) THEN
1199
1200        DO ind=1,ndomain
1201          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1202         
1203          rval2d=>field(ind)%rval2d
1204       
1205          req=>message%request(ind)
1206          DO isend=1,req%nsend
1207            send=>req%send(isend)
1208            value=>send%value
1209
1210           
1211            IF (send%rank/=mpi_rank) THEN
1212              ireq=send%ireq
1213
1214              buffer_r=>message%buffers(ireq)%r
1215              offset=send%offset
1216                                           
1217              DO n=1,send%size
1218                buffer_r(offset+n)=rval2d(value(n))
1219              ENDDO
1220
1221              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1222                !$OMP CRITICAL           
1223                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1224                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1225                !$OMP END CRITICAL
1226              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1227                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               &
1228                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1229              ENDIF
1230             
1231            ENDIF
1232          ENDDO
1233        ENDDO
1234       
1235        DO ind=1,ndomain
1236          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1237          rval2d=>field(ind)%rval2d
1238          req=>message%request(ind)       
1239
1240          DO irecv=1,req%nrecv
1241            recv=>req%recv(irecv)
1242
1243            IF (recv%rank==mpi_rank) THEN
1244
1245              value=>recv%value
1246              src_value => recv%src_value
1247              src_rval2d=>field(recv%domain)%rval2d
1248              sgn=>recv%sign
1249
1250              DO n=1,recv%size
1251                rval2d(value(n))=src_rval2d(src_value(n))*sgn(n)
1252              ENDDO
1253               
1254                   
1255            ELSE
1256           
1257              ireq=recv%ireq
1258              buffer_r=>message%buffers(ireq)%r
1259              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1260               !$OMP CRITICAL           
1261                CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,               &
1262                  recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1263               !$OMP END CRITICAL
1264              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1265                 CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,              &
1266                   recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1267              ENDIF
1268           
1269            ENDIF
1270          ENDDO
1271       
1272        ENDDO
1273       
1274     
1275      ELSE  IF (field(1)%ndim==3) THEN
1276        max_req=0
1277        DO ind=1,ndomain
1278          req=>message%request(ind)
1279          IF (req%nsend>max_req) max_req=req%nsend
1280        ENDDO
1281             
1282        DO ind=1,ndomain
1283          IF (.NOT. assigned_domain(ind) ) CYCLE
1284
1285          dim3=size(field(ind)%rval3d,2)
1286          CALL distrib_level(1,dim3, lbegin,lend)
1287
1288          rval3d=>field(ind)%rval3d
1289          req=>message%request(ind)
1290 
1291          DO isend=1,req%nsend
1292            send=>req%send(isend)
1293            value=>send%value
1294
1295            IF (send%rank/=mpi_rank) THEN
1296              ireq=send%ireq
1297              buffer_r=>message%buffers(ireq)%r
1298
1299              offset=send%offset*dim3 + (lbegin-1)*send%size
1300         
1301              CALL trace_start("copy_to_buffer")
1302
1303              DO d3=lbegin,lend
1304                DO n=1,send%size
1305                  buffer_r(n+offset)=rval3d(value(n),d3)
1306                ENDDO
1307                offset=offset+send%size
1308              ENDDO
1309              CALL trace_end("copy_to_buffer")
1310
1311              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1312                !$OMP BARRIER
1313              ENDIF
1314             
1315              IF (is_omp_level_master) THEN
1316                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1317                  !$OMP CRITICAL   
1318                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1319                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1320                  !$OMP END CRITICAL
1321                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1322                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        &
1323                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1324                ENDIF
1325              ENDIF
1326            ELSE
1327              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1328                !$OMP BARRIER
1329              ENDIF
1330            ENDIF
1331          ENDDO
1332
1333          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1334            DO isend=req%nsend+1,max_req
1335              !$OMP BARRIER
1336            ENDDO
1337          ENDIF
1338
1339        ENDDO
1340         
1341        DO ind=1,ndomain
1342          IF (.NOT. assigned_domain(ind) ) CYCLE
1343          dim3=size(field(ind)%rval3d,2)
1344          CALL distrib_level(1,dim3, lbegin,lend)
1345          rval3d=>field(ind)%rval3d
1346          req=>message%request(ind)
1347
1348          DO irecv=1,req%nrecv
1349            recv=>req%recv(irecv)
1350
1351            IF (recv%rank==mpi_rank) THEN
1352              value=>recv%value
1353              src_value => recv%src_value
1354              src_rval3d=>field(recv%domain)%rval3d
1355              sgn=>recv%sign
1356
1357              CALL trace_start("copy_data")
1358
1359              DO d3=lbegin,lend
1360                DO n=1,recv%size
1361                  rval3d(value(n),d3)=src_rval3d(src_value(n),d3)*sgn(n)
1362                ENDDO
1363              ENDDO
1364              CALL trace_end("copy_data")
1365
1366            ELSE
1367              ireq=recv%ireq
1368              buffer_r=>message%buffers(ireq)%r
1369 
1370              IF (is_omp_level_master) THEN
1371                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1372                  !$OMP CRITICAL
1373                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1374                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1375                  !$OMP END CRITICAL
1376                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1377                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        &
1378                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1379                ENDIF
1380              ENDIF
1381            ENDIF 
1382          ENDDO
1383       
1384        ENDDO
1385
1386      ELSE  IF (field(1)%ndim==4) THEN
1387
1388        max_req=0
1389        DO ind=1,ndomain
1390          req=>message%request(ind)
1391          IF (req%nsend>max_req) max_req=req%nsend
1392        ENDDO
1393   
1394        DO ind=1,ndomain
1395          IF (.NOT. assigned_domain(ind) ) CYCLE
1396
1397          dim3=size(field(ind)%rval4d,2)
1398          CALL distrib_level(1,dim3, lbegin,lend)
1399          dim4=size(field(ind)%rval4d,3)
1400          rval4d=>field(ind)%rval4d
1401          req=>message%request(ind)
1402
1403          DO isend=1,req%nsend
1404            send=>req%send(isend)
1405            value=>send%value
1406
1407            IF (send%rank/=mpi_rank) THEN
1408
1409              ireq=send%ireq
1410              buffer_r=>message%buffers(ireq)%r
1411
1412              CALL trace_start("copy_to_buffer")
1413              DO d4=1,dim4
1414                offset=send%offset*dim3*dim4 + dim3*send%size*(d4-1) +               &
1415                  (lbegin-1)*send%size
1416
1417                DO d3=lbegin,lend
1418                  DO n=1,send%size
1419                    buffer_r(n+offset)=rval4d(value(n),d3,d4)
1420                  ENDDO
1421                  offset=offset+send%size
1422                ENDDO
1423              ENDDO
1424              CALL trace_end("copy_to_buffer")
1425
1426              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1427                !$OMP BARRIER
1428              ENDIF
1429
1430              IF (is_omp_level_master) THEN
1431                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1432                  !$OMP CRITICAL
1433                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1434                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1435                  !$OMP END CRITICAL
1436                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1437                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   &
1438                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)
1439                ENDIF
1440              ENDIF
1441            ELSE
1442              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1443                !$OMP BARRIER
1444              ENDIF
1445            ENDIF
1446          ENDDO
1447         
1448          IF (mpi_threading_mode==MPI_THREAD_SERIALIZED.OR.mpi_threading_mode==MPI_THREAD_MULTIPLE .AND. omp_level_size>1) THEN
1449            DO isend=req%nsend+1,max_req
1450              !$OMP BARRIER
1451            ENDDO
1452          ENDIF
1453
1454        ENDDO
1455       
1456        DO ind=1,ndomain
1457          IF (.NOT. assigned_domain(ind) ) CYCLE
1458         
1459          dim3=size(field(ind)%rval4d,2)
1460          CALL distrib_level(1,dim3, lbegin,lend)
1461          dim4=size(field(ind)%rval4d,3)
1462          rval4d=>field(ind)%rval4d
1463          req=>message%request(ind)
1464
1465          DO irecv=1,req%nrecv
1466            recv=>req%recv(irecv)
1467            IF (recv%rank==mpi_rank) THEN
1468              value=>recv%value
1469              src_value => recv%src_value
1470              src_rval4d=>field(recv%domain)%rval4d
1471              sgn=>recv%sign
1472
1473              CALL trace_start("copy_data")
1474              DO d4=1,dim4
1475                DO d3=lbegin,lend
1476                  DO n=1,recv%size
1477                    rval4d(value(n),d3,d4)=src_rval4d(src_value(n),d3,d4)*sgn(n)
1478                  ENDDO
1479                ENDDO
1480              ENDDO
1481               
1482              CALL trace_end("copy_data")
1483                   
1484            ELSE
1485
1486              ireq=recv%ireq
1487              buffer_r=>message%buffers(ireq)%r
1488              IF (is_omp_level_master) THEN
1489                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN
1490                 !$OMP CRITICAL           
1491                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1492                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1493                  !$OMP END CRITICAL
1494                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN
1495                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   &
1496                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
1497                ENDIF
1498              ENDIF
1499            ENDIF
1500          ENDDO
1501        ENDDO
1502
1503      ENDIF     
1504
1505      IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN
1506!$OMP BARRIER
1507!$OMP MASTER       
1508
1509        DO ireq=1,message%nreq_send
1510          buffer_r=>message%buffers(ireq)%r
1511          msize=message%buffers(ireq)%size
1512          rank=message%buffers(ireq)%rank
1513          CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1514            comm_icosa, message%mpi_req(ireq),ierr)
1515        ENDDO
1516
1517        DO ireq=message%nreq_send+1,message%nreq_send+message%nreq_recv
1518          buffer_r=>message%buffers(ireq)%r
1519          msize=message%buffers(ireq)%size
1520          rank=message%buffers(ireq)%rank
1521          CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          &
1522            comm_icosa, message%mpi_req(ireq),ierr)
1523        ENDDO
1524
1525!$OMP END MASTER
1526      ENDIF             
1527    ENDIF
1528   
1529!$OMP BARRIER
1530!    CALL trace_end("send_message_mpi")
1531
1532    CALL exit_profile(id_mpi)
1533   
1534  END SUBROUTINE send_message_mpi
1535 
1536  SUBROUTINE test_message_mpi(message)
1537  IMPLICIT NONE
1538    TYPE(t_message) :: message
1539   
1540    INTEGER :: ierr
1541
1542!$OMP MASTER
1543    IF (message%pending .AND. .NOT. message%completed) CALL MPI_TESTALL(message%nreq,&
1544      message%mpi_req,message%completed,message%status,ierr)
1545!$OMP END MASTER
1546  END SUBROUTINE  test_message_mpi
1547 
1548   
1549  SUBROUTINE wait_message_mpi(message)
1550  USE profiling_mod
1551  USE field_mod
1552  USE domain_mod
1553  USE mpi_mod
1554  USE mpipara
1555  USE omp_para
1556  USE trace
1557  USE grid_param
1558  IMPLICIT NONE
1559    TYPE(t_message) :: message
1560
1561    TYPE(t_field),POINTER :: field(:)
1562    REAL(rstd),POINTER :: rval2d(:) 
1563    REAL(rstd),POINTER :: rval3d(:,:) 
1564    REAL(rstd),POINTER :: rval4d(:,:,:) 
1565    REAL(rstd),POINTER :: buffer_r(:) 
1566    INTEGER,POINTER :: value(:) 
1567    INTEGER,POINTER :: sgn(:) 
1568    TYPE(ARRAY),POINTER :: recv,send 
1569    TYPE(t_request),POINTER :: req
1570    INTEGER, ALLOCATABLE :: mpi_req(:)
1571    INTEGER, ALLOCATABLE :: status(:,:)
1572    INTEGER :: irecv,isend
1573    INTEGER :: ireq,nreq
1574    INTEGER :: ind,n,l,m,i
1575    INTEGER :: dim3,dim4,d3,d4,lbegin,lend
1576    INTEGER :: offset
1577
1578    IF(grid_type == grid_unst) RETURN
1579
1580    message%open=.FALSE.
1581    IF (.NOT. message%pending) RETURN
1582
1583    CALL enter_profile(id_mpi)
1584
1585!    CALL trace_start("wait_message_mpi")
1586
1587    field=>message%field
1588    nreq=message%nreq
1589   
1590    IF (field(1)%data_type==type_real) THEN
1591      IF (field(1)%ndim==2) THEN
1592
1593!$OMP MASTER
1594        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1595          message%status,ierr)
1596!$OMP END MASTER
1597!$OMP BARRIER
1598       
1599        DO ind=1,ndomain
1600          IF (.NOT. assigned_domain(ind) .OR. .NOT. is_omp_level_master) CYCLE
1601         
1602          rval2d=>field(ind)%rval2d
1603          req=>message%request(ind)
1604          DO irecv=1,req%nrecv
1605            recv=>req%recv(irecv)
1606            IF (recv%rank/=mpi_rank) THEN
1607              ireq=recv%ireq
1608              buffer_r=>message%buffers(ireq)%r
1609              value=>recv%value
1610              sgn=>recv%sign
1611              offset=recv%offset
1612              DO n=1,recv%size
1613                rval2d(value(n))=buffer_r(n+offset)*sgn(n) 
1614              ENDDO
1615
1616            ENDIF
1617          ENDDO
1618       
1619        ENDDO
1620     
1621     
1622      ELSE  IF (field(1)%ndim==3) THEN
1623
1624!$OMP MASTER
1625        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1626          message%status,ierr)
1627!$OMP END MASTER
1628!$OMP BARRIER
1629
1630       
1631        DO ind=1,ndomain
1632          IF (.NOT. assigned_domain(ind) ) CYCLE
1633
1634          rval3d=>field(ind)%rval3d
1635          req=>message%request(ind)
1636          DO irecv=1,req%nrecv
1637            recv=>req%recv(irecv)
1638            IF (recv%rank/=mpi_rank) THEN
1639              ireq=recv%ireq
1640              buffer_r=>message%buffers(ireq)%r
1641              value=>recv%value
1642              sgn=>recv%sign
1643             
1644              dim3=size(rval3d,2)
1645   
1646              CALL distrib_level(1,dim3, lbegin,lend)
1647              offset=recv%offset*dim3 + (lbegin-1)*recv%size
1648              CALL trace_start("copy_from_buffer")
1649             
1650              IF (req%vector) THEN
1651                DO d3=lbegin,lend
1652                  DO n=1,recv%size
1653                    rval3d(value(n),d3)=buffer_r(n+offset)*sgn(n) 
1654                  ENDDO
1655                  offset=offset+recv%size
1656                ENDDO
1657              ELSE
1658                DO d3=lbegin,lend
1659                  DO n=1,recv%size
1660                    rval3d(value(n),d3)=buffer_r(n+offset) 
1661                  ENDDO
1662                  offset=offset+recv%size
1663                ENDDO
1664              ENDIF
1665               
1666              CALL trace_end("copy_from_buffer")
1667            ENDIF
1668          ENDDO
1669       
1670        ENDDO
1671
1672      ELSE  IF (field(1)%ndim==4) THEN
1673!$OMP MASTER
1674        IF (.NOT. message%completed) CALL MPI_WAITALL(nreq,message%mpi_req,          &
1675          message%status,ierr)
1676!$OMP END MASTER
1677!$OMP BARRIER
1678
1679               
1680        DO ind=1,ndomain
1681          IF (.NOT. assigned_domain(ind) ) CYCLE
1682
1683          rval4d=>field(ind)%rval4d
1684          req=>message%request(ind)
1685          DO irecv=1,req%nrecv
1686            recv=>req%recv(irecv)
1687            IF (recv%rank/=mpi_rank) THEN
1688              ireq=recv%ireq
1689              buffer_r=>message%buffers(ireq)%r
1690              value=>recv%value
1691              sgn=>recv%sign
1692
1693              dim3=size(rval4d,2)
1694              CALL distrib_level(1,dim3, lbegin,lend)
1695              dim4=size(rval4d,3)
1696              CALL trace_start("copy_from_buffer")
1697              DO d4=1,dim4
1698                offset=recv%offset*dim3*dim4 + dim3*recv%size*(d4-1) +               &
1699                  (lbegin-1)*recv%size
1700                DO d3=lbegin,lend
1701                  DO n=1,recv%size
1702                    rval4d(value(n),d3,d4)=buffer_r(n+offset)*sgn(n) 
1703                  ENDDO
1704                  offset=offset+recv%size
1705                ENDDO
1706              ENDDO
1707              CALL trace_end("copy_from_buffer")
1708            ENDIF
1709          ENDDO
1710       
1711        ENDDO
1712     
1713      ENDIF     
1714     
1715    ENDIF
1716
1717!$OMP MASTER
1718    message%pending=.FALSE.
1719!$OMP END MASTER
1720
1721!    CALL trace_end("wait_message_mpi")
1722!$OMP BARRIER
1723   
1724    CALL exit_profile(id_mpi)
1725
1726  END SUBROUTINE wait_message_mpi
1727
1728  SUBROUTINE transfert_request_mpi(field,request)
1729  USE field_mod
1730  USE grid_param
1731  IMPLICIT NONE
1732    TYPE(t_field),POINTER :: field(:)
1733    TYPE(t_request),POINTER :: request(:)
1734
1735    TYPE(t_message),SAVE :: message
1736   
1737    IF(grid_type == grid_unst) RETURN
1738
1739    CALL init_message_mpi(field,request, message)
1740    CALL transfert_message_mpi(field,message)
1741    CALL finalize_message_mpi(field,message)
1742   
1743  END SUBROUTINE transfert_request_mpi
1744 
1745   
1746   
1747  SUBROUTINE transfert_request_seq(field,request)
1748  USE field_mod
1749  USE domain_mod
1750  IMPLICIT NONE
1751    TYPE(t_field),POINTER :: field(:)
1752    TYPE(t_request),POINTER :: request(:)
1753    REAL(rstd),POINTER :: rval2d(:) 
1754    REAL(rstd),POINTER :: rval3d(:,:) 
1755    REAL(rstd),POINTER :: rval4d(:,:,:) 
1756    INTEGER :: ind
1757    TYPE(t_request),POINTER :: req
1758    INTEGER :: n
1759    REAL(rstd) :: var1,var2
1760   
1761    DO ind=1,ndomain
1762      req=>request(ind)
1763      rval2d=>field(ind)%rval2d
1764      rval3d=>field(ind)%rval3d
1765      rval4d=>field(ind)%rval4d
1766     
1767      IF (field(ind)%data_type==type_real) THEN
1768        IF (field(ind)%ndim==2) THEN
1769          DO n=1,req%size
1770            rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*&
1771              req%target_sign(n)
1772          ENDDO
1773        ELSE IF (field(ind)%ndim==3) THEN
1774          DO n=1,req%size
1775            rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*&
1776              req%target_sign(n)
1777          ENDDO
1778        ELSE IF (field(ind)%ndim==4) THEN
1779          DO n=1,req%size
1780            rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*&
1781              req%target_sign(n)
1782          ENDDO
1783        ENDIF
1784      ENDIF       
1785
1786    ENDDO
1787   
1788  END SUBROUTINE transfert_request_seq
1789 
1790 
1791  SUBROUTINE gather_field(field_loc,field_glo)
1792  USE field_mod
1793  USE domain_mod
1794  USE mpi_mod
1795  USE mpipara
1796  IMPLICIT NONE
1797    TYPE(t_field),POINTER :: field_loc(:)
1798    TYPE(t_field),POINTER :: field_glo(:)
1799    INTEGER, ALLOCATABLE :: mpi_req(:)
1800    INTEGER, ALLOCATABLE :: status(:,:)
1801    INTEGER :: ireq,nreq
1802    INTEGER :: ind_glo,ind_loc   
1803 
1804    IF (.NOT. using_mpi) THEN
1805   
1806      DO ind_loc=1,ndomain
1807        IF (field_loc(ind_loc)%ndim==2) field_glo(ind_loc)%rval2d=field_loc(ind_loc)%rval2d
1808        IF (field_loc(ind_loc)%ndim==3) field_glo(ind_loc)%rval3d=field_loc(ind_loc)%rval3d
1809        IF (field_loc(ind_loc)%ndim==4) field_glo(ind_loc)%rval4d=field_loc(ind_loc)%rval4d
1810      ENDDO
1811   
1812    ELSE
1813         
1814      nreq=ndomain
1815      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
1816      ALLOCATE(mpi_req(nreq))
1817      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
1818   
1819   
1820      ireq=0
1821      IF (mpi_rank==0) THEN
1822        DO ind_glo=1,ndomain_glo
1823          ireq=ireq+1
1824
1825          IF (field_glo(ind_glo)%ndim==2) THEN
1826            CALL MPI_IRECV(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
1827                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1828   
1829          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1830            CALL MPI_IRECV(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
1831                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1832
1833          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1834            CALL MPI_IRECV(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
1835                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1836          ENDIF
1837         
1838        ENDDO
1839      ENDIF
1840 
1841      DO ind_loc=1,ndomain
1842        ireq=ireq+1
1843
1844        IF (field_loc(ind_loc)%ndim==2) THEN
1845          CALL MPI_ISEND(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
1846                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1847        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
1848          CALL MPI_ISEND(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
1849                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1850        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
1851          CALL MPI_ISEND(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
1852                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1853        ENDIF
1854     
1855      ENDDO
1856   
1857      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
1858
1859    ENDIF
1860       
1861  END SUBROUTINE gather_field
1862
1863  SUBROUTINE bcast_field(field_glo)
1864  USE field_mod
1865  USE domain_mod
1866  USE mpi_mod
1867  USE mpipara
1868  IMPLICIT NONE
1869    TYPE(t_field),POINTER :: field_glo(:)
1870    INTEGER :: ind_glo   
1871 
1872    IF (.NOT. using_mpi) THEN
1873   
1874! nothing to do
1875   
1876    ELSE
1877         
1878      DO ind_glo=1,ndomain_glo
1879
1880          IF (field_glo(ind_glo)%ndim==2) THEN
1881            CALL MPI_BCAST(field_glo(ind_glo)%rval2d, size(field_glo(ind_glo)%rval2d) , MPI_REAL8, 0, comm_icosa, ierr)
1882          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1883            CALL MPI_BCAST(field_glo(ind_glo)%rval3d, size(field_glo(ind_glo)%rval3d) , MPI_REAL8, 0, comm_icosa, ierr)
1884          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1885            CALL MPI_BCAST(field_glo(ind_glo)%rval4d, size(field_glo(ind_glo)%rval4d) , MPI_REAL8, 0, comm_icosa, ierr)
1886          ENDIF
1887         
1888        ENDDO
1889      ENDIF
1890       
1891  END SUBROUTINE bcast_field
1892
1893  SUBROUTINE scatter_field(field_glo,field_loc)
1894  USE field_mod
1895  USE domain_mod
1896  USE mpi_mod
1897  USE mpipara
1898  IMPLICIT NONE
1899    TYPE(t_field),POINTER :: field_glo(:)
1900    TYPE(t_field),POINTER :: field_loc(:)
1901    INTEGER, ALLOCATABLE :: mpi_req(:)
1902    INTEGER, ALLOCATABLE :: status(:,:)
1903    INTEGER :: ireq,nreq
1904    INTEGER :: ind_glo,ind_loc   
1905 
1906    IF (.NOT. using_mpi) THEN
1907   
1908      DO ind_loc=1,ndomain
1909        IF (field_loc(ind_loc)%ndim==2) field_loc(ind_loc)%rval2d=field_glo(ind_loc)%rval2d
1910        IF (field_loc(ind_loc)%ndim==3) field_loc(ind_loc)%rval3d=field_glo(ind_loc)%rval3d
1911        IF (field_loc(ind_loc)%ndim==4) field_loc(ind_loc)%rval4d=field_glo(ind_loc)%rval4d
1912      ENDDO
1913   
1914    ELSE
1915         
1916      nreq=ndomain
1917      IF (mpi_rank==0) nreq=nreq+ndomain_glo 
1918      ALLOCATE(mpi_req(nreq))
1919      ALLOCATE(status(MPI_STATUS_SIZE,nreq))
1920   
1921   
1922      ireq=0
1923      IF (mpi_rank==0) THEN
1924        DO ind_glo=1,ndomain_glo
1925          ireq=ireq+1
1926
1927          IF (field_glo(ind_glo)%ndim==2) THEN
1928            CALL MPI_ISEND(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 ,   &
1929                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1930   
1931          ELSE IF (field_glo(ind_glo)%ndim==3) THEN
1932            CALL MPI_ISEND(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 ,   &
1933                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1934
1935          ELSE IF (field_glo(ind_glo)%ndim==4) THEN
1936            CALL MPI_ISEND(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 ,   &
1937                         domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr)
1938          ENDIF
1939         
1940        ENDDO
1941      ENDIF
1942 
1943      DO ind_loc=1,ndomain
1944        ireq=ireq+1
1945
1946        IF (field_loc(ind_loc)%ndim==2) THEN
1947          CALL MPI_IRECV(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 ,   &
1948                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1949        ELSE IF (field_loc(ind_loc)%ndim==3) THEN
1950          CALL MPI_IRECV(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 ,   &
1951                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1952        ELSE IF (field_loc(ind_loc)%ndim==4) THEN
1953          CALL MPI_IRECV(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 ,   &
1954                         0, ind_loc, comm_icosa, mpi_req(ireq), ierr)
1955        ENDIF
1956     
1957      ENDDO
1958   
1959      CALL MPI_WAITALL(nreq,mpi_req,status,ierr)
1960
1961    ENDIF
1962       
1963  END SUBROUTINE scatter_field
1964 
1965  SUBROUTINE trace_in
1966  USE trace
1967  IMPLICIT NONE
1968 
1969    CALL trace_start("transfert_buffer")
1970  END SUBROUTINE trace_in             
1971
1972  SUBROUTINE trace_out
1973  USE trace
1974  IMPLICIT NONE
1975 
1976    CALL trace_end("transfert_buffer")
1977  END SUBROUTINE trace_out             
1978
1979
1980
1981
1982!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1983!! Definition des Broadcast --> 4D   !!
1984!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1985
1986!! -- Les chaine de charactï¿œre -- !!
1987
1988  SUBROUTINE bcast_mpi_c(var1)
1989  IMPLICIT NONE
1990    CHARACTER(LEN=*),INTENT(INOUT) :: Var1
1991   
1992    CALL bcast_mpi_cgen(Var1,len(Var1))
1993
1994  END SUBROUTINE bcast_mpi_c
1995
1996!! -- Les entiers -- !!
1997 
1998  SUBROUTINE bcast_mpi_i(var)
1999  USE mpipara
2000  IMPLICIT NONE
2001    INTEGER,INTENT(INOUT) :: Var
2002   
2003    INTEGER               :: var_tmp(1)
2004   
2005    IF (is_mpi_master) var_tmp(1)=var
2006    CALL bcast_mpi_igen(Var_tmp,1)
2007    var=var_tmp(1)
2008   
2009  END SUBROUTINE bcast_mpi_i
2010
2011  SUBROUTINE bcast_mpi_i1(var)
2012  IMPLICIT NONE
2013    INTEGER,INTENT(INOUT) :: Var(:)
2014
2015    CALL bcast_mpi_igen(Var,size(Var))
2016   
2017  END SUBROUTINE bcast_mpi_i1
2018
2019  SUBROUTINE bcast_mpi_i2(var)
2020  IMPLICIT NONE
2021    INTEGER,INTENT(INOUT) :: Var(:,:)
2022   
2023    CALL bcast_mpi_igen(Var,size(Var))
2024 
2025  END SUBROUTINE bcast_mpi_i2
2026
2027  SUBROUTINE bcast_mpi_i3(var)
2028  IMPLICIT NONE
2029    INTEGER,INTENT(INOUT) :: Var(:,:,:)
2030   
2031    CALL bcast_mpi_igen(Var,size(Var))
2032
2033  END SUBROUTINE bcast_mpi_i3
2034
2035  SUBROUTINE bcast_mpi_i4(var)
2036  IMPLICIT NONE
2037    INTEGER,INTENT(INOUT) :: Var(:,:,:,:)
2038   
2039    CALL bcast_mpi_igen(Var,size(Var))
2040
2041  END SUBROUTINE bcast_mpi_i4
2042
2043
2044!! -- Les reels -- !!
2045
2046  SUBROUTINE bcast_mpi_r(var)
2047  USE mpipara
2048  IMPLICIT NONE
2049    REAL,INTENT(INOUT) :: Var
2050    REAL               :: var_tmp(1)
2051   
2052    IF (is_mpi_master) var_tmp(1)=var
2053    CALL bcast_mpi_rgen(Var_tmp,1)
2054    var=var_tmp(1)   
2055
2056  END SUBROUTINE bcast_mpi_r
2057
2058  SUBROUTINE bcast_mpi_r1(var)
2059  IMPLICIT NONE
2060    REAL,INTENT(INOUT) :: Var(:)
2061   
2062    CALL bcast_mpi_rgen(Var,size(Var))
2063
2064  END SUBROUTINE bcast_mpi_r1
2065
2066  SUBROUTINE bcast_mpi_r2(var)
2067  IMPLICIT NONE
2068    REAL,INTENT(INOUT) :: Var(:,:)
2069   
2070    CALL bcast_mpi_rgen(Var,size(Var))
2071
2072  END SUBROUTINE bcast_mpi_r2
2073
2074  SUBROUTINE bcast_mpi_r3(var)
2075  IMPLICIT NONE
2076    REAL,INTENT(INOUT) :: Var(:,:,:)
2077   
2078    CALL bcast_mpi_rgen(Var,size(Var))
2079
2080  END SUBROUTINE bcast_mpi_r3
2081
2082  SUBROUTINE bcast_mpi_r4(var)
2083  IMPLICIT NONE
2084    REAL,INTENT(INOUT) :: Var(:,:,:,:)
2085   
2086    CALL bcast_mpi_rgen(Var,size(Var))
2087
2088  END SUBROUTINE bcast_mpi_r4
2089 
2090!! -- Les booleans -- !!
2091
2092  SUBROUTINE bcast_mpi_l(var)
2093  USE mpipara
2094  IMPLICIT NONE
2095    LOGICAL,INTENT(INOUT) :: Var
2096    LOGICAL               :: var_tmp(1)
2097   
2098    IF (is_mpi_master) var_tmp(1)=var
2099    CALL bcast_mpi_lgen(Var_tmp,1)
2100    var=var_tmp(1)   
2101
2102  END SUBROUTINE bcast_mpi_l
2103
2104  SUBROUTINE bcast_mpi_l1(var)
2105  IMPLICIT NONE
2106    LOGICAL,INTENT(INOUT) :: Var(:)
2107   
2108    CALL bcast_mpi_lgen(Var,size(Var))
2109
2110  END SUBROUTINE bcast_mpi_l1
2111
2112  SUBROUTINE bcast_mpi_l2(var)
2113  IMPLICIT NONE
2114    LOGICAL,INTENT(INOUT) :: Var(:,:)
2115   
2116    CALL bcast_mpi_lgen(Var,size(Var))
2117
2118  END SUBROUTINE bcast_mpi_l2
2119
2120  SUBROUTINE bcast_mpi_l3(var)
2121  IMPLICIT NONE
2122    LOGICAL,INTENT(INOUT) :: Var(:,:,:)
2123   
2124    CALL bcast_mpi_lgen(Var,size(Var))
2125
2126  END SUBROUTINE bcast_mpi_l3
2127
2128  SUBROUTINE bcast_mpi_l4(var)
2129  IMPLICIT NONE
2130    LOGICAL,INTENT(INOUT) :: Var(:,:,:,:)
2131   
2132    CALL bcast_mpi_lgen(Var,size(Var))
2133
2134  END SUBROUTINE bcast_mpi_l4
2135 
2136!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2137!! DEFINITION DES FONCTIONS DE TRANSFERT GENERIQUES !
2138!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2139
2140  SUBROUTINE bcast_mpi_cgen(var,nb)
2141    USE mpi_mod
2142    USE mpipara
2143    IMPLICIT NONE
2144   
2145    CHARACTER(LEN=*),INTENT(INOUT) :: Var
2146    INTEGER,INTENT(IN) :: nb
2147
2148    IF (.NOT. using_mpi) RETURN
2149   
2150    CALL MPI_BCAST(Var,nb,MPI_CHARACTER,mpi_master,comm_icosa,ierr)
2151       
2152  END SUBROUTINE bcast_mpi_cgen
2153
2154
2155     
2156  SUBROUTINE bcast_mpi_igen(var,nb)
2157    USE mpi_mod
2158    USE mpipara
2159    IMPLICIT NONE
2160   
2161    INTEGER,INTENT(IN) :: nb
2162    INTEGER,DIMENSION(nb),INTENT(INOUT) :: Var
2163
2164    IF (.NOT. using_mpi) RETURN
2165
2166    CALL MPI_BCAST(Var,nb,MPI_INTEGER,mpi_master,comm_icosa,ierr)
2167       
2168  END SUBROUTINE bcast_mpi_igen
2169
2170
2171
2172 
2173  SUBROUTINE bcast_mpi_rgen(var,nb)
2174    USE mpi_mod
2175    USE mpipara
2176    IMPLICIT NONE
2177   
2178    INTEGER,INTENT(IN) :: nb
2179    REAL,DIMENSION(nb),INTENT(INOUT) :: Var
2180
2181    IF (.NOT. using_mpi) RETURN
2182
2183    CALL MPI_BCAST(Var,nb,MPI_REAL8,mpi_master,comm_icosa,ierr)
2184   
2185  END SUBROUTINE bcast_mpi_rgen
2186 
2187
2188
2189
2190  SUBROUTINE bcast_mpi_lgen(var,nb)
2191    USE mpi_mod
2192    USE mpipara
2193    IMPLICIT NONE
2194   
2195    INTEGER,INTENT(IN) :: nb
2196    LOGICAL,DIMENSION(nb),INTENT(INOUT) :: Var
2197
2198    IF (.NOT. using_mpi) RETURN
2199
2200    CALL MPI_BCAST(Var,nb,MPI_LOGICAL,mpi_master,comm_icosa,ierr)
2201
2202  END SUBROUTINE bcast_mpi_lgen
2203 
2204   
2205END MODULE transfert_mpi_mod
2206     
2207       
2208       
2209       
2210     
Note: See TracBrowser for help on using the repository browser.