source: LMDZ4/branches/LMDZ4-dev/libf/dyn3dpar/parallel.F90 @ 1200

Last change on this file since 1200 was 1089, checked in by yann meurdesoif, 16 years ago

correction Othman Bouizi : Initialisation de is_using_omp pour la dynamique

YM

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