source: LMDZ4/trunk/libf/dyn3dpar/mod_hallo.F90 @ 785

Last change on this file since 785 was 764, checked in by Laurent Fairhead, 17 years ago

Merge entre la version V3_conv et le HEAD
YM, JG, LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 17.7 KB
RevLine 
[630]1module mod_Hallo
2USE parallel
3implicit none
4!  include 'mpif.h'
[764]5  integer, parameter :: MaxRequest=200
[630]6  integer, parameter :: MaxProc=80
7  integer, parameter :: MaxBufferSize=1024*1024*16
8  integer, parameter :: ListSize=1000
9 
10  integer,save       :: MaxBufferSize_Used
11 
12    real,save,pointer,dimension(:) :: Buffer
13!   pointer (Pbuffer,Buffer)
14!   real,dimension(:) :: Buffer
15   integer,dimension(Listsize) :: Buffer_Pos
16   integer :: Index_Pos
17   
18  type Hallo
19    real, dimension(:,:),pointer :: Field
20    integer :: offset
21    integer :: size
22    integer :: NbLevel
23    integer :: Stride
24  end type Hallo
25 
26  type request_SR
27    integer :: NbRequest=0
28!    real,dimension(:),allocatable :: Buffer
29!    real,dimension(:),pointer :: Buffer
30    integer :: Pos
31    integer :: Index
32    type(Hallo),dimension(MaxRequest) :: Hallo
33    integer :: MSG_Request
34  end type request_SR
35
36  type request
37    type(request_SR),dimension(0:MaxProc-1) :: RequestSend
38    type(request_SR),dimension(0:MaxProc-1) :: RequestRecv
39    integer :: tag=1
40  end type request
41 
42   
43  contains
44
45  subroutine Init_mod_hallo(MPI_Buffer)
46    implicit none
47    real,dimension(:),target :: MPI_BUFFER   
48    integer i
49    Index_Pos=1
50    Buffer_Pos(Index_Pos)=1
51   
52    Buffer=>MPI_Buffer
53    MaxBufferSize_Used=0
54  end subroutine init_mod_hallo
55     
56  subroutine allocate_buffer(Size,Index,Pos)
57  implicit none
58    integer :: Size
59    integer :: Index
60    integer :: Pos
61
62    if (Buffer_pos(Index_pos)+Size>MaxBufferSize_Used) MaxBufferSize_Used=Buffer_pos(Index_pos)+Size 
63    if (Buffer_pos(Index_pos)+Size>MaxBufferSize) then
64      print *,'STOP :: La taille de MaxBufferSize dans mod_hallo.F90 est trop petite !!!!'
65      stop
66    endif
67   
68    if (Index_pos>=ListSize) then
69      print *,'STOP :: La taille de ListSize dans mod_hallo.F90 est trop petite !!!!'
70      stop
71    endif
72     
73    Pos=Buffer_Pos(Index_Pos)
74    Buffer_Pos(Index_pos+1)=Buffer_Pos(Index_Pos)+Size
75    Index_Pos=Index_Pos+1
76    Index=Index_Pos
77   
78  end subroutine allocate_buffer
79     
80  subroutine deallocate_buffer(Index)
81  implicit none
82    integer :: Index
83   
84    Buffer_Pos(Index)=-1
85   
86    do while (Buffer_Pos(Index_Pos)==-1 .and. Index_Pos>1)
87      Index_Pos=Index_Pos-1
88    end do
89
90  end subroutine deallocate_buffer 
91 
92  subroutine SetTag(a_request,tag)
93  implicit none
94    type(request):: a_request
95    integer :: tag
96   
97    a_request%tag=tag
98  end subroutine SetTag
99 
100 
101  subroutine Init_Hallo(Field,Stride,NbLevel,offset,size,NewHallo)
102    integer :: Stride
103    integer :: NbLevel
104    integer :: size
105    integer :: offset
106    real, dimension(Stride,NbLevel),target :: Field
107    type(Hallo) :: NewHallo
108   
109    NewHallo%Field=>Field
110    NewHallo%Stride=Stride
111    NewHallo%NbLevel=NbLevel
112    NewHallo%size=size
113    NewHallo%offset=offset
114   
115   
116  end subroutine Init_Hallo
117 
118  subroutine Register_SendField(Field,ij,ll,offset,size,target,a_request)
119  implicit none
120
121#include "dimensions90.h"
122#include "paramet90.h"   
123    include 'mpif.h'
124   
125      INTEGER :: ij,ll,offset,size,target
126      REAL, dimension(ij,ll) :: Field
127      type(request),target :: a_request
128      type(request_SR),pointer :: Ptr_request
129
130      Ptr_Request=>a_request%RequestSend(target)
131      Ptr_Request%NbRequest=Ptr_Request%NbRequest+1
132      if (Ptr_Request%NbRequest>=MaxRequest) then
133        print *,'STOP :: La taille de MaxRequest dans mod_hallo.F90 est trop petite !!!!'
134        stop
135      endif     
136      call Init_Hallo(Field,ij,ll,offset,size,Ptr_request%Hallo(Ptr_Request%NbRequest))
137     
138   end subroutine Register_SendField     
139     
140  subroutine Register_RecvField(Field,ij,ll,offset,size,target,a_request)
141  implicit none
142
143#include "dimensions90.h"
144#include "paramet90.h"   
145    include 'mpif.h'
146   
147      INTEGER :: ij,ll,offset,size,target
148      REAL, dimension(ij,ll) :: Field
149      type(request),target :: a_request
150      type(request_SR),pointer :: Ptr_request
151
152      Ptr_Request=>a_request%RequestRecv(target)
153      Ptr_Request%NbRequest=Ptr_Request%NbRequest+1
154     
155      if (Ptr_Request%NbRequest>=MaxRequest) then
156        print *,'STOP :: La taille de MaxRequest dans mod_hallo.F90 est trop petite !!!!'
157        stop
158      endif   
159           
160      call Init_Hallo(Field,ij,ll,offset,size,Ptr_request%Hallo(Ptr_Request%NbRequest))
161
162     
163   end subroutine Register_RecvField     
164 
165  subroutine Register_SwapField(FieldS,FieldR,ij,ll,jj_Nb_New,a_request)
166 
167      implicit none
168#include "dimensions90.h"
169#include "paramet90.h"   
170    include 'mpif.h'
171   
172    INTEGER :: ij,ll
173    REAL, dimension(ij,ll) :: FieldS
174    REAL, dimension(ij,ll) :: FieldR
175    type(request) :: a_request
176    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
177    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
178   
179    integer ::i,jje,jjb
180   
181    jj_begin_New(0)=1
182    jj_End_New(0)=jj_Nb_New(0)
183    do i=1,MPI_Size-1
184      jj_begin_New(i)=jj_end_New(i-1)+1
185      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
186    enddo
187   
188    do i=0,MPI_Size-1
189      if (i /= MPI_Rank) then
190        jjb=max(jj_begin_new(i),jj_begin)
191        jje=min(jj_end_new(i),jj_end)
192       
193        if (ij==ip1jm .and. jje==jjp1) jje=jjm
194       
195        if (jje >= jjb) then
196          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
197        endif
198       
199        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
200        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
201       
202        if (ij==ip1jm .and. jje==jjp1) jje=jjm
203       
204        if (jje >= jjb) then
205          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
206        endif
207       
208      endif
209    enddo
210   
211  end subroutine Register_SwapField   
212 
213 
214    subroutine Register_SwapFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down,a_request)
215 
216      implicit none
217#include "dimensions90.h"
218#include "paramet90.h"   
219    include 'mpif.h'
220   
221    INTEGER :: ij,ll,Up,Down
222    REAL, dimension(ij,ll) :: FieldS
223    REAL, dimension(ij,ll) :: FieldR
224    type(request) :: a_request
225    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
226    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
227   
228    integer ::i,jje,jjb
229   
230    jj_begin_New(0)=1
231    jj_End_New(0)=jj_Nb_New(0)
232    do i=1,MPI_Size-1
233      jj_begin_New(i)=jj_end_New(i-1)+1
234      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
235    enddo
236   
237    do i=0,MPI_Size-1
238      jj_begin_New(i)=max(1,jj_begin_New(i)-Up)
239      jj_end_New(i)=min(jjp1,jj_end_new(i)+Down)
240    enddo
241   
242    do i=0,MPI_Size-1
243      if (i /= MPI_Rank) then
244        jjb=max(jj_begin_new(i),jj_begin)
245        jje=min(jj_end_new(i),jj_end)
246       
247        if (ij==ip1jm .and. jje==jjp1) jje=jjm
248       
249        if (jje >= jjb) then
250          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
251        endif
252       
253        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
254        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
255       
256        if (ij==ip1jm .and. jje==jjp1) jje=jjm
257       
258        if (jje >= jjb) then
259          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
260        endif
261       
262      endif
263    enddo
264   
265  end subroutine Register_SwapFieldHallo
266 
267  subroutine Register_Hallo(Field,ij,ll,RUp,Rdown,SUp,SDown,a_request)
268 
269      implicit none
270#include "dimensions90.h"
271#include "paramet90.h"   
272    include 'mpif.h'
273   
274      INTEGER :: ij,ll
275      REAL, dimension(ij,ll) :: Field
276      INTEGER :: Sup,Sdown,rup,rdown
277      type(request) :: a_request
278      type(Hallo),pointer :: PtrHallo
279      LOGICAL :: SendUp,SendDown
280      LOGICAL :: RecvUp,RecvDown
281   
282 
283      SendUp=.TRUE.
284      SendDown=.TRUE.
285      RecvUp=.TRUE.
286      RecvDown=.TRUE.
287       
288      IF (pole_nord) THEN
289        SendUp=.FALSE.
290        RecvUp=.FALSE.
291      ENDIF
292 
293      IF (pole_sud) THEN
294        SendDown=.FALSE.
295        RecvDown=.FALSE.
296      ENDIF
297     
298      if (Sup.eq.0) then
299        SendUp=.FALSE.
300       endif
301     
302      if (Sdown.eq.0) then
303        SendDown=.FALSE.
304      endif
305
306      if (Rup.eq.0) then
307        RecvUp=.FALSE.
308      endif
309     
310      if (Rdown.eq.0) then
311        RecvDown=.FALSE.
312      endif
313     
314      IF (SendUp) THEN
315        call Register_SendField(Field,ij,ll,jj_begin,SUp,MPI_Rank-1,a_request)
316      ENDIF
317 
318      IF (SendDown) THEN
319        call Register_SendField(Field,ij,ll,jj_end-SDown+1,SDown,MPI_Rank+1,a_request)
320      ENDIF
321   
322 
323      IF (RecvUp) THEN
324        call Register_RecvField(Field,ij,ll,jj_begin-Rup,RUp,MPI_Rank-1,a_request)
325      ENDIF
326 
327      IF (RecvDown) THEN
328        call Register_RecvField(Field,ij,ll,jj_end+1,RDown,MPI_Rank+1,a_request)
329      ENDIF
330 
331    end subroutine Register_Hallo
332   
333    subroutine SendRequest(a_Request)
334      implicit none
335
336#include "dimensions90.h"
337#include "paramet90.h"
338      include 'mpif.h'
339
340      type(request),target :: a_request
341      type(request_SR),pointer :: Req
342      type(Hallo),pointer :: PtrHallo
343      integer :: SizeBuffer
344      integer :: i,rank,l,ij,Pos,ierr
345      integer :: offset
346!      real,dimension(:),pointer :: Buffer
347      real,dimension(:,:),pointer :: Field
348      integer :: Nb
349       
350      do rank=0,MPI_SIZE-1
351     
352        Req=>a_Request%RequestSend(rank)
353       
354        SizeBuffer=0
355        do i=1,Req%NbRequest
356          PtrHallo=>Req%Hallo(i)
357          SizeBuffer=SizeBuffer+PtrHallo%size*PtrHallo%NbLevel*iip1
358        enddo
359     
360        if (SizeBuffer>0) then
361       
362!          allocate(Req%Buffer(SizeBuffer))
363          call allocate_buffer(SizeBuffer,Req%Index,Req%pos)
364
365          Pos=Req%Pos
366!          Buffer=>req%Buffer
367          do i=1,Req%NbRequest
368            PtrHallo=>Req%Hallo(i)
369            offset=(PtrHallo%offset-1)*iip1+1
370            Nb=iip1*PtrHallo%size-1
371            Field=>PtrHallo%Field
372           
373            do l=1,PtrHallo%NbLevel
374!cdir NODEP
375              do ij=0,Nb
376                Buffer(Pos+ij)=Field(Offset+ij,l)
377              enddo
378!              Buffer(Pos:Pos+Nb)=Field(offset:offset+Nb,l)
379             
380              Pos=Pos+Nb+1
381            enddo
382           
383          enddo
384   
385!         print *, 'process',MPI_RANK,'ISSEND: requette ',a_request%tag,'au process',rank,'de taille',SizeBuffer
386!         call MPI_ISSEND(Req%Buffer,SizeBuffer,MPI_REAL8,rank,a_request%tag,     &
[764]387!                         COMM_LMDZ,Req%MSG_Request,ierr)
[630]388         call MPI_ISSEND(Buffer(req%Pos),SizeBuffer,MPI_REAL8,rank,a_request%tag,     &
[764]389                         COMM_LMDZ,Req%MSG_Request,ierr)
[630]390
391        endif
392
393    enddo
394   
395           
396      do rank=0,MPI_SIZE-1
397         
398          Req=>a_Request%RequestRecv(rank)
399          SizeBuffer=0
400         
401          do i=1,Req%NbRequest
402            PtrHallo=>Req%Hallo(i)
403            SizeBuffer=SizeBuffer+PtrHallo%size*PtrHallo%NbLevel*iip1
404          enddo
405       
406          if (SizeBuffer>0) then
407!            allocate(Req%Buffer(SizeBuffer))
408             call allocate_buffer(SizeBuffer,Req%Index,Req%Pos)
409!            print *, 'process',MPI_RANK,'IRECV: requette ',a_request%tag,'au process',rank,'de taille',SizeBuffer
410           
411!           call MPI_IRECV(Req%Buffer,SizeBuffer,MPI_REAL8,rank,a_request%tag,     &
[764]412!                           COMM_LMDZ,Req%MSG_Request,ierr)
[630]413            call MPI_IRECV(Buffer(Req%Pos),SizeBuffer,MPI_REAL8,rank,a_request%tag,     &
[764]414                           COMM_LMDZ,Req%MSG_Request,ierr)
[630]415
416          endif
417     
418      enddo
419                       
420   end subroutine SendRequest
421   
422   subroutine WaitRequest(a_Request)
423   implicit none
424   
425#include "dimensions90.h"
426#include "paramet90.h"
427      include 'mpif.h'   
428     
429      type(request),target :: a_request
430      type(request_SR),pointer :: Req
431      type(Hallo),pointer :: PtrHallo
[764]432      integer, dimension(2*mpi_size) :: TabRequest
433      integer, dimension(MPI_STATUS_SIZE,2*mpi_size) :: TabStatus
[630]434      integer :: NbRequest
435      integer :: i,rank,pos,ij,l,ierr
436      integer :: offset
437      integer :: Nb
438     
439     
440      NbRequest=0
441      do rank=0,MPI_SIZE-1
442        Req=>a_request%RequestSend(rank)
443        if (Req%NbRequest>0) then
444          NbRequest=NbRequest+1
445          TabRequest(NbRequest)=Req%MSG_Request
446        endif
447      enddo
448     
449      do rank=0,MPI_SIZE-1
450        Req=>a_request%RequestRecv(rank)
451        if (Req%NbRequest>0) then
452          NbRequest=NbRequest+1
453          TabRequest(NbRequest)=Req%MSG_Request
454        endif
455      enddo
456     
457      if (NbRequest>0) call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
458     
459      do rank=0,MPI_Size-1
460        Req=>a_request%RequestRecv(rank)
461        if (Req%NbRequest>0) then
462          Pos=Req%Pos
463          do i=1,Req%NbRequest
464            PtrHallo=>Req%Hallo(i)
465            offset=(PtrHallo%offset-1)*iip1+1
466            Nb=iip1*PtrHallo%size-1
467           
468            do l=1,PtrHallo%NbLevel
469!cdir NODEP
470              do ij=0,Nb
471                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
472              enddo
473!              PtrHallo%Field(offset:offset+Nb,l)=Buffer(Pos:Pos+Nb)
474!             do ij=offset,offset+iip1*PtrHallo%size-1
475!                PtrHallo%Field(ij,l)=Buffer(Pos)
476!                Pos=Pos+1
477!              enddo
478              Pos=Pos+Nb+1
479            enddo
480           
481          enddo
482        endif
483      enddo
484     
485      do rank=0,MPI_SIZE-1
486        Req=>a_request%RequestSend(rank)
487        if (Req%NbRequest>0) then
488          call deallocate_buffer(Req%Index)
489          Req%NbRequest=0
490        endif
491      enddo
492             
493      do rank=0,MPI_SIZE-1
494        Req=>a_request%RequestRecv(rank)
495        if (Req%NbRequest>0) then
496          call deallocate_buffer(Req%Index)
497          Req%NbRequest=0
498        endif
499      enddo
500     
501      a_request%tag=1
502    end subroutine WaitRequest
503     
504   subroutine WaitSendRequest(a_Request)
505   implicit none
506   
507#include "dimensions90.h"
508#include "paramet90.h"
509      include 'mpif.h'   
510     
511      type(request),target :: a_request
512      type(request_SR),pointer :: Req
513      type(Hallo),pointer :: PtrHallo
[764]514      integer, dimension(mpi_size) :: TabRequest
515      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
[630]516      integer :: NbRequest
517      integer :: i,rank,pos,ij,l,ierr
518      integer :: offset
519     
520     
521      NbRequest=0
522      do rank=0,MPI_SIZE-1
523        Req=>a_request%RequestSend(rank)
524        if (Req%NbRequest>0) then
525          NbRequest=NbRequest+1
526          TabRequest(NbRequest)=Req%MSG_Request
527        endif
528      enddo
529     
530
531      if (NbRequest>0) call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
532     
533     
534      do rank=0,MPI_SIZE-1
535        Req=>a_request%RequestSend(rank)
536        if (Req%NbRequest>0) then
537          call deallocate_buffer(Req%Index)
538          Req%NbRequest=0
539        endif
540      enddo
541             
542      a_request%tag=1
543    end subroutine WaitSendRequest
544   
545   subroutine WaitRecvRequest(a_Request)
546   implicit none
547   
548#include "dimensions90.h"
549#include "paramet90.h"
550      include 'mpif.h'   
551     
552      type(request),target :: a_request
553      type(request_SR),pointer :: Req
554      type(Hallo),pointer :: PtrHallo
[764]555      integer, dimension(mpi_size) :: TabRequest
556      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
[630]557      integer :: NbRequest
558      integer :: i,rank,pos,ij,l,ierr
559      integer :: offset,Nb
560     
561     
562      NbRequest=0
563     
564      do rank=0,MPI_SIZE-1
565        Req=>a_request%RequestRecv(rank)
566        if (Req%NbRequest>0) then
567          NbRequest=NbRequest+1
568          TabRequest(NbRequest)=Req%MSG_Request
569        endif
570      enddo
571     
572     
573      if (NbRequest>0) call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
574     
575      do rank=0,MPI_Size-1
576        Req=>a_request%RequestRecv(rank)
577        if (Req%NbRequest>0) then
578          Pos=Req%Pos
579          do i=1,Req%NbRequest
580            PtrHallo=>Req%Hallo(i)
581            offset=(PtrHallo%offset-1)*iip1+1
582            Nb=iip1*PtrHallo%size-1
583           
584            do l=1,PtrHallo%NbLevel
585!cdir NODEP
586              do ij=0,Nb
587                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
588              enddo
589                 Pos=Pos+Nb+1
590            enddo
591          enddo
592        endif
593      enddo
594     
595           
596      do rank=0,MPI_SIZE-1
597        Req=>a_request%RequestRecv(rank)
598        if (Req%NbRequest>0) then
599          call deallocate_buffer(Req%Index)
600          Req%NbRequest=0
601        endif
602      enddo
603     
604      a_request%tag=1
605    end subroutine WaitRecvRequest
606   
607   
608   
609    subroutine CopyField(FieldS,FieldR,ij,ll,jj_Nb_New)
610 
611      implicit none
612#include "dimensions90.h"
613#include "paramet90.h"   
614    include 'mpif.h'
615   
616    INTEGER :: ij,ll
617    REAL, dimension(ij,ll) :: FieldS
618    REAL, dimension(ij,ll) :: FieldR
619    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
620    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
621   
622    integer ::i,jje,jjb,ijb,ije
623   
624    jj_begin_New(0)=1
625    jj_End_New(0)=jj_Nb_New(0)
626    do i=1,MPI_Size-1
627      jj_begin_New(i)=jj_end_New(i-1)+1
628      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
629    enddo
630   
631    jjb=max(jj_begin,jj_begin_new(MPI_Rank))
632    jje=min(jj_end,jj_end_new(MPI_Rank))
633    if (ij==ip1jm) jje=min(jje,jjm)
634
635    if (jje >= jjb) then
636      ijb=(jjb-1)*iip1+1
637      ije=jje*iip1
638      FieldR(ijb:ije,1:ll)=FieldS(ijb:ije,1:ll)
639    endif
640
641  end subroutine CopyField   
642
643  subroutine CopyFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down)
644 
645      implicit none
646#include "dimensions90.h"
647#include "paramet90.h"   
648    include 'mpif.h'
649   
650    INTEGER :: ij,ll,Up,Down
651    REAL, dimension(ij,ll) :: FieldS
652    REAL, dimension(ij,ll) :: FieldR
653    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
654    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
655
656    integer ::i,jje,jjb,ijb,ije
657
658     
659    jj_begin_New(0)=1
660    jj_End_New(0)=jj_Nb_New(0)
661    do i=1,MPI_Size-1
662      jj_begin_New(i)=jj_end_New(i-1)+1
663      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
664    enddo
665
666       
667    jjb=max(jj_begin,jj_begin_new(MPI_Rank)-Up)
668    jje=min(jj_end,jj_end_new(MPI_Rank)+Down)
669    if (ij==ip1jm) jje=min(jje,jjm)
670   
671   
672    if (jje >= jjb) then
673      ijb=(jjb-1)*iip1+1
674      ije=jje*iip1
675      FieldR(ijb:ije,1:ll)=FieldS(ijb:ije,1:ll)
676    endif
677   end subroutine CopyFieldHallo       
678         
679end module mod_Hallo
680   
Note: See TracBrowser for help on using the repository browser.