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

Last change on this file since 803 was 711, checked in by ymipsl, 6 years ago

Adding halo transfer for scalar field ad vorticity point.
=> Use request "req_z1_scal" with transfer function

YM

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