source: LMDZ5/branches/LMDZ_tree_FC/libf/dyn3dmem/bands.F90 @ 2924

Last change on this file since 2924 was 2771, checked in by Laurent Fairhead, 8 years ago

Modifications to allow MPI domain partition on 2 latitude bands
rather than 3
LF

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 13.3 KB
Line 
1!
2! $Id$
3!
4  module Bands
5  USE parallel_lmdz
6    integer, parameter :: bands_caldyn=1
7    integer, parameter :: bands_vanleer=2
8    integer, parameter :: bands_dissip=3
9   
10    INTEGER,dimension(:),allocatable :: jj_Nb_Caldyn
11    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer
12    INTEGER,dimension(:),allocatable :: jj_Nb_vanleer2
13    INTEGER,dimension(:),allocatable :: jj_Nb_dissip
14    INTEGER,dimension(:),allocatable :: jj_Nb_physic
15    INTEGER,dimension(:),allocatable :: jj_Nb_physic_bis
16   
17    TYPE(distrib),SAVE,TARGET :: distrib_Caldyn
18    TYPE(distrib),SAVE,TARGET :: distrib_vanleer
19    TYPE(distrib),SAVE,TARGET :: distrib_vanleer2
20    TYPE(distrib),SAVE,TARGET :: distrib_dissip
21    TYPE(distrib),SAVE,TARGET :: distrib_physic
22    TYPE(distrib),SAVE,TARGET :: distrib_physic_bis
23
24    INTEGER,dimension(:),allocatable :: distrib_phys
25 
26  contains
27 
28  subroutine AllocateBands
29    USE parallel_lmdz
30    implicit none
31   
32    allocate(jj_Nb_Caldyn(0:MPI_Size-1))
33    allocate(jj_Nb_vanleer(0:MPI_Size-1))
34    allocate(jj_Nb_vanleer2(0:MPI_Size-1))
35    allocate(jj_Nb_dissip(0:MPI_Size-1))
36    allocate(jj_Nb_physic(0:MPI_Size-1))
37    allocate(jj_Nb_physic_bis(0:MPI_Size-1))
38    allocate(distrib_phys(0:MPI_Size-1))
39 
40  end subroutine AllocateBands
41 
42  subroutine Read_distrib
43    USE parallel_lmdz
44    implicit none
45
46    include "dimensions.h"
47      integer :: i,j
48      character (len=4) :: siim,sjjm,sllm,sproc
49      character (len=255) :: filename
50      integer :: unit_number=10
51      integer :: ierr
52   
53      call AllocateBands
54      write(siim,'(i3)') iim
55      write(sjjm,'(i3)') jjm
56      write(sllm,'(i3)') llm
57      write(sproc,'(i3)') mpi_size
58      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
59                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
60       
61      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='old',FORM='formatted',IOSTAT=ierr)
62     
63      if (ierr==0) then
64     
65         do i=0,mpi_size-1
66          read (unit_number,*) j,jj_nb_caldyn(i)
67        enddo
68     
69        do i=0,mpi_size-1
70          read (unit_number,*) j,jj_nb_vanleer(i)
71        enddo
72     
73        do i=0,mpi_size-1
74          read (unit_number,*) j,jj_nb_dissip(i)
75        enddo
76     
77        do i=0,mpi_size-1
78          read (unit_number,*) j,distrib_phys(i)
79        enddo
80       
81        CLOSE(unit_number) 
82 
83      else
84        do i=0,mpi_size-1
85          jj_nb_caldyn(i)=(jjm+1)/mpi_size
86          if (i<MOD(jjm+1,mpi_size)) jj_nb_caldyn(i)=jj_nb_caldyn(i)+1
87        enddo
88     
89        jj_nb_vanleer(:)=jj_nb_caldyn(:)
90        jj_nb_dissip(:)=jj_nb_caldyn(:)
91       
92        do i=0,mpi_size-1
93          distrib_phys(i)=(iim*(jjm-1)+2)/mpi_size
94          IF (i<MOD(iim*(jjm-1)+2,mpi_size)) distrib_phys(i)=distrib_phys(i)+1
95        enddo
96      endif
97     
98!      distrib_phys(:)=jj_nb_caldyn(:)*iim
99!      distrib_phys(0) = distrib_phys(0) - (iim-1)
100!      distrib_phys(mpi_size-1) = distrib_phys(mpi_size-1) - (iim-1)
101     
102   end subroutine Read_distrib
103   
104   
105   SUBROUTINE  Set_Bands
106     USE parallel_lmdz
107     IMPLICIT NONE
108     INCLUDE 'dimensions.h'   
109     INTEGER :: i, ij
110     INTEGER :: jj_para_begin(0:mpi_size-1)
111     INTEGER :: jj_para_end(0:mpi_size-1)
112       
113      do i=0,mpi_size-1
114         jj_nb_vanleer2(i)=(jjm+1)/mpi_size
115         if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1
116      enddo
117         
118      jj_para_begin(0)=1
119      ij=distrib_phys(0)+iim-1
120      jj_para_end(0)=((ij-1)/iim)+1
121     
122      DO i=1,mpi_Size-1
123        ij=ij+1
124        jj_para_begin(i)=((ij-1)/iim)+1
125        ij=ij+distrib_phys(i)-1
126        jj_para_end(i)=((ij-1)/iim)+1
127      ENDDO
128 
129       do i=0,MPI_Size-1
130        jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1
131        if (i/=0) then
132          if (jj_para_begin(i)==jj_para_end(i-1)) then
133            jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1
134          endif
135        endif
136      enddo
137     
138      do i=0,MPI_Size-1
139        jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1
140        if (i/=0) then
141          if (jj_para_begin(i)==jj_para_end(i-1)) then
142            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
143          else
144            jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1
145            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
146          endif
147        endif
148      enddo
149
150      CALL create_distrib(jj_Nb_Caldyn,distrib_caldyn)
151      CALL create_distrib(jj_Nb_vanleer,distrib_vanleer)
152      CALL create_distrib(jj_Nb_vanleer2,distrib_vanleer2)
153      CALL create_distrib(jj_Nb_dissip,distrib_dissip)
154      CALL create_distrib(jj_Nb_physic,distrib_physic)
155      CALL create_distrib(jj_Nb_physic_bis,distrib_physic_bis)
156     
157      distrib_physic_bis%jjb_u=distrib_physic%jjb_u
158      distrib_physic_bis%jje_u=distrib_physic%jje_u
159      distrib_physic_bis%jjnb_u=distrib_physic%jjnb_u
160
161      distrib_physic_bis%ijb_u=distrib_physic%ijb_u
162      distrib_physic_bis%ije_u=distrib_physic%ije_u
163      distrib_physic_bis%ijnb_u=distrib_physic%ijnb_u
164
165      distrib_physic_bis%jjb_v=distrib_physic%jjb_v
166      distrib_physic_bis%jje_v=distrib_physic%jje_v
167      distrib_physic_bis%jjnb_v=distrib_physic%jjnb_v
168
169      distrib_physic_bis%ijb_v=distrib_physic%ijb_v
170      distrib_physic_bis%ije_v=distrib_physic%ije_v
171      distrib_physic_bis%ijnb_v=distrib_physic%ijnb_v
172     
173    end subroutine Set_Bands
174
175
176    subroutine AdjustBands_caldyn(new_dist)
177      use times
178      USE parallel_lmdz
179      implicit none
180      TYPE(distrib),INTENT(INOUT) :: new_dist
181
182      real :: minvalue,maxvalue
183      integer :: min_proc,max_proc
184      integer :: i,j
185      real,allocatable,dimension(:) :: value
186      integer,allocatable,dimension(:) :: index
187      real :: tmpvalue
188      integer :: tmpindex
189     
190      allocate(value(0:mpi_size-1))
191      allocate(index(0:mpi_size-1))
192       
193 
194      call allgather_timer_average
195
196      do i=0,mpi_size-1
197        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
198        index(i)=i
199      enddo
200     
201      do i=0,mpi_size-2
202        do j=i+1,mpi_size-1
203          if (value(i)>value(j)) then
204            tmpvalue=value(i)
205            value(i)=value(j)
206            value(j)=tmpvalue
207           
208            tmpindex=index(i)
209            index(i)=index(j)
210            index(j)=tmpindex
211           endif
212         enddo
213      enddo
214     
215      maxvalue=value(mpi_size-1)
216      max_proc=index(mpi_size-1)           
217           
218      do i=0,mpi_size-2
219        minvalue=value(i)
220        min_proc=index(i)
221        if (jj_nb_caldyn(max_proc)>2) then
222          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
223             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
224             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
225             exit
226           else
227             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
228                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
229               jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
230               jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
231               exit
232             endif
233           endif
234         endif
235      enddo
236     
237      deallocate(value)
238      deallocate(index)
239      CALL create_distrib(jj_nb_caldyn,new_dist)
240       
241    end subroutine AdjustBands_caldyn
242   
243    subroutine AdjustBands_vanleer(new_dist)
244      use times
245      USE parallel_lmdz
246      implicit none
247      TYPE(distrib),INTENT(INOUT) :: new_dist
248
249      real :: minvalue,maxvalue
250      integer :: min_proc,max_proc
251      integer :: i,j
252      real,allocatable,dimension(:) :: value
253      integer,allocatable,dimension(:) :: index
254      real :: tmpvalue
255      integer :: tmpindex
256     
257      allocate(value(0:mpi_size-1))
258      allocate(index(0:mpi_size-1))
259       
260 
261      call allgather_timer_average
262
263      do i=0,mpi_size-1
264        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
265        index(i)=i
266      enddo
267     
268      do i=0,mpi_size-2
269        do j=i+1,mpi_size-1
270          if (value(i)>value(j)) then
271            tmpvalue=value(i)
272            value(i)=value(j)
273            value(j)=tmpvalue
274           
275            tmpindex=index(i)
276            index(i)=index(j)
277            index(j)=tmpindex
278           endif
279         enddo
280      enddo
281     
282      maxvalue=value(mpi_size-1)
283      max_proc=index(mpi_size-1)           
284           
285      do i=0,mpi_size-2
286        minvalue=value(i)
287        min_proc=index(i)
288
289        if (jj_nb_vanleer(max_proc)>2) then
290          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
291             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
292             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
293             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
294             exit
295           else
296             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
297               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
298               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
299               exit
300             endif
301           endif
302         endif
303      enddo
304     
305      deallocate(value)
306      deallocate(index)
307 
308      CALL create_distrib(jj_nb_vanleer,new_dist)
309         
310    end subroutine AdjustBands_vanleer
311
312    subroutine AdjustBands_dissip(new_dist)
313      use times
314      USE parallel_lmdz
315      implicit none
316      TYPE(distrib),INTENT(INOUT) :: new_dist
317     
318      real :: minvalue,maxvalue
319      integer :: min_proc,max_proc
320      integer :: i,j
321      real,allocatable,dimension(:) :: value
322      integer,allocatable,dimension(:) :: index
323      real :: tmpvalue
324      integer :: tmpindex
325     
326      allocate(value(0:mpi_size-1))
327      allocate(index(0:mpi_size-1))
328       
329 
330      call allgather_timer_average
331
332      do i=0,mpi_size-1
333        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
334        index(i)=i
335      enddo
336     
337      do i=0,mpi_size-2
338        do j=i+1,mpi_size-1
339          if (value(i)>value(j)) then
340            tmpvalue=value(i)
341            value(i)=value(j)
342            value(j)=tmpvalue
343           
344            tmpindex=index(i)
345            index(i)=index(j)
346            index(j)=tmpindex
347           endif
348         enddo
349      enddo
350     
351      maxvalue=value(mpi_size-1)
352      max_proc=index(mpi_size-1)           
353           
354      do i=0,mpi_size-2
355        minvalue=value(i)
356        min_proc=index(i)
357
358        if (jj_nb_dissip(max_proc)>3) then
359          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
360             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
361             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
362             exit
363           else
364             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
365                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
366               jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
367               jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
368               exit
369             endif
370           endif
371         endif
372      enddo
373     
374      deallocate(value)
375      deallocate(index)
376 
377      CALL create_distrib(jj_nb_dissip,new_dist)
378         
379    end subroutine AdjustBands_dissip
380
381    subroutine AdjustBands_physic
382      use times
383#ifdef CPP_PHYS
384! Ehouarn: what follows is only related to // physics
385      USE mod_phys_lmdz_para, only : klon_mpi_para_nb
386#endif
387      USE parallel_lmdz
388      implicit none
389
390      integer :: i,Index
391      real,allocatable,dimension(:) :: value
392      integer,allocatable,dimension(:) :: Inc
393      real :: medium
394      integer :: NbTot,sgn
395     
396      allocate(value(0:mpi_size-1))
397      allocate(Inc(0:mpi_size-1))
398       
399 
400      call allgather_timer_average
401     
402      medium=0
403      do i=0,mpi_size-1
404        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
405        medium=medium+value(i)
406      enddo   
407     
408      medium=medium/mpi_size     
409      NbTot=0
410#ifdef CPP_PHYS
411      do i=0,mpi_size-1
412        Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i))
413        NbTot=NbTot+Inc(i) 
414      enddo
415     
416      if (NbTot>=0) then
417        Sgn=1
418      else
419        Sgn=-1
420        NbTot=-NbTot
421      endif
422     
423      Index=0
424      do i=1,NbTot
425        Inc(Index)=Inc(Index)-Sgn
426        Index=Index+1
427        if (Index>mpi_size-1) Index=0
428      enddo
429     
430      do i=0,mpi_size-1
431        distrib_phys(i)=klon_mpi_para_nb(i)+inc(i)
432      enddo
433#endif 
434         
435    end subroutine AdjustBands_physic
436
437    subroutine WriteBands
438    USE parallel_lmdz
439    implicit none
440    include "dimensions.h"
441
442      integer :: i,j
443      character (len=4) :: siim,sjjm,sllm,sproc
444      character (len=255) :: filename
445      integer :: unit_number=10
446      integer :: ierr
447 
448      write(siim,'(i3)') iim
449      write(sjjm,'(i3)') jjm
450      write(sllm,'(i3)') llm
451      write(sproc,'(i3)') mpi_size
452
453      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
454                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
455     
456      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
457     
458      if (ierr==0) then
459       
460!       write (unit_number,*) '*** Bandes caldyn ***'
461        do i=0,mpi_size-1
462          write (unit_number,*) i,jj_nb_caldyn(i)
463        enddo
464       
465!       write (unit_number,*) '*** Bandes vanleer ***'
466        do i=0,mpi_size-1
467          write (unit_number,*) i,jj_nb_vanleer(i)
468        enddo
469       
470!        write (unit_number,*) '*** Bandes dissip ***'
471        do i=0,mpi_size-1
472          write (unit_number,*) i,jj_nb_dissip(i)
473        enddo
474       
475        do i=0,mpi_size-1
476          write (unit_number,*) i,distrib_phys(i)
477        enddo
478       
479        CLOSE(unit_number)   
480      else
481        print *,'probleme lors de l ecriture des bandes'
482      endif
483       
484    end subroutine WriteBands
485 
486  end module Bands
487 
488 
489
Note: See TracBrowser for help on using the repository browser.