[1279] | 1 | ! |
---|
| 2 | ! $Id: bands.F90 1907 2013-11-26 13:10:46Z fairhead $ |
---|
| 3 | ! |
---|
[630] | 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 |
---|
[774] | 16 | INTEGER,dimension(:),allocatable :: distrib_phys |
---|
[630] | 17 | |
---|
| 18 | contains |
---|
| 19 | |
---|
| 20 | subroutine AllocateBands |
---|
[1823] | 21 | USE parallel_lmdz |
---|
[630] | 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)) |
---|
[774] | 30 | allocate(distrib_phys(0:MPI_Size-1)) |
---|
[630] | 31 | |
---|
| 32 | end subroutine AllocateBands |
---|
| 33 | |
---|
[774] | 34 | subroutine Read_distrib |
---|
[1823] | 35 | USE parallel_lmdz |
---|
[630] | 36 | implicit none |
---|
| 37 | |
---|
[792] | 38 | include "dimensions.h" |
---|
[630] | 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 |
---|
[774] | 70 | read (unit_number,*) j,distrib_phys(i) |
---|
[630] | 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(:) |
---|
[774] | 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 |
---|
[630] | 91 | |
---|
[774] | 92 | |
---|
| 93 | SUBROUTINE Set_Bands |
---|
[1823] | 94 | USE parallel_lmdz |
---|
[1615] | 95 | #ifdef CPP_PHYS |
---|
| 96 | ! Ehouarn: what follows is only related to // physics |
---|
[774] | 97 | USE mod_phys_lmdz_para, ONLY : jj_para_begin,jj_para_end |
---|
[1279] | 98 | #endif |
---|
[774] | 99 | IMPLICIT NONE |
---|
[792] | 100 | INCLUDE 'dimensions.h' |
---|
[774] | 101 | INTEGER :: i |
---|
| 102 | |
---|
[630] | 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 | |
---|
[1615] | 108 | #ifdef CPP_PHYS |
---|
[630] | 109 | do i=0,MPI_Size-1 |
---|
[774] | 110 | jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
[630] | 111 | if (i/=0) then |
---|
[774] | 112 | if (jj_para_begin(i)==jj_para_end(i-1)) then |
---|
[630] | 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 |
---|
[774] | 119 | jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
[630] | 120 | if (i/=0) then |
---|
[774] | 121 | if (jj_para_begin(i)==jj_para_end(i-1)) then |
---|
[630] | 122 | jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1 |
---|
[764] | 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 |
---|
[630] | 127 | endif |
---|
| 128 | enddo |
---|
[1279] | 129 | #endif |
---|
[630] | 130 | |
---|
[774] | 131 | end subroutine Set_Bands |
---|
[630] | 132 | |
---|
[774] | 133 | |
---|
[630] | 134 | subroutine AdjustBands_caldyn |
---|
| 135 | use times |
---|
[1823] | 136 | USE parallel_lmdz |
---|
[630] | 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 | |
---|
| 151 | call allgather_timer_average |
---|
| 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 | |
---|
| 199 | subroutine AdjustBands_vanleer |
---|
| 200 | use times |
---|
[1823] | 201 | USE parallel_lmdz |
---|
[630] | 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 | |
---|
| 216 | call allgather_timer_average |
---|
| 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 | |
---|
| 265 | subroutine AdjustBands_dissip |
---|
| 266 | use times |
---|
[1823] | 267 | USE parallel_lmdz |
---|
[630] | 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 | |
---|
| 282 | call allgather_timer_average |
---|
| 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 | |
---|
| 331 | subroutine AdjustBands_physic |
---|
| 332 | use times |
---|
[1615] | 333 | #ifdef CPP_PHYS |
---|
| 334 | ! Ehouarn: what follows is only related to // physics |
---|
[774] | 335 | USE mod_phys_lmdz_para, only : klon_mpi_para_nb |
---|
[1279] | 336 | #endif |
---|
[1823] | 337 | USE parallel_lmdz |
---|
[630] | 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 | |
---|
| 350 | call allgather_timer_average |
---|
| 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 |
---|
[1615] | 360 | #ifdef CPP_PHYS |
---|
[630] | 361 | do i=0,mpi_size-1 |
---|
[774] | 362 | Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i)) |
---|
[630] | 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 |
---|
[774] | 381 | distrib_phys(i)=klon_mpi_para_nb(i)+inc(i) |
---|
[630] | 382 | enddo |
---|
[1279] | 383 | #endif |
---|
[630] | 384 | end subroutine AdjustBands_physic |
---|
| 385 | |
---|
| 386 | subroutine WriteBands |
---|
[1823] | 387 | USE parallel_lmdz |
---|
[630] | 388 | implicit none |
---|
[792] | 389 | include "dimensions.h" |
---|
[630] | 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 |
---|
[774] | 425 | write (unit_number,*) i,distrib_phys(i) |
---|
[630] | 426 | enddo |
---|
| 427 | |
---|
| 428 | CLOSE(unit_number) |
---|
| 429 | else |
---|
[764] | 430 | print *,'probleme lors de l ecriture des bandes' |
---|
[630] | 431 | endif |
---|
| 432 | |
---|
| 433 | end subroutine WriteBands |
---|
| 434 | |
---|
| 435 | end module Bands |
---|
| 436 | |
---|
| 437 | |
---|