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

Last change on this file since 5106 was 5105, checked in by abarral, 4 months ago

Replace 1DUTILS.h by module lmdz_1dutils.f90
Replace 1DConv.h by module lmdz_old_1dconv.f90 (it's only used by old_* files)
Convert *.F to *.f90
Fix gradsdef.h formatting
Remove unnecessary "RETURN" at the end of functions/subroutines

  • 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.5 KB
Line 
1
2! $Id$
3
4  MODULE parallel_lmdz
5  USE mod_const_mpi
6  USE lmdz_mpi, ONLY: using_mpi
7      use IOIPSL
8    INTEGER,PARAMETER :: halo_max=3
9   
10    LOGICAL,SAVE :: using_omp ! .TRUE. if using OpenMP
11    LOGICAL,SAVE :: is_master ! .TRUE. if the core is both MPI & OpenMP master
12!$OMP THREADPRIVATE(is_master)
13   
14    integer, save :: mpi_size
15    integer, save :: mpi_rank
16    integer, save :: jj_begin
17    integer, save :: jj_end
18    integer, save :: jj_nb
19    integer, save :: ij_begin
20    integer, save :: ij_end
21    logical, save :: pole_nord
22    logical, save :: pole_sud
23
24    integer,save  :: jjb_u
25    integer,save  :: jje_u
26    integer,save  :: jjnb_u
27    integer,save  :: jjb_v
28    integer,save  :: jje_v
29    integer,save  :: jjnb_v   
30
31    integer,save  :: ijb_u
32    integer,save  :: ije_u
33    integer,save  :: ijnb_u   
34   
35    integer,save  :: ijb_v
36    integer,save  :: ije_v
37    integer,save  :: ijnb_v   
38     
39   
40    integer, allocatable, save, dimension(:) :: jj_begin_para
41    integer, allocatable, save, dimension(:) :: jj_end_para
42    integer, allocatable, save, dimension(:) :: jj_nb_para
43    integer, save :: OMP_CHUNK
44    integer, save :: omp_rank
45    integer, save :: omp_size 
46!$OMP THREADPRIVATE(omp_rank)
47
48    TYPE distrib
49      integer :: jj_begin
50      integer :: jj_end
51      integer :: jj_nb
52      integer :: ij_begin
53      integer :: ij_end
54
55      integer  :: jjb_u
56      integer  :: jje_u
57      integer  :: jjnb_u
58      integer  :: jjb_v
59      integer  :: jje_v
60      integer  :: jjnb_v   
61 
62      integer  :: ijb_u
63      integer  :: ije_u
64      integer  :: ijnb_u   
65   
66      integer  :: ijb_v
67      integer  :: ije_v
68      integer  :: ijnb_v   
69     
70   
71      integer, pointer :: jj_begin_para(:) => NULL()
72      integer, pointer :: jj_end_para(:) => NULL()
73      integer, pointer :: jj_nb_para(:) => NULL()
74    END TYPE distrib 
75   
76    INTERFACE ASSIGNMENT (=)
77      MODULE PROCEDURE copy_distrib
78    END INTERFACE
79    TYPE(distrib),SAVE :: current_dist
80   
81 contains
82 
83    SUBROUTINE init_parallel
84    USE vampir
85    USE lmdz_mpi
86    implicit none
87      INCLUDE "dimensions.h"
88      INCLUDE "paramet.h"
89      INCLUDE "iniprint.h"
90
91      integer :: ierr
92      integer :: i,j
93      integer :: type_size
94      integer, dimension(3) :: blocklen,type
95      integer :: comp_id
96      character(len=4)  :: num
97      character(len=20) :: filename
98 
99#ifdef CPP_OMP   
100      INTEGER :: OMP_GET_NUM_THREADS
101      EXTERNAL OMP_GET_NUM_THREADS
102      INTEGER :: OMP_GET_THREAD_NUM
103      EXTERNAL OMP_GET_THREAD_NUM
104#endif 
105
106#ifdef CPP_OMP
107       using_OMP=.TRUE.
108#else
109       using_OMP=.FALSE.
110#endif
111     
112      CALL InitVampir
113     
114      IF (using_mpi) THEN
115        CALL MPI_COMM_SIZE(COMM_LMDZ,mpi_size,ierr)
116        CALL MPI_COMM_RANK(COMM_LMDZ,mpi_rank,ierr)
117      ELSE
118        mpi_size=1
119        mpi_rank=0
120      ENDIF
121
122
123! Open text output file with mpi_rank in suffix of file name
124      IF (lunout /= 5 .and. lunout /= 6) THEN
125         WRITE(num,'(I4.4)') mpi_rank
126         filename='lmdz.out_'//num
127         IF (mpi_rank /= 0) THEN
128            OPEN(UNIT=lunout,FILE=TRIM(filename),ACTION='write', &
129               STATUS='unknown',FORM='formatted',IOSTAT=ierr)
130         ENDIF
131      ENDIF
132
133     
134      allocate(jj_begin_para(0:mpi_size-1))
135      allocate(jj_end_para(0:mpi_size-1))
136      allocate(jj_nb_para(0:mpi_size-1))
137     
138      do i=0,mpi_size-1
139        jj_nb_para(i)=(jjm+1)/mpi_size
140        if ( i < MOD((jjm+1),mpi_size) ) jj_nb_para(i)=jj_nb_para(i)+1
141       
142        if (jj_nb_para(i) <= 1 ) then
143         
144         write(lunout,*)"Arret : le nombre de bande de lattitude par process est trop faible (<2)."
145         write(lunout,*)" ---> diminuez le nombre de CPU ou augmentez la taille en lattitude"
146         
147          IF (using_mpi) CALL MPI_ABORT(COMM_LMDZ,-1, ierr)
148
149        endif
150       
151      enddo
152     
153!      jj_nb_para(0)=11
154!      jj_nb_para(1)=25
155!      jj_nb_para(2)=25
156!      jj_nb_para(3)=12     
157
158      j=1
159     
160      do i=0,mpi_size-1
161       
162        jj_begin_para(i)=j
163        jj_end_para(i)=j+jj_Nb_para(i)-1
164        j=j+jj_Nb_para(i)
165     
166      enddo
167     
168      jj_begin = jj_begin_para(mpi_rank)
169      jj_end   = jj_end_para(mpi_rank)
170      jj_nb    = jj_nb_para(mpi_rank)
171     
172      ij_begin=(jj_begin-1)*iip1+1
173      ij_end=jj_end*iip1
174     
175      if (mpi_rank==0) then
176        pole_nord=.TRUE.
177      else
178        pole_nord=.FALSE.
179      endif
180     
181      if (mpi_rank==mpi_size-1) then
182        pole_sud=.TRUE.
183      else
184        pole_sud=.FALSE.
185      endif
186       
187      write(lunout,*)"init_parallel: jj_begin",jj_begin
188      write(lunout,*)"init_parallel: jj_end",jj_end
189      write(lunout,*)"init_parallel: ij_begin",ij_begin
190      write(lunout,*)"init_parallel: ij_end",ij_end
191      jjb_u=MAX(jj_begin-halo_max,1)
192      jje_u=MIN(jj_end+halo_max,jjp1)
193      jjnb_u=jje_u-jjb_u+1
194
195      jjb_v=MAX(jj_begin-halo_max,1)
196      jje_v=MIN(jj_end+halo_max,jjm)
197      jjnb_v=jje_v-jjb_v+1
198
199      ijb_u=MAX(ij_begin-halo_max*iip1,1)
200      ije_u=MIN(ij_end+halo_max*iip1,ip1jmp1)
201      ijnb_u=ije_u-ijb_u+1
202
203      ijb_v=MAX(ij_begin-halo_max*iip1,1)
204      ije_v=MIN(ij_end+halo_max*iip1,ip1jm)
205      ijnb_v=ije_v-ijb_v+1
206     
207!$OMP PARALLEL
208
209#ifdef CPP_OMP
210!$OMP MASTER
211        omp_size=OMP_GET_NUM_THREADS()
212!$OMP END MASTER
213!$OMP BARRIER
214        omp_rank=OMP_GET_THREAD_NUM() 
215
216!Config  Key  = omp_chunk
217!Config  Desc = taille des blocs openmp
218!Config  Def  = 1
219!Config  Help = defini la taille des packets d'itération openmp
220!Config         distribue a chaque tache lors de l'entree dans une
221!Config         boucle parallelisee
222
223!$OMP MASTER
224      omp_chunk=(llm+1)/omp_size
225      IF (MOD(llm+1,omp_size)/=0) omp_chunk=omp_chunk+1
226      CALL getin('omp_chunk',omp_chunk)
227!$OMP END MASTER
228!$OMP BARRIER       
229#else   
230        omp_size=1
231        omp_rank=0
232#endif
233!$OMP END PARALLEL         
234      CALL create_distrib(jj_nb_para,current_dist)
235     
236      IF ((mpi_rank==0).and.(omp_rank==0)) THEN
237        is_master=.TRUE.
238      ELSE
239        is_master=.FALSE.
240      ENDIF
241     
242    END SUBROUTINE  init_parallel
243
244    SUBROUTINE create_distrib(jj_nb_new,d)
245    IMPLICIT NONE
246      INCLUDE "dimensions.h"
247      INCLUDE "paramet.h"
248     
249      INTEGER,INTENT(IN) :: jj_Nb_New(0:MPI_Size-1)
250      TYPE(distrib),INTENT(INOUT) :: d
251      INTEGER :: i 
252 
253      IF (.NOT. ASSOCIATED(d%jj_nb_para)) ALLOCATE(d%jj_nb_para(0:MPI_Size-1))
254      IF (.NOT. ASSOCIATED(d%jj_begin_para)) ALLOCATE(d%jj_begin_para(0:MPI_Size-1))
255      IF (.NOT. ASSOCIATED(d%jj_end_para)) ALLOCATE(d%jj_end_para(0:MPI_Size-1))
256     
257      d%jj_Nb_Para=jj_Nb_New
258     
259      d%jj_begin_para(0)=1
260      d%jj_end_para(0)=d%jj_Nb_Para(0)
261     
262      do i=1,mpi_size-1
263       
264        d%jj_begin_para(i)=d%jj_end_para(i-1)+1
265        d%jj_end_para(i)=d%jj_begin_para(i)+d%jj_Nb_para(i)-1
266     
267      enddo
268     
269      d%jj_begin = d%jj_begin_para(mpi_rank)
270      d%jj_end   = d%jj_end_para(mpi_rank)
271      d%jj_nb    = d%jj_nb_para(mpi_rank)
272     
273      d%ij_begin=(d%jj_begin-1)*iip1+1
274      d%ij_end=d%jj_end*iip1
275
276      d%jjb_u=MAX(d%jj_begin-halo_max,1)
277      d%jje_u=MIN(d%jj_end+halo_max,jjp1)
278      d%jjnb_u=d%jje_u-d%jjb_u+1
279
280      d%jjb_v=MAX(d%jj_begin-halo_max,1)
281      d%jje_v=MIN(d%jj_end+halo_max,jjm)
282      d%jjnb_v=d%jje_v-d%jjb_v+1
283
284      d%ijb_u=MAX(d%ij_begin-halo_max*iip1,1)
285      d%ije_u=MIN(d%ij_end+halo_max*iip1,ip1jmp1)
286      d%ijnb_u=d%ije_u-d%ijb_u+1
287
288      d%ijb_v=MAX(d%ij_begin-halo_max*iip1,1)
289      d%ije_v=MIN(d%ij_end+halo_max*iip1,ip1jm)
290      d%ijnb_v=d%ije_v-d%ijb_v+1     
291
292    END SUBROUTINE create_distrib
293
294     
295    SUBROUTINE Set_Distrib(d)
296    IMPLICIT NONE
297
298    INCLUDE "dimensions.h"
299    INCLUDE "paramet.h"
300    TYPE(distrib),INTENT(IN) :: d
301
302      jj_begin = d%jj_begin
303      jj_end = d%jj_end
304      jj_nb = d%jj_nb
305      ij_begin = d%ij_begin
306      ij_end = d%ij_end
307
308      jjb_u = d%jjb_u
309      jje_u = d%jje_u
310      jjnb_u = d%jjnb_u
311      jjb_v = d%jjb_v
312      jje_v = d%jje_v
313      jjnb_v = d%jjnb_v
314 
315      ijb_u = d%ijb_u
316      ije_u = d%ije_u
317      ijnb_u = d%ijnb_u
318   
319      ijb_v = d%ijb_v
320      ije_v = d%ije_v
321      ijnb_v = d%ijnb_v
322     
323   
324      jj_begin_para(:) = d%jj_begin_para(:)
325      jj_end_para(:) = d%jj_end_para(:)
326      jj_nb_para(:) = d%jj_nb_para(:)
327      current_dist=d
328
329    END SUBROUTINE Set_Distrib
330
331    SUBROUTINE copy_distrib(dist,new_dist)
332    IMPLICIT NONE
333
334    INCLUDE "dimensions.h"
335    INCLUDE "paramet.h"
336    TYPE(distrib),INTENT(INOUT) :: dist
337    TYPE(distrib),INTENT(IN) :: new_dist
338
339     dist%jj_begin = new_dist%jj_begin
340     dist%jj_end = new_dist%jj_end
341     dist%jj_nb = new_dist%jj_nb
342     dist%ij_begin = new_dist%ij_begin
343     dist%ij_end = new_dist%ij_end
344
345     dist%jjb_u = new_dist%jjb_u
346     dist%jje_u = new_dist%jje_u
347     dist%jjnb_u = new_dist%jjnb_u
348     dist%jjb_v = new_dist%jjb_v
349     dist%jje_v = new_dist%jje_v
350     dist%jjnb_v = new_dist%jjnb_v
351   
352     dist%ijb_u = new_dist%ijb_u
353     dist%ije_u = new_dist%ije_u
354     dist%ijnb_u = new_dist%ijnb_u
355     
356     dist%ijb_v = new_dist%ijb_v
357     dist%ije_v = new_dist%ije_v
358     dist%ijnb_v = new_dist%ijnb_v
359         
360     
361     dist%jj_begin_para(:) = new_dist%jj_begin_para(:)
362     dist%jj_end_para(:) = new_dist%jj_end_para(:)
363     dist%jj_nb_para(:) = new_dist%jj_nb_para(:)
364 
365    END SUBROUTINE copy_distrib
366   
367   
368    SUBROUTINE get_current_distrib(d)
369    IMPLICIT NONE
370
371    INCLUDE "dimensions.h"
372    INCLUDE "paramet.h"
373    TYPE(distrib),INTENT(OUT) :: d
374
375     d=current_dist
376
377    END SUBROUTINE get_current_distrib
378   
379    SUBROUTINE Finalize_parallel
380    USE lmdz_mpi
381    ! ug Pour les sorties XIOS
382        USE wxios
383
384#ifdef CPP_COUPLE
385! Use of Oasis-MCT coupler
386#if defined CPP_OMCT
387    use mod_prism
388#else
389    use mod_prism_proto
390#endif
391! Ehouarn: surface_data module is in 'phylmd' ...
392      use surface_data, ONLY: type_ocean
393      implicit none
394#else
395      implicit none
396! without the surface_data module, we declare (and set) a dummy 'type_ocean'
397      character(len=6),parameter :: type_ocean="dummy"
398#endif
399
400      include "dimensions.h"
401      include "paramet.h"
402
403      integer :: ierr
404      integer :: i
405
406      if (allocated(jj_begin_para)) deallocate(jj_begin_para)
407      if (allocated(jj_end_para))   deallocate(jj_end_para)
408      if (allocated(jj_nb_para))    deallocate(jj_nb_para)
409
410      if (type_ocean == 'couple') then
411#ifdef CPP_COUPLE
412        IF (using_xios) THEN
413          !Fermeture propre de XIOS
414          CALL wxios_close()
415          CALL prism_terminate_proto(ierr)
416          IF (ierr .ne. PRISM_Ok) THEN
417            CALL abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
418          ENDIF
419        ELSE
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        ENDIF
425#else
426        CALL abort_gcm('Finalize_parallel','type_ocean = couple but CPP_COUPLE not defined',1)
427#endif
428      else
429        IF (using_xios) THEN
430          !Fermeture propre de XIOS
431          CALL wxios_close()
432        ENDIF
433        IF (using_mpi) CALL MPI_FINALIZE(ierr)
434      end if
435     
436    END SUBROUTINE  Finalize_parallel
437       
438    SUBROUTINE Pack_Data(Field,ij,ll,row,Buffer)
439    implicit none
440
441      INCLUDE "dimensions.h"
442      INCLUDE "paramet.h"
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
464      INCLUDE "dimensions.h"
465      INCLUDE "paramet.h"
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
487    USE lmdz_mpi
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)
499    USE lmdz_mpi
500    USE Vampir
501    implicit none
502      INCLUDE "dimensions.h"
503      INCLUDE "paramet.h"   
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
513
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==0) then
541          SendDown=.FALSE.
542          RecvUp=.FALSE.
543        endif
544     
545        if (down==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)
558          CALL MPI_ISEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
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)
571          CALL MPI_ISEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
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
612     
613    END SUBROUTINE  exchange_Hallo
614   
615
616    SUBROUTINE Gather_Field(Field,ij,ll,rank)
617    USE lmdz_mpi
618    implicit none
619    INCLUDE "dimensions.h"
620    INCLUDE "paramet.h"
621    INCLUDE "iniprint.h"
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 
640        CALL abort_gcm("parallel_lmdz","erreur dans Gather_Field",1)
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
654              CALL abort_gcm("parallel_lmdz","erreur dans Gather_Field",1)
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)
697    USE lmdz_mpi
698    implicit none
699    INCLUDE "dimensions.h"
700    INCLUDE "paramet.h"   
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)
715    USE lmdz_mpi
716    implicit none
717    INCLUDE "dimensions.h"
718    INCLUDE "paramet.h"   
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)
735!    USE lmdz_mpi
736!    implicit none
737!      INCLUDE "dimensions.h"
738!      INCLUDE "paramet.h"   
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.