source: LMDZ6/trunk/libf/dyn3dmem/mod_hallo.f90 @ 5322

Last change on this file since 5322 was 5285, checked in by abarral, 3 weeks ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

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