source: LMDZ5/trunk/libf/dyn3dmem/mod_hallo.F90 @ 5435

Last change on this file since 5435 was 2620, checked in by fhourdin, 8 years ago

Correction pour compenser un bug dans la bibliothèque MPI sur ada.

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