My Project
 All Classes Files Functions Variables Macros
bands.F90
Go to the documentation of this file.
1 !
2 ! $Id: bands.F90 1615 2012-02-10 15:42:26Z emillour $
3 !
4  module bands
5 
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  INTEGER,dimension(:),allocatable :: distrib_phys
17 
18  contains
19 
20  subroutine allocatebands
21  use parallel
22  implicit none
23 
24  allocate(jj_nb_caldyn(0:mpi_size-1))
25  allocate(jj_nb_vanleer(0:mpi_size-1))
26  allocate(jj_nb_vanleer2(0:mpi_size-1))
27  allocate(jj_nb_dissip(0:mpi_size-1))
28  allocate(jj_nb_physic(0:mpi_size-1))
29  allocate(jj_nb_physic_bis(0:mpi_size-1))
30  allocate(distrib_phys(0:mpi_size-1))
31 
32  end subroutine allocatebands
33 
34  subroutine read_distrib
35  use parallel
36  implicit none
37 
38  include "dimensions.h"
39  integer :: i,j
40  character (len=4) :: siim,sjjm,sllm,sproc
41  character (len=255) :: filename
42  integer :: unit_number=10
43  integer :: ierr
44 
45  call allocatebands
46  write(siim,'(i3)') iim
47  write(sjjm,'(i3)') jjm
48  write(sllm,'(i3)') llm
49  write(sproc,'(i3)') mpi_size
50  filename='Bands_'//trim(adjustl(siim))//'x'//trim(adjustl(sjjm))//'x'//trim(adjustl(sllm))//'_' &
51  //trim(adjustl(sproc))//'prc.dat'
52 
53  OPEN(unit=unit_number,file=trim(filename),status='old',form='formatted',iostat=ierr)
54 
55  if (ierr==0) then
56 
57  do i=0,mpi_size-1
58  read (unit_number,*) j,jj_nb_caldyn(i)
59  enddo
60 
61  do i=0,mpi_size-1
62  read (unit_number,*) j,jj_nb_vanleer(i)
63  enddo
64 
65  do i=0,mpi_size-1
66  read (unit_number,*) j,jj_nb_dissip(i)
67  enddo
68 
69  do i=0,mpi_size-1
70  read (unit_number,*) j,distrib_phys(i)
71  enddo
72 
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  do i=0,mpi_size-1
85  distrib_phys(i)=(iim*(jjm-1)+2)/mpi_size
86  IF (i<mod(iim*(jjm-1)+2,mpi_size)) distrib_phys(i)=distrib_phys(i)+1
87  enddo
88  endif
89 
90  end subroutine read_distrib
91 
92 
93  SUBROUTINE set_bands
94  USE parallel
95 #ifdef CPP_PHYS
96 ! Ehouarn: what follows is only related to // physics
97  USE mod_phys_lmdz_para, ONLY : jj_para_begin,jj_para_end
98 #endif
99  IMPLICIT NONE
100  include 'dimensions.h'
101  INTEGER :: i
102 
103  do i=0,mpi_size-1
104  jj_nb_vanleer2(i)=(jjm+1)/mpi_size
105  if (i<mod(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1
106  enddo
107 
108 #ifdef CPP_PHYS
109  do i=0,mpi_size-1
110  jj_nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1
111  if (i/=0) then
112  if (jj_para_begin(i)==jj_para_end(i-1)) then
113  jj_nb_physic(i-1)=jj_nb_physic(i-1)-1
114  endif
115  endif
116  enddo
117 
118  do i=0,mpi_size-1
119  jj_nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1
120  if (i/=0) then
121  if (jj_para_begin(i)==jj_para_end(i-1)) then
122  jj_nb_physic_bis(i)=jj_nb_physic_bis(i)-1
123  else
124  jj_nb_physic_bis(i-1)=jj_nb_physic_bis(i-1)+1
125  jj_nb_physic_bis(i)=jj_nb_physic_bis(i)-1
126  endif
127  endif
128  enddo
129 #endif
130 
131  end subroutine set_bands
132 
133 
135  use times
136  use parallel
137  implicit none
138 
139  real :: minvalue,maxvalue
140  integer :: min_proc,max_proc
141  integer :: i,j
142  real,allocatable,dimension(:) :: value
143  integer,allocatable,dimension(:) :: index
144  real :: tmpvalue
145  integer :: tmpindex
146 
147  allocate(value(0:mpi_size-1))
148  allocate(index(0:mpi_size-1))
149 
150 
152 
153  do i=0,mpi_size-1
154  value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i)
155  index(i)=i
156  enddo
157 
158  do i=0,mpi_size-2
159  do j=i+1,mpi_size-1
160  if (value(i)>value(j)) then
161  tmpvalue=value(i)
162  value(i)=value(j)
163  value(j)=tmpvalue
164 
165  tmpindex=index(i)
166  index(i)=index(j)
167  index(j)=tmpindex
168  endif
169  enddo
170  enddo
171 
172  maxvalue=value(mpi_size-1)
173  max_proc=index(mpi_size-1)
174 
175  do i=0,mpi_size-2
176  minvalue=value(i)
177  min_proc=index(i)
178  if (jj_nb_caldyn(max_proc)>3) then
179  if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then
180  jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
181  jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
182  exit
183  else
184  if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) &
185  -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then
186  jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1
187  jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1
188  exit
189  endif
190  endif
191  endif
192  enddo
193 
194  deallocate(value)
195  deallocate(index)
196 
197  end subroutine adjustbands_caldyn
198 
200  use times
201  use parallel
202  implicit none
203 
204  real :: minvalue,maxvalue
205  integer :: min_proc,max_proc
206  integer :: i,j
207  real,allocatable,dimension(:) :: value
208  integer,allocatable,dimension(:) :: index
209  real :: tmpvalue
210  integer :: tmpindex
211 
212  allocate(value(0:mpi_size-1))
213  allocate(index(0:mpi_size-1))
214 
215 
217 
218  do i=0,mpi_size-1
219  value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i)
220  index(i)=i
221  enddo
222 
223  do i=0,mpi_size-2
224  do j=i+1,mpi_size-1
225  if (value(i)>value(j)) then
226  tmpvalue=value(i)
227  value(i)=value(j)
228  value(j)=tmpvalue
229 
230  tmpindex=index(i)
231  index(i)=index(j)
232  index(j)=tmpindex
233  endif
234  enddo
235  enddo
236 
237  maxvalue=value(mpi_size-1)
238  max_proc=index(mpi_size-1)
239 
240  do i=0,mpi_size-2
241  minvalue=value(i)
242  min_proc=index(i)
243 
244  if (jj_nb_vanleer(max_proc)>3) then
245  if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. &
246  timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then
247  jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
248  jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
249  exit
250  else
251  if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then
252  jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1
253  jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1
254  exit
255  endif
256  endif
257  endif
258  enddo
259 
260  deallocate(value)
261  deallocate(index)
262 
263  end subroutine adjustbands_vanleer
264 
266  use times
267  use parallel
268  implicit none
269 
270  real :: minvalue,maxvalue
271  integer :: min_proc,max_proc
272  integer :: i,j
273  real,allocatable,dimension(:) :: value
274  integer,allocatable,dimension(:) :: index
275  real :: tmpvalue
276  integer :: tmpindex
277 
278  allocate(value(0:mpi_size-1))
279  allocate(index(0:mpi_size-1))
280 
281 
283 
284  do i=0,mpi_size-1
285  value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i)
286  index(i)=i
287  enddo
288 
289  do i=0,mpi_size-2
290  do j=i+1,mpi_size-1
291  if (value(i)>value(j)) then
292  tmpvalue=value(i)
293  value(i)=value(j)
294  value(j)=tmpvalue
295 
296  tmpindex=index(i)
297  index(i)=index(j)
298  index(j)=tmpindex
299  endif
300  enddo
301  enddo
302 
303  maxvalue=value(mpi_size-1)
304  max_proc=index(mpi_size-1)
305 
306  do i=0,mpi_size-2
307  minvalue=value(i)
308  min_proc=index(i)
309 
310  if (jj_nb_dissip(max_proc)>3) then
311  if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then
312  jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
313  jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
314  exit
315  else
316  if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) &
317  - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then
318  jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1
319  jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1
320  exit
321  endif
322  endif
323  endif
324  enddo
325 
326  deallocate(value)
327  deallocate(index)
328 
329  end subroutine adjustbands_dissip
330 
332  use times
333 #ifdef CPP_PHYS
334 ! Ehouarn: what follows is only related to // physics
335  USE mod_phys_lmdz_para, only : klon_mpi_para_nb
336 #endif
337  USE parallel
338  implicit none
339 
340  integer :: i,index
341  real,allocatable,dimension(:) :: value
342  integer,allocatable,dimension(:) :: inc
343  real :: medium
344  integer :: nbtot,sgn
345 
346  allocate(value(0:mpi_size-1))
347  allocate(inc(0:mpi_size-1))
348 
349 
351 
352  medium=0
353  do i=0,mpi_size-1
354  value(i)=timer_average(jj_nb_physic(i),timer_physic,i)
355  medium=medium+value(i)
356  enddo
357 
358  medium=medium/mpi_size
359  nbtot=0
360 #ifdef CPP_PHYS
361  do i=0,mpi_size-1
362  inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i))
363  nbtot=nbtot+inc(i)
364  enddo
365 
366  if (nbtot>=0) then
367  sgn=1
368  else
369  sgn=-1
370  nbtot=-nbtot
371  endif
372 
373  index=0
374  do i=1,nbtot
375  inc(index)=inc(index)-sgn
376  index=index+1
377  if (index>mpi_size-1) index=0
378  enddo
379 
380  do i=0,mpi_size-1
381  distrib_phys(i)=klon_mpi_para_nb(i)+inc(i)
382  enddo
383 #endif
384  end subroutine adjustbands_physic
385 
386  subroutine writebands
387  USE parallel
388  implicit none
389  include "dimensions.h"
390 
391  integer :: i,j
392  character (len=4) :: siim,sjjm,sllm,sproc
393  character (len=255) :: filename
394  integer :: unit_number=10
395  integer :: ierr
396 
397  write(siim,'(i3)') iim
398  write(sjjm,'(i3)') jjm
399  write(sllm,'(i3)') llm
400  write(sproc,'(i3)') mpi_size
401 
402  filename='Bands_'//trim(adjustl(siim))//'x'//trim(adjustl(sjjm))//'x'//trim(adjustl(sllm))//'_' &
403  //trim(adjustl(sproc))//'prc.dat'
404 
405  OPEN(unit=unit_number,file=trim(filename),status='replace',form='formatted',iostat=ierr)
406 
407  if (ierr==0) then
408 
409 ! write (unit_number,*) '*** Bandes caldyn ***'
410  do i=0,mpi_size-1
411  write (unit_number,*) i,jj_nb_caldyn(i)
412  enddo
413 
414 ! write (unit_number,*) '*** Bandes vanleer ***'
415  do i=0,mpi_size-1
416  write (unit_number,*) i,jj_nb_vanleer(i)
417  enddo
418 
419 ! write (unit_number,*) '*** Bandes dissip ***'
420  do i=0,mpi_size-1
421  write (unit_number,*) i,jj_nb_dissip(i)
422  enddo
423 
424  do i=0,mpi_size-1
425  write (unit_number,*) i,distrib_phys(i)
426  enddo
427 
428  CLOSE(unit_number)
429  else
430  print *,'probleme lors de l ecriture des bandes'
431  endif
432 
433  end subroutine writebands
434 
435  end module bands
436 
437