source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/dyn3dpar/parallel_lmdz.F90 @ 5409

Last change on this file since 5409 was 2239, checked in by Ehouarn Millour, 10 years ago

Reorganizing physics/dynamics interface:

  • what is related to dynamics-physics interface is now in a seperate directory: dynlmdz_phy* for physics in phy*
  • 1d model and related dependencies (including a couple from "dynamics", set up as symbolic links) is now in subdirectory "dyn1d" of phy*.
  • "bibio" directory is now "misc" and should only contain autonomous utilities.
  • "cosp" is now a subdirectory of phylmd.

EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 15.8 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#ifdef CPP_IOIPSL
7      use IOIPSL, only: getin
8#else
9! if not using IOIPSL, we still need to use (a local version of) getin
10      use ioipsl_getincom, only: getin
11#endif   
12   
13    LOGICAL,SAVE :: using_mpi=.TRUE.
14    LOGICAL,SAVE :: using_omp
15   
16    integer, save :: mpi_size
17    integer, save :: mpi_rank
18    integer, save :: jj_begin
19    integer, save :: jj_end
20    integer, save :: jj_nb
21    integer, save :: ij_begin
22    integer, save :: ij_end
23    logical, save :: pole_nord
24    logical, save :: pole_sud
25   
26    integer, allocatable, save, dimension(:) :: jj_begin_para
27    integer, allocatable, save, dimension(:) :: jj_end_para
28    integer, allocatable, save, dimension(:) :: jj_nb_para
29    integer, save :: OMP_CHUNK
30    integer, save :: omp_rank
31    integer, save :: omp_size 
32!$OMP THREADPRIVATE(omp_rank)
33
34! Ehouarn: add "dummy variables" (which are in dyn3d_mem/parallel_lmdz.F90)
35! so that calfis_loc compiles even if using dyn3dpar
36    integer,save  :: jjb_u
37    integer,save  :: jje_u
38    integer,save  :: jjnb_u
39    integer,save  :: jjb_v
40    integer,save  :: jje_v
41    integer,save  :: jjnb_v   
42
43    integer,save  :: ijb_u
44    integer,save  :: ije_u
45    integer,save  :: ijnb_u   
46   
47    integer,save  :: ijb_v
48    integer,save  :: ije_v
49    integer,save  :: ijnb_v   
50
51 contains
52 
53    subroutine init_parallel
54    USE vampir
55    implicit none
56#ifdef CPP_MPI
57      include 'mpif.h'
58#endif
59#include "dimensions.h"
60#include "paramet.h"
61#include "iniprint.h"
62
63      integer :: ierr
64      integer :: i,j
65      integer :: type_size
66      integer, dimension(3) :: blocklen,type
67      integer :: comp_id
68      character(len=4)  :: num
69      character(len=20) :: filename
70 
71#ifdef CPP_OMP   
72      INTEGER :: OMP_GET_NUM_THREADS
73      EXTERNAL OMP_GET_NUM_THREADS
74      INTEGER :: OMP_GET_THREAD_NUM
75      EXTERNAL OMP_GET_THREAD_NUM
76#endif 
77
78#ifdef CPP_MPI
79       using_mpi=.TRUE.
80#else
81       using_mpi=.FALSE.
82#endif
83     
84
85#ifdef CPP_OMP
86       using_OMP=.TRUE.
87#else
88       using_OMP=.FALSE.
89#endif
90     
91      call InitVampir
92     
93      IF (using_mpi) THEN
94#ifdef CPP_MPI
95        call MPI_COMM_SIZE(COMM_LMDZ,mpi_size,ierr)
96        call MPI_COMM_RANK(COMM_LMDZ,mpi_rank,ierr)
97#endif
98      ELSE
99        mpi_size=1
100        mpi_rank=0
101      ENDIF
102
103
104! Open text output file with mpi_rank in suffix of file name
105      IF (lunout /= 5 .and. lunout /= 6) THEN
106         WRITE(num,'(I4.4)') mpi_rank
107         filename='lmdz.out_'//num
108         IF (mpi_rank .NE. 0) THEN
109            OPEN(UNIT=lunout,FILE=TRIM(filename),ACTION='write', &
110               STATUS='unknown',FORM='formatted',IOSTAT=ierr)
111         ENDIF
112      ENDIF
113
114     
115      allocate(jj_begin_para(0:mpi_size-1))
116      allocate(jj_end_para(0:mpi_size-1))
117      allocate(jj_nb_para(0:mpi_size-1))
118     
119      do i=0,mpi_size-1
120        jj_nb_para(i)=(jjm+1)/mpi_size
121        if ( i < MOD((jjm+1),mpi_size) ) jj_nb_para(i)=jj_nb_para(i)+1
122       
123        if (jj_nb_para(i) <= 2 ) then
124         
125         write(lunout,*)"Arret : le nombre de bande de lattitude par process est trop faible (<2)."
126         write(lunout,*)" ---> diminuez le nombre de CPU ou augmentez la taille en lattitude"
127         
128#ifdef CPP_MPI
129          IF (using_mpi) call MPI_ABORT(COMM_LMDZ,-1, ierr)
130#endif         
131        endif
132       
133      enddo
134     
135!      jj_nb_para(0)=11
136!      jj_nb_para(1)=25
137!      jj_nb_para(2)=25
138!      jj_nb_para(3)=12     
139
140      j=1
141     
142      do i=0,mpi_size-1
143       
144        jj_begin_para(i)=j
145        jj_end_para(i)=j+jj_Nb_para(i)-1
146        j=j+jj_Nb_para(i)
147     
148      enddo
149     
150      jj_begin = jj_begin_para(mpi_rank)
151      jj_end   = jj_end_para(mpi_rank)
152      jj_nb    = jj_nb_para(mpi_rank)
153     
154      ij_begin=(jj_begin-1)*iip1+1
155      ij_end=jj_end*iip1
156     
157      if (mpi_rank.eq.0) then
158        pole_nord=.TRUE.
159      else
160        pole_nord=.FALSE.
161      endif
162     
163      if (mpi_rank.eq.mpi_size-1) then
164        pole_sud=.TRUE.
165      else
166        pole_sud=.FALSE.
167      endif
168       
169      write(lunout,*)"init_parallel: jj_begin",jj_begin
170      write(lunout,*)"init_parallel: jj_end",jj_end
171      write(lunout,*)"init_parallel: ij_begin",ij_begin
172      write(lunout,*)"init_parallel: ij_end",ij_end
173
174!$OMP PARALLEL
175
176#ifdef CPP_OMP
177!$OMP MASTER
178        omp_size=OMP_GET_NUM_THREADS()
179!$OMP END MASTER
180!$OMP BARRIER
181        omp_rank=OMP_GET_THREAD_NUM()   
182
183!Config  Key  = omp_chunk
184!Config  Desc = taille des blocs openmp
185!Config  Def  = 1
186!Config  Help = defini la taille des packets d'it�ration openmp
187!Config         distribue a chaque tache lors de l'entree dans une
188!Config         boucle parallelisee
189
190!$OMP MASTER
191      omp_chunk=(llm+1)/omp_size
192      IF (MOD(llm+1,omp_size)/=0) omp_chunk=omp_chunk+1
193      CALL getin('omp_chunk',omp_chunk)
194!$OMP END MASTER
195!$OMP BARRIER       
196#else   
197        omp_size=1
198        omp_rank=0
199#endif
200!$OMP END PARALLEL         
201   
202    end subroutine init_parallel
203
204   
205    subroutine SetDistrib(jj_Nb_New)
206    implicit none
207
208#include "dimensions.h"
209#include "paramet.h"
210
211      INTEGER,dimension(0:MPI_Size-1) :: jj_Nb_New
212      INTEGER :: i 
213 
214      jj_Nb_Para=jj_Nb_New
215     
216      jj_begin_para(0)=1
217      jj_end_para(0)=jj_Nb_Para(0)
218     
219      do i=1,mpi_size-1
220       
221        jj_begin_para(i)=jj_end_para(i-1)+1
222        jj_end_para(i)=jj_begin_para(i)+jj_Nb_para(i)-1
223     
224      enddo
225     
226      jj_begin = jj_begin_para(mpi_rank)
227      jj_end   = jj_end_para(mpi_rank)
228      jj_nb    = jj_nb_para(mpi_rank)
229     
230      ij_begin=(jj_begin-1)*iip1+1
231      ij_end=jj_end*iip1
232
233    end subroutine SetDistrib
234
235
236
237   
238    subroutine Finalize_parallel
239#ifdef CPP_XIOS
240    ! ug Pour les sorties XIOS
241        USE wxios
242#endif
243#ifdef CPP_COUPLE
244! Use of Oasis-MCT coupler
245#if defined CPP_OMCT
246    use mod_prism
247#else
248    use mod_prism_proto
249#endif
250! Ehouarn: surface_data module is in 'phylmd' ...
251      use surface_data, only : type_ocean
252      implicit none
253#else
254      implicit none
255! without the surface_data module, we declare (and set) a dummy 'type_ocean'
256      character(len=6),parameter :: type_ocean="dummy"
257#endif
258! #endif of #ifdef CPP_EARTH
259
260      include "dimensions.h"
261      include "paramet.h"
262#ifdef CPP_MPI
263      include 'mpif.h'
264#endif     
265
266      integer :: ierr
267      integer :: i
268
269      if (allocated(jj_begin_para)) deallocate(jj_begin_para)
270      if (allocated(jj_end_para))   deallocate(jj_end_para)
271      if (allocated(jj_nb_para))    deallocate(jj_nb_para)
272
273      if (type_ocean == 'couple') then
274#ifdef CPP_XIOS
275    !Fermeture propre de XIOS
276      CALL wxios_close()
277#else
278#ifdef CPP_COUPLE
279         call prism_terminate_proto(ierr)
280         IF (ierr .ne. PRISM_Ok) THEN
281            call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
282         endif
283#endif
284#endif
285      else
286#ifdef CPP_XIOS
287    !Fermeture propre de XIOS
288      CALL wxios_close()
289#endif
290#ifdef CPP_MPI
291         IF (using_mpi) call MPI_FINALIZE(ierr)
292#endif
293      end if
294     
295    end subroutine Finalize_parallel
296       
297    subroutine Pack_Data(Field,ij,ll,row,Buffer)
298    implicit none
299
300#include "dimensions.h"
301#include "paramet.h"
302
303      integer, intent(in) :: ij,ll,row
304      real,dimension(ij,ll),intent(in) ::Field
305      real,dimension(ll*iip1*row), intent(out) :: Buffer
306           
307      integer :: Pos
308      integer :: i,l
309     
310      Pos=0
311      do l=1,ll
312        do i=1,row*iip1
313          Pos=Pos+1
314          Buffer(Pos)=Field(i,l)
315        enddo
316      enddo
317     
318    end subroutine Pack_data
319     
320    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
321    implicit none
322
323#include "dimensions.h"
324#include "paramet.h"
325
326      integer, intent(in) :: ij,ll,row
327      real,dimension(ij,ll),intent(out) ::Field
328      real,dimension(ll*iip1*row), intent(in) :: Buffer
329           
330      integer :: Pos
331      integer :: i,l
332     
333      Pos=0
334     
335      do l=1,ll
336        do i=1,row*iip1
337          Pos=Pos+1
338          Field(i,l)=Buffer(Pos)
339        enddo
340      enddo
341     
342    end subroutine UnPack_data
343
344   
345    SUBROUTINE barrier
346    IMPLICIT NONE
347#ifdef CPP_MPI
348    INCLUDE 'mpif.h'
349#endif
350    INTEGER :: ierr
351   
352!$OMP CRITICAL (MPI)     
353#ifdef CPP_MPI
354      IF (using_mpi) CALL MPI_Barrier(COMM_LMDZ,ierr)
355#endif
356!$OMP END CRITICAL (MPI)
357   
358    END SUBROUTINE barrier
359       
360     
361    subroutine exchange_hallo(Field,ij,ll,up,down)
362    USE Vampir
363    implicit none
364#include "dimensions.h"
365#include "paramet.h"   
366#ifdef CPP_MPI
367    include 'mpif.h'
368#endif   
369      INTEGER :: ij,ll
370      REAL, dimension(ij,ll) :: Field
371      INTEGER :: up,down
372     
373      INTEGER :: ierr
374      LOGICAL :: SendUp,SendDown
375      LOGICAL :: RecvUp,RecvDown
376      INTEGER, DIMENSION(4) :: Request
377#ifdef CPP_MPI
378      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
379#else
380      INTEGER, DIMENSION(1,4) :: Status
381#endif
382      INTEGER :: NbRequest
383      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
384      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
385      INTEGER :: Buffer_size     
386
387      IF (using_mpi) THEN
388
389        CALL barrier
390     
391        call VTb(VThallo)
392     
393        SendUp=.TRUE.
394        SendDown=.TRUE.
395        RecvUp=.TRUE.
396        RecvDown=.TRUE.
397         
398        IF (pole_nord) THEN
399          SendUp=.FALSE.
400          RecvUp=.FALSE.
401        ENDIF
402   
403        IF (pole_sud) THEN
404          SendDown=.FALSE.
405          RecvDown=.FALSE.
406        ENDIF
407       
408        if (up.eq.0) then
409          SendDown=.FALSE.
410          RecvUp=.FALSE.
411        endif
412     
413        if (down.eq.0) then
414          SendUp=.FALSE.
415          RecvDown=.FALSE.
416        endif
417     
418        NbRequest=0
419 
420        IF (SendUp) THEN
421          NbRequest=NbRequest+1
422          buffer_size=down*iip1*ll
423          allocate(Buffer_Send_up(Buffer_size))
424          call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
425!$OMP CRITICAL (MPI)
426#ifdef CPP_MPI
427          call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
428                          COMM_LMDZ,Request(NbRequest),ierr)
429#endif
430!$OMP END CRITICAL (MPI)
431        ENDIF
432 
433        IF (SendDown) THEN
434          NbRequest=NbRequest+1
435           
436          buffer_size=up*iip1*ll
437          allocate(Buffer_Send_down(Buffer_size))
438          call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
439       
440!$OMP CRITICAL (MPI)
441#ifdef CPP_MPI
442          call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
443                          COMM_LMDZ,Request(NbRequest),ierr)
444#endif
445!$OMP END CRITICAL (MPI)
446        ENDIF
447   
448 
449        IF (RecvUp) THEN
450          NbRequest=NbRequest+1
451          buffer_size=up*iip1*ll
452          allocate(Buffer_recv_up(Buffer_size))
453             
454!$OMP CRITICAL (MPI)
455#ifdef CPP_MPI
456          call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
457                          COMM_LMDZ,Request(NbRequest),ierr)
458#endif
459!$OMP END CRITICAL (MPI)
460     
461       
462        ENDIF
463 
464        IF (RecvDown) THEN
465          NbRequest=NbRequest+1
466          buffer_size=down*iip1*ll
467          allocate(Buffer_recv_down(Buffer_size))
468       
469!$OMP CRITICAL (MPI)
470#ifdef CPP_MPI
471          call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
472                          COMM_LMDZ,Request(NbRequest),ierr)
473#endif
474!$OMP END CRITICAL (MPI)
475       
476        ENDIF
477 
478#ifdef CPP_MPI
479        if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
480#endif
481        IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
482        IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
483
484        call VTe(VThallo)
485        call barrier
486     
487      ENDIF  ! using_mpi
488     
489      RETURN
490     
491    end subroutine exchange_Hallo
492   
493
494    subroutine Gather_Field(Field,ij,ll,rank)
495    implicit none
496#include "dimensions.h"
497#include "paramet.h"
498#include "iniprint.h"
499#ifdef CPP_MPI
500    include 'mpif.h'
501#endif   
502      INTEGER :: ij,ll,rank
503      REAL, dimension(ij,ll) :: Field
504      REAL, dimension(:),allocatable :: Buffer_send   
505      REAL, dimension(:),allocatable :: Buffer_Recv
506      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
507      INTEGER :: ierr
508      INTEGER ::i
509     
510      IF (using_mpi) THEN
511
512        if (ij==ip1jmp1) then
513           allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
514           call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
515        else if (ij==ip1jm) then
516           allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
517           call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
518        else
519           write(lunout,*)ij 
520        stop 'erreur dans Gather_Field'
521        endif
522       
523        if (MPI_Rank==rank) then
524          allocate(Buffer_Recv(ij*ll))
525
526!CDIR NOVECTOR
527          do i=0,MPI_Size-1
528             
529            if (ij==ip1jmp1) then
530              Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
531            else if (ij==ip1jm) then
532              Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
533            else
534              stop 'erreur dans Gather_Field'
535            endif
536                   
537            if (i==0) then
538              displ(i)=0
539            else
540              displ(i)=displ(i-1)+Recv_count(i-1)
541            endif
542           
543          enddo
544         
545        else
546          ! Ehouarn: When in debug mode, ifort complains (for call MPI_GATHERV
547          !          below) about Buffer_Recv() being not allocated.
548          !          So make a dummy allocation.
549          allocate(Buffer_Recv(1))
550        endif ! of if (MPI_Rank==rank)
551 
552!$OMP CRITICAL (MPI)
553#ifdef CPP_MPI
554        call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
555                          Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
556#endif
557!$OMP END CRITICAL (MPI)
558     
559        if (MPI_Rank==rank) then                 
560     
561          if (ij==ip1jmp1) then
562            do i=0,MPI_Size-1
563              call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
564                               jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
565            enddo
566          else if (ij==ip1jm) then
567            do i=0,MPI_Size-1
568               call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
569                               min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
570            enddo
571          endif
572        endif
573      ENDIF ! using_mpi
574     
575    end subroutine Gather_Field
576
577
578    subroutine AllGather_Field(Field,ij,ll)
579    implicit none
580#include "dimensions.h"
581#include "paramet.h"   
582#ifdef CPP_MPI
583    include 'mpif.h'
584#endif   
585      INTEGER :: ij,ll
586      REAL, dimension(ij,ll) :: Field
587      INTEGER :: ierr
588     
589      IF (using_mpi) THEN
590        call Gather_Field(Field,ij,ll,0)
591!$OMP CRITICAL (MPI)
592#ifdef CPP_MPI
593      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
594#endif
595!$OMP END CRITICAL (MPI)
596      ENDIF
597     
598    end subroutine AllGather_Field
599   
600   subroutine Broadcast_Field(Field,ij,ll,rank)
601    implicit none
602#include "dimensions.h"
603#include "paramet.h"   
604#ifdef CPP_MPI
605    include 'mpif.h'
606#endif   
607      INTEGER :: ij,ll
608      REAL, dimension(ij,ll) :: Field
609      INTEGER :: rank
610      INTEGER :: ierr
611     
612      IF (using_mpi) THEN
613     
614!$OMP CRITICAL (MPI)
615#ifdef CPP_MPI
616      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
617#endif
618!$OMP END CRITICAL (MPI)
619     
620      ENDIF
621    end subroutine Broadcast_Field
622       
623   
624!  Subroutine verif_hallo(Field,ij,ll,up,down)
625!    implicit none
626!#include "dimensions.h"
627!#include "paramet.h"   
628!    include 'mpif.h'
629!   
630!      INTEGER :: ij,ll
631!      REAL, dimension(ij,ll) :: Field
632!      INTEGER :: up,down
633!     
634!      REAL,dimension(ij,ll): NewField
635!     
636!      NewField=0
637!     
638!      ijb=ij_begin
639!      ije=ij_end
640!      if (pole_nord)
641!      NewField(ij_be       
642
643  end MODULE parallel_lmdz
Note: See TracBrowser for help on using the repository browser.