source: codes/icosagcm/trunk/src/parallel/transfert_mpi.f90 @ 922

Last change on this file since 922 was 901, checked in by adurocher, 5 years ago

trunk : Fixed compilation with --std=f2008 with gfortran

Added dynamico_abort() to replace non standard ABORT() intrinsic
Other modifications to respect the fortran standard

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