! $Id$ module Bands USE parallel_lmdz integer, parameter :: bands_caldyn=1 integer, parameter :: bands_vanleer=2 integer, parameter :: bands_dissip=3 INTEGER,dimension(:),allocatable :: jj_Nb_Caldyn INTEGER,dimension(:),allocatable :: jj_Nb_vanleer INTEGER,dimension(:),allocatable :: jj_Nb_vanleer2 INTEGER,dimension(:),allocatable :: jj_Nb_dissip INTEGER,dimension(:),allocatable :: jj_Nb_physic INTEGER,dimension(:),allocatable :: jj_Nb_physic_bis TYPE(distrib),SAVE,TARGET :: distrib_Caldyn TYPE(distrib),SAVE,TARGET :: distrib_vanleer TYPE(distrib),SAVE,TARGET :: distrib_vanleer2 TYPE(distrib),SAVE,TARGET :: distrib_dissip TYPE(distrib),SAVE,TARGET :: distrib_physic TYPE(distrib),SAVE,TARGET :: distrib_physic_bis INTEGER,dimension(:),allocatable :: distrib_phys contains SUBROUTINE AllocateBands USE parallel_lmdz IMPLICIT NONE allocate(jj_Nb_Caldyn(0:MPI_Size-1)) allocate(jj_Nb_vanleer(0:MPI_Size-1)) allocate(jj_Nb_vanleer2(0:MPI_Size-1)) allocate(jj_Nb_dissip(0:MPI_Size-1)) allocate(jj_Nb_physic(0:MPI_Size-1)) allocate(jj_Nb_physic_bis(0:MPI_Size-1)) allocate(distrib_phys(0:MPI_Size-1)) END SUBROUTINE AllocateBands SUBROUTINE Read_distrib USE parallel_lmdz IMPLICIT NONE include "dimensions.h" INTEGER :: i,j character (len=4) :: siim,sjjm,sllm,sproc character (len=255) :: filename INTEGER :: unit_number=10 INTEGER :: ierr CALL AllocateBands WRITE(siim,'(i3)') iim WRITE(sjjm,'(i3)') jjm WRITE(sllm,'(i3)') llm WRITE(sproc,'(i3)') mpi_size filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_' & //TRIM(ADJUSTL(sproc))//'prc.dat' OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='old',FORM='formatted',IOSTAT=ierr) if (ierr==0) THEN do i=0,mpi_size-1 read (unit_number,*) j,jj_nb_caldyn(i) enddo do i=0,mpi_size-1 read (unit_number,*) j,jj_nb_vanleer(i) enddo do i=0,mpi_size-1 read (unit_number,*) j,jj_nb_dissip(i) enddo do i=0,mpi_size-1 read (unit_number,*) j,distrib_phys(i) enddo CLOSE(unit_number) else do i=0,mpi_size-1 jj_nb_caldyn(i)=(jjm+1)/mpi_size if (ivalue(j)) THEN tmpvalue=value(i) value(i)=value(j) value(j)=tmpvalue tmpindex=index(i) index(i)=index(j) index(j)=tmpindex endif enddo enddo maxvalue=value(mpi_size-1) max_proc=index(mpi_size-1) do i=0,mpi_size-2 minvalue=value(i) min_proc=index(i) if (jj_nb_caldyn(max_proc)>2) THEN if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) THEN jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 exit else if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) & -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) THEN jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 exit endif endif endif enddo deallocate(value) deallocate(index) CALL create_distrib(jj_nb_caldyn,new_dist) END SUBROUTINE AdjustBands_caldyn SUBROUTINE AdjustBands_vanleer(new_dist) use times USE parallel_lmdz IMPLICIT NONE TYPE(distrib),INTENT(INOUT) :: new_dist REAL :: minvalue,maxvalue INTEGER :: min_proc,max_proc INTEGER :: i,j real,allocatable,dimension(:) :: value integer,allocatable,dimension(:) :: index REAL :: tmpvalue INTEGER :: tmpindex allocate(value(0:mpi_size-1)) allocate(index(0:mpi_size-1)) CALL allgather_timer_average do i=0,mpi_size-1 value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i) index(i)=i enddo do i=0,mpi_size-2 do j=i+1,mpi_size-1 if (value(i)>value(j)) THEN tmpvalue=value(i) value(i)=value(j) value(j)=tmpvalue tmpindex=index(i) index(i)=index(j) index(j)=tmpindex endif enddo enddo maxvalue=value(mpi_size-1) max_proc=index(mpi_size-1) do i=0,mpi_size-2 minvalue=value(i) min_proc=index(i) if (jj_nb_vanleer(max_proc)>2) THEN if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. & timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) THEN jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 exit else if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) THEN jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 exit endif endif endif enddo deallocate(value) deallocate(index) CALL create_distrib(jj_nb_vanleer,new_dist) END SUBROUTINE AdjustBands_vanleer SUBROUTINE AdjustBands_dissip(new_dist) use times USE parallel_lmdz IMPLICIT NONE TYPE(distrib),INTENT(INOUT) :: new_dist REAL :: minvalue,maxvalue INTEGER :: min_proc,max_proc INTEGER :: i,j real,allocatable,dimension(:) :: value integer,allocatable,dimension(:) :: index REAL :: tmpvalue INTEGER :: tmpindex allocate(value(0:mpi_size-1)) allocate(index(0:mpi_size-1)) CALL allgather_timer_average do i=0,mpi_size-1 value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i) index(i)=i enddo do i=0,mpi_size-2 do j=i+1,mpi_size-1 if (value(i)>value(j)) THEN tmpvalue=value(i) value(i)=value(j) value(j)=tmpvalue tmpindex=index(i) index(i)=index(j) index(j)=tmpindex endif enddo enddo maxvalue=value(mpi_size-1) max_proc=index(mpi_size-1) do i=0,mpi_size-2 minvalue=value(i) min_proc=index(i) if (jj_nb_dissip(max_proc)>3) THEN if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) THEN jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 exit else if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) & - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) THEN jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 exit endif endif endif enddo deallocate(value) deallocate(index) CALL create_distrib(jj_nb_dissip,new_dist) END SUBROUTINE AdjustBands_dissip SUBROUTINE AdjustBands_physic use times ! Ehouarn: what follows is only related to // physics USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_PHYS USE lmdz_phys_para, ONLY: klon_mpi_para_nb USE parallel_lmdz IMPLICIT NONE INTEGER :: i,Index real,allocatable,dimension(:) :: value integer,allocatable,dimension(:) :: Inc REAL :: medium INTEGER :: NbTot,sgn allocate(value(0:mpi_size-1)) allocate(Inc(0:mpi_size-1)) CALL allgather_timer_average medium=0 do i=0,mpi_size-1 value(i)=timer_average(jj_nb_physic(i),timer_physic,i) medium=medium+value(i) enddo medium=medium/mpi_size NbTot=0 IF (CPPKEY_PHYS) THEN do i=0,mpi_size-1 Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i)) NbTot=NbTot+Inc(i) enddo if (NbTot>=0) THEN Sgn=1 else Sgn=-1 NbTot=-NbTot endif Index=0 do i=1,NbTot Inc(Index)=Inc(Index)-Sgn Index=Index+1 if (Index>mpi_size-1) Index=0 enddo do i=0,mpi_size-1 distrib_phys(i)=klon_mpi_para_nb(i)+inc(i) enddo END IF END SUBROUTINE AdjustBands_physic SUBROUTINE WriteBands USE parallel_lmdz IMPLICIT NONE include "dimensions.h" INTEGER :: i,j character (len=4) :: siim,sjjm,sllm,sproc character (len=255) :: filename INTEGER :: unit_number=10 INTEGER :: ierr WRITE(siim,'(i3)') iim WRITE(sjjm,'(i3)') jjm WRITE(sllm,'(i3)') llm WRITE(sproc,'(i3)') mpi_size filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_' & //TRIM(ADJUSTL(sproc))//'prc.dat' OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr) if (ierr==0) THEN ! write (unit_number,*) '*** Bandes caldyn ***' do i=0,mpi_size-1 write (unit_number,*) i,jj_nb_caldyn(i) enddo ! write (unit_number,*) '*** Bandes vanleer ***' do i=0,mpi_size-1 write (unit_number,*) i,jj_nb_vanleer(i) enddo ! write (unit_number,*) '*** Bandes dissip ***' do i=0,mpi_size-1 write (unit_number,*) i,jj_nb_dissip(i) enddo do i=0,mpi_size-1 write (unit_number,*) i,distrib_phys(i) enddo CLOSE(unit_number) else print *,'probleme lors de l ecriture des bandes' endif END SUBROUTINE WriteBands end module Bands