source: LMDZ5/trunk/libf/dyn3dpar/parallel_lmdz.F90 @ 1823

Last change on this file since 1823 was 1823, checked in by Ehouarn Millour, 11 years ago

Remplacement de parallel.F90 (en conflit avec orchidée) par parallel_lmdz.F90.
UG
.........................................
Renaming parallel.F90 (conflicting with orchidée) into parallel_lmdz.F90.
UG

File size: 14.4 KB
Line 
1!
2! $Id: parallel.F90 1810 2013-07-24 08:06:39Z emillour $
3!
4  MODULE parallel_lmdz
5  USE mod_const_mpi
6   
7    LOGICAL,SAVE :: using_mpi=.TRUE.
8    LOGICAL,SAVE :: using_omp
9   
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
23    integer, save :: OMP_CHUNK
24    integer, save :: omp_rank
25    integer, save :: omp_size 
26!$OMP THREADPRIVATE(omp_rank)
27
28 contains
29 
30    subroutine init_parallel
31    USE vampir
32    implicit none
33#ifdef CPP_MPI
34      include 'mpif.h'
35#endif
36#include "dimensions.h"
37#include "paramet.h"
38#include "iniprint.h"
39
40      integer :: ierr
41      integer :: i,j
42      integer :: type_size
43      integer, dimension(3) :: blocklen,type
44      integer :: comp_id
45      character(len=4)  :: num
46      character(len=20) :: filename
47 
48#ifdef CPP_OMP   
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 
54
55#ifdef CPP_MPI
56       using_mpi=.TRUE.
57#else
58       using_mpi=.FALSE.
59#endif
60     
61
62#ifdef CPP_OMP
63       using_OMP=.TRUE.
64#else
65       using_OMP=.FALSE.
66#endif
67     
68      call InitVampir
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
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
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         
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"
104         
105#ifdef CPP_MPI
106          IF (using_mpi) call MPI_ABORT(COMM_LMDZ,-1, ierr)
107#endif         
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       
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
150
151!$OMP PARALLEL
152
153#ifdef CPP_OMP
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         
163   
164    end subroutine init_parallel
165
166   
167    subroutine SetDistrib(jj_Nb_New)
168    implicit none
169
170#include "dimensions.h"
171#include "paramet.h"
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
201#ifdef CPP_COUPLE
202    use mod_prism_proto
203! Ehouarn: surface_data module is in 'phylmd' ...
204      use surface_data, only : type_ocean
205      implicit none
206#else
207      implicit none
208! without the surface_data module, we declare (and set) a dummy 'type_ocean'
209      character(len=6),parameter :: type_ocean="dummy"
210#endif
211! #endif of #ifdef CPP_EARTH
212
213      include "dimensions.h"
214      include "paramet.h"
215#ifdef CPP_MPI
216      include 'mpif.h'
217#endif     
218
219      integer :: ierr
220      integer :: i
221
222      if (allocated(jj_begin_para)) deallocate(jj_begin_para)
223      if (allocated(jj_end_para))   deallocate(jj_end_para)
224      if (allocated(jj_nb_para))    deallocate(jj_nb_para)
225
226      if (type_ocean == 'couple') then
227#ifdef CPP_COUPLE
228         call prism_terminate_proto(ierr)
229         IF (ierr .ne. PRISM_Ok) THEN
230            call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
231         endif
232#endif
233      else
234#ifdef CPP_MPI
235         IF (using_mpi) call MPI_FINALIZE(ierr)
236#endif
237      end if
238     
239    end subroutine Finalize_parallel
240       
241    subroutine Pack_Data(Field,ij,ll,row,Buffer)
242    implicit none
243
244#include "dimensions.h"
245#include "paramet.h"
246
247      integer, intent(in) :: ij,ll,row
248      real,dimension(ij,ll),intent(in) ::Field
249      real,dimension(ll*iip1*row), intent(out) :: Buffer
250           
251      integer :: Pos
252      integer :: i,l
253     
254      Pos=0
255      do l=1,ll
256        do i=1,row*iip1
257          Pos=Pos+1
258          Buffer(Pos)=Field(i,l)
259        enddo
260      enddo
261     
262    end subroutine Pack_data
263     
264    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
265    implicit none
266
267#include "dimensions.h"
268#include "paramet.h"
269
270      integer, intent(in) :: ij,ll,row
271      real,dimension(ij,ll),intent(out) ::Field
272      real,dimension(ll*iip1*row), intent(in) :: Buffer
273           
274      integer :: Pos
275      integer :: i,l
276     
277      Pos=0
278     
279      do l=1,ll
280        do i=1,row*iip1
281          Pos=Pos+1
282          Field(i,l)=Buffer(Pos)
283        enddo
284      enddo
285     
286    end subroutine UnPack_data
287
288   
289    SUBROUTINE barrier
290    IMPLICIT NONE
291#ifdef CPP_MPI
292    INCLUDE 'mpif.h'
293#endif
294    INTEGER :: ierr
295   
296!$OMP CRITICAL (MPI)     
297#ifdef CPP_MPI
298      IF (using_mpi) CALL MPI_Barrier(COMM_LMDZ,ierr)
299#endif
300!$OMP END CRITICAL (MPI)
301   
302    END SUBROUTINE barrier
303       
304     
305    subroutine exchange_hallo(Field,ij,ll,up,down)
306    USE Vampir
307    implicit none
308#include "dimensions.h"
309#include "paramet.h"   
310#ifdef CPP_MPI
311    include 'mpif.h'
312#endif   
313      INTEGER :: ij,ll
314      REAL, dimension(ij,ll) :: Field
315      INTEGER :: up,down
316     
317      INTEGER :: ierr
318      LOGICAL :: SendUp,SendDown
319      LOGICAL :: RecvUp,RecvDown
320      INTEGER, DIMENSION(4) :: Request
321#ifdef CPP_MPI
322      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
323#else
324      INTEGER, DIMENSION(1,4) :: Status
325#endif
326      INTEGER :: NbRequest
327      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
328      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
329      INTEGER :: Buffer_size     
330
331      IF (using_mpi) THEN
332
333        CALL barrier
334     
335        call VTb(VThallo)
336     
337        SendUp=.TRUE.
338        SendDown=.TRUE.
339        RecvUp=.TRUE.
340        RecvDown=.TRUE.
341         
342        IF (pole_nord) THEN
343          SendUp=.FALSE.
344          RecvUp=.FALSE.
345        ENDIF
346   
347        IF (pole_sud) THEN
348          SendDown=.FALSE.
349          RecvDown=.FALSE.
350        ENDIF
351       
352        if (up.eq.0) then
353          SendDown=.FALSE.
354          RecvUp=.FALSE.
355        endif
356     
357        if (down.eq.0) then
358          SendUp=.FALSE.
359          RecvDown=.FALSE.
360        endif
361     
362        NbRequest=0
363 
364        IF (SendUp) THEN
365          NbRequest=NbRequest+1
366          buffer_size=down*iip1*ll
367          allocate(Buffer_Send_up(Buffer_size))
368          call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
369!$OMP CRITICAL (MPI)
370#ifdef CPP_MPI
371          call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
372                          COMM_LMDZ,Request(NbRequest),ierr)
373#endif
374!$OMP END CRITICAL (MPI)
375        ENDIF
376 
377        IF (SendDown) THEN
378          NbRequest=NbRequest+1
379           
380          buffer_size=up*iip1*ll
381          allocate(Buffer_Send_down(Buffer_size))
382          call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
383       
384!$OMP CRITICAL (MPI)
385#ifdef CPP_MPI
386          call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
387                          COMM_LMDZ,Request(NbRequest),ierr)
388#endif
389!$OMP END CRITICAL (MPI)
390        ENDIF
391   
392 
393        IF (RecvUp) THEN
394          NbRequest=NbRequest+1
395          buffer_size=up*iip1*ll
396          allocate(Buffer_recv_up(Buffer_size))
397             
398!$OMP CRITICAL (MPI)
399#ifdef CPP_MPI
400          call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
401                          COMM_LMDZ,Request(NbRequest),ierr)
402#endif
403!$OMP END CRITICAL (MPI)
404     
405       
406        ENDIF
407 
408        IF (RecvDown) THEN
409          NbRequest=NbRequest+1
410          buffer_size=down*iip1*ll
411          allocate(Buffer_recv_down(Buffer_size))
412       
413!$OMP CRITICAL (MPI)
414#ifdef CPP_MPI
415          call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
416                          COMM_LMDZ,Request(NbRequest),ierr)
417#endif
418!$OMP END CRITICAL (MPI)
419       
420        ENDIF
421 
422#ifdef CPP_MPI
423        if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
424#endif
425        IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
426        IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
427
428        call VTe(VThallo)
429        call barrier
430     
431      ENDIF  ! using_mpi
432     
433      RETURN
434     
435    end subroutine exchange_Hallo
436   
437
438    subroutine Gather_Field(Field,ij,ll,rank)
439    implicit none
440#include "dimensions.h"
441#include "paramet.h"
442#include "iniprint.h"
443#ifdef CPP_MPI
444    include 'mpif.h'
445#endif   
446      INTEGER :: ij,ll,rank
447      REAL, dimension(ij,ll) :: Field
448      REAL, dimension(:),allocatable :: Buffer_send   
449      REAL, dimension(:),allocatable :: Buffer_Recv
450      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
451      INTEGER :: ierr
452      INTEGER ::i
453     
454      IF (using_mpi) THEN
455
456        if (ij==ip1jmp1) then
457           allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
458           call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
459        else if (ij==ip1jm) then
460           allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
461           call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
462        else
463           write(lunout,*)ij 
464        stop 'erreur dans Gather_Field'
465        endif
466       
467        if (MPI_Rank==rank) then
468          allocate(Buffer_Recv(ij*ll))
469
470!CDIR NOVECTOR
471          do i=0,MPI_Size-1
472             
473            if (ij==ip1jmp1) then
474              Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
475            else if (ij==ip1jm) then
476              Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
477            else
478              stop 'erreur dans Gather_Field'
479            endif
480                   
481            if (i==0) then
482              displ(i)=0
483            else
484              displ(i)=displ(i-1)+Recv_count(i-1)
485            endif
486           
487          enddo
488         
489        else
490          ! Ehouarn: When in debug mode, ifort complains (for call MPI_GATHERV
491          !          below) about Buffer_Recv() being not allocated.
492          !          So make a dummy allocation.
493          allocate(Buffer_Recv(1))
494        endif ! of if (MPI_Rank==rank)
495 
496!$OMP CRITICAL (MPI)
497#ifdef CPP_MPI
498        call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
499                          Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
500#endif
501!$OMP END CRITICAL (MPI)
502     
503        if (MPI_Rank==rank) then                 
504     
505          if (ij==ip1jmp1) then
506            do i=0,MPI_Size-1
507              call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
508                               jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
509            enddo
510          else if (ij==ip1jm) then
511            do i=0,MPI_Size-1
512               call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
513                               min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
514            enddo
515          endif
516        endif
517      ENDIF ! using_mpi
518     
519    end subroutine Gather_Field
520
521
522    subroutine AllGather_Field(Field,ij,ll)
523    implicit none
524#include "dimensions.h"
525#include "paramet.h"   
526#ifdef CPP_MPI
527    include 'mpif.h'
528#endif   
529      INTEGER :: ij,ll
530      REAL, dimension(ij,ll) :: Field
531      INTEGER :: ierr
532     
533      IF (using_mpi) THEN
534        call Gather_Field(Field,ij,ll,0)
535!$OMP CRITICAL (MPI)
536#ifdef CPP_MPI
537      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
538#endif
539!$OMP END CRITICAL (MPI)
540      ENDIF
541     
542    end subroutine AllGather_Field
543   
544   subroutine Broadcast_Field(Field,ij,ll,rank)
545    implicit none
546#include "dimensions.h"
547#include "paramet.h"   
548#ifdef CPP_MPI
549    include 'mpif.h'
550#endif   
551      INTEGER :: ij,ll
552      REAL, dimension(ij,ll) :: Field
553      INTEGER :: rank
554      INTEGER :: ierr
555     
556      IF (using_mpi) THEN
557     
558!$OMP CRITICAL (MPI)
559#ifdef CPP_MPI
560      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
561#endif
562!$OMP END CRITICAL (MPI)
563     
564      ENDIF
565    end subroutine Broadcast_Field
566       
567   
568!  Subroutine verif_hallo(Field,ij,ll,up,down)
569!    implicit none
570!#include "dimensions.h"
571!#include "paramet.h"   
572!    include 'mpif.h'
573!   
574!      INTEGER :: ij,ll
575!      REAL, dimension(ij,ll) :: Field
576!      INTEGER :: up,down
577!     
578!      REAL,dimension(ij,ll): NewField
579!     
580!      NewField=0
581!     
582!      ijb=ij_begin
583!      ije=ij_end
584!      if (pole_nord)
585!      NewField(ij_be       
586
587  end MODULE parallel_lmdz
Note: See TracBrowser for help on using the repository browser.