LMDZ
bands.F90
Go to the documentation of this file.
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
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
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 
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
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
143  else
146  endif
147  endif
148  enddo
149 
156 
159  distrib_physic_bis%jjnb_u=distrib_physic%jjnb_u
160 
163  distrib_physic_bis%ijnb_u=distrib_physic%ijnb_u
164 
167  distrib_physic_bis%jjnb_v=distrib_physic%jjnb_v
168 
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 
195 
196  do i=0,mpi_size-1
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)>3) 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 
262 
263  do i=0,mpi_size-1
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)>3) 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 
331 
332  do i=0,mpi_size-1
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 
401 
402  medium=0
403  do i=0,mpi_size-1
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
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 
Definition: bands.F90:4
type(distrib), target, save distrib_vanleer2
Definition: bands.F90:19
integer, dimension(:), allocatable distrib_phys
Definition: bands.F90:24
integer, parameter timer_physic
Definition: times.F90:10
integer, dimension(:), allocatable jj_nb_caldyn
Definition: bands.F90:10
subroutine adjustbands_dissip(new_dist)
Definition: bands.F90:313
type(distrib), target, save distrib_dissip
Definition: bands.F90:20
subroutine writebands
Definition: bands.F90:438
integer, save mpi_size
real, dimension(:,:,:), allocatable timer_average
Definition: times.F90:19
!$Id Turb_fcg_gcssold get_uvd hqturb_gcssold endif!large scale llm day day1 day day1 *dt_toga endif!time annee_ref dt_toga u_toga vq_toga w_prof vq_prof llm day day1 day day1 *dt_dice endif!time annee_ref dt_dice swup_dice vg_dice omega_dice tg_prof vg_profd w_profd omega_profd!do llm!print llm l llm
integer, dimension(:), allocatable jj_nb_dissip
Definition: bands.F90:13
integer, parameter timer_caldyn
Definition: times.F90:7
subroutine allgather_timer_average
Definition: times.F90:173
type(distrib), target, save distrib_vanleer
Definition: bands.F90:18
integer, dimension(:,:,:), allocatable timer_iteration
Definition: times.F90:18
subroutine create_distrib(jj_nb_new, d)
integer, parameter timer_dissip
Definition: times.F90:9
real, dimension(:,:,:), allocatable timer_delta
Definition: times.F90:20
integer, dimension(:), allocatable jj_nb_physic
Definition: bands.F90:14
type(distrib), target, save distrib_physic
Definition: bands.F90:21
subroutine adjustbands_vanleer(new_dist)
Definition: bands.F90:244
Definition: times.F90:1
subroutine allocatebands
Definition: bands.F90:29
integer, dimension(:), allocatable jj_nb_vanleer2
Definition: bands.F90:12
subroutine set_bands
Definition: bands.F90:106
integer, parameter bands_vanleer
Definition: bands.F90:7
c c zjulian c cym CALL iim cym klev iim
Definition: ini_bilKP_ave.h:24
integer, dimension(:), allocatable jj_nb_physic_bis
Definition: bands.F90:15
subroutine adjustbands_caldyn(new_dist)
Definition: bands.F90:177
integer, parameter bands_dissip
Definition: bands.F90:8
subroutine read_distrib
Definition: bands.F90:43
type(distrib), target, save distrib_physic_bis
Definition: bands.F90:22
integer, parameter timer_vanleer
Definition: times.F90:8
integer, dimension(:), allocatable jj_nb_vanleer
Definition: bands.F90:11
type(distrib), target, save distrib_caldyn
Definition: bands.F90:17
integer, parameter bands_caldyn
Definition: bands.F90:6
!$Header!integer nvarmx s s unit
Definition: gradsdef.h:20
subroutine adjustbands_physic
Definition: bands.F90:382