source: LMDZ4/trunk/libf/dyn3dpar/parallel.F90 @ 1134

Last change on this file since 1134 was 1000, checked in by Laurent Fairhead, 16 years ago
  • Modifs sur le parallelisme: masquage dans la physique
  • Inclusion strato
  • mise en coherence etat0
  • le mode offline fonctionne maintenant en parallele,
  • les fichiers de la dynamiques sont correctement sortis et peuvent etre reconstruit avec rebuild
  • la version parallele de la dynamique peut s'executer sans MPI (sur 1 proc)
  • L'OPENMP fonctionne maintenant sans la parallelisation MPI.

YM
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 13.1 KB
Line 
1  module parallel
2  USE mod_const_mpi
3   
4    LOGICAL,SAVE :: using_mpi
5    LOGICAL,SAVE :: using_omp
6   
7    integer, save :: mpi_size
8    integer, save :: mpi_rank
9    integer, save :: jj_begin
10    integer, save :: jj_end
11    integer, save :: jj_nb
12    integer, save :: ij_begin
13    integer, save :: ij_end
14    logical, save :: pole_nord
15    logical, save :: pole_sud
16   
17    integer, allocatable, save, dimension(:) :: jj_begin_para
18    integer, allocatable, save, dimension(:) :: jj_end_para
19    integer, allocatable, save, dimension(:) :: jj_nb_para
20    integer, save :: OMP_CHUNK
21    integer, save :: omp_rank
22    integer, save :: omp_size 
23!$OMP THREADPRIVATE(omp_rank)
24
25 contains
26 
27    subroutine init_parallel
28    USE vampir
29    implicit none
30#ifdef CPP_MPI
31      include 'mpif.h'
32#endif
33#include "dimensions.h"
34#include "paramet.h"
35   
36      integer :: ierr
37      integer :: i,j
38      integer :: type_size
39      integer, dimension(3) :: blocklen,type
40      integer :: comp_id
41
42#ifdef CPP_OMP   
43      INTEGER :: OMP_GET_NUM_THREADS
44      EXTERNAL OMP_GET_NUM_THREADS
45      INTEGER :: OMP_GET_THREAD_NUM
46      EXTERNAL OMP_GET_THREAD_NUM
47#endif 
48
49#ifdef CPP_MPI
50       using_mpi=.TRUE.
51#else
52       using_mpi=.FALSE.
53#endif
54     
55      call InitVampir
56     
57      IF (using_mpi) THEN
58#ifdef CPP_MPI
59        call MPI_COMM_SIZE(COMM_LMDZ,mpi_size,ierr)
60        call MPI_COMM_RANK(COMM_LMDZ,mpi_rank,ierr)
61#endif
62      ELSE
63        mpi_size=1
64        mpi_rank=0
65      ENDIF
66 
67     
68      allocate(jj_begin_para(0:mpi_size-1))
69      allocate(jj_end_para(0:mpi_size-1))
70      allocate(jj_nb_para(0:mpi_size-1))
71     
72      do i=0,mpi_size-1
73        jj_nb_para(i)=(jjm+1)/mpi_size
74        if ( i < MOD((jjm+1),mpi_size) ) jj_nb_para(i)=jj_nb_para(i)+1
75       
76        if (jj_nb_para(i) <= 2 ) then
77         
78         print *,"Arret : le nombre de bande de lattitude par process est trop faible (<2)."
79         print *," ---> diminuez le nombre de CPU ou augmentez la taille en lattitude"
80         
81#ifdef CPP_MPI
82          IF (using_mpi) call MPI_ABORT(COMM_LMDZ,-1, ierr)
83#endif         
84        endif
85       
86      enddo
87     
88!      jj_nb_para(0)=11
89!      jj_nb_para(1)=25
90!      jj_nb_para(2)=25
91!      jj_nb_para(3)=12     
92
93      j=1
94     
95      do i=0,mpi_size-1
96       
97        jj_begin_para(i)=j
98        jj_end_para(i)=j+jj_Nb_para(i)-1
99        j=j+jj_Nb_para(i)
100     
101      enddo
102     
103      jj_begin = jj_begin_para(mpi_rank)
104      jj_end   = jj_end_para(mpi_rank)
105      jj_nb    = jj_nb_para(mpi_rank)
106     
107      ij_begin=(jj_begin-1)*iip1+1
108      ij_end=jj_end*iip1
109     
110      if (mpi_rank.eq.0) then
111        pole_nord=.TRUE.
112      else
113        pole_nord=.FALSE.
114      endif
115     
116      if (mpi_rank.eq.mpi_size-1) then
117        pole_sud=.TRUE.
118      else
119        pole_sud=.FALSE.
120      endif
121       
122      print *,"jj_begin",jj_begin
123      print *,"jj_end",jj_end
124      print *,"ij_begin",ij_begin
125      print *,"ij_end",ij_end
126
127!$OMP PARALLEL
128
129#ifdef CPP_OMP
130!$OMP MASTER
131        omp_size=OMP_GET_NUM_THREADS()
132!$OMP END MASTER
133        omp_rank=OMP_GET_THREAD_NUM()   
134#else   
135        omp_size=1
136        omp_rank=0
137#endif
138!$OMP END PARALLEL         
139   
140    end subroutine init_parallel
141
142   
143    subroutine SetDistrib(jj_Nb_New)
144    implicit none
145
146#include "dimensions.h"
147#include "paramet.h"
148
149      INTEGER,dimension(0:MPI_Size-1) :: jj_Nb_New
150      INTEGER :: i 
151 
152      jj_Nb_Para=jj_Nb_New
153     
154      jj_begin_para(0)=1
155      jj_end_para(0)=jj_Nb_Para(0)
156     
157      do i=1,mpi_size-1
158       
159        jj_begin_para(i)=jj_end_para(i-1)+1
160        jj_end_para(i)=jj_begin_para(i)+jj_Nb_para(i)-1
161     
162      enddo
163     
164      jj_begin = jj_begin_para(mpi_rank)
165      jj_end   = jj_end_para(mpi_rank)
166      jj_nb    = jj_nb_para(mpi_rank)
167     
168      ij_begin=(jj_begin-1)*iip1+1
169      ij_end=jj_end*iip1
170
171    end subroutine SetDistrib
172
173
174
175   
176    subroutine Finalize_parallel
177#ifdef CPP_COUPLE
178    use mod_prism_proto
179#endif
180      use surface_data, only : type_ocean
181      implicit none
182
183      include "dimensions.h"
184      include "paramet.h"
185#ifdef CPP_MPI
186      include 'mpif.h'
187#endif     
188
189      integer :: ierr
190      integer :: i
191      deallocate(jj_begin_para)
192      deallocate(jj_end_para)
193      deallocate(jj_nb_para)
194
195      if (type_ocean == 'couple') then
196#ifdef CPP_COUPLE
197         call prism_terminate_proto(ierr)
198         IF (ierr .ne. PRISM_Ok) THEN
199            call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1)
200         endif
201#endif
202      else
203#ifdef CPP_MPI
204         IF (using_mpi) call MPI_FINALIZE(ierr)
205#endif
206      end if
207     
208    end subroutine Finalize_parallel
209       
210    subroutine Pack_Data(Field,ij,ll,row,Buffer)
211    implicit none
212
213#include "dimensions.h"
214#include "paramet.h"
215
216      integer, intent(in) :: ij,ll,row
217      real,dimension(ij,ll),intent(in) ::Field
218      real,dimension(ll*iip1*row), intent(out) :: Buffer
219           
220      integer :: Pos
221      integer :: i,l
222     
223      Pos=0
224      do l=1,ll
225        do i=1,row*iip1
226          Pos=Pos+1
227          Buffer(Pos)=Field(i,l)
228        enddo
229      enddo
230     
231    end subroutine Pack_data
232     
233    subroutine Unpack_Data(Field,ij,ll,row,Buffer)
234    implicit none
235
236#include "dimensions.h"
237#include "paramet.h"
238
239      integer, intent(in) :: ij,ll,row
240      real,dimension(ij,ll),intent(out) ::Field
241      real,dimension(ll*iip1*row), intent(in) :: Buffer
242           
243      integer :: Pos
244      integer :: i,l
245     
246      Pos=0
247     
248      do l=1,ll
249        do i=1,row*iip1
250          Pos=Pos+1
251          Field(i,l)=Buffer(Pos)
252        enddo
253      enddo
254     
255    end subroutine UnPack_data
256
257   
258    SUBROUTINE barrier
259    IMPLICIT NONE
260#ifdef CPP_MPI
261    INCLUDE 'mpif.h'
262#endif
263    INTEGER :: ierr
264   
265!$OMP CRITICAL (MPI)     
266#ifdef CPP_MPI
267      IF (using_mpi) CALL MPI_Barrier(COMM_LMDZ,ierr)
268#endif
269!$OMP END CRITICAL (MPI)
270   
271    END SUBROUTINE barrier
272       
273     
274    subroutine exchange_hallo(Field,ij,ll,up,down)
275    USE Vampir
276    implicit none
277#include "dimensions.h"
278#include "paramet.h"   
279#ifdef CPP_MPI
280    include 'mpif.h'
281#endif   
282      INTEGER :: ij,ll
283      REAL, dimension(ij,ll) :: Field
284      INTEGER :: up,down
285     
286      INTEGER :: ierr
287      LOGICAL :: SendUp,SendDown
288      LOGICAL :: RecvUp,RecvDown
289      INTEGER, DIMENSION(4) :: Request
290#ifdef CPP_MPI
291      INTEGER, DIMENSION(MPI_STATUS_SIZE,4) :: Status
292#else
293      INTEGER, DIMENSION(1,4) :: Status
294#endif
295      INTEGER :: NbRequest
296      REAL, dimension(:),allocatable :: Buffer_Send_up,Buffer_Send_down
297      REAL, dimension(:),allocatable :: Buffer_Recv_up,Buffer_Recv_down
298      INTEGER :: Buffer_size     
299
300      IF (using_mpi) THEN
301
302        CALL barrier
303     
304        call VTb(VThallo)
305     
306        SendUp=.TRUE.
307        SendDown=.TRUE.
308        RecvUp=.TRUE.
309        RecvDown=.TRUE.
310         
311        IF (pole_nord) THEN
312          SendUp=.FALSE.
313          RecvUp=.FALSE.
314        ENDIF
315   
316        IF (pole_sud) THEN
317          SendDown=.FALSE.
318          RecvDown=.FALSE.
319        ENDIF
320       
321        if (up.eq.0) then
322          SendDown=.FALSE.
323          RecvUp=.FALSE.
324        endif
325     
326        if (down.eq.0) then
327          SendUp=.FALSE.
328          RecvDown=.FALSE.
329        endif
330     
331        NbRequest=0
332 
333        IF (SendUp) THEN
334          NbRequest=NbRequest+1
335          buffer_size=down*iip1*ll
336          allocate(Buffer_Send_up(Buffer_size))
337          call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up)
338!$OMP CRITICAL (MPI)
339#ifdef CPP_MPI
340          call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1,     &
341                          COMM_LMDZ,Request(NbRequest),ierr)
342#endif
343!$OMP END CRITICAL (MPI)
344        ENDIF
345 
346        IF (SendDown) THEN
347          NbRequest=NbRequest+1
348           
349          buffer_size=up*iip1*ll
350          allocate(Buffer_Send_down(Buffer_size))
351          call PACK_Data(Field(ij_end+1-up*iip1,1),ij,ll,up,Buffer_send_down)
352       
353!$OMP CRITICAL (MPI)
354#ifdef CPP_MPI
355          call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1,     &
356                          COMM_LMDZ,Request(NbRequest),ierr)
357#endif
358!$OMP END CRITICAL (MPI)
359        ENDIF
360   
361 
362        IF (RecvUp) THEN
363          NbRequest=NbRequest+1
364          buffer_size=up*iip1*ll
365          allocate(Buffer_recv_up(Buffer_size))
366             
367!$OMP CRITICAL (MPI)
368#ifdef CPP_MPI
369          call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1,  &
370                          COMM_LMDZ,Request(NbRequest),ierr)
371#endif
372!$OMP END CRITICAL (MPI)
373     
374       
375        ENDIF
376 
377        IF (RecvDown) THEN
378          NbRequest=NbRequest+1
379          buffer_size=down*iip1*ll
380          allocate(Buffer_recv_down(Buffer_size))
381       
382!$OMP CRITICAL (MPI)
383#ifdef CPP_MPI
384          call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1,     &
385                          COMM_LMDZ,Request(NbRequest),ierr)
386#endif
387!$OMP END CRITICAL (MPI)
388       
389        ENDIF
390 
391#ifdef CPP_MPI
392        if (NbRequest > 0) call MPI_WAITALL(NbRequest,Request,Status,ierr)
393#endif
394        IF (RecvUp)  call Unpack_Data(Field(ij_begin-up*iip1,1),ij,ll,up,Buffer_Recv_up)
395        IF (RecvDown) call Unpack_Data(Field(ij_end+1,1),ij,ll,down,Buffer_Recv_down) 
396
397        call VTe(VThallo)
398        call barrier
399     
400      ENDIF  ! using_mpi
401     
402      RETURN
403     
404    end subroutine exchange_Hallo
405   
406
407    subroutine Gather_Field(Field,ij,ll,rank)
408    implicit none
409#include "dimensions.h"
410#include "paramet.h"   
411#ifdef CPP_MPI
412    include 'mpif.h'
413#endif   
414      INTEGER :: ij,ll,rank
415      REAL, dimension(ij,ll) :: Field
416      REAL, dimension(:),allocatable :: Buffer_send   
417      REAL, dimension(:),allocatable :: Buffer_Recv
418      INTEGER, dimension(0:MPI_Size-1) :: Recv_count, displ
419      INTEGER :: ierr
420      INTEGER ::i
421     
422      IF (using_mpi) THEN
423
424        if (ij==ip1jmp1) then
425           allocate(Buffer_send(iip1*ll*(jj_end-jj_begin+1)))
426           call Pack_Data(Field(ij_begin,1),ij,ll,jj_end-jj_begin+1,Buffer_send)
427        else if (ij==ip1jm) then
428           allocate(Buffer_send(iip1*ll*(min(jj_end,jjm)-jj_begin+1)))
429           call Pack_Data(Field(ij_begin,1),ij,ll,min(jj_end,jjm)-jj_begin+1,Buffer_send)
430        else
431           print *,ij 
432        stop 'erreur dans Gather_Field'
433        endif
434       
435        if (MPI_Rank==rank) then
436          allocate(Buffer_Recv(ij*ll))
437
438!CDIR NOVECTOR
439          do i=0,MPI_Size-1
440             
441            if (ij==ip1jmp1) then
442              Recv_count(i)=(jj_end_para(i)-jj_begin_para(i)+1)*ll*iip1
443            else if (ij==ip1jm) then
444              Recv_count(i)=(min(jj_end_para(i),jjm)-jj_begin_para(i)+1)*ll*iip1
445            else
446              stop 'erreur dans Gather_Field'
447            endif
448                   
449            if (i==0) then
450              displ(i)=0
451            else
452              displ(i)=displ(i-1)+Recv_count(i-1)
453            endif
454           
455          enddo
456         
457        endif
458 
459!$OMP CRITICAL (MPI)
460#ifdef CPP_MPI
461        call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8,   &
462                          Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr)
463#endif
464!$OMP END CRITICAL (MPI)
465     
466        if (MPI_Rank==rank) then                 
467     
468          if (ij==ip1jmp1) then
469            do i=0,MPI_Size-1
470              call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                 &
471                               jj_end_para(i)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
472            enddo
473          else if (ij==ip1jm) then
474            do i=0,MPI_Size-1
475               call Unpack_Data(Field((jj_begin_para(i)-1)*iip1+1,1),ij,ll,                       &
476                               min(jj_end_para(i),jjm)-jj_begin_para(i)+1,Buffer_Recv(displ(i)+1))
477            enddo
478          endif
479        endif
480      ENDIF ! using_mpi
481     
482    end subroutine Gather_Field
483
484
485    subroutine AllGather_Field(Field,ij,ll)
486    implicit none
487#include "dimensions.h"
488#include "paramet.h"   
489#ifdef CPP_MPI
490    include 'mpif.h'
491#endif   
492      INTEGER :: ij,ll
493      REAL, dimension(ij,ll) :: Field
494      INTEGER :: ierr
495     
496      IF (using_mpi) THEN
497        call Gather_Field(Field,ij,ll,0)
498!$OMP CRITICAL (MPI)
499#ifdef CPP_MPI
500      call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
501#endif
502!$OMP END CRITICAL (MPI)
503      ENDIF
504     
505    end subroutine AllGather_Field
506   
507   subroutine Broadcast_Field(Field,ij,ll,rank)
508    implicit none
509#include "dimensions.h"
510#include "paramet.h"   
511#ifdef CPP_MPI
512    include 'mpif.h'
513#endif   
514      INTEGER :: ij,ll
515      REAL, dimension(ij,ll) :: Field
516      INTEGER :: rank
517      INTEGER :: ierr
518     
519      IF (using_mpi) THEN
520     
521!$OMP CRITICAL (MPI)
522#ifdef CPP_MPI
523      call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr)
524#endif
525!$OMP END CRITICAL (MPI)
526     
527      ENDIF
528    end subroutine Broadcast_Field
529       
530   
531    /* 
532  Subroutine verif_hallo(Field,ij,ll,up,down)
533    implicit none
534#include "dimensions.h"
535#include "paramet.h"   
536    include 'mpif.h'
537   
538      INTEGER :: ij,ll
539      REAL, dimension(ij,ll) :: Field
540      INTEGER :: up,down
541     
542      REAL,dimension(ij,ll): NewField
543     
544      NewField=0
545     
546      ijb=ij_begin
547      ije=ij_end
548      if (pole_nord)
549      NewField(ij_be       
550*/
551  end module parallel
Note: See TracBrowser for help on using the repository browser.