source: LMDZ6/branches/blowing_snow/libf/dyn3dmem/parallel_lmdz.F90 @ 5453

Last change on this file since 5453 was 4469, checked in by Laurent Fairhead, 22 months ago

Replaced STOP instructions by calls to abort_gcm for a cleaner exit

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