source: LMDZ6/trunk/libf/dyn3dmem/bands.f90

Last change on this file was 5271, checked in by abarral, 28 hours ago

Move dimensions.h into a module
Nb: doesn't compile yet

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