source: LMDZ6/branches/Amaury_dev/libf/dyn3dmem/parallel_lmdz.F90 @ 5099

Last change on this file since 5099 was 5099, checked in by abarral, 2 months ago

Replace most uses of CPP_DUST by the corresponding logical defined in lmdz_cppkeys_wrapper.F90
Convert several files from .F to .f90 to allow Dust to compile w/o rrtm/ecrad
Create lmdz_yoerad.f90
(lint) Remove "!" on otherwise empty line

  • 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.6 KB
Line 
1
2! $Id$
3
4  MODULE parallel_lmdz
5  USE mod_const_mpi
6  USE lmdz_mpi, ONLY : using_mpi
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   
13    INTEGER,PARAMETER :: halo_max=3
14   
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)
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
90    USE lmdz_mpi
91    implicit none
92      INCLUDE "dimensions.h"
93      INCLUDE "paramet.h"
94      INCLUDE "iniprint.h"
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 /= 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       
147        if (jj_nb_para(i) <= 1 ) then
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)
153
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==0) then
181        pole_nord=.TRUE.
182      else
183        pole_nord=.FALSE.
184      endif
185     
186      if (mpi_rank==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
218!$OMP BARRIER
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
227
228!$OMP MASTER
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)
232!$OMP END MASTER
233!$OMP BARRIER       
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     
241      IF ((mpi_rank==0).and.(omp_rank==0)) THEN
242        is_master=.true.
243      ELSE
244        is_master=.false.
245      ENDIF
246     
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
385    USE lmdz_mpi
386    ! ug Pour les sorties XIOS
387        USE wxios
388
389#ifdef CPP_COUPLE
390! Use of Oasis-MCT coupler
391#if defined CPP_OMCT
392    use mod_prism
393#else
394    use mod_prism_proto
395#endif
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
405      include "dimensions.h"
406      include "paramet.h"
407
408      integer :: ierr
409      integer :: i
410
411      if (allocated(jj_begin_para)) deallocate(jj_begin_para)
412      if (allocated(jj_end_para))   deallocate(jj_end_para)
413      if (allocated(jj_nb_para))    deallocate(jj_nb_para)
414
415      if (type_ocean == 'couple') then
416#ifdef CPP_COUPLE
417        IF (using_xios) THEN
418          !Fermeture propre de XIOS
419          CALL wxios_close()
420          CALL prism_terminate_proto(ierr)
421          IF (ierr .ne. PRISM_Ok) THEN
422            CALL abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
423          ENDIF
424        ELSE
425           call prism_terminate_proto(ierr)
426           IF (ierr .ne. PRISM_Ok) THEN
427              call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
428           endif
429        ENDIF
430#else
431        call abort_gcm('Finalize_parallel','type_ocean = couple but CPP_COUPLE not defined',1)
432#endif
433      else
434        IF (using_xios) THEN
435          !Fermeture propre de XIOS
436          CALL wxios_close()
437        ENDIF
438        IF (using_mpi) call MPI_FINALIZE(ierr)
439      end if
440     
441    end subroutine Finalize_parallel
442       
443    subroutine Pack_Data(Field,ij,ll,row,Buffer)
444    implicit none
445
446      INCLUDE "dimensions.h"
447      INCLUDE "paramet.h"
448
449      integer, intent(in) :: ij,ll,row
450      real,dimension(ij,ll),intent(in) ::Field
451      real,dimension(ll*iip1*row), intent(out) :: Buffer
452           
453      integer :: Pos
454      integer :: i,l
455     
456      Pos=0
457      do l=1,ll
458        do i=1,row*iip1
459          Pos=Pos+1
460          Buffer(Pos)=Field(i,l)
461        enddo
462      enddo
463     
464    end subroutine Pack_data
465     
466    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
467    implicit none
468
469      INCLUDE "dimensions.h"
470      INCLUDE "paramet.h"
471
472      integer, intent(in) :: ij,ll,row
473      real,dimension(ij,ll),intent(out) ::Field
474      real,dimension(ll*iip1*row), intent(in) :: Buffer
475           
476      integer :: Pos
477      integer :: i,l
478     
479      Pos=0
480     
481      do l=1,ll
482        do i=1,row*iip1
483          Pos=Pos+1
484          Field(i,l)=Buffer(Pos)
485        enddo
486      enddo
487     
488    end subroutine UnPack_data
489
490   
491    SUBROUTINE barrier
492    USE lmdz_mpi
493    IMPLICIT NONE
494    INTEGER :: ierr
495   
496!$OMP CRITICAL (MPI)     
497      IF (using_mpi) CALL MPI_Barrier(COMM_LMDZ,ierr)
498!$OMP END CRITICAL (MPI)
499   
500    END SUBROUTINE barrier
501       
502     
503    subroutine exchange_hallo(Field,ij,ll,up,down)
504    USE lmdz_mpi
505    USE Vampir
506    implicit none
507      INCLUDE "dimensions.h"
508      INCLUDE "paramet.h"   
509      INTEGER :: ij,ll
510      REAL, dimension(ij,ll) :: Field
511      INTEGER :: up,down
512     
513      INTEGER :: ierr
514      LOGICAL :: SendUp,SendDown
515      LOGICAL :: RecvUp,RecvDown
516      INTEGER, DIMENSION(4) :: Request
517      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
518
519      INTEGER :: NbRequest
520      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
521      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
522      INTEGER :: Buffer_size     
523
524      IF (using_mpi) THEN
525
526        CALL barrier
527     
528        call VTb(VThallo)
529     
530        SendUp=.TRUE.
531        SendDown=.TRUE.
532        RecvUp=.TRUE.
533        RecvDown=.TRUE.
534         
535        IF (pole_nord) THEN
536          SendUp=.FALSE.
537          RecvUp=.FALSE.
538        ENDIF
539   
540        IF (pole_sud) THEN
541          SendDown=.FALSE.
542          RecvDown=.FALSE.
543        ENDIF
544       
545        if (up==0) then
546          SendDown=.FALSE.
547          RecvUp=.FALSE.
548        endif
549     
550        if (down==0) then
551          SendUp=.FALSE.
552          RecvDown=.FALSE.
553        endif
554     
555        NbRequest=0
556 
557        IF (SendUp) THEN
558          NbRequest=NbRequest+1
559          buffer_size=down*iip1*ll
560          allocate(Buffer_Send_up(Buffer_size))
561          call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
562!$OMP CRITICAL (MPI)
563          call MPI_ISEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
564                          COMM_LMDZ,Request(NbRequest),ierr)
565!$OMP END CRITICAL (MPI)
566        ENDIF
567 
568        IF (SendDown) THEN
569          NbRequest=NbRequest+1
570           
571          buffer_size=up*iip1*ll
572          allocate(Buffer_Send_down(Buffer_size))
573          call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
574       
575!$OMP CRITICAL (MPI)
576          call MPI_ISEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
577                          COMM_LMDZ,Request(NbRequest),ierr)
578!$OMP END CRITICAL (MPI)
579        ENDIF
580   
581 
582        IF (RecvUp) THEN
583          NbRequest=NbRequest+1
584          buffer_size=up*iip1*ll
585          allocate(Buffer_recv_up(Buffer_size))
586             
587!$OMP CRITICAL (MPI)
588          call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
589                          COMM_LMDZ,Request(NbRequest),ierr)
590!$OMP END CRITICAL (MPI)
591     
592       
593        ENDIF
594 
595        IF (RecvDown) THEN
596          NbRequest=NbRequest+1
597          buffer_size=down*iip1*ll
598          allocate(Buffer_recv_down(Buffer_size))
599       
600!$OMP CRITICAL (MPI)
601          call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
602                          COMM_LMDZ,Request(NbRequest),ierr)
603!$OMP END CRITICAL (MPI)
604       
605        ENDIF
606 
607        if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
608        IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
609        IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
610
611        call VTe(VThallo)
612        call barrier
613     
614      ENDIF  ! using_mpi
615     
616      RETURN
617     
618    end subroutine exchange_Hallo
619   
620
621    subroutine Gather_Field(Field,ij,ll,rank)
622    USE lmdz_mpi
623    implicit none
624    INCLUDE "dimensions.h"
625    INCLUDE "paramet.h"
626    INCLUDE "iniprint.h"
627      INTEGER :: ij,ll,rank
628      REAL, dimension(ij,ll) :: Field
629      REAL, dimension(:),allocatable :: Buffer_send   
630      REAL, dimension(:),allocatable :: Buffer_Recv
631      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
632      INTEGER :: ierr
633      INTEGER ::i
634     
635      IF (using_mpi) THEN
636
637        if (ij==ip1jmp1) then
638           allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
639           call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
640        else if (ij==ip1jm) then
641           allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
642           call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
643        else
644           write(lunout,*)ij 
645        CALL abort_gcm("parallel_lmdz","erreur dans Gather_Field",1)
646        endif
647       
648        if (MPI_Rank==rank) then
649          allocate(Buffer_Recv(ij*ll))
650
651!CDIR NOVECTOR
652          do i=0,MPI_Size-1
653             
654            if (ij==ip1jmp1) then
655              Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
656            else if (ij==ip1jm) then
657              Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
658            else
659              CALL abort_gcm("parallel_lmdz","erreur dans Gather_Field",1)
660            endif
661                   
662            if (i==0) then
663              displ(i)=0
664            else
665              displ(i)=displ(i-1)+Recv_count(i-1)
666            endif
667           
668          enddo
669         
670        else
671          ! Ehouarn: When in debug mode, ifort complains (for call MPI_GATHERV
672          !          below) about Buffer_Recv() being not allocated.
673          !          So make a dummy allocation.
674          allocate(Buffer_Recv(1))
675        endif ! of if (MPI_Rank==rank)
676 
677!$OMP CRITICAL (MPI)
678        call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
679                          Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
680!$OMP END CRITICAL (MPI)
681     
682        if (MPI_Rank==rank) then                 
683     
684          if (ij==ip1jmp1) then
685            do i=0,MPI_Size-1
686              call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
687                               jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
688            enddo
689          else if (ij==ip1jm) then
690            do i=0,MPI_Size-1
691               call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
692                               min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
693            enddo
694          endif
695        endif
696      ENDIF ! using_mpi
697     
698    end subroutine Gather_Field
699
700
701    subroutine AllGather_Field(Field,ij,ll)
702    USE lmdz_mpi
703    implicit none
704    INCLUDE "dimensions.h"
705    INCLUDE "paramet.h"   
706      INTEGER :: ij,ll
707      REAL, dimension(ij,ll) :: Field
708      INTEGER :: ierr
709     
710      IF (using_mpi) THEN
711        call Gather_Field(Field,ij,ll,0)
712!$OMP CRITICAL (MPI)
713      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
714!$OMP END CRITICAL (MPI)
715      ENDIF
716     
717    end subroutine AllGather_Field
718   
719   subroutine Broadcast_Field(Field,ij,ll,rank)
720    USE lmdz_mpi
721    implicit none
722    INCLUDE "dimensions.h"
723    INCLUDE "paramet.h"   
724      INTEGER :: ij,ll
725      REAL, dimension(ij,ll) :: Field
726      INTEGER :: rank
727      INTEGER :: ierr
728     
729      IF (using_mpi) THEN
730     
731!$OMP CRITICAL (MPI)
732      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
733!$OMP END CRITICAL (MPI)
734     
735      ENDIF
736    end subroutine Broadcast_Field
737       
738   
739!  Subroutine verif_hallo(Field,ij,ll,up,down)
740!    USE lmdz_mpi
741!    implicit none
742!      INCLUDE "dimensions.h"
743!      INCLUDE "paramet.h"   
744
745!      INTEGER :: ij,ll
746!      REAL, dimension(ij,ll) :: Field
747!      INTEGER :: up,down
748
749!      REAL,dimension(ij,ll): NewField
750
751!      NewField=0
752
753!      ijb=ij_begin
754!      ije=ij_end
755!      if (pole_nord)
756!      NewField(ij_be       
757
758  end MODULE parallel_lmdz
Note: See TracBrowser for help on using the repository browser.