source: LMDZ6/branches/LMDZ_ECRad/libf/dyn3dmem/parallel_lmdz.F90 @ 5452

Last change on this file since 5452 was 4727, checked in by idelkadi, 15 months ago

Merged trunk changes -r4488:4726 LMDZ_ECRad branch

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