source: LMDZ4/trunk/libf/dyn3dpar/bands.F90 @ 1187

Last change on this file since 1187 was 792, checked in by Laurent Fairhead, 17 years ago

Modifications suite a la transformation des fichiers include pour
qu'ils soient compatibles a la fois au format fixe et au format libre
Un bon nombre de fichier *.inc du coup disparaissent
LF

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