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

Last change on this file since 807 was 806, checked in by lsce, 17 years ago

ACo + YM : correction bug communicateur pour le couple

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