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

Last change on this file since 825 was 811, checked in by lsce, 17 years ago

ACo + YM : correction bug

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