source: LMDZ5/trunk/libf/dyn3dmem/bands.F90 @ 1632

Last change on this file since 1632 was 1632, checked in by Laurent Fairhead, 12 years ago

Import initial du répertoire dyn3dmem

Attention! ceci n'est qu'une version préliminaire du code "basse mémoire":
le code contenu dans ce répertoire est basé sur la r1320 et a donc besoin
d'être mis à jour par rapport à la dynamique parallèle d'aujourd'hui.
Ce code est toutefois mis à disposition pour circonvenir à des problèmes
de mémoire que certaines configurations du modèle pourraient rencontrer.
Dans l'état, il compile et tourne sur vargas et au CCRT


Initial import of dyn3dmem

Warning! this is just a preliminary version of the memory light code:
it is based on r1320 of the code and thus needs to be updated before
it can replace the present dyn3dpar code. It is nevertheless put at your
disposal to circumvent some memory problems some LMDZ configurations may
encounter. In its present state, it will compile and run on vargas and CCRT

File size: 13.3 KB
Line 
1!
2! $Id: bands.F90 1279 2009-12-10 09:02:56Z fairhead $
3!
4  module Bands
5  USE parallel
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
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
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
107#ifdef CPP_EARTH
108! Ehouarn: what follows is only related to // physics; for now only for Earth
109     USE mod_phys_lmdz_para, ONLY : jj_para_begin,jj_para_end
110#endif
111     IMPLICIT NONE
112     INCLUDE 'dimensions.h'   
113     INTEGER :: i
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#ifdef CPP_EARTH
121! Ehouarn: what follows is only related to // physics; for now only for Earth         
122      do i=0,MPI_Size-1
123        jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1
124        if (i/=0) then
125          if (jj_para_begin(i)==jj_para_end(i-1)) then
126            jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1
127          endif
128        endif
129      enddo
130     
131      do i=0,MPI_Size-1
132        jj_Nb_physic_bis(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_bis(i)=jj_Nb_physic_bis(i)-1
136          else
137            jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1
138            jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1
139          endif
140        endif
141      enddo
142#endif     
143      CALL create_distrib(jj_Nb_Caldyn,distrib_caldyn)
144      CALL create_distrib(jj_Nb_vanleer,distrib_vanleer)
145      CALL create_distrib(jj_Nb_vanleer2,distrib_vanleer2)
146      CALL create_distrib(jj_Nb_dissip,distrib_dissip)
147      CALL create_distrib(jj_Nb_physic,distrib_physic)
148      CALL create_distrib(jj_Nb_physic_bis,distrib_physic_bis)
149     
150      distrib_physic_bis%jjb_u=distrib_physic%jjb_u
151      distrib_physic_bis%jje_u=distrib_physic%jje_u
152      distrib_physic_bis%jjnb_u=distrib_physic%jjnb_u
153
154      distrib_physic_bis%ijb_u=distrib_physic%ijb_u
155      distrib_physic_bis%ije_u=distrib_physic%ije_u
156      distrib_physic_bis%ijnb_u=distrib_physic%ijnb_u
157
158      distrib_physic_bis%jjb_v=distrib_physic%jjb_v
159      distrib_physic_bis%jje_v=distrib_physic%jje_v
160      distrib_physic_bis%jjnb_v=distrib_physic%jjnb_v
161
162      distrib_physic_bis%ijb_v=distrib_physic%ijb_v
163      distrib_physic_bis%ije_v=distrib_physic%ije_v
164      distrib_physic_bis%ijnb_v=distrib_physic%ijnb_v
165     
166    end subroutine Set_Bands
167
168
169    subroutine AdjustBands_caldyn(new_dist)
170      use times
171      use parallel
172      implicit none
173      TYPE(distrib),INTENT(INOUT) :: new_dist
174
175      real :: minvalue,maxvalue
176      integer :: min_proc,max_proc
177      integer :: i,j
178      real,allocatable,dimension(:) :: value
179      integer,allocatable,dimension(:) :: index
180      real :: tmpvalue
181      integer :: tmpindex
182     
183      allocate(value(0:mpi_size-1))
184      allocate(index(0:mpi_size-1))
185       
186 
187      call allgather_timer_average
188
189      do i=0,mpi_size-1
190        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
191        index(i)=i
192      enddo
193     
194      do i=0,mpi_size-2
195        do j=i+1,mpi_size-1
196          if (value(i)>value(j)) then
197            tmpvalue=value(i)
198            value(i)=value(j)
199            value(j)=tmpvalue
200           
201            tmpindex=index(i)
202            index(i)=index(j)
203            index(j)=tmpindex
204           endif
205         enddo
206      enddo
207     
208      maxvalue=value(mpi_size-1)
209      max_proc=index(mpi_size-1)           
210           
211      do i=0,mpi_size-2
212        minvalue=value(i)
213        min_proc=index(i)
214        if (jj_nb_caldyn(max_proc)>3) then
215          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
216             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
217             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
218             exit
219           else
220             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
221                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
222               jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
223               jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
224               exit
225             endif
226           endif
227         endif
228      enddo
229     
230      deallocate(value)
231      deallocate(index)
232      CALL create_distrib(jj_nb_caldyn,new_dist)
233       
234    end subroutine AdjustBands_caldyn
235   
236    subroutine AdjustBands_vanleer(new_dist)
237      use times
238      use parallel
239      implicit none
240      TYPE(distrib),INTENT(INOUT) :: new_dist
241
242      real :: minvalue,maxvalue
243      integer :: min_proc,max_proc
244      integer :: i,j
245      real,allocatable,dimension(:) :: value
246      integer,allocatable,dimension(:) :: index
247      real :: tmpvalue
248      integer :: tmpindex
249     
250      allocate(value(0:mpi_size-1))
251      allocate(index(0:mpi_size-1))
252       
253 
254      call allgather_timer_average
255
256      do i=0,mpi_size-1
257        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
258        index(i)=i
259      enddo
260     
261      do i=0,mpi_size-2
262        do j=i+1,mpi_size-1
263          if (value(i)>value(j)) then
264            tmpvalue=value(i)
265            value(i)=value(j)
266            value(j)=tmpvalue
267           
268            tmpindex=index(i)
269            index(i)=index(j)
270            index(j)=tmpindex
271           endif
272         enddo
273      enddo
274     
275      maxvalue=value(mpi_size-1)
276      max_proc=index(mpi_size-1)           
277           
278      do i=0,mpi_size-2
279        minvalue=value(i)
280        min_proc=index(i)
281
282        if (jj_nb_vanleer(max_proc)>3) then
283          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
284             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
285             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
286             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
287             exit
288           else
289             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
290               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
291               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
292               exit
293             endif
294           endif
295         endif
296      enddo
297     
298      deallocate(value)
299      deallocate(index)
300 
301      CALL create_distrib(jj_nb_vanleer,new_dist)
302         
303    end subroutine AdjustBands_vanleer
304
305    subroutine AdjustBands_dissip(new_dist)
306      use times
307      use parallel
308      implicit none
309      TYPE(distrib),INTENT(INOUT) :: new_dist
310     
311      real :: minvalue,maxvalue
312      integer :: min_proc,max_proc
313      integer :: i,j
314      real,allocatable,dimension(:) :: value
315      integer,allocatable,dimension(:) :: index
316      real :: tmpvalue
317      integer :: tmpindex
318     
319      allocate(value(0:mpi_size-1))
320      allocate(index(0:mpi_size-1))
321       
322 
323      call allgather_timer_average
324
325      do i=0,mpi_size-1
326        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
327        index(i)=i
328      enddo
329     
330      do i=0,mpi_size-2
331        do j=i+1,mpi_size-1
332          if (value(i)>value(j)) then
333            tmpvalue=value(i)
334            value(i)=value(j)
335            value(j)=tmpvalue
336           
337            tmpindex=index(i)
338            index(i)=index(j)
339            index(j)=tmpindex
340           endif
341         enddo
342      enddo
343     
344      maxvalue=value(mpi_size-1)
345      max_proc=index(mpi_size-1)           
346           
347      do i=0,mpi_size-2
348        minvalue=value(i)
349        min_proc=index(i)
350
351        if (jj_nb_dissip(max_proc)>3) then
352          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
353             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
354             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
355             exit
356           else
357             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
358                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
359               jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
360               jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
361               exit
362             endif
363           endif
364         endif
365      enddo
366     
367      deallocate(value)
368      deallocate(index)
369 
370      CALL create_distrib(jj_nb_dissip,new_dist)
371         
372    end subroutine AdjustBands_dissip
373
374    subroutine AdjustBands_physic
375      use times
376#ifdef CPP_EARTH
377! Ehouarn: what follows is only related to // physics; for now only for Earth
378      USE mod_phys_lmdz_para, only : klon_mpi_para_nb
379#endif
380      USE parallel
381      implicit none
382
383      integer :: i,Index
384      real,allocatable,dimension(:) :: value
385      integer,allocatable,dimension(:) :: Inc
386      real :: medium
387      integer :: NbTot,sgn
388     
389      allocate(value(0:mpi_size-1))
390      allocate(Inc(0:mpi_size-1))
391       
392 
393      call allgather_timer_average
394     
395      medium=0
396      do i=0,mpi_size-1
397        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
398        medium=medium+value(i)
399      enddo   
400     
401      medium=medium/mpi_size     
402      NbTot=0
403#ifdef CPP_EARTH
404! Ehouarn: what follows is only related to // physics; for now only for Earth
405      do i=0,mpi_size-1
406        Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i))
407        NbTot=NbTot+Inc(i) 
408      enddo
409     
410      if (NbTot>=0) then
411        Sgn=1
412      else
413        Sgn=-1
414        NbTot=-NbTot
415      endif
416     
417      Index=0
418      do i=1,NbTot
419        Inc(Index)=Inc(Index)-Sgn
420        Index=Index+1
421        if (Index>mpi_size-1) Index=0
422      enddo
423     
424      do i=0,mpi_size-1
425        distrib_phys(i)=klon_mpi_para_nb(i)+inc(i)
426      enddo
427#endif 
428         
429    end subroutine AdjustBands_physic
430
431    subroutine WriteBands
432    USE parallel
433    implicit none
434    include "dimensions.h"
435
436      integer :: i,j
437      character (len=4) :: siim,sjjm,sllm,sproc
438      character (len=255) :: filename
439      integer :: unit_number=10
440      integer :: ierr
441 
442      write(siim,'(i3)') iim
443      write(sjjm,'(i3)') jjm
444      write(sllm,'(i3)') llm
445      write(sproc,'(i3)') mpi_size
446
447      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
448                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
449     
450      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
451     
452      if (ierr==0) then
453       
454!       write (unit_number,*) '*** Bandes caldyn ***'
455        do i=0,mpi_size-1
456          write (unit_number,*) i,jj_nb_caldyn(i)
457        enddo
458       
459!       write (unit_number,*) '*** Bandes vanleer ***'
460        do i=0,mpi_size-1
461          write (unit_number,*) i,jj_nb_vanleer(i)
462        enddo
463       
464!        write (unit_number,*) '*** Bandes dissip ***'
465        do i=0,mpi_size-1
466          write (unit_number,*) i,jj_nb_dissip(i)
467        enddo
468       
469        do i=0,mpi_size-1
470          write (unit_number,*) i,distrib_phys(i)
471        enddo
472       
473        CLOSE(unit_number)   
474      else
475        print *,'probleme lors de l ecriture des bandes'
476      endif
477       
478    end subroutine WriteBands
479 
480  end module Bands
481 
482 
Note: See TracBrowser for help on using the repository browser.