source: LMDZ4/branches/V3_test/libf/dyn3dpar/bands.F90 @ 4462

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

Nettoyage version parallele
YM/LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.8 KB
Line 
1  module Bands
2 
3    integer, parameter :: bands_caldyn=1
4    integer, parameter :: bands_vanleer=2
5    integer, parameter :: bands_dissip=3
6   
7    INTEGER,dimension(:),allocatable :: jj_Nb_Caldyn
8    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer
9    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer2
10    INTEGER,dimension(:),allocatable :: jj_Nb_dissip
11    INTEGER,dimension(:),allocatable :: jj_Nb_physic
12    INTEGER,dimension(:),allocatable :: jj_Nb_physic_bis
13    INTEGER,dimension(:),allocatable :: New_klon_para_nb
14 
15  contains
16 
17  subroutine AllocateBands
18    use parallel
19    implicit none
20   
21    allocate(jj_Nb_Caldyn(0:MPI_Size-1))
22    allocate(jj_Nb_vanleer(0:MPI_Size-1))
23    allocate(jj_Nb_vanleer2(0:MPI_Size-1))
24    allocate(jj_Nb_dissip(0:MPI_Size-1))
25    allocate(jj_Nb_physic(0:MPI_Size-1))
26    allocate(jj_Nb_physic_bis(0:MPI_Size-1))
27    allocate(New_klon_para_nb(0:MPI_Size-1))
28 
29  end subroutine AllocateBands
30 
31  subroutine InitBands
32    use parallel
33    use dimphy
34    implicit none
35
36#include "dimensions90.h"
37!#include "paramet90.h"   
38      integer :: i,j
39      character (len=4) :: siim,sjjm,sllm,sproc
40      character (len=255) :: filename
41      integer :: unit_number=10
42      integer :: ierr
43   
44      call AllocateBands
45      write(siim,'(i3)') iim
46      write(sjjm,'(i3)') jjm
47      write(sllm,'(i3)') llm
48      write(sproc,'(i3)') mpi_size
49      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
50                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
51       
52      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='old',FORM='formatted',IOSTAT=ierr)
53     
54      if (ierr==0) then
55     
56         do i=0,mpi_size-1
57          read (unit_number,*) j,jj_nb_caldyn(i)
58        enddo
59     
60        do i=0,mpi_size-1
61          read (unit_number,*) j,jj_nb_vanleer(i)
62        enddo
63     
64        do i=0,mpi_size-1
65          read (unit_number,*) j,jj_nb_dissip(i)
66        enddo
67     
68        do i=0,mpi_size-1
69          read (unit_number,*) j,klon_para_nb(i)
70        enddo
71       
72        call InitDimphy
73        CLOSE(unit_number) 
74 
75      else
76        do i=0,mpi_size-1
77          jj_nb_caldyn(i)=(jjm+1)/mpi_size
78          if (i<MOD(jjm+1,mpi_size)) jj_nb_caldyn(i)=jj_nb_caldyn(i)+1
79        enddo
80     
81        jj_nb_vanleer(:)=jj_nb_caldyn(:)
82        jj_nb_dissip(:)=jj_nb_caldyn(:)
83   
84      endif
85   
86      do i=0,mpi_size-1
87         jj_nb_vanleer2(i)=(jjm+1)/mpi_size
88         if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1
89      enddo
90         
91      do i=0,MPI_Size-1
92        jj_Nb_physic(i)=jjphy_para_end(i)-jjphy_para_begin(i)+1
93        if (i/=0) then
94          if (jjphy_para_begin(i)==jjphy_para_end(i-1)) then
95            jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1
96          endif
97        endif
98      enddo
99     
100      do i=0,MPI_Size-1
101        jj_Nb_physic_bis(i)=jjphy_para_end(i)-jjphy_para_begin(i)+1
102        if (i/=0) then
103          if (jjphy_para_begin(i)==jjphy_para_end(i-1)) then
104            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
105          else
106            jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1
107            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
108          endif
109        endif
110      enddo
111     
112      do i=0,MPI_Size-1
113        New_klon_para_nb(i)=klon_para_nb(i)
114      enddo
115       
116    end subroutine InitBands
117
118    subroutine AdjustBands_caldyn
119      use times
120      implicit none
121
122      real :: minvalue,maxvalue
123      integer :: min_proc,max_proc
124      integer :: i,j
125      real,allocatable,dimension(:) :: value
126      integer,allocatable,dimension(:) :: index
127      real :: tmpvalue
128      integer :: tmpindex
129     
130      allocate(value(0:mpi_size-1))
131      allocate(index(0:mpi_size-1))
132       
133 
134      call allgather_timer_average
135
136      do i=0,mpi_size-1
137        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
138        index(i)=i
139      enddo
140     
141      do i=0,mpi_size-2
142        do j=i+1,mpi_size-1
143          if (value(i)>value(j)) then
144            tmpvalue=value(i)
145            value(i)=value(j)
146            value(j)=tmpvalue
147           
148            tmpindex=index(i)
149            index(i)=index(j)
150            index(j)=tmpindex
151           endif
152         enddo
153      enddo
154     
155      maxvalue=value(mpi_size-1)
156      max_proc=index(mpi_size-1)           
157           
158      do i=0,mpi_size-2
159        minvalue=value(i)
160        min_proc=index(i)
161        if (jj_nb_caldyn(max_proc)>3) then
162          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
163             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
164             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
165             exit
166           else
167             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
168                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
169               jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
170               jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
171               exit
172             endif
173           endif
174         endif
175      enddo
176     
177      deallocate(value)
178      deallocate(index)
179         
180    end subroutine AdjustBands_caldyn
181   
182    subroutine AdjustBands_vanleer
183      use times
184      implicit none
185
186      real :: minvalue,maxvalue
187      integer :: min_proc,max_proc
188      integer :: i,j
189      real,allocatable,dimension(:) :: value
190      integer,allocatable,dimension(:) :: index
191      real :: tmpvalue
192      integer :: tmpindex
193     
194      allocate(value(0:mpi_size-1))
195      allocate(index(0:mpi_size-1))
196       
197 
198      call allgather_timer_average
199
200      do i=0,mpi_size-1
201        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
202        index(i)=i
203      enddo
204     
205      do i=0,mpi_size-2
206        do j=i+1,mpi_size-1
207          if (value(i)>value(j)) then
208            tmpvalue=value(i)
209            value(i)=value(j)
210            value(j)=tmpvalue
211           
212            tmpindex=index(i)
213            index(i)=index(j)
214            index(j)=tmpindex
215           endif
216         enddo
217      enddo
218     
219      maxvalue=value(mpi_size-1)
220      max_proc=index(mpi_size-1)           
221           
222      do i=0,mpi_size-2
223        minvalue=value(i)
224        min_proc=index(i)
225
226        if (jj_nb_vanleer(max_proc)>3) then
227          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
228             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
229             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
230             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
231             exit
232           else
233             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
234               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
235               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
236               exit
237             endif
238           endif
239         endif
240      enddo
241     
242      deallocate(value)
243      deallocate(index)
244         
245    end subroutine AdjustBands_vanleer
246
247    subroutine AdjustBands_dissip
248      use times
249      implicit none
250
251      real :: minvalue,maxvalue
252      integer :: min_proc,max_proc
253      integer :: i,j
254      real,allocatable,dimension(:) :: value
255      integer,allocatable,dimension(:) :: index
256      real :: tmpvalue
257      integer :: tmpindex
258     
259      allocate(value(0:mpi_size-1))
260      allocate(index(0:mpi_size-1))
261       
262 
263      call allgather_timer_average
264
265      do i=0,mpi_size-1
266        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
267        index(i)=i
268      enddo
269     
270      do i=0,mpi_size-2
271        do j=i+1,mpi_size-1
272          if (value(i)>value(j)) then
273            tmpvalue=value(i)
274            value(i)=value(j)
275            value(j)=tmpvalue
276           
277            tmpindex=index(i)
278            index(i)=index(j)
279            index(j)=tmpindex
280           endif
281         enddo
282      enddo
283     
284      maxvalue=value(mpi_size-1)
285      max_proc=index(mpi_size-1)           
286           
287      do i=0,mpi_size-2
288        minvalue=value(i)
289        min_proc=index(i)
290
291        if (jj_nb_dissip(max_proc)>3) then
292          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
293             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
294             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
295             exit
296           else
297             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
298                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
299               jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
300               jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
301               exit
302             endif
303           endif
304         endif
305      enddo
306     
307      deallocate(value)
308      deallocate(index)
309         
310    end subroutine AdjustBands_dissip
311
312    subroutine AdjustBands_physic
313      use times
314      use dimphy
315      implicit none
316
317      integer :: i,Index
318      real,allocatable,dimension(:) :: value
319      integer,allocatable,dimension(:) :: Inc
320      real :: medium
321      integer :: NbTot,sgn
322     
323      allocate(value(0:mpi_size-1))
324      allocate(Inc(0:mpi_size-1))
325       
326 
327      call allgather_timer_average
328     
329      medium=0
330      do i=0,mpi_size-1
331        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
332        medium=medium+value(i)
333      enddo   
334     
335      medium=medium/mpi_size     
336      NbTot=0
337      do i=0,mpi_size-1
338        Inc(i)=nint(klon_para_nb(i)*(medium-value(i))/value(i))
339        NbTot=NbTot+Inc(i) 
340      enddo
341     
342      if (NbTot>=0) then
343        Sgn=1
344      else
345        Sgn=-1
346        NbTot=-NbTot
347      endif
348     
349      Index=0
350      do i=1,NbTot
351        Inc(Index)=Inc(Index)-Sgn
352        Index=Index+1
353        if (Index>mpi_size-1) Index=0
354      enddo
355     
356      do i=0,mpi_size-1
357        New_klon_para_nb(i)=klon_para_nb(i)+inc(i)
358      enddo
359     
360    end subroutine AdjustBands_physic
361
362    subroutine WriteBands
363    USE parallel
364    implicit none
365#include "dimensions90.h"
366
367      integer :: i,j
368      character (len=4) :: siim,sjjm,sllm,sproc
369      character (len=255) :: filename
370      integer :: unit_number=10
371      integer :: ierr
372 
373      write(siim,'(i3)') iim
374      write(sjjm,'(i3)') jjm
375      write(sllm,'(i3)') llm
376      write(sproc,'(i3)') mpi_size
377
378      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
379                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
380     
381      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
382     
383      if (ierr==0) then
384       
385!       write (unit_number,*) '*** Bandes caldyn ***'
386        do i=0,mpi_size-1
387          write (unit_number,*) i,jj_nb_caldyn(i)
388        enddo
389       
390!       write (unit_number,*) '*** Bandes vanleer ***'
391        do i=0,mpi_size-1
392          write (unit_number,*) i,jj_nb_vanleer(i)
393        enddo
394       
395!        write (unit_number,*) '*** Bandes dissip ***'
396        do i=0,mpi_size-1
397          write (unit_number,*) i,jj_nb_dissip(i)
398        enddo
399       
400        do i=0,mpi_size-1
401          write (unit_number,*) i,New_klon_para_nb(i)
402        enddo
403       
404        CLOSE(unit_number)   
405      else
406        print *,'probleme lors de l ecriture des bandes'
407      endif
408       
409    end subroutine WriteBands
410 
411  end module Bands
412 
413 
Note: See TracBrowser for help on using the repository browser.