source: LMDZ4/branches/LMDZ4_par_0/libf/dyn3dpar/bands.F90 @ 5379

Last change on this file since 5379 was 630, checked in by Laurent Fairhead, 20 years ago

Import d'une version parallele de la dynamique YM
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.7 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          endif
106        endif
107      enddo
108     
109      do i=0,MPI_Size-1
110        New_klon_para_nb(i)=klon_para_nb(i)
111      enddo
112       
113    end subroutine InitBands
114
115    subroutine AdjustBands_caldyn
116      use times
117      implicit none
118
119      real :: minvalue,maxvalue
120      integer :: min_proc,max_proc
121      integer :: i,j
122      real,allocatable,dimension(:) :: value
123      integer,allocatable,dimension(:) :: index
124      real :: tmpvalue
125      integer :: tmpindex
126     
127      allocate(value(0:mpi_size-1))
128      allocate(index(0:mpi_size-1))
129       
130 
131      call allgather_timer_average
132
133      do i=0,mpi_size-1
134        value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
135        index(i)=i
136      enddo
137     
138      do i=0,mpi_size-2
139        do j=i+1,mpi_size-1
140          if (value(i)>value(j)) then
141            tmpvalue=value(i)
142            value(i)=value(j)
143            value(j)=tmpvalue
144           
145            tmpindex=index(i)
146            index(i)=index(j)
147            index(j)=tmpindex
148           endif
149         enddo
150      enddo
151     
152      maxvalue=value(mpi_size-1)
153      max_proc=index(mpi_size-1)           
154           
155      do i=0,mpi_size-2
156        minvalue=value(i)
157        min_proc=index(i)
158        if (jj_nb_caldyn(max_proc)>3) then
159          if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
160             jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
161             jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
162             exit
163           else
164             if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)                 &
165                -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
166               jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
167               jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
168               exit
169             endif
170           endif
171         endif
172      enddo
173     
174      deallocate(value)
175      deallocate(index)
176         
177    end subroutine AdjustBands_caldyn
178   
179    subroutine AdjustBands_vanleer
180      use times
181      implicit none
182
183      real :: minvalue,maxvalue
184      integer :: min_proc,max_proc
185      integer :: i,j
186      real,allocatable,dimension(:) :: value
187      integer,allocatable,dimension(:) :: index
188      real :: tmpvalue
189      integer :: tmpindex
190     
191      allocate(value(0:mpi_size-1))
192      allocate(index(0:mpi_size-1))
193       
194 
195      call allgather_timer_average
196
197      do i=0,mpi_size-1
198        value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
199        index(i)=i
200      enddo
201     
202      do i=0,mpi_size-2
203        do j=i+1,mpi_size-1
204          if (value(i)>value(j)) then
205            tmpvalue=value(i)
206            value(i)=value(j)
207            value(j)=tmpvalue
208           
209            tmpindex=index(i)
210            index(i)=index(j)
211            index(j)=tmpindex
212           endif
213         enddo
214      enddo
215     
216      maxvalue=value(mpi_size-1)
217      max_proc=index(mpi_size-1)           
218           
219      do i=0,mpi_size-2
220        minvalue=value(i)
221        min_proc=index(i)
222
223        if (jj_nb_vanleer(max_proc)>3) then
224          if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
225             timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
226             jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
227             jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
228             exit
229           else
230             if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
231               jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
232               jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
233               exit
234             endif
235           endif
236         endif
237      enddo
238     
239      deallocate(value)
240      deallocate(index)
241         
242    end subroutine AdjustBands_vanleer
243
244    subroutine AdjustBands_dissip
245      use times
246      implicit none
247
248      real :: minvalue,maxvalue
249      integer :: min_proc,max_proc
250      integer :: i,j
251      real,allocatable,dimension(:) :: value
252      integer,allocatable,dimension(:) :: index
253      real :: tmpvalue
254      integer :: tmpindex
255     
256      allocate(value(0:mpi_size-1))
257      allocate(index(0:mpi_size-1))
258       
259 
260      call allgather_timer_average
261
262      do i=0,mpi_size-1
263        value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
264        index(i)=i
265      enddo
266     
267      do i=0,mpi_size-2
268        do j=i+1,mpi_size-1
269          if (value(i)>value(j)) then
270            tmpvalue=value(i)
271            value(i)=value(j)
272            value(j)=tmpvalue
273           
274            tmpindex=index(i)
275            index(i)=index(j)
276            index(j)=tmpindex
277           endif
278         enddo
279      enddo
280     
281      maxvalue=value(mpi_size-1)
282      max_proc=index(mpi_size-1)           
283           
284      do i=0,mpi_size-2
285        minvalue=value(i)
286        min_proc=index(i)
287
288        if (jj_nb_dissip(max_proc)>3) then
289          if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
290             jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
291             jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
292             exit
293           else
294             if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)         &
295                - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
296               jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
297               jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
298               exit
299             endif
300           endif
301         endif
302      enddo
303     
304      deallocate(value)
305      deallocate(index)
306         
307    end subroutine AdjustBands_dissip
308
309    subroutine AdjustBands_physic
310      use times
311      use dimphy
312      implicit none
313
314      integer :: i,Index
315      real,allocatable,dimension(:) :: value
316      integer,allocatable,dimension(:) :: Inc
317      real :: medium
318      integer :: NbTot,sgn
319     
320      allocate(value(0:mpi_size-1))
321      allocate(Inc(0:mpi_size-1))
322       
323 
324      call allgather_timer_average
325     
326      medium=0
327      do i=0,mpi_size-1
328        value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
329        medium=medium+value(i)
330      enddo   
331     
332      medium=medium/mpi_size     
333      NbTot=0
334      do i=0,mpi_size-1
335        Inc(i)=nint(klon_para_nb(i)*(medium-value(i))/value(i))
336        NbTot=NbTot+Inc(i) 
337      enddo
338     
339      if (NbTot>=0) then
340        Sgn=1
341      else
342        Sgn=-1
343        NbTot=-NbTot
344      endif
345     
346      Index=0
347      do i=1,NbTot
348        Inc(Index)=Inc(Index)-Sgn
349        Index=Index+1
350        if (Index>mpi_size-1) Index=0
351      enddo
352     
353      do i=0,mpi_size-1
354        New_klon_para_nb(i)=klon_para_nb(i)+inc(i)
355      enddo
356     
357    end subroutine AdjustBands_physic
358
359    subroutine WriteBands
360    USE parallel
361    implicit none
362#include "dimensions90.h"
363
364      integer :: i,j
365      character (len=4) :: siim,sjjm,sllm,sproc
366      character (len=255) :: filename
367      integer :: unit_number=10
368      integer :: ierr
369 
370      write(siim,'(i3)') iim
371      write(sjjm,'(i3)') jjm
372      write(sllm,'(i3)') llm
373      write(sproc,'(i3)') mpi_size
374
375      filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_'  &
376                        //TRIM(ADJUSTL(sproc))//'prc.dat'   
377     
378      OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr)
379     
380      if (ierr==0) then
381       
382!       write (unit_number,*) '*** Bandes caldyn ***'
383        do i=0,mpi_size-1
384          write (unit_number,*) i,jj_nb_caldyn(i)
385        enddo
386       
387!       write (unit_number,*) '*** Bandes vanleer ***'
388        do i=0,mpi_size-1
389          write (unit_number,*) i,jj_nb_vanleer(i)
390        enddo
391       
392!        write (unit_number,*) '*** Bandes dissip ***'
393        do i=0,mpi_size-1
394          write (unit_number,*) i,jj_nb_dissip(i)
395        enddo
396       
397        do i=0,mpi_size-1
398          write (unit_number,*) i,New_klon_para_nb(i)
399        enddo
400       
401        CLOSE(unit_number)   
402      else
403        print *,'problème lors de l écriture des bandes'
404      endif
405       
406    end subroutine WriteBands
407 
408  end module Bands
409 
410 
Note: See TracBrowser for help on using the repository browser.