source: LMDZ4/tags/LMDZ4_V3_3/libf/dyn3dpar/parallel.F90

Last change on this file was 919, checked in by (none), 16 years ago

This commit was manufactured by cvs2svn to create tag 'LMDZ4_V3_3'.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.9 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      use surface_data, only : ocean
140      implicit none
141
142      include "dimensions.h"
143      include "paramet.h"
144      integer :: ierr
145      integer :: i
146      include 'mpif.h'
147     
148      deallocate(jj_begin_para)
149      deallocate(jj_end_para)
150      deallocate(jj_nb_para)
151
152      if (ocean == 'couple') then
153#ifdef CPP_COUPLE
154         call prism_terminate_proto(ierr)
155         IF (ierr .ne. PRISM_Ok) THEN
156            call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
157         endif
158#endif
159      else
160         call MPI_FINALIZE(ierr)
161      end if
162     
163    end subroutine Finalize_parallel
164       
165    subroutine Pack_Data(Field,ij,ll,row,Buffer)
166    implicit none
167
168#include "dimensions.h"
169#include "paramet.h"
170
171      integer, intent(in) :: ij,ll,row
172      real,dimension(ij,ll),intent(in) ::Field
173      real,dimension(ll*iip1*row), intent(out) :: Buffer
174           
175      integer :: Pos
176      integer :: i,l
177     
178      Pos=0
179      do l=1,ll
180        do i=1,row*iip1
181          Pos=Pos+1
182          Buffer(Pos)=Field(i,l)
183        enddo
184      enddo
185     
186    end subroutine Pack_data
187     
188    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
189    implicit none
190
191#include "dimensions.h"
192#include "paramet.h"
193
194      integer, intent(in) :: ij,ll,row
195      real,dimension(ij,ll),intent(out) ::Field
196      real,dimension(ll*iip1*row), intent(in) :: Buffer
197           
198      integer :: Pos
199      integer :: i,l
200     
201      Pos=0
202     
203      do l=1,ll
204        do i=1,row*iip1
205          Pos=Pos+1
206          Field(i,l)=Buffer(Pos)
207        enddo
208      enddo
209     
210    end subroutine UnPack_data
211     
212    subroutine exchange_hallo(Field,ij,ll,up,down)
213    USE Vampir
214    implicit none
215#include "dimensions.h"
216#include "paramet.h"   
217    include 'mpif.h'
218   
219      INTEGER :: ij,ll
220      REAL, dimension(ij,ll) :: Field
221      INTEGER :: up,down
222     
223      INTEGER :: ierr
224      LOGICAL :: SendUp,SendDown
225      LOGICAL :: RecvUp,RecvDown
226      INTEGER, DIMENSION(4) :: Request
227      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
228      INTEGER :: NbRequest
229      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
230      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
231      INTEGER :: Buffer_size     
232     
233      call MPI_Barrier(COMM_LMDZ,ierr)
234      call VTb(VThallo)
235     
236      SendUp=.TRUE.
237      SendDown=.TRUE.
238      RecvUp=.TRUE.
239      RecvDown=.TRUE.
240       
241      IF (pole_nord) THEN
242        SendUp=.FALSE.
243        RecvUp=.FALSE.
244      ENDIF
245 
246      IF (pole_sud) THEN
247        SendDown=.FALSE.
248        RecvDown=.FALSE.
249      ENDIF
250     
251      if (up.eq.0) then
252        SendDown=.FALSE.
253        RecvUp=.FALSE.
254      endif
255     
256      if (down.eq.0) then
257        SendUp=.FALSE.
258        RecvDown=.FALSE.
259      endif
260     
261      NbRequest=0
262 
263      IF (SendUp) THEN
264        NbRequest=NbRequest+1
265        buffer_size=down*iip1*ll
266        allocate(Buffer_Send_up(Buffer_size))
267        call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
268        call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
269                        COMM_LMDZ,Request(NbRequest),ierr)
270      ENDIF
271 
272      IF (SendDown) THEN
273        NbRequest=NbRequest+1
274       
275        buffer_size=up*iip1*ll
276        allocate(Buffer_Send_down(Buffer_size))
277        call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
278       
279        call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
280                        COMM_LMDZ,Request(NbRequest),ierr)
281      ENDIF
282   
283 
284      IF (RecvUp) THEN
285        NbRequest=NbRequest+1
286        buffer_size=up*iip1*ll
287        allocate(Buffer_recv_up(Buffer_size))
288             
289        call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
290                        COMM_LMDZ,Request(NbRequest),ierr)
291     
292       
293      ENDIF
294 
295      IF (RecvDown) THEN
296        NbRequest=NbRequest+1
297        buffer_size=down*iip1*ll
298        allocate(Buffer_recv_down(Buffer_size))
299       
300        call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
301                        COMM_LMDZ,Request(NbRequest),ierr)
302     
303       
304      ENDIF
305 
306      if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
307      IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
308      IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
309
310      call VTe(VThallo)
311      call MPI_Barrier(COMM_LMDZ,ierr)
312      RETURN
313     
314    end subroutine exchange_Hallo
315   
316   
317    subroutine Gather_Field(Field,ij,ll,rank)
318    implicit none
319#include "dimensions.h"
320#include "paramet.h"   
321    include 'mpif.h'
322   
323      INTEGER :: ij,ll,rank
324      REAL, dimension(ij,ll) :: Field
325      REAL, dimension(:),allocatable :: Buffer_send   
326      REAL, dimension(:),allocatable :: Buffer_Recv
327      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
328      INTEGER :: ierr
329      INTEGER ::i
330     
331      if (ij==ip1jmp1) then
332         allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
333         call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
334      else if (ij==ip1jm) then
335         allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
336         call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
337      else
338         print *,ij
339         stop 'erreur dans Gather_Field'
340      endif
341     
342      if (MPI_Rank==rank) then
343        allocate(Buffer_Recv(ij*ll))
344        do i=0,MPI_Size-1
345         
346          if (ij==ip1jmp1) then
347            Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
348          else if (ij==ip1jm) then
349            Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
350          else
351            stop 'erreur dans Gather_Field'
352          endif
353         
354          if (i==0) then
355            displ(i)=0
356          else
357            displ(i)=displ(i-1)+Recv_count(i-1)
358          endif
359         
360        enddo
361       
362      endif
363     
364      call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
365                        Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
366     
367      if (MPI_Rank==rank) then                 
368     
369        if (ij==ip1jmp1) then
370          do i=0,MPI_Size-1
371            call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
372                             jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
373          enddo
374        else if (ij==ip1jm) then
375          do i=0,MPI_Size-1
376             call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
377                             min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
378          enddo
379        endif
380     
381      endif
382     
383    end subroutine Gather_Field
384   
385    subroutine AllGather_Field(Field,ij,ll)
386    implicit none
387#include "dimensions.h"
388#include "paramet.h"   
389    include 'mpif.h'
390   
391      INTEGER :: ij,ll
392      REAL, dimension(ij,ll) :: Field
393      INTEGER :: ierr
394     
395      call Gather_Field(Field,ij,ll,0)
396      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
397     
398    end subroutine AllGather_Field
399   
400   subroutine Broadcast_Field(Field,ij,ll,rank)
401    implicit none
402#include "dimensions.h"
403#include "paramet.h"   
404    include 'mpif.h'
405   
406      INTEGER :: ij,ll
407      REAL, dimension(ij,ll) :: Field
408      INTEGER :: rank
409      INTEGER :: ierr
410     
411      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
412     
413    end subroutine Broadcast_Field
414       
415   
416    /* 
417  Subroutine verif_hallo(Field,ij,ll,up,down)
418    implicit none
419#include "dimensions.h"
420#include "paramet.h"   
421    include 'mpif.h'
422   
423      INTEGER :: ij,ll
424      REAL, dimension(ij,ll) :: Field
425      INTEGER :: up,down
426     
427      REAL,dimension(ij,ll): NewField
428     
429      NewField=0
430     
431      ijb=ij_begin
432      ije=ij_end
433      if (pole_nord)
434      NewField(ij_be       
435*/
436  end module parallel
Note: See TracBrowser for help on using the repository browser.