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

Last change on this file since 1821 was 1810, checked in by Ehouarn Millour, 11 years ago

Updating makelmdz and create_makelmdz for :

1 parallelism
2 compilation with various versions of orchidee
3 compilation in 1D mode
4 some cleaning

Also some updates in dyn3dmem:
1 allocate_field_mod.f90 renamed allocate_field_mod.F90
2 module dimensions renamed dimensions_mod
3 module allocate_field renamed allocate_field

FH

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