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

Last change on this file since 995 was 995, checked in by lsce, 16 years ago
  • Modifications liées au calcul des nouveau sous-fractions
  • Nettoyage de ocean slab : il reste uniquement la version avec glace de mer forcé
  • Nouveaux variables pour distiguer la version et type d'ocean : type_ocean=force/slab/couple, version_ocean=opa8/nemo pour couplé ou version_ocean=sicOBS pour slab

JG

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