source: LMDZ4/branches/V3_test/libf/dyn3dpar/parallel.F90 @ 709

Last change on this file since 709 was 709, checked in by Laurent Fairhead, 18 years ago

Nouvelles versions de la dynamique YM
LF

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