source: LMDZ5/trunk/libf/parallel.F90 @ 1630

Last change on this file since 1630 was 1630, checked in by Laurent Fairhead, 12 years ago

Importation initiale du répertoire dyn3dmem


Initial import of dyn3dmem directory

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