source: LMDZ5/branches/LMDZ5-DOFOCO/libf/dyn3dmem/mod_hallo.F90 @ 3575

Last change on this file since 3575 was 1632, checked in by Laurent Fairhead, 12 years ago

Import initial du répertoire dyn3dmem

Attention! ceci n'est qu'une version préliminaire du code "basse mémoire":
le code contenu dans ce répertoire est basé sur la r1320 et a donc besoin
d'être mis à jour par rapport à la dynamique parallèle d'aujourd'hui.
Ce code est toutefois mis à disposition pour circonvenir à des problèmes
de mémoire que certaines configurations du modèle pourraient rencontrer.
Dans l'état, il compile et tourne sur vargas et au CCRT


Initial import of dyn3dmem

Warning! this is just a preliminary version of the memory light code:
it is based on r1320 of the code and thus needs to be updated before
it can replace the present dyn3dpar code. It is nevertheless put at your
disposal to circumvent some memory problems some LMDZ configurations may
encounter. In its present state, it will compile and run on vargas and CCRT

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