source: LMDZ4/branches/pre_V3/libf/phylmd/dimphy.F90 @ 4249

Last change on this file since 4249 was 642, checked in by Laurent Fairhead, 20 years ago

Pour une compilation correcte quand on est en monocpu
LF

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