source: LMDZ4/branches/LMDZ4_AR5/libf/dyn3dpar/parallel.F90 @ 1476

Last change on this file since 1476 was 1476, checked in by jghattas, 13 years ago

Bug correction for ce0l in dyn3dpar :

  • added missing MPI_FINALIZE
  • removed STOP before end of program


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