source: LMDZ4/trunk/libf/dyn3dpar/parallel.F90 @ 792

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

Modifications suite a la transformation des fichiers include pour
qu'ils soient compatibles a la fois au format fixe et au format libre
Un bon nombre de fichier *.inc du coup disparaissent
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 11.1 KB
RevLine 
[630]1  module parallel
2 
3    integer, save :: mpi_size
4    integer, save :: mpi_rank
[764]5    integer, save :: COMM_LMDZ
[630]6    integer, save :: jj_begin
7    integer, save :: jj_end
8    integer, save :: jj_nb
9    integer, save :: ij_begin
10    integer, save :: ij_end
11    logical, save :: pole_nord
12    logical, save :: pole_sud
13   
14    integer, allocatable, save, dimension(:) :: jj_begin_para
15    integer, allocatable, save, dimension(:) :: jj_end_para
16    integer, allocatable, save, dimension(:) :: jj_nb_para
[764]17    integer, save :: OMP_CHUNK
[630]18   
19 contains
20 
21    subroutine init_parallel
22    USE vampir
[764]23#ifdef CPP_COUPLE
24    USE mod_prism_proto
25#endif
[630]26    implicit none
27   
28      integer :: ierr
29      integer :: i,j
30      integer :: type_size
31      integer, dimension(3) :: blocklen,type
[764]32      integer :: comp_id
[630]33     
34     
35      include 'mpif.h'
[792]36#include "dimensions.h"
37#include "paramet.h"
[764]38
39#ifdef CPP_COUPLE
40       call prism_init_comp_proto (comp_id, 'lmdz.x', ierr)
41       call prism_get_localcomm_proto(COMM_LMDZ,ierr)
42#else     
[630]43      call MPI_INIT(ierr)
[764]44      COMM_LMDZ=MPI_COMM_WORLD
45#endif
[630]46      call InitVampir
[764]47      call MPI_COMM_SIZE(COMM_LMDZ,mpi_size,ierr)
48      call MPI_COMM_RANK(COMM_LMDZ,mpi_rank,ierr)
[630]49 
50     
51      allocate(jj_begin_para(0:mpi_size-1))
52      allocate(jj_end_para(0:mpi_size-1))
53      allocate(jj_nb_para(0:mpi_size-1))
54     
55      do i=0,mpi_size-1
56        jj_nb_para(i)=(jjm+1)/mpi_size
57        if ( i < MOD((jjm+1),mpi_size) ) jj_nb_para(i)=jj_nb_para(i)+1
58       
59        if (jj_nb_para(i) <= 2 ) then
60         
61         print *,"Arret : le nombre de bande de lattitude par process est trop faible (<2)."
62          print *," ---> diminuez le nombre de CPU ou augmentez la taille en lattitude"
63         
[764]64          call MPI_ABORT(COMM_LMDZ,-1, ierr)
[630]65         
66        endif
67       
68      enddo
69     
70!      jj_nb_para(0)=11
71!      jj_nb_para(1)=25
72!      jj_nb_para(2)=25
73!      jj_nb_para(3)=12     
74
75      j=1
76     
77      do i=0,mpi_size-1
78       
79        jj_begin_para(i)=j
80        jj_end_para(i)=j+jj_Nb_para(i)-1
81        j=j+jj_Nb_para(i)
82     
83      enddo
84     
85      jj_begin = jj_begin_para(mpi_rank)
86      jj_end   = jj_end_para(mpi_rank)
87      jj_nb    = jj_nb_para(mpi_rank)
88     
89      ij_begin=(jj_begin-1)*iip1+1
90      ij_end=jj_end*iip1
91     
92      if (mpi_rank.eq.0) then
93        pole_nord=.TRUE.
94      else
95        pole_nord=.FALSE.
96      endif
97     
98      if (mpi_rank.eq.mpi_size-1) then
99        pole_sud=.TRUE.
100      else
101        pole_sud=.FALSE.
102      endif
103       
104      print *,"jj_begin",jj_begin
105      print *,"jj_end",jj_end
106      print *,"ij_begin",ij_begin
107      print *,"ij_end",ij_end
108   
109   
110    end subroutine init_parallel
111
112   
113    subroutine SetDistrib(jj_Nb_New)
114    implicit none
115
[792]116#include "dimensions.h"
117#include "paramet.h"
[630]118
119      INTEGER,dimension(0:MPI_Size-1) :: jj_Nb_New
120      INTEGER :: i 
121 
122      jj_Nb_Para=jj_Nb_New
123     
124      jj_begin_para(0)=1
125      jj_end_para(0)=jj_Nb_Para(0)
126     
127      do i=1,mpi_size-1
128       
129        jj_begin_para(i)=jj_end_para(i-1)+1
130        jj_end_para(i)=jj_begin_para(i)+jj_Nb_para(i)-1
131     
132      enddo
133     
134      jj_begin = jj_begin_para(mpi_rank)
135      jj_end   = jj_end_para(mpi_rank)
136      jj_nb    = jj_nb_para(mpi_rank)
137     
138      ij_begin=(jj_begin-1)*iip1+1
139      ij_end=jj_end*iip1
140
141    end subroutine SetDistrib
142
143
144
145   
146    subroutine Finalize_parallel
[764]147#ifdef CPP_COUPLE
148    use mod_prism_proto
149#endif
[630]150    implicit none
151
[792]152#include "dimensions.h"
153#include "paramet.h"
[630]154      integer :: ierr
155      integer :: i
156      include 'mpif.h'
157     
158      deallocate(jj_begin_para)
159      deallocate(jj_end_para)
160      deallocate(jj_nb_para)
[764]161
162#ifdef CPP_COUPLE
163     call prism_terminate_proto(ierr)
164     IF (ierr .ne. PRISM_Ok) THEN
165       call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
166     endif
167#else         
[630]168      call MPI_FINALIZE(ierr)
[764]169#endif
[630]170     
171    end subroutine Finalize_parallel
[764]172       
[630]173    subroutine Pack_Data(Field,ij,ll,row,Buffer)
174    implicit none
175
[792]176#include "dimensions.h"
177#include "paramet.h"
[630]178
179      integer, intent(in) :: ij,ll,row
180      real,dimension(ij,ll),intent(in) ::Field
181      real,dimension(ll*iip1*row), intent(out) :: Buffer
182           
183      integer :: Pos
184      integer :: i,l
185     
186      Pos=0
187      do l=1,ll
188        do i=1,row*iip1
189          Pos=Pos+1
190          Buffer(Pos)=Field(i,l)
191        enddo
192      enddo
193     
194    end subroutine Pack_data
195     
196    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
197    implicit none
198
[792]199#include "dimensions.h"
200#include "paramet.h"
[630]201
202      integer, intent(in) :: ij,ll,row
203      real,dimension(ij,ll),intent(out) ::Field
204      real,dimension(ll*iip1*row), intent(in) :: Buffer
205           
206      integer :: Pos
207      integer :: i,l
208     
209      Pos=0
210     
211      do l=1,ll
212        do i=1,row*iip1
213          Pos=Pos+1
214          Field(i,l)=Buffer(Pos)
215        enddo
216      enddo
217     
218    end subroutine UnPack_data
219     
220    subroutine exchange_hallo(Field,ij,ll,up,down)
221    USE Vampir
222    implicit none
[792]223#include "dimensions.h"
224#include "paramet.h"   
[630]225    include 'mpif.h'
226   
227      INTEGER :: ij,ll
228      REAL, dimension(ij,ll) :: Field
229      INTEGER :: up,down
230     
231      INTEGER :: ierr
232      LOGICAL :: SendUp,SendDown
233      LOGICAL :: RecvUp,RecvDown
234      INTEGER, DIMENSION(4) :: Request
235      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
236      INTEGER :: NbRequest
237      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
238      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
239      INTEGER :: Buffer_size     
240     
[764]241      call MPI_Barrier(COMM_LMDZ,ierr)
[630]242      call VTb(VThallo)
243     
244      SendUp=.TRUE.
245      SendDown=.TRUE.
246      RecvUp=.TRUE.
247      RecvDown=.TRUE.
248       
249      IF (pole_nord) THEN
250        SendUp=.FALSE.
251        RecvUp=.FALSE.
252      ENDIF
253 
254      IF (pole_sud) THEN
255        SendDown=.FALSE.
256        RecvDown=.FALSE.
257      ENDIF
258     
259      if (up.eq.0) then
260        SendDown=.FALSE.
261        RecvUp=.FALSE.
262      endif
263     
264      if (down.eq.0) then
265        SendUp=.FALSE.
266        RecvDown=.FALSE.
267      endif
268     
269      NbRequest=0
270 
271      IF (SendUp) THEN
272        NbRequest=NbRequest+1
273        buffer_size=down*iip1*ll
274        allocate(Buffer_Send_up(Buffer_size))
275        call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
276        call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
[764]277                        COMM_LMDZ,Request(NbRequest),ierr)
[630]278      ENDIF
279 
280      IF (SendDown) THEN
281        NbRequest=NbRequest+1
282       
283        buffer_size=up*iip1*ll
284        allocate(Buffer_Send_down(Buffer_size))
285        call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
286       
287        call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
[764]288                        COMM_LMDZ,Request(NbRequest),ierr)
[630]289      ENDIF
290   
291 
292      IF (RecvUp) THEN
293        NbRequest=NbRequest+1
294        buffer_size=up*iip1*ll
295        allocate(Buffer_recv_up(Buffer_size))
296             
297        call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
[764]298                        COMM_LMDZ,Request(NbRequest),ierr)
[630]299     
300       
301      ENDIF
302 
303      IF (RecvDown) THEN
304        NbRequest=NbRequest+1
305        buffer_size=down*iip1*ll
306        allocate(Buffer_recv_down(Buffer_size))
307       
308        call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
[764]309                        COMM_LMDZ,Request(NbRequest),ierr)
[630]310     
311       
312      ENDIF
313 
314      if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
315      IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
316      IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
317
318      call VTe(VThallo)
[764]319      call MPI_Barrier(COMM_LMDZ,ierr)
[630]320      RETURN
321     
322    end subroutine exchange_Hallo
323   
324   
325    subroutine Gather_Field(Field,ij,ll,rank)
326    implicit none
[792]327#include "dimensions.h"
328#include "paramet.h"   
[630]329    include 'mpif.h'
330   
331      INTEGER :: ij,ll,rank
332      REAL, dimension(ij,ll) :: Field
333      REAL, dimension(:),allocatable :: Buffer_send   
334      REAL, dimension(:),allocatable :: Buffer_Recv
335      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
336      INTEGER :: ierr
337      INTEGER ::i
338     
339      if (ij==ip1jmp1) then
340         allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
341         call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
342      else if (ij==ip1jm) then
343         allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
344         call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
345      else
346         print *,ij
347         stop 'erreur dans Gather_Field'
348      endif
349     
350      if (MPI_Rank==rank) then
351        allocate(Buffer_Recv(ij*ll))
352        do i=0,MPI_Size-1
353         
354          if (ij==ip1jmp1) then
355            Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
356          else if (ij==ip1jm) then
357            Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
358          else
359            stop 'erreur dans Gather_Field'
360          endif
361         
362          if (i==0) then
363            displ(i)=0
364          else
365            displ(i)=displ(i-1)+Recv_count(i-1)
366          endif
367         
368        enddo
369       
370      endif
371     
372      call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
[764]373                        Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
[630]374     
375      if (MPI_Rank==rank) then                 
376     
377        if (ij==ip1jmp1) then
378          do i=0,MPI_Size-1
379            call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
380                             jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
381          enddo
382        else if (ij==ip1jm) then
383          do i=0,MPI_Size-1
384             call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
385                             min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
386          enddo
387        endif
388     
389      endif
390     
391    end subroutine Gather_Field
392   
393    subroutine AllGather_Field(Field,ij,ll)
394    implicit none
[792]395#include "dimensions.h"
396#include "paramet.h"   
[630]397    include 'mpif.h'
398   
399      INTEGER :: ij,ll
400      REAL, dimension(ij,ll) :: Field
401      INTEGER :: ierr
402     
403      call Gather_Field(Field,ij,ll,0)
[764]404      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
[630]405     
406    end subroutine AllGather_Field
407   
408   subroutine Broadcast_Field(Field,ij,ll,rank)
409    implicit none
[792]410#include "dimensions.h"
411#include "paramet.h"   
[630]412    include 'mpif.h'
413   
414      INTEGER :: ij,ll
415      REAL, dimension(ij,ll) :: Field
416      INTEGER :: rank
417      INTEGER :: ierr
418     
[764]419      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
[630]420     
421    end subroutine Broadcast_Field
422       
423   
424    /* 
425  Subroutine verif_hallo(Field,ij,ll,up,down)
426    implicit none
[792]427#include "dimensions.h"
428#include "paramet.h"   
[630]429    include 'mpif.h'
430   
431      INTEGER :: ij,ll
432      REAL, dimension(ij,ll) :: Field
433      INTEGER :: up,down
434     
435      REAL,dimension(ij,ll): NewField
436     
437      NewField=0
438     
439      ijb=ij_begin
440      ije=ij_end
441      if (pole_nord)
442      NewField(ij_be       
443*/
444  end module parallel
Note: See TracBrowser for help on using the repository browser.