source: LMDZ4/branches/V3_test/libf/phylmd/mod_phys_mpi.F90 @ 1358

Last change on this file since 1358 was 718, checked in by Laurent Fairhead, 18 years ago

Corrections bugs divers YM
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 7.5 KB
Line 
1module mod_phys_mpi
2 
3  logical, save :: monocpu
4!  logical, save :: parallel
5   
6  INTEGER,SAVE :: iiphy_begin
7  INTEGER,SAVE :: iiphy_end
8  INTEGER,SAVE :: jjphy_begin
9  INTEGER,SAVE :: jjphy_end
10  INTEGER,SAVE :: jjphy_nb
11  INTEGER,SAVE :: ijphy_begin
12  INTEGER,SAVE :: ijphy_end
13  INTEGER,SAVE :: ijphy_nb
14  INTEGER,SAVE :: klon_begin
15  INTEGER,SAVE :: klon_end
16  INTEGER,SAVE :: klon_mpi
17 
18  INTEGER,SAVE,allocatable,dimension(:) :: jjphy_para_nb
19  INTEGER,SAVE,allocatable,dimension(:) :: jjphy_para_begin
20  INTEGER,SAVE,allocatable,dimension(:) :: jjphy_para_end
21  INTEGER,SAVE,dimension(:),allocatable :: Liste_i
22  INTEGER,SAVE,dimension(:),allocatable :: Liste_j
23  INTEGER,SAVE,dimension(:),allocatable :: klon_para_nb
24  INTEGER,SAVE,dimension(:),allocatable :: klon_para_begin
25  INTEGER,SAVE,dimension(:),allocatable :: klon_para_end
26
27 
28  INTEGER,save :: phy_rank
29  INTEGER,save :: phy_size
30 
31
32  INTEGER,save :: klev, klevp1,klevm1,klon,kfdia,kidia
33  INTEGER,save :: klon2
34  INTEGER,save :: nbtr
35  INTEGER,save :: kdlon, kflev
36!$OMP THREADPRIVATE(klon,kfdia,kidia,kdlon)
37
38 
39contains
40 
41   subroutine InitDimphy
42
43#ifdef CPP_PARA
44   use parallel
45#endif
46
47   implicit none
48#include "dimensions90.h"
49#include "paramet90.h"
50    integer :: rank
51    integer :: i,j,Pos,Index
52    logical,save :: First=.true.
53#ifndef CPP_PARA
54   logical, save :: pole_nord = .true.
55   logical, save :: pole_sud = .true.
56#endif
57
58#ifdef CPP_PARA     
59      monocpu=.false.
60#else
61      monocpu = .true.
62#endif
63     
64      if (monocpu) then
65       
66        phy_rank=0
67        phy_size=1
68       
69#ifdef CPP_PARA     
70      else
71 
72        phy_rank=mpi_rank
73        phy_size=mpi_size
74#endif
75     
76      endif
77
78      if (First) then
79        allocate(jjphy_para_nb(0:phy_size-1))
80        allocate(jjphy_para_begin(0:phy_size-1))
81        allocate(jjphy_para_end(0:phy_size-1))
82        allocate(klon_para_nb(0:phy_size-1))
83        allocate(klon_para_begin(0:phy_size-1))
84        allocate(klon_para_end(0:phy_size-1))
85      endif
86       
87      klev=llm
88      klevp1=klev+1
89      klevm1=klev-1
90      klon2=iim*(jjm-1)+2
91     
92     
93!        do i=0,phy_size-1
94!          klon_para_nb(i)=(klon2/phy_size)
95!          if (i<MOD(klon2,phy_size)) klon_para_nb(i)=klon_para_nb(i)+1
96!          if (i==0) then
97!            klon_para_begin(i)=1
98!          else
99!            klon_para_begin(i)=klon_para_end(i-1)+1
100!          endif
101!          klon_para_end(i)=klon_para_begin(i)+klon_para_nb(i)-1
102!        enddo
103     
104      if (First) then
105        do i=0,phy_size-1
106          klon_para_nb(i)=(klon2/phy_size)
107          if (i<MOD(klon2,phy_size)) klon_para_nb(i)=klon_para_nb(i)+1
108        enddo
109      endif
110       
111      do i=0,phy_size-1
112        if (i==0) then
113          klon_para_begin(i)=1
114        else
115          klon_para_begin(i)=klon_para_end(i-1)+1
116        endif
117        klon_para_end(i)=klon_para_begin(i)+klon_para_nb(i)-1
118      enddo
119       
120      klon=klon_para_nb(phy_rank)
121      klon_mpi=klon_para_nb(phy_rank)
122      klon_begin=klon_para_begin(phy_rank)
123      klon_end=klon_para_end(phy_rank)
124       
125      nbtr=nqmx-2+1/(nqmx-1)
126
127      kflev=klev
128   
129      if (.not.first) then
130        deallocate(liste_i,liste_j)
131      endif
132     
133     
134      allocate(Liste_i(klon))
135      allocate(Liste_j(klon))
136   
137      Index=1
138      if (Pole_Nord) then
139        Liste_i(1)=1
140        Liste_j(1)=1
141        iiphy_begin=1
142        jjphy_begin=1
143        ijphy_begin=1
144        Index=Index+1
145      endif
146   
147      Pos=2
148      do j=2,jjm
149        do i=1,iim
150     
151          if (Pos==klon_begin) then
152            iiphy_begin=i
153            jjphy_begin=j
154            ijphy_begin=j*iip1+i
155          endif
156       
157          if (Pos==klon_end) then
158            iiphy_end=i
159            jjphy_end=j
160            ijphy_end=j*iip1+i
161          endif
162       
163          if (Pos>=klon_begin .AND. Pos<=klon_end) then
164            Liste_i(Index)=i
165            Liste_j(Index)=j
166            Index=Index+1
167          endif
168          Pos=Pos+1
169        enddo
170      enddo
171   
172      if (pole_sud) then
173        Liste_i(Index)=1
174        Liste_j(Index)=jjp1
175        iiphy_end=1
176        jjphy_end=jjp1
177        ijphy_end=jjp1*iip1+1
178      endif
179   
180      ijphy_nb=ijphy_end-ijphy_begin+1
181      jjphy_nb=jjphy_end-jjphy_begin+1
182       
183           
184      Pos=2
185      rank=0
186      jjphy_para_begin(rank)=1
187   
188      do j=2,jjm
189        do i=1,iim
190          if (Pos==klon_para_begin(rank)) then
191            jjphy_para_begin(rank)=j
192          endif
193          if (Pos==klon_para_end(rank)) then
194            jjphy_para_end(rank)=j
195            rank=rank+1
196          endif
197          Pos=Pos+1
198        enddo
199      enddo
200      jjphy_para_end(rank)=jjp1
201     
202      First=.false.   
203   end subroutine InitDimphy
204   
205
206   subroutine GatherField(Fields,Fieldr,ll)
207#ifdef CPP_PARA
208      USE parallel, ONLY : COMM_LMDZ
209#endif
210     implicit none
211#ifdef CPP_PARA     
212     include 'mpif.h'
213#endif
214     INTEGER :: ll
215     REAL, dimension(klon_mpi,ll) :: Fields
216     REAL, dimension(klon2,ll) :: Fieldr     
217     REAL, dimension(klon2*ll) :: Field_tmp     
218     INTEGER, dimension(0:phy_size-1) :: displs
219     INTEGER, dimension(0:phy_size-1) :: sendcounts
220     INTEGER :: l,Pos,rank,ierr
221     INTEGER :: klon_b,klon_e, Nb
222     
223     Pos=1     
224     do rank=0,phy_size-1
225       klon_b=klon_para_begin(rank)
226       klon_e=klon_para_end(rank)
227       Nb=klon_para_nb(rank)
228       displs(rank)=Pos-1
229       sendcounts(rank)=Nb*ll
230       Pos=Pos+Nb*ll
231     enddo
232         
233      if (monocpu) then
234        Fieldr(:,:)=Fields(:,:)
235      else
236#ifdef CPP_PARA
237        call MPI_Gatherv(Fields,klon_mpi*ll,MPI_REAL8,Field_tmp,sendcounts,  &
238                         displs,MPI_REAL8,0,COMM_LMDZ,ierr)
239        Pos=1     
240        do rank=0,phy_size-1
241          klon_b=klon_para_begin(rank)
242          klon_e=klon_para_end(rank)
243          Nb=klon_para_nb(rank)
244          do l=1,ll
245           Fieldr(klon_b:klon_e,l)=Field_tmp(Pos:Pos+Nb-1)
246           Pos=Pos+Nb
247          enddo
248        enddo
249#endif
250     endif
251     
252   end subroutine GatherField
253
254   subroutine AllGatherField(Fields,Fieldr,ll)
255#ifdef CPP_PARA
256      USE parallel, ONLY : COMM_LMDZ
257#endif
258     implicit none
259#ifdef CPP_PARA     
260     include 'mpif.h'
261#endif   
262     INTEGER :: ll
263     REAL, dimension(klon_mpi,ll) :: Fields
264     REAL, dimension(klon2,ll) :: Fieldr     
265     
266     INTEGER :: ierr
267     
268      IF (monocpu) THEN
269        Fieldr(:,:)=Fields(:,:)
270      ELSE
271       call GatherField(Fields,Fieldr,ll)
272#ifdef CPP_PARA             
273       call MPI_BCAST(Fieldr,klon2*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
274#endif
275      ENDIF 
276   end subroutine AllGatherField
277   
278           
279   subroutine ScatterField(Fields,Fieldr,ll)
280#ifdef CPP_PARA
281      USE parallel, ONLY : COMM_LMDZ
282#endif
283     implicit none
284#ifdef CPP_PARA     
285     include 'mpif.h'
286#endif
287     INTEGER :: ll
288     REAL, dimension(klon2,ll) :: Fields
289     REAL, dimension(klon_mpi,ll) :: Fieldr     
290     REAL, dimension(klon2*ll) :: Field_tmp     
291     INTEGER, dimension(0:phy_size-1) :: displs
292     INTEGER, dimension(0:phy_size-1) :: sendcounts
293     INTEGER :: l,Pos,rank,ierr
294     INTEGER :: klon_b,klon_e, Nb
295
296     Pos=1   
297     do rank=0,phy_size-1
298       klon_b=klon_para_begin(rank)
299       klon_e=klon_para_end(rank)
300       Nb=klon_para_nb(rank)
301       displs(rank)=Pos-1
302       sendcounts(rank)=Nb*ll
303       do l=1,ll
304         Field_tmp(Pos:Pos+Nb-1)=Fields(klon_b:klon_e,l)
305         Pos=Pos+Nb
306       enddo
307     enddo
308           
309     if (monocpu) then
310       Fieldr(:,:)=Fields(:,:)
311     else
312#ifdef CPP_PARA     
313       call MPI_Scatterv(Field_tmp,sendcounts,displs,MPI_REAL8,         &
314                         Fieldr,klon_mpi*ll,MPI_REAL8,0,COMM_LMDZ,ierr)
315#endif
316     endif
317     
318       
319   end subroutine ScatterField
320   
321   
322end module mod_phys_mpi
Note: See TracBrowser for help on using the repository browser.