source: LMDZ5/trunk/libf/dyn3dmem/parallel_lmdz.F90 @ 1823

Last change on this file since 1823 was 1823, checked in by Ehouarn Millour, 11 years ago

Remplacement de parallel.F90 (en conflit avec orchidée) par parallel_lmdz.F90.
UG
.........................................
Renaming parallel.F90 (conflicting with orchidée) into parallel_lmdz.F90.
UG

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