source: LMDZ4/branches/LMDZ4_par_0/libf/phylmd/dimphy.F90 @ 1566

Last change on this file since 1566 was 633, checked in by (none), 20 years ago

This commit was manufactured by cvs2svn to create branch 'LMDZ4_par_0'.

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