source: LMDZ5/branches/LF-private/libf/dyn3dmem/mod_hallo.F90 @ 5440

Last change on this file since 5440 was 1847, checked in by yann meurdesoif, 11 years ago

Improved memory management of buffer for MPi request

YM

File size: 44.9 KB
Line 
1module mod_Hallo
2USE parallel_lmdz
3implicit none
4  logical,save :: use_mpi_alloc
5  integer, parameter :: MaxProc=512
6  integer, parameter :: DefaultMaxBufferSize=1024*1024*100
7  integer, SAVE :: MaxBufferSize=0
8  integer, parameter :: ListSize=1000
9 
10  integer,save       :: MaxBufferSize_Used
11!$OMP THREADPRIVATE( MaxBufferSize_Used)
12
13   real,save,pointer,dimension(:) :: Buffer
14!$OMP THREADPRIVATE(Buffer)
15
16   integer,save,dimension(Listsize) :: Buffer_Pos
17   integer,save :: Index_Pos
18!$OMP THREADPRIVATE(Buffer_Pos,Index_pos)
19   
20  type Hallo
21    real, dimension(:,:),pointer :: Field
22    integer :: offset
23    integer :: size
24    integer :: NbLevel
25    integer :: Stride
26  end type Hallo
27 
28  type request_SR
29    integer :: NbRequest=0
30    integer :: NbRequestMax=0
31    integer :: BufferSize
32    integer :: Pos
33    integer :: Index
34    type(Hallo), POINTER :: Hallo(:)
35    integer :: MSG_Request
36  end type request_SR
37
38  type request
39    type(request_SR),dimension(0:MaxProc-1) :: RequestSend
40    type(request_SR),dimension(0:MaxProc-1) :: RequestRecv
41    integer :: tag=1
42  end type request
43 
44   TYPE(distrib),SAVE :: distrib_gather
45
46
47  INTERFACE Register_SwapField_u
48    MODULE PROCEDURE Register_SwapField1d_u,Register_SwapField2d_u1d,Register_SwapField3d_u
49  END INTERFACE Register_SwapField_u
50
51  INTERFACE Register_SwapField_v
52    MODULE PROCEDURE Register_SwapField1d_v,Register_SwapField2d_v1d,Register_SwapField3d_v
53  END INTERFACE Register_SwapField_v
54
55  INTERFACE Register_SwapField2d_u
56    MODULE PROCEDURE Register_SwapField1d_u2d,Register_SwapField2d_u2d,Register_SwapField3d_u2d
57  END INTERFACE Register_SwapField2d_u
58
59  INTERFACE Register_SwapField2d_v
60    MODULE PROCEDURE Register_SwapField1d_v2d,Register_SwapField2d_v2d,Register_SwapField3d_v2d
61  END INTERFACE Register_SwapField2d_v
62
63  contains
64
65  subroutine Init_mod_hallo
66  USE dimensions_mod
67  USE IOIPSL
68    implicit none
69    integer :: jj_nb_gather(0:mpi_size-1)
70   
71    Index_Pos=1
72    Buffer_Pos(Index_Pos)=1
73    MaxBufferSize_Used=0
74!$OMP MASTER     
75    MaxBufferSize=DefaultMaxBufferSize
76    CALL getin("mpi_buffer_size",MaxBufferSize)
77!$OMP END MASTER
78!$OMP BARRIER
79   
80    IF (use_mpi_alloc .AND. using_mpi) THEN
81      CALL create_global_mpi_buffer
82    ELSE
83      CALL create_standard_mpi_buffer
84    ENDIF
85     
86!$OMP MASTER     
87     jj_nb_gather(:)=0
88     jj_nb_gather(0)=jjp1
89     
90     CALL create_distrib(jj_nb_gather,distrib_gather)
91!$OMP END MASTER
92!$OMP BARRIER
93
94  end subroutine init_mod_hallo
95
96  SUBROUTINE create_standard_mpi_buffer
97  IMPLICIT NONE
98   
99    ALLOCATE(Buffer(MaxBufferSize))
100   
101  END SUBROUTINE create_standard_mpi_buffer
102 
103  SUBROUTINE create_global_mpi_buffer
104  IMPLICIT NONE
105#ifdef CPP_MPI
106  INCLUDE 'mpif.h'
107#endif 
108    POINTER (Pbuffer,MPI_Buffer(MaxBufferSize))
109    REAL :: MPI_Buffer
110#ifdef CPP_MPI
111    INTEGER(KIND=MPI_ADDRESS_KIND) :: BS
112#else
113    INTEGER(KIND=8) :: BS
114#endif
115    INTEGER :: i,ierr
116
117!  Allocation du buffer MPI
118      Bs=8*MaxBufferSize
119!$OMP CRITICAL (MPI)
120#ifdef CPP_MPI
121      CALL MPI_ALLOC_MEM(BS,MPI_INFO_NULL,Pbuffer,ierr)
122#endif
123!$OMP END CRITICAL (MPI)
124      DO i=1,MaxBufferSize
125        MPI_Buffer(i)=i
126      ENDDO
127     
128      CALL  Associate_buffer(MPI_Buffer)
129     
130  CONTAINS
131     
132     SUBROUTINE Associate_buffer(MPI_Buffer)
133     IMPLICIT NONE
134       REAL,DIMENSION(:),target :: MPI_Buffer 
135
136         Buffer=>MPI_Buffer
137 
138      END SUBROUTINE  Associate_buffer
139                                     
140  END SUBROUTINE create_global_mpi_buffer
141 
142     
143  subroutine allocate_buffer(Size,Index,Pos)
144  implicit none
145    integer :: Size
146    integer :: Index
147    integer :: Pos
148
149    if (Buffer_pos(Index_pos)+Size>MaxBufferSize_Used) MaxBufferSize_Used=Buffer_pos(Index_pos)+Size 
150    if (Buffer_pos(Index_pos)+Size>MaxBufferSize) then
151      print *,'STOP :: La taille de MaxBufferSize dans mod_hallo.F90 est trop petite !!!!'
152      stop
153    endif
154   
155    if (Index_pos>=ListSize) then
156      print *,'STOP :: La taille de ListSize dans mod_hallo.F90 est trop petite !!!!'
157      stop
158    endif
159     
160    Pos=Buffer_Pos(Index_Pos)
161    Buffer_Pos(Index_pos+1)=Buffer_Pos(Index_Pos)+Size
162    Index_Pos=Index_Pos+1
163    Index=Index_Pos
164   
165  end subroutine allocate_buffer
166     
167  subroutine deallocate_buffer(Index)
168  implicit none
169    integer :: Index
170   
171    Buffer_Pos(Index)=-1
172   
173    do while (Buffer_Pos(Index_Pos)==-1 .and. Index_Pos>1)
174      Index_Pos=Index_Pos-1
175    end do
176
177  end subroutine deallocate_buffer 
178 
179  subroutine SetTag(a_request,tag)
180  implicit none
181    type(request):: a_request
182    integer :: tag
183   
184    a_request%tag=tag
185  end subroutine SetTag
186 
187 
188  subroutine New_Hallo(Field,Stride,NbLevel,offset,size,Ptr_request)
189    integer :: Stride
190    integer :: NbLevel
191    integer :: size
192    integer :: offset
193    real, dimension(Stride,NbLevel),target :: Field
194    type(request_SR),pointer :: Ptr_request
195    type(Hallo),POINTER :: NewHallos(:),HalloSwitch(:), NewHallo
196   
197    Ptr_Request%NbRequest=Ptr_Request%NbRequest+1
198    IF(Ptr_Request%NbRequestMax==0) THEN
199       Ptr_Request%NbRequestMax=10
200       ALLOCATE(Ptr_Request%Hallo(Ptr_Request%NbRequestMax))
201    ELSE IF ( Ptr_Request%NbRequest > Ptr_Request%NbRequestMax) THEN
202      Ptr_Request%NbRequestMax=INT(Ptr_Request%NbRequestMax*1.2)
203      ALLOCATE(NewHallos(Ptr_Request%NbRequestMax))
204      NewHallos(1:Ptr_Request%NbRequest-1)=Ptr_Request%hallo(1:Ptr_Request%NbRequest-1)
205      HalloSwitch=>Ptr_Request%hallo
206      Ptr_Request%hallo=>NewHallos
207      DEALLOCATE(HalloSwitch)
208    ENDIF
209   
210    NewHallo=>Ptr_Request%hallo(Ptr_Request%NbRequest)
211         
212    NewHallo%Field=>Field
213    NewHallo%Stride=Stride
214    NewHallo%NbLevel=NbLevel
215    NewHallo%size=size
216    NewHallo%offset=offset
217   
218  end subroutine New_Hallo
219 
220  subroutine Register_SendField(Field,ij,ll,offset,size,target,a_request)
221  USE dimensions_mod
222  implicit none
223
224   
225      INTEGER :: ij,ll,offset,size,target
226      REAL, dimension(ij,ll) :: Field
227      type(request),target :: a_request
228      type(request_SR),pointer :: Ptr_request
229
230      Ptr_Request=>a_request%RequestSend(target)
231
232      call New_Hallo(Field,ij,ll,offset,size,Ptr_request)
233     
234   end subroutine Register_SendField     
235     
236  subroutine Register_RecvField(Field,ij,ll,offset,size,target,a_request)
237  USE dimensions_mod
238  implicit none
239
240   
241      INTEGER :: ij,ll,offset,size,target
242      REAL, dimension(ij,ll) :: Field
243      type(request),target :: a_request
244      type(request_SR),pointer :: Ptr_request
245
246      Ptr_Request=>a_request%RequestRecv(target)
247           
248      call New_Hallo(Field,ij,ll,offset,size,Ptr_request)
249
250     
251   end subroutine Register_RecvField     
252 
253  subroutine Register_SwapField(FieldS,FieldR,ij,ll,jj_Nb_New,a_request)
254  USE dimensions_mod
255      implicit none
256
257   
258    INTEGER :: ij,ll
259    REAL, dimension(ij,ll) :: FieldS
260    REAL, dimension(ij,ll) :: FieldR
261    type(request) :: a_request
262    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
263    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
264   
265    integer ::i,jje,jjb
266   
267    jj_begin_New(0)=1
268    jj_End_New(0)=jj_Nb_New(0)
269    do i=1,MPI_Size-1
270      jj_begin_New(i)=jj_end_New(i-1)+1
271      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
272    enddo
273   
274    do i=0,MPI_Size-1
275      if (i /= MPI_Rank) then
276        jjb=max(jj_begin_new(i),jj_begin)
277        jje=min(jj_end_new(i),jj_end)
278       
279        if (ij==ip1jm .and. jje==jjp1) jje=jjm
280       
281        if (jje >= jjb) then
282          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
283        endif
284       
285        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
286        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
287       
288        if (ij==ip1jm .and. jje==jjp1) jje=jjm
289       
290        if (jje >= jjb) then
291          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
292        endif
293       
294      endif
295    enddo
296   
297  end subroutine Register_SwapField   
298 
299
300 
301  subroutine Register_SwapFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down,a_request)
302  USE dimensions_mod
303 
304      implicit none
305   
306    INTEGER :: ij,ll,Up,Down
307    REAL, dimension(ij,ll) :: FieldS
308    REAL, dimension(ij,ll) :: FieldR
309    type(request) :: a_request
310    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
311    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
312   
313    integer ::i,jje,jjb
314   
315    jj_begin_New(0)=1
316    jj_End_New(0)=jj_Nb_New(0)
317    do i=1,MPI_Size-1
318      jj_begin_New(i)=jj_end_New(i-1)+1
319      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
320    enddo
321   
322    do i=0,MPI_Size-1
323      jj_begin_New(i)=max(1,jj_begin_New(i)-Up)
324      jj_end_New(i)=min(jjp1,jj_end_new(i)+Down)
325    enddo
326   
327    do i=0,MPI_Size-1
328      if (i /= MPI_Rank) then
329        jjb=max(jj_begin_new(i),jj_begin)
330        jje=min(jj_end_new(i),jj_end)
331       
332        if (ij==ip1jm .and. jje==jjp1) jje=jjm
333       
334        if (jje >= jjb) then
335          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
336        endif
337       
338        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
339        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
340       
341        if (ij==ip1jm .and. jje==jjp1) jje=jjm
342       
343        if (jje >= jjb) then
344          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
345        endif
346       
347      endif
348    enddo
349   
350  end subroutine Register_SwapFieldHallo
351
352
353
354  SUBROUTINE Register_SwapField1d_u(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
355  USE parallel_lmdz
356  USE dimensions_mod
357      IMPLICIT NONE
358   
359    REAL, DIMENSION(:),INTENT(IN)     :: FieldS
360    REAL, DIMENSION(:),INTENT(OUT)    :: FieldR
361    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
362    TYPE(distrib),INTENT(IN)          :: new_dist
363    INTEGER,OPTIONAL,INTENT(IN)       :: up
364    INTEGER,OPTIONAL,INTENT(IN)       :: down     
365    TYPE(request),INTENT(INOUT)         :: a_request
366
367    INTEGER                           :: halo_up
368    INTEGER                           :: halo_down
369   
370   
371    halo_up=0
372    halo_down=0
373    IF (PRESENT(up))   halo_up=up
374    IF (PRESENT(down)) halo_down=down
375
376    IF (PRESENT(old_dist)) THEN
377      CALL  Register_SwapField_gen_u(FieldS,FieldR,1,old_dist,new_dist,halo_up,halo_down,a_request)
378    ELSE
379      CALL  Register_SwapField_gen_u(FieldS,FieldR,1,current_dist,new_dist,halo_up,halo_down,a_request)
380    ENDIF
381       
382  END SUBROUTINE  Register_SwapField1d_u
383
384
385  SUBROUTINE Register_SwapField2d_u1d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
386  USE parallel_lmdz
387  USE dimensions_mod
388    IMPLICIT NONE
389   
390    REAL, DIMENSION(:,:),INTENT(IN)     :: FieldS
391    REAL, DIMENSION(:,:),INTENT(OUT)    :: FieldR
392    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
393    TYPE(distrib),INTENT(IN)          :: new_dist
394    INTEGER,OPTIONAL,INTENT(IN)       :: up
395    INTEGER,OPTIONAL,INTENT(IN)       :: down     
396    TYPE(request),INTENT(INOUT)         :: a_request
397
398    INTEGER                           :: halo_up
399    INTEGER                           :: halo_down
400    INTEGER                           :: ll
401       
402   
403    halo_up=0
404    halo_down=0
405    IF (PRESENT(up))   halo_up=up
406    IF (PRESENT(down)) halo_down=down
407   
408    ll=size(FieldS,2)
409   
410    IF (PRESENT(old_dist)) THEN
411      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
412    ELSE
413      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
414    ENDIF
415   
416  END SUBROUTINE  Register_SwapField2d_u1d
417   
418
419  SUBROUTINE Register_SwapField3d_u(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
420  USE parallel_lmdz
421  USE dimensions_mod
422      IMPLICIT NONE
423   
424    REAL, DIMENSION(:,:,:),INTENT(IN)     :: FieldS
425    REAL, DIMENSION(:,:,:),INTENT(OUT)    :: FieldR
426    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
427    TYPE(distrib),INTENT(IN)          :: new_dist
428    INTEGER,OPTIONAL,INTENT(IN)       :: up
429    INTEGER,OPTIONAL,INTENT(IN)       :: down     
430    TYPE(request),INTENT(INOUT)         :: a_request
431
432    INTEGER                           :: halo_up
433    INTEGER                           :: halo_down
434    INTEGER                           :: ll
435       
436   
437    halo_up=0
438    halo_down=0
439    IF (PRESENT(up))   halo_up=up
440    IF (PRESENT(down)) halo_down=down
441   
442    ll=size(FieldS,2)*size(FieldS,3)
443   
444    IF (PRESENT(old_dist)) THEN
445      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
446    ELSE
447      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
448    ENDIF
449   
450  END SUBROUTINE  Register_SwapField3d_u
451 
452
453
454 SUBROUTINE Register_SwapField1d_u2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
455  USE parallel_lmdz
456  USE dimensions_mod
457
458      IMPLICIT NONE
459   
460    REAL, DIMENSION(:,:),INTENT(IN)     :: FieldS
461    REAL, DIMENSION(:,:),INTENT(OUT)    :: FieldR
462    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
463    TYPE(distrib),OPTIONAL,INTENT(IN)          :: new_dist !LF
464    INTEGER,OPTIONAL,INTENT(IN)       :: up
465    INTEGER,OPTIONAL,INTENT(IN)       :: down     
466    TYPE(request),INTENT(INOUT)         :: a_request
467
468    INTEGER                           :: halo_up
469    INTEGER                           :: halo_down
470   
471   
472    halo_up=0
473    halo_down=0
474    IF (PRESENT(up))   halo_up=up
475    IF (PRESENT(down)) halo_down=down
476
477    IF (PRESENT(old_dist)) THEN
478      CALL  Register_SwapField_gen_u(FieldS,FieldR,1,old_dist,new_dist,halo_up,halo_down,a_request)
479    ELSE
480      CALL  Register_SwapField_gen_u(FieldS,FieldR,1,current_dist,new_dist,halo_up,halo_down,a_request)
481    ENDIF
482       
483  END SUBROUTINE  Register_SwapField1d_u2d
484
485
486  SUBROUTINE Register_SwapField2d_u2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
487  USE parallel_lmdz
488  USE dimensions_mod
489
490      IMPLICIT NONE
491   
492    REAL, DIMENSION(:,:,:),INTENT(IN)     :: FieldS
493    REAL, DIMENSION(:,:,:),INTENT(OUT)    :: FieldR
494    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
495    TYPE(distrib),INTENT(IN)          :: new_dist
496    INTEGER,OPTIONAL,INTENT(IN)       :: up
497    INTEGER,OPTIONAL,INTENT(IN)       :: down     
498    TYPE(request),INTENT(INOUT)         :: a_request
499
500    INTEGER                           :: halo_up
501    INTEGER                           :: halo_down
502    INTEGER                           :: ll
503       
504   
505    halo_up=0
506    halo_down=0
507    IF (PRESENT(up))   halo_up=up
508    IF (PRESENT(down)) halo_down=down
509   
510    ll=size(FieldS,3)
511   
512    IF (PRESENT(old_dist)) THEN
513      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
514    ELSE
515      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
516    ENDIF
517   
518  END SUBROUTINE  Register_SwapField2d_u2d
519   
520
521  SUBROUTINE Register_SwapField3d_u2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
522  USE parallel_lmdz
523  USE dimensions_mod
524      IMPLICIT NONE
525   
526    REAL, DIMENSION(:,:,:,:),INTENT(IN)     :: FieldS
527    REAL, DIMENSION(:,:,:,:),INTENT(OUT)    :: FieldR
528    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
529    TYPE(distrib),INTENT(IN)          :: new_dist
530    INTEGER,OPTIONAL,INTENT(IN)       :: up
531    INTEGER,OPTIONAL,INTENT(IN)       :: down     
532    TYPE(request),INTENT(INOUT)         :: a_request
533
534    INTEGER                           :: halo_up
535    INTEGER                           :: halo_down
536    INTEGER                           :: ll
537       
538   
539    halo_up=0
540    halo_down=0
541    IF (PRESENT(up))   halo_up=up
542    IF (PRESENT(down)) halo_down=down
543   
544    ll=size(FieldS,3)*size(FieldS,4)
545   
546    IF (PRESENT(old_dist)) THEN
547      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
548    ELSE
549      CALL  Register_SwapField_gen_u(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
550    ENDIF
551   
552  END SUBROUTINE  Register_SwapField3d_u2d
553
554
555
556
557
558
559
560  SUBROUTINE Register_SwapField1d_v(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
561  USE parallel_lmdz
562  USE dimensions_mod
563      IMPLICIT NONE
564   
565    REAL, DIMENSION(:),INTENT(IN)     :: FieldS
566    REAL, DIMENSION(:),INTENT(OUT)    :: FieldR
567    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
568    TYPE(distrib),INTENT(IN)          :: new_dist
569    INTEGER,OPTIONAL,INTENT(IN)       :: up
570    INTEGER,OPTIONAL,INTENT(IN)       :: down     
571    TYPE(request),INTENT(INOUT)         :: a_request
572
573    INTEGER                           :: halo_up
574    INTEGER                           :: halo_down
575   
576   
577    halo_up=0
578    halo_down=0
579    IF (PRESENT(up))   halo_up=up
580    IF (PRESENT(down)) halo_down=down
581
582    IF (PRESENT(old_dist)) THEN
583      CALL  Register_SwapField_gen_v(FieldS,FieldR,1,old_dist,new_dist,halo_up,halo_down,a_request)
584    ELSE
585      CALL  Register_SwapField_gen_v(FieldS,FieldR,1,current_dist,new_dist,halo_up,halo_down,a_request)
586    ENDIF
587       
588  END SUBROUTINE  Register_SwapField1d_v
589
590
591  SUBROUTINE Register_SwapField2d_v1d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
592  USE parallel_lmdz
593  USE dimensions_mod
594      IMPLICIT NONE
595   
596    REAL, DIMENSION(:,:),INTENT(IN)     :: FieldS
597    REAL, DIMENSION(:,:),INTENT(OUT)    :: FieldR
598    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
599    TYPE(distrib),INTENT(IN)          :: new_dist
600    INTEGER,OPTIONAL,INTENT(IN)       :: up
601    INTEGER,OPTIONAL,INTENT(IN)       :: down     
602    TYPE(request),INTENT(INOUT)         :: a_request
603
604    INTEGER                           :: halo_up
605    INTEGER                           :: halo_down
606    INTEGER                           :: ll
607       
608   
609    halo_up=0
610    halo_down=0
611    IF (PRESENT(up))   halo_up=up
612    IF (PRESENT(down)) halo_down=down
613   
614    ll=size(FieldS,2)
615   
616    IF (PRESENT(old_dist)) THEN
617      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
618    ELSE
619      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
620    ENDIF
621   
622  END SUBROUTINE  Register_SwapField2d_v1d
623   
624
625  SUBROUTINE Register_SwapField3d_v(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
626  USE parallel_lmdz
627  USE dimensions_mod
628      IMPLICIT NONE
629   
630    REAL, DIMENSION(:,:,:),INTENT(IN)     :: FieldS
631    REAL, DIMENSION(:,:,:),INTENT(OUT)    :: FieldR
632    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
633    TYPE(distrib),INTENT(IN)          :: new_dist
634    INTEGER,OPTIONAL,INTENT(IN)       :: up
635    INTEGER,OPTIONAL,INTENT(IN)       :: down     
636    TYPE(request),INTENT(INOUT)         :: a_request
637
638    INTEGER                           :: halo_up
639    INTEGER                           :: halo_down
640    INTEGER                           :: ll
641       
642   
643    halo_up=0
644    halo_down=0
645    IF (PRESENT(up))   halo_up=up
646    IF (PRESENT(down)) halo_down=down
647   
648    ll=size(FieldS,2)*size(FieldS,3)
649   
650    IF (PRESENT(old_dist)) THEN
651      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
652    ELSE
653      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
654    ENDIF
655   
656  END SUBROUTINE  Register_SwapField3d_v
657
658
659
660
661  SUBROUTINE Register_SwapField1d_v2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
662  USE parallel_lmdz
663  USE dimensions_mod
664      IMPLICIT NONE
665   
666    REAL, DIMENSION(:,:),INTENT(IN)     :: FieldS
667    REAL, DIMENSION(:,:),INTENT(OUT)    :: FieldR
668    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
669    TYPE(distrib),OPTIONAL,INTENT(IN)          :: new_dist !LF
670    INTEGER,OPTIONAL,INTENT(IN)       :: up
671    INTEGER,OPTIONAL,INTENT(IN)       :: down     
672    TYPE(request),INTENT(INOUT)         :: a_request
673
674    INTEGER                           :: halo_up
675    INTEGER                           :: halo_down
676   
677   
678    halo_up=0
679    halo_down=0
680    IF (PRESENT(up))   halo_up=up
681    IF (PRESENT(down)) halo_down=down
682
683    IF (PRESENT(old_dist)) THEN
684      CALL  Register_SwapField_gen_v(FieldS,FieldR,1,old_dist,new_dist,halo_up,halo_down,a_request)
685    ELSE
686      CALL  Register_SwapField_gen_v(FieldS,FieldR,1,current_dist,new_dist,halo_up,halo_down,a_request)
687    ENDIF
688       
689  END SUBROUTINE  Register_SwapField1d_v2d
690
691
692  SUBROUTINE Register_SwapField2d_v2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
693  USE parallel_lmdz
694  USE dimensions_mod
695      IMPLICIT NONE
696   
697    REAL, DIMENSION(:,:,:),INTENT(IN)     :: FieldS
698    REAL, DIMENSION(:,:,:),INTENT(OUT)    :: FieldR
699    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
700    TYPE(distrib),INTENT(IN)          :: new_dist
701    INTEGER,OPTIONAL,INTENT(IN)       :: up
702    INTEGER,OPTIONAL,INTENT(IN)       :: down     
703    TYPE(request),INTENT(INOUT)         :: a_request
704
705    INTEGER                           :: halo_up
706    INTEGER                           :: halo_down
707    INTEGER                           :: ll
708       
709   
710    halo_up=0
711    halo_down=0
712    IF (PRESENT(up))   halo_up=up
713    IF (PRESENT(down)) halo_down=down
714   
715    ll=size(FieldS,3)
716   
717    IF (PRESENT(old_dist)) THEN
718      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
719    ELSE
720      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
721    ENDIF
722   
723  END SUBROUTINE  Register_SwapField2d_v2d
724   
725
726  SUBROUTINE Register_SwapField3d_v2d(FieldS,FieldR,new_dist,a_request,old_dist,up,down)
727  USE parallel_lmdz
728  USE dimensions_mod
729      IMPLICIT NONE
730   
731    REAL, DIMENSION(:,:,:,:),INTENT(IN)     :: FieldS
732    REAL, DIMENSION(:,:,:,:),INTENT(OUT)    :: FieldR
733    TYPE(distrib),OPTIONAL,INTENT(IN) :: old_dist
734    TYPE(distrib),INTENT(IN)          :: new_dist
735    INTEGER,OPTIONAL,INTENT(IN)       :: up
736    INTEGER,OPTIONAL,INTENT(IN)       :: down     
737    TYPE(request),INTENT(INOUT)         :: a_request
738
739    INTEGER                           :: halo_up
740    INTEGER                           :: halo_down
741    INTEGER                           :: ll
742       
743   
744    halo_up=0
745    halo_down=0
746    IF (PRESENT(up))   halo_up=up
747    IF (PRESENT(down)) halo_down=down
748   
749    ll=size(FieldS,3)*size(FieldS,4)
750   
751    IF (PRESENT(old_dist)) THEN
752      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,old_dist,new_dist,halo_up,halo_down,a_request)
753    ELSE
754      CALL  Register_SwapField_gen_v(FieldS,FieldR,ll,current_dist,new_dist,halo_up,halo_down,a_request)
755    ENDIF
756   
757  END SUBROUTINE  Register_SwapField3d_v2d
758 
759 
760
761  SUBROUTINE Register_SwapField_gen_u(FieldS,FieldR,ll,old_dist,new_dist,Up,Down,a_request)
762  USE parallel_lmdz
763  USE dimensions_mod
764      IMPLICIT NONE
765   
766    INTEGER :: ll,Up,Down
767    TYPE(distrib)  :: old_dist
768    TYPE(distrib)  :: new_dist
769    REAL, DIMENSION(old_dist%ijb_u:old_dist%ije_u,ll) :: FieldS
770    REAL, DIMENSION(new_dist%ijb_u:new_dist%ije_u,ll) :: FieldR
771    TYPE(request) :: a_request
772    INTEGER,DIMENSION(0:MPI_Size-1) :: jj_Nb_New   
773    INTEGER,DIMENSION(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
774   
775    INTEGER ::i,l,jje,jjb,ijb,ije
776   
777    DO i=0,MPI_Size-1
778      jj_begin_New(i)=max(1,new_dist%jj_begin_para(i)-Up)
779      jj_end_New(i)=min(jjp1,new_dist%jj_end_para(i)+Down)
780    ENDDO
781   
782    DO i=0,MPI_Size-1
783      IF (i /= MPI_Rank) THEN
784        jjb=max(jj_begin_new(i),old_dist%jj_begin)
785        jje=min(jj_end_new(i),old_dist%jj_end)
786       
787        IF (jje >= jjb) THEN
788          CALL Register_SendField(FieldS,old_dist%ijnb_u,ll,jjb-old_dist%jjb_u+1,jje-jjb+1,i,a_request)
789        ENDIF
790       
791        jjb=max(jj_begin_new(MPI_Rank),old_dist%jj_begin_Para(i))
792        jje=min(jj_end_new(MPI_Rank),old_dist%jj_end_Para(i))
793       
794        IF (jje >= jjb) THEN
795          CALL Register_RecvField(FieldR,new_dist%ijnb_u,ll,jjb-new_dist%jjb_u+1,jje-jjb+1,i,a_request)
796        ENDIF
797      ELSE
798        jjb=max(jj_begin_new(i),old_dist%jj_begin)
799        jje=min(jj_end_new(i),old_dist%jj_end)
800        ijb=(jjb-1)*iip1+1
801        ije=jje*iip1
802!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
803        DO l=1,ll
804          FieldR(ijb:ije,l)=FieldS(ijb:ije,l)             
805        ENDDO
806!$OMP END DO NOWAIT
807      ENDIF
808    ENDDO
809   
810  END SUBROUTINE Register_SwapField_gen_u
811
812
813
814  SUBROUTINE Register_SwapField_gen_v(FieldS,FieldR,ll,old_dist,new_dist,Up,Down,a_request)
815  USE parallel_lmdz
816  USE dimensions_mod
817    IMPLICIT NONE
818   
819    INTEGER :: ll,Up,Down
820    TYPE(distrib)  :: old_dist
821    TYPE(distrib)  :: new_dist
822    REAL, DIMENSION(old_dist%ijb_v:old_dist%ije_v,ll) :: FieldS
823    REAL, DIMENSION(new_dist%ijb_v:new_dist%ije_v,ll) :: FieldR
824    TYPE(request) :: a_request
825    INTEGER,DIMENSION(0:MPI_Size-1) :: jj_Nb_New   
826    INTEGER,DIMENSION(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
827   
828    INTEGER ::i,l,jje,jjb,ijb,ije
829   
830    DO i=0,MPI_Size-1
831      jj_begin_New(i)=max(1,new_dist%jj_begin_para(i)-Up)
832      jj_end_New(i)=min(jjp1,new_dist%jj_end_para(i)+Down)
833    ENDDO
834   
835    DO i=0,MPI_Size-1
836      IF (i /= MPI_Rank) THEN
837        jjb=max(jj_begin_new(i),old_dist%jj_begin)
838        jje=min(jj_end_new(i),old_dist%jj_end)
839
840        IF (jje==jjp1) jje=jjm       
841
842        IF (jje >= jjb) THEN
843          CALL Register_SendField(FieldS,old_dist%ijnb_v,ll,jjb-old_dist%jjb_v+1,jje-jjb+1,i,a_request)
844        ENDIF
845       
846        jjb=max(jj_begin_new(MPI_Rank),old_dist%jj_begin_Para(i))
847        jje=min(jj_end_new(MPI_Rank),old_dist%jj_end_Para(i))
848
849        IF (jje==jjp1) jje=jjm
850       
851        IF (jje >= jjb) THEN
852          CALL Register_RecvField(FieldR,new_dist%ijnb_v,ll,jjb-new_dist%jjb_v+1,jje-jjb+1,i,a_request)
853        ENDIF
854      ELSE
855        jjb=max(jj_begin_new(i),old_dist%jj_begin)
856        jje=min(jj_end_new(i),old_dist%jj_end)
857        IF (jje==jjp1) jje=jjm
858        ijb=(jjb-1)*iip1+1
859        ije=jje*iip1
860!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
861        DO l=1,ll
862          FieldR(ijb:ije,l)=FieldS(ijb:ije,l)
863        ENDDO             
864!$OMP END DO NOWAIT
865      ENDIF
866    ENDDO
867   
868  END SUBROUTINE Register_SwapField_gen_v
869
870
871 
872
873 
874  subroutine Register_Hallo(Field,ij,ll,RUp,Rdown,SUp,SDown,a_request)
875  USE dimensions_mod
876      implicit none
877
878#ifdef CPP_MPI
879    include 'mpif.h'
880#endif   
881      INTEGER :: ij,ll
882      REAL, dimension(ij,ll) :: Field
883      INTEGER :: Sup,Sdown,rup,rdown
884      type(request) :: a_request
885      type(Hallo),pointer :: PtrHallo
886      LOGICAL :: SendUp,SendDown
887      LOGICAL :: RecvUp,RecvDown
888   
889 
890      SendUp=.TRUE.
891      SendDown=.TRUE.
892      RecvUp=.TRUE.
893      RecvDown=.TRUE.
894       
895      IF (pole_nord) THEN
896        SendUp=.FALSE.
897        RecvUp=.FALSE.
898      ENDIF
899 
900      IF (pole_sud) THEN
901        SendDown=.FALSE.
902        RecvDown=.FALSE.
903      ENDIF
904     
905      if (Sup.eq.0) then
906        SendUp=.FALSE.
907       endif
908     
909      if (Sdown.eq.0) then
910        SendDown=.FALSE.
911      endif
912
913      if (Rup.eq.0) then
914        RecvUp=.FALSE.
915      endif
916     
917      if (Rdown.eq.0) then
918        RecvDown=.FALSE.
919      endif
920     
921      IF (SendUp) THEN
922        call Register_SendField(Field,ij,ll,jj_begin,SUp,MPI_Rank-1,a_request)
923      ENDIF
924 
925      IF (SendDown) THEN
926        call Register_SendField(Field,ij,ll,jj_end-SDown+1,SDown,MPI_Rank+1,a_request)
927      ENDIF
928   
929 
930      IF (RecvUp) THEN
931        call Register_RecvField(Field,ij,ll,jj_begin-Rup,RUp,MPI_Rank-1,a_request)
932      ENDIF
933 
934      IF (RecvDown) THEN
935        call Register_RecvField(Field,ij,ll,jj_end+1,RDown,MPI_Rank+1,a_request)
936      ENDIF
937 
938    end subroutine Register_Hallo
939
940
941  subroutine Register_Hallo_u(Field,ll,RUp,Rdown,SUp,SDown,a_request)
942  USE dimensions_mod
943      implicit none
944#ifdef CPP_MPI
945    include 'mpif.h'
946#endif   
947      INTEGER :: ll
948      REAL, dimension(ijb_u:ije_u,ll) :: Field
949      INTEGER :: Sup,Sdown,rup,rdown
950      type(request) :: a_request
951      type(Hallo),pointer :: PtrHallo
952      LOGICAL :: SendUp,SendDown
953      LOGICAL :: RecvUp,RecvDown
954   
955 
956      SendUp=.TRUE.
957      SendDown=.TRUE.
958      RecvUp=.TRUE.
959      RecvDown=.TRUE.
960       
961      IF (pole_nord) THEN
962        SendUp=.FALSE.
963        RecvUp=.FALSE.
964      ENDIF
965 
966      IF (pole_sud) THEN
967        SendDown=.FALSE.
968        RecvDown=.FALSE.
969      ENDIF
970     
971      if (Sup.eq.0) then
972        SendUp=.FALSE.
973       endif
974     
975      if (Sdown.eq.0) then
976        SendDown=.FALSE.
977      endif
978
979      if (Rup.eq.0) then
980        RecvUp=.FALSE.
981      endif
982     
983      if (Rdown.eq.0) then
984        RecvDown=.FALSE.
985      endif
986     
987      IF (SendUp) THEN
988        call Register_SendField(Field,ijnb_u,ll,jj_begin-jjb_u+1,SUp,MPI_Rank-1,a_request)
989      ENDIF
990 
991      IF (SendDown) THEN
992        call Register_SendField(Field,ijnb_u,ll,jj_end-SDown+1-jjb_u+1,SDown,MPI_Rank+1,a_request)
993      ENDIF
994   
995 
996      IF (RecvUp) THEN
997        call Register_RecvField(Field,ijnb_u,ll,jj_begin-Rup-jjb_u+1,RUp,MPI_Rank-1,a_request)
998      ENDIF
999 
1000      IF (RecvDown) THEN
1001        call Register_RecvField(Field,ijnb_u,ll,jj_end+1-jjb_u+1,RDown,MPI_Rank+1,a_request)
1002      ENDIF
1003 
1004    end subroutine Register_Hallo_u
1005
1006  subroutine Register_Hallo_v(Field,ll,RUp,Rdown,SUp,SDown,a_request)
1007  USE dimensions_mod
1008      implicit none
1009#ifdef CPP_MPI
1010    include 'mpif.h'
1011#endif   
1012      INTEGER :: ll
1013      REAL, dimension(ijb_v:ije_v,ll) :: Field
1014      INTEGER :: Sup,Sdown,rup,rdown
1015      type(request) :: a_request
1016      type(Hallo),pointer :: PtrHallo
1017      LOGICAL :: SendUp,SendDown
1018      LOGICAL :: RecvUp,RecvDown
1019   
1020 
1021      SendUp=.TRUE.
1022      SendDown=.TRUE.
1023      RecvUp=.TRUE.
1024      RecvDown=.TRUE.
1025       
1026      IF (pole_nord) THEN
1027        SendUp=.FALSE.
1028        RecvUp=.FALSE.
1029      ENDIF
1030 
1031      IF (pole_sud) THEN
1032        SendDown=.FALSE.
1033        RecvDown=.FALSE.
1034      ENDIF
1035     
1036      if (Sup.eq.0) then
1037        SendUp=.FALSE.
1038       endif
1039     
1040      if (Sdown.eq.0) then
1041        SendDown=.FALSE.
1042      endif
1043
1044      if (Rup.eq.0) then
1045        RecvUp=.FALSE.
1046      endif
1047     
1048      if (Rdown.eq.0) then
1049        RecvDown=.FALSE.
1050      endif
1051     
1052      IF (SendUp) THEN
1053        call Register_SendField(Field,ijnb_v,ll,jj_begin-jjb_v+1,SUp,MPI_Rank-1,a_request)
1054      ENDIF
1055 
1056      IF (SendDown) THEN
1057        call Register_SendField(Field,ijnb_v,ll,jj_end-SDown+1-jjb_v+1,SDown,MPI_Rank+1,a_request)
1058      ENDIF
1059   
1060 
1061      IF (RecvUp) THEN
1062        call Register_RecvField(Field,ijnb_v,ll,jj_begin-Rup-jjb_v+1,RUp,MPI_Rank-1,a_request)
1063      ENDIF
1064 
1065      IF (RecvDown) THEN
1066        call Register_RecvField(Field,ijnb_v,ll,jj_end+1-jjb_v+1,RDown,MPI_Rank+1,a_request)
1067      ENDIF
1068 
1069    end subroutine Register_Hallo_v
1070   
1071    subroutine SendRequest(a_Request)
1072    USE dimensions_mod
1073      implicit none
1074
1075#ifdef CPP_MPI
1076      include 'mpif.h'
1077#endif
1078
1079      type(request),target :: a_request
1080      type(request_SR),pointer :: Req
1081      type(Hallo),pointer :: PtrHallo
1082      integer :: SizeBuffer
1083      integer :: i,rank,l,ij,Pos,ierr
1084      integer :: offset
1085      real,dimension(:,:),pointer :: Field
1086      integer :: Nb
1087       
1088      do rank=0,MPI_SIZE-1
1089     
1090        Req=>a_Request%RequestSend(rank)
1091       
1092        SizeBuffer=0
1093        do i=1,Req%NbRequest
1094          PtrHallo=>Req%Hallo(i)
1095!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1096          DO l=1,PtrHallo%NbLevel
1097            SizeBuffer=SizeBuffer+PtrHallo%size*iip1
1098          ENDDO
1099!$OMP ENDDO NOWAIT         
1100        enddo
1101     
1102         Req%BufferSize=SizeBuffer
1103         if (Req%NbRequest>0) then
1104       
1105          call allocate_buffer(SizeBuffer,Req%Index,Req%pos)
1106
1107          Pos=Req%Pos
1108          do i=1,Req%NbRequest
1109            PtrHallo=>Req%Hallo(i)
1110            offset=(PtrHallo%offset-1)*iip1+1
1111            Nb=iip1*PtrHallo%size-1
1112            Field=>PtrHallo%Field
1113
1114!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
1115            do l=1,PtrHallo%NbLevel
1116!cdir NODEP
1117              do ij=0,Nb
1118                Buffer(Pos+ij)=Field(Offset+ij,l)
1119              enddo
1120             
1121              Pos=Pos+Nb+1
1122            enddo
1123!$OMP END DO NOWAIT           
1124          enddo
1125   
1126         if (SizeBuffer>0) then
1127!$OMP CRITICAL (MPI)
1128         
1129#ifdef CPP_MPI
1130         call MPI_ISSEND(Buffer(req%Pos),SizeBuffer,MPI_REAL_LMDZ,rank,a_request%tag+1000*omp_rank,     &
1131                         COMM_LMDZ,Req%MSG_Request,ierr)
1132#endif
1133         IF (.NOT.using_mpi) THEN
1134           PRINT *,'Erreur, echange MPI en mode sequentiel !!!'
1135           STOP
1136         ENDIF
1137!         PRINT *,"-------------------------------------------------------------------"
1138!         PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->"
1139!         PRINT *,"Requete envoye au proc :",rank,"tag :",a_request%tag+1000*omp_rank
1140!         PRINT *,"Taille du message :",SizeBuffer,"requete no :",Req%MSG_Request
1141!         PRINT *,"-------------------------------------------------------------------"
1142!$OMP END CRITICAL (MPI)
1143        endif
1144       endif
1145    enddo
1146   
1147           
1148      do rank=0,MPI_SIZE-1
1149         
1150          Req=>a_Request%RequestRecv(rank)
1151          SizeBuffer=0
1152         
1153          do i=1,Req%NbRequest
1154            PtrHallo=>Req%Hallo(i)
1155
1156!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1157            DO l=1,PtrHallo%NbLevel
1158              SizeBuffer=SizeBuffer+PtrHallo%size*iip1
1159            ENDDO
1160!$OMP ENDDO NOWAIT         
1161          enddo
1162         
1163          Req%BufferSize=SizeBuffer
1164         
1165          if (Req%NbRequest>0) then
1166          call allocate_buffer(SizeBuffer,Req%Index,Req%Pos)
1167   
1168          if (SizeBuffer>0) then
1169
1170!$OMP CRITICAL (MPI)
1171
1172#ifdef CPP_MPI
1173             call MPI_IRECV(Buffer(Req%Pos),SizeBuffer,MPI_REAL_LMDZ,rank,a_request%tag+1000*omp_rank,     &
1174                           COMM_LMDZ,Req%MSG_Request,ierr)
1175#endif             
1176             IF (.NOT.using_mpi) THEN
1177               PRINT *,'Erreur, echange MPI en mode sequentiel !!!'
1178               STOP
1179             ENDIF
1180
1181!         PRINT *,"-------------------------------------------------------------------"
1182!         PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->"
1183!         PRINT *,"Requete en attente du proc :",rank,"tag :",a_request%tag+1000*omp_rank
1184!         PRINT *,"Taille du message :",SizeBuffer,"requete no :",Req%MSG_Request
1185!         PRINT *,"-------------------------------------------------------------------"
1186
1187!$OMP END CRITICAL (MPI)
1188          endif
1189        endif
1190     
1191      enddo
1192                       
1193   end subroutine SendRequest
1194   
1195   subroutine WaitRequest(a_Request)
1196   USE dimensions_mod
1197   implicit none
1198   
1199#ifdef CPP_MPI
1200      include 'mpif.h'   
1201#endif
1202     
1203      type(request),target :: a_request
1204      type(request_SR),pointer :: Req
1205      type(Hallo),pointer :: PtrHallo
1206      integer, dimension(2*mpi_size) :: TabRequest
1207#ifdef CPP_MPI
1208      integer, dimension(MPI_STATUS_SIZE,2*mpi_size) :: TabStatus
1209#else
1210      integer, dimension(1,2*mpi_size) :: TabStatus
1211#endif
1212      integer :: NbRequest
1213      integer :: i,rank,pos,ij,l,ierr
1214      integer :: offset
1215      integer :: Nb
1216     
1217     
1218      NbRequest=0
1219      do rank=0,MPI_SIZE-1
1220        Req=>a_request%RequestSend(rank)
1221        if (Req%NbRequest>0 .AND. Req%BufferSize > 0) then
1222          NbRequest=NbRequest+1
1223          TabRequest(NbRequest)=Req%MSG_Request
1224        endif
1225      enddo
1226     
1227      do rank=0,MPI_SIZE-1
1228        Req=>a_request%RequestRecv(rank)
1229        if (Req%NbRequest>0 .AND. Req%BufferSize > 0 ) then
1230          NbRequest=NbRequest+1
1231          TabRequest(NbRequest)=Req%MSG_Request
1232        endif
1233      enddo
1234     
1235      if (NbRequest>0) then
1236!$OMP CRITICAL (MPI)
1237!        PRINT *,"-------------------------------------------------------------------"
1238!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
1239!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
1240#ifdef CPP_MPI
1241        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
1242#endif
1243!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
1244!        PRINT *,"-------------------------------------------------------------------"
1245!$OMP END CRITICAL (MPI)
1246      endif
1247      do rank=0,MPI_Size-1
1248        Req=>a_request%RequestRecv(rank)
1249        if (Req%NbRequest>0) then
1250          Pos=Req%Pos
1251          do i=1,Req%NbRequest
1252            PtrHallo=>Req%Hallo(i)
1253            offset=(PtrHallo%offset-1)*iip1+1
1254            Nb=iip1*PtrHallo%size-1
1255
1256!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
1257            do l=1,PtrHallo%NbLevel
1258!cdir NODEP
1259              do ij=0,Nb
1260                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
1261              enddo
1262
1263              Pos=Pos+Nb+1
1264            enddo
1265!$OMP ENDDO NOWAIT         
1266          enddo
1267        endif
1268      enddo
1269     
1270      do rank=0,MPI_SIZE-1
1271        Req=>a_request%RequestSend(rank)
1272        if (Req%NbRequest>0) then
1273          call deallocate_buffer(Req%Index)
1274          Req%NbRequest=0
1275        endif
1276      enddo
1277             
1278      do rank=0,MPI_SIZE-1
1279        Req=>a_request%RequestRecv(rank)
1280        if (Req%NbRequest>0) then
1281          call deallocate_buffer(Req%Index)
1282          Req%NbRequest=0
1283        endif
1284      enddo
1285     
1286      a_request%tag=1
1287    end subroutine WaitRequest
1288     
1289   subroutine WaitSendRequest(a_Request)
1290   USE dimensions_mod
1291   implicit none
1292   
1293#ifdef CPP_MPI
1294      include 'mpif.h'   
1295#endif     
1296      type(request),target :: a_request
1297      type(request_SR),pointer :: Req
1298      type(Hallo),pointer :: PtrHallo
1299      integer, dimension(mpi_size) :: TabRequest
1300#ifdef CPP_MPI
1301      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
1302#else
1303      integer, dimension(1,mpi_size) :: TabStatus
1304#endif
1305      integer :: NbRequest
1306      integer :: i,rank,pos,ij,l,ierr
1307      integer :: offset
1308     
1309     
1310      NbRequest=0
1311      do rank=0,MPI_SIZE-1
1312        Req=>a_request%RequestSend(rank)
1313        if (Req%NbRequest>0) then
1314          NbRequest=NbRequest+1
1315          TabRequest(NbRequest)=Req%MSG_Request
1316        endif
1317      enddo
1318     
1319
1320      if (NbRequest>0 .AND. Req%BufferSize > 0 ) THEN
1321!$OMP CRITICAL (MPI)     
1322!        PRINT *,"-------------------------------------------------------------------"
1323!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
1324!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
1325#ifdef CPP_MPI
1326        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
1327#endif
1328!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
1329!        PRINT *,"-------------------------------------------------------------------"
1330
1331!$OMP END CRITICAL (MPI)
1332      endif     
1333     
1334      do rank=0,MPI_SIZE-1
1335        Req=>a_request%RequestSend(rank)
1336        if (Req%NbRequest>0) then
1337          call deallocate_buffer(Req%Index)
1338          Req%NbRequest=0
1339        endif
1340      enddo
1341             
1342      a_request%tag=1
1343    end subroutine WaitSendRequest
1344   
1345   subroutine WaitRecvRequest(a_Request)
1346   USE dimensions_mod
1347   implicit none
1348   
1349#ifdef CPP_MPI
1350      include 'mpif.h'   
1351#endif
1352     
1353      type(request),target :: a_request
1354      type(request_SR),pointer :: Req
1355      type(Hallo),pointer :: PtrHallo
1356      integer, dimension(mpi_size) :: TabRequest
1357#ifdef CPP_MPI
1358      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
1359#else
1360      integer, dimension(1,mpi_size) :: TabStatus
1361#endif
1362      integer :: NbRequest
1363      integer :: i,rank,pos,ij,l,ierr
1364      integer :: offset,Nb
1365     
1366     
1367      NbRequest=0
1368     
1369      do rank=0,MPI_SIZE-1
1370        Req=>a_request%RequestRecv(rank)
1371        if (Req%NbRequest>0 .AND. Req%BufferSize > 0 ) then
1372          NbRequest=NbRequest+1
1373          TabRequest(NbRequest)=Req%MSG_Request
1374        endif
1375      enddo
1376     
1377     
1378      if (NbRequest>0) then
1379!$OMP CRITICAL (MPI)     
1380!        PRINT *,"-------------------------------------------------------------------"
1381!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
1382!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
1383#ifdef CPP_MPI
1384        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
1385#endif
1386!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
1387!        PRINT *,"-------------------------------------------------------------------"
1388!$OMP END CRITICAL (MPI)     
1389      endif
1390     
1391      do rank=0,MPI_Size-1
1392        Req=>a_request%RequestRecv(rank)
1393        if (Req%NbRequest>0) then
1394          Pos=Req%Pos
1395          do i=1,Req%NbRequest
1396            PtrHallo=>Req%Hallo(i)
1397            offset=(PtrHallo%offset-1)*iip1+1
1398            Nb=iip1*PtrHallo%size-1
1399!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
1400            do l=1,PtrHallo%NbLevel
1401!cdir NODEP
1402              do ij=0,Nb
1403                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
1404              enddo
1405                 Pos=Pos+Nb+1
1406            enddo
1407!$OMP END DO NOWAIT
1408          enddo
1409        endif
1410      enddo
1411     
1412           
1413      do rank=0,MPI_SIZE-1
1414        Req=>a_request%RequestRecv(rank)
1415        if (Req%NbRequest>0) then
1416          call deallocate_buffer(Req%Index)
1417          Req%NbRequest=0
1418        endif
1419      enddo
1420     
1421      a_request%tag=1
1422    end subroutine WaitRecvRequest
1423   
1424   
1425   
1426    subroutine CopyField(FieldS,FieldR,ij,ll,jj_Nb_New)
1427    USE dimensions_mod
1428 
1429      implicit none
1430   
1431    INTEGER :: ij,ll,l
1432    REAL, dimension(ij,ll) :: FieldS
1433    REAL, dimension(ij,ll) :: FieldR
1434    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
1435    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
1436   
1437    integer ::i,jje,jjb,ijb,ije
1438   
1439    jj_begin_New(0)=1
1440    jj_End_New(0)=jj_Nb_New(0)
1441    do i=1,MPI_Size-1
1442      jj_begin_New(i)=jj_end_New(i-1)+1
1443      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
1444    enddo
1445   
1446    jjb=max(jj_begin,jj_begin_new(MPI_Rank))
1447    jje=min(jj_end,jj_end_new(MPI_Rank))
1448    if (ij==ip1jm) jje=min(jje,jjm)
1449
1450    if (jje >= jjb) then
1451      ijb=(jjb-1)*iip1+1
1452      ije=jje*iip1
1453
1454!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1455      do l=1,ll
1456        FieldR(ijb:ije,l)=FieldS(ijb:ije,l)
1457      enddo
1458!$OMP ENDDO NOWAIT
1459    endif
1460
1461
1462  end subroutine CopyField   
1463
1464  subroutine CopyFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down)
1465  USE dimensions_mod
1466 
1467      implicit none
1468   
1469    INTEGER :: ij,ll,Up,Down
1470    REAL, dimension(ij,ll) :: FieldS
1471    REAL, dimension(ij,ll) :: FieldR
1472    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
1473    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
1474
1475    integer ::i,jje,jjb,ijb,ije,l
1476
1477     
1478    jj_begin_New(0)=1
1479    jj_End_New(0)=jj_Nb_New(0)
1480    do i=1,MPI_Size-1
1481      jj_begin_New(i)=jj_end_New(i-1)+1
1482      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
1483    enddo
1484
1485       
1486    jjb=max(jj_begin,jj_begin_new(MPI_Rank)-Up)
1487    jje=min(jj_end,jj_end_new(MPI_Rank)+Down)
1488    if (ij==ip1jm) jje=min(jje,jjm)
1489   
1490   
1491    if (jje >= jjb) then
1492      ijb=(jjb-1)*iip1+1
1493      ije=jje*iip1
1494
1495!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1496      do l=1,ll
1497        FieldR(ijb:ije,l)=FieldS(ijb:ije,l)
1498      enddo
1499!$OMP ENDDO NOWAIT
1500
1501    endif
1502   end subroutine CopyFieldHallo       
1503
1504   subroutine Gather_field_u(field_loc,field_glo,ll)
1505   USE dimensions_mod
1506   implicit none
1507     integer :: ll
1508     real :: field_loc(ijb_u:ije_u,ll)
1509     real :: field_glo(ip1jmp1,ll)
1510     type(request) :: request_gather
1511     integer       :: l
1512
1513
1514!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1515     DO l=1,ll
1516       field_glo(ij_begin:ij_end,l)=field_loc(ij_begin:ij_end,l)
1517     ENDDO
1518     
1519     call register_SwapField(field_glo,field_glo,ip1jmp1,ll,distrib_gather%jj_nb_para,request_gather)
1520     call SendRequest(request_gather)
1521!$OMP BARRIER
1522     call WaitRequest(request_gather)       
1523!$OMP BARRIER
1524
1525    end subroutine Gather_field_u
1526       
1527   subroutine Gather_field_v(field_loc,field_glo,ll)
1528   USE dimensions_mod
1529   implicit none
1530     integer :: ll
1531     real :: field_loc(ijb_v:ije_v,ll)
1532     real :: field_glo(ip1jm,ll)
1533     type(request) :: request_gather
1534     integer :: ijb,ije
1535     integer       :: l
1536     
1537   
1538     ijb=ij_begin
1539     ije=ij_end
1540     if (pole_sud) ije=ij_end-iip1
1541       
1542!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1543     DO l=1,ll
1544       field_glo(ijb:ije,l)=field_loc(ijb:ije,l)
1545     ENDDO
1546     
1547     call register_SwapField(field_glo,field_glo,ip1jm,ll,distrib_gather%jj_nb_para,request_gather)
1548     call SendRequest(request_gather)
1549!$OMP BARRIER
1550     call WaitRequest(request_gather)       
1551!$OMP BARRIER
1552
1553    end subroutine Gather_field_v
1554     
1555   subroutine Scatter_field_u(field_glo,field_loc,ll)
1556   USE dimensions_mod
1557   implicit none
1558     integer :: ll
1559     real :: field_glo(ip1jmp1,ll)
1560     real :: field_loc(ijb_u:ije_u,ll)
1561     type(request) :: request_gather
1562     TYPE(distrib) :: distrib_swap
1563     integer       :: l
1564     
1565!$OMP BARRIER
1566!$OMP MASTER     
1567     call get_current_distrib(distrib_swap)
1568     call set_Distrib(distrib_gather)
1569!$OMP END MASTER
1570!$OMP BARRIER
1571 
1572     call register_SwapField(field_glo,field_glo,ip1jmp1,ll,distrib_swap%jj_nb_para,request_gather)
1573     call SendRequest(request_gather)
1574!$OMP BARRIER
1575     call WaitRequest(request_gather)       
1576!$OMP BARRIER
1577!$OMP MASTER     
1578     call set_Distrib(distrib_swap)
1579!$OMP END MASTER
1580!$OMP BARRIER
1581
1582!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1583       DO l=1,ll
1584         field_loc(ij_begin:ij_end,l)=field_glo(ij_begin:ij_end,l)
1585       ENDDO
1586
1587    end subroutine Scatter_field_u
1588
1589   subroutine Scatter_field_v(field_glo,field_loc,ll)
1590   USE dimensions_mod
1591   implicit none
1592     integer :: ll
1593     real :: field_glo(ip1jmp1,ll)
1594     real :: field_loc(ijb_v:ije_v,ll)
1595     type(request) :: request_gather
1596     TYPE(distrib) :: distrib_swap
1597     integer       :: ijb,ije,l
1598     
1599
1600!$OMP BARRIER
1601!$OMP MASTER     
1602     call get_current_distrib(distrib_swap)
1603     call set_Distrib(distrib_gather)
1604!$OMP END MASTER
1605!$OMP BARRIER
1606     call register_SwapField(field_glo,field_glo,ip1jm,ll,distrib_swap%jj_nb_para,request_gather)
1607     call SendRequest(request_gather)
1608!$OMP BARRIER
1609     call WaitRequest(request_gather)       
1610!$OMP BARRIER
1611!$OMP MASTER
1612     call set_Distrib(distrib_swap)
1613!$OMP END MASTER
1614!$OMP BARRIER
1615     ijb=ij_begin
1616     ije=ij_end
1617     if (pole_sud) ije=ij_end-iip1
1618     
1619!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1620       DO l=1,ll
1621         field_loc(ijb:ije,l)=field_glo(ijb:ije,l)
1622       ENDDO
1623
1624    end subroutine Scatter_field_v
1625             
1626end module mod_Hallo
1627   
Note: See TracBrowser for help on using the repository browser.