source: LMDZ5/trunk/libf/dyn3dpar/parallel.F90 @ 1638

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