[1632] | 1 | ! |
---|
[1707] | 2 | ! $Id$ |
---|
[1632] | 3 | ! |
---|
| 4 | module Bands |
---|
[1864] | 5 | USE parallel_lmdz |
---|
[1632] | 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 |
---|
[1864] | 29 | USE parallel_lmdz |
---|
[1632] | 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 |
---|
[1864] | 43 | USE parallel_lmdz |
---|
[1632] | 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 | |
---|
| 89 | jj_nb_vanleer(:)=jj_nb_caldyn(:) |
---|
| 90 | jj_nb_dissip(:)=jj_nb_caldyn(:) |
---|
| 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 |
---|
[1864] | 106 | USE parallel_lmdz |
---|
[1707] | 107 | #ifdef CPP_PHYS |
---|
| 108 | ! Ehouarn: what follows is only related to // physics |
---|
[1632] | 109 | USE mod_phys_lmdz_para, ONLY : jj_para_begin,jj_para_end |
---|
| 110 | #endif |
---|
| 111 | IMPLICIT NONE |
---|
| 112 | INCLUDE 'dimensions.h' |
---|
| 113 | INTEGER :: i |
---|
| 114 | |
---|
| 115 | do i=0,mpi_size-1 |
---|
| 116 | jj_nb_vanleer2(i)=(jjm+1)/mpi_size |
---|
| 117 | if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1 |
---|
| 118 | enddo |
---|
| 119 | |
---|
[1707] | 120 | #ifdef CPP_PHYS |
---|
[1632] | 121 | do i=0,MPI_Size-1 |
---|
| 122 | jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
| 123 | if (i/=0) then |
---|
| 124 | if (jj_para_begin(i)==jj_para_end(i-1)) then |
---|
| 125 | jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1 |
---|
| 126 | endif |
---|
| 127 | endif |
---|
| 128 | enddo |
---|
| 129 | |
---|
| 130 | do i=0,MPI_Size-1 |
---|
| 131 | jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
| 132 | if (i/=0) then |
---|
| 133 | if (jj_para_begin(i)==jj_para_end(i-1)) then |
---|
| 134 | jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1 |
---|
| 135 | else |
---|
| 136 | jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1 |
---|
| 137 | jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1 |
---|
| 138 | endif |
---|
| 139 | endif |
---|
| 140 | enddo |
---|
| 141 | #endif |
---|
| 142 | CALL create_distrib(jj_Nb_Caldyn,distrib_caldyn) |
---|
| 143 | CALL create_distrib(jj_Nb_vanleer,distrib_vanleer) |
---|
| 144 | CALL create_distrib(jj_Nb_vanleer2,distrib_vanleer2) |
---|
| 145 | CALL create_distrib(jj_Nb_dissip,distrib_dissip) |
---|
| 146 | CALL create_distrib(jj_Nb_physic,distrib_physic) |
---|
| 147 | CALL create_distrib(jj_Nb_physic_bis,distrib_physic_bis) |
---|
| 148 | |
---|
| 149 | distrib_physic_bis%jjb_u=distrib_physic%jjb_u |
---|
| 150 | distrib_physic_bis%jje_u=distrib_physic%jje_u |
---|
| 151 | distrib_physic_bis%jjnb_u=distrib_physic%jjnb_u |
---|
| 152 | |
---|
| 153 | distrib_physic_bis%ijb_u=distrib_physic%ijb_u |
---|
| 154 | distrib_physic_bis%ije_u=distrib_physic%ije_u |
---|
| 155 | distrib_physic_bis%ijnb_u=distrib_physic%ijnb_u |
---|
| 156 | |
---|
| 157 | distrib_physic_bis%jjb_v=distrib_physic%jjb_v |
---|
| 158 | distrib_physic_bis%jje_v=distrib_physic%jje_v |
---|
| 159 | distrib_physic_bis%jjnb_v=distrib_physic%jjnb_v |
---|
| 160 | |
---|
| 161 | distrib_physic_bis%ijb_v=distrib_physic%ijb_v |
---|
| 162 | distrib_physic_bis%ije_v=distrib_physic%ije_v |
---|
| 163 | distrib_physic_bis%ijnb_v=distrib_physic%ijnb_v |
---|
| 164 | |
---|
| 165 | end subroutine Set_Bands |
---|
| 166 | |
---|
| 167 | |
---|
| 168 | subroutine AdjustBands_caldyn(new_dist) |
---|
| 169 | use times |
---|
[1864] | 170 | USE parallel_lmdz |
---|
[1632] | 171 | implicit none |
---|
| 172 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
| 173 | |
---|
| 174 | real :: minvalue,maxvalue |
---|
| 175 | integer :: min_proc,max_proc |
---|
| 176 | integer :: i,j |
---|
| 177 | real,allocatable,dimension(:) :: value |
---|
| 178 | integer,allocatable,dimension(:) :: index |
---|
| 179 | real :: tmpvalue |
---|
| 180 | integer :: tmpindex |
---|
| 181 | |
---|
| 182 | allocate(value(0:mpi_size-1)) |
---|
| 183 | allocate(index(0:mpi_size-1)) |
---|
| 184 | |
---|
| 185 | |
---|
| 186 | call allgather_timer_average |
---|
| 187 | |
---|
| 188 | do i=0,mpi_size-1 |
---|
| 189 | value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i) |
---|
| 190 | index(i)=i |
---|
| 191 | enddo |
---|
| 192 | |
---|
| 193 | do i=0,mpi_size-2 |
---|
| 194 | do j=i+1,mpi_size-1 |
---|
| 195 | if (value(i)>value(j)) then |
---|
| 196 | tmpvalue=value(i) |
---|
| 197 | value(i)=value(j) |
---|
| 198 | value(j)=tmpvalue |
---|
| 199 | |
---|
| 200 | tmpindex=index(i) |
---|
| 201 | index(i)=index(j) |
---|
| 202 | index(j)=tmpindex |
---|
| 203 | endif |
---|
| 204 | enddo |
---|
| 205 | enddo |
---|
| 206 | |
---|
| 207 | maxvalue=value(mpi_size-1) |
---|
| 208 | max_proc=index(mpi_size-1) |
---|
| 209 | |
---|
| 210 | do i=0,mpi_size-2 |
---|
| 211 | minvalue=value(i) |
---|
| 212 | min_proc=index(i) |
---|
| 213 | if (jj_nb_caldyn(max_proc)>3) then |
---|
| 214 | if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) then |
---|
| 215 | jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 |
---|
| 216 | jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 |
---|
| 217 | exit |
---|
| 218 | else |
---|
| 219 | if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) & |
---|
| 220 | -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) then |
---|
| 221 | jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 |
---|
| 222 | jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 |
---|
| 223 | exit |
---|
| 224 | endif |
---|
| 225 | endif |
---|
| 226 | endif |
---|
| 227 | enddo |
---|
| 228 | |
---|
| 229 | deallocate(value) |
---|
| 230 | deallocate(index) |
---|
| 231 | CALL create_distrib(jj_nb_caldyn,new_dist) |
---|
| 232 | |
---|
| 233 | end subroutine AdjustBands_caldyn |
---|
| 234 | |
---|
| 235 | subroutine AdjustBands_vanleer(new_dist) |
---|
| 236 | use times |
---|
[1864] | 237 | USE parallel_lmdz |
---|
[1632] | 238 | implicit none |
---|
| 239 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
| 240 | |
---|
| 241 | real :: minvalue,maxvalue |
---|
| 242 | integer :: min_proc,max_proc |
---|
| 243 | integer :: i,j |
---|
| 244 | real,allocatable,dimension(:) :: value |
---|
| 245 | integer,allocatable,dimension(:) :: index |
---|
| 246 | real :: tmpvalue |
---|
| 247 | integer :: tmpindex |
---|
| 248 | |
---|
| 249 | allocate(value(0:mpi_size-1)) |
---|
| 250 | allocate(index(0:mpi_size-1)) |
---|
| 251 | |
---|
| 252 | |
---|
| 253 | call allgather_timer_average |
---|
| 254 | |
---|
| 255 | do i=0,mpi_size-1 |
---|
| 256 | value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i) |
---|
| 257 | index(i)=i |
---|
| 258 | enddo |
---|
| 259 | |
---|
| 260 | do i=0,mpi_size-2 |
---|
| 261 | do j=i+1,mpi_size-1 |
---|
| 262 | if (value(i)>value(j)) then |
---|
| 263 | tmpvalue=value(i) |
---|
| 264 | value(i)=value(j) |
---|
| 265 | value(j)=tmpvalue |
---|
| 266 | |
---|
| 267 | tmpindex=index(i) |
---|
| 268 | index(i)=index(j) |
---|
| 269 | index(j)=tmpindex |
---|
| 270 | endif |
---|
| 271 | enddo |
---|
| 272 | enddo |
---|
| 273 | |
---|
| 274 | maxvalue=value(mpi_size-1) |
---|
| 275 | max_proc=index(mpi_size-1) |
---|
| 276 | |
---|
| 277 | do i=0,mpi_size-2 |
---|
| 278 | minvalue=value(i) |
---|
| 279 | min_proc=index(i) |
---|
| 280 | |
---|
| 281 | if (jj_nb_vanleer(max_proc)>3) then |
---|
| 282 | if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. & |
---|
| 283 | timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) then |
---|
| 284 | jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 |
---|
| 285 | jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 |
---|
| 286 | exit |
---|
| 287 | else |
---|
| 288 | if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) then |
---|
| 289 | jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 |
---|
| 290 | jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 |
---|
| 291 | exit |
---|
| 292 | endif |
---|
| 293 | endif |
---|
| 294 | endif |
---|
| 295 | enddo |
---|
| 296 | |
---|
| 297 | deallocate(value) |
---|
| 298 | deallocate(index) |
---|
| 299 | |
---|
| 300 | CALL create_distrib(jj_nb_vanleer,new_dist) |
---|
| 301 | |
---|
| 302 | end subroutine AdjustBands_vanleer |
---|
| 303 | |
---|
| 304 | subroutine AdjustBands_dissip(new_dist) |
---|
| 305 | use times |
---|
[1864] | 306 | USE parallel_lmdz |
---|
[1632] | 307 | implicit none |
---|
| 308 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
| 309 | |
---|
| 310 | real :: minvalue,maxvalue |
---|
| 311 | integer :: min_proc,max_proc |
---|
| 312 | integer :: i,j |
---|
| 313 | real,allocatable,dimension(:) :: value |
---|
| 314 | integer,allocatable,dimension(:) :: index |
---|
| 315 | real :: tmpvalue |
---|
| 316 | integer :: tmpindex |
---|
| 317 | |
---|
| 318 | allocate(value(0:mpi_size-1)) |
---|
| 319 | allocate(index(0:mpi_size-1)) |
---|
| 320 | |
---|
| 321 | |
---|
| 322 | call allgather_timer_average |
---|
| 323 | |
---|
| 324 | do i=0,mpi_size-1 |
---|
| 325 | value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i) |
---|
| 326 | index(i)=i |
---|
| 327 | enddo |
---|
| 328 | |
---|
| 329 | do i=0,mpi_size-2 |
---|
| 330 | do j=i+1,mpi_size-1 |
---|
| 331 | if (value(i)>value(j)) then |
---|
| 332 | tmpvalue=value(i) |
---|
| 333 | value(i)=value(j) |
---|
| 334 | value(j)=tmpvalue |
---|
| 335 | |
---|
| 336 | tmpindex=index(i) |
---|
| 337 | index(i)=index(j) |
---|
| 338 | index(j)=tmpindex |
---|
| 339 | endif |
---|
| 340 | enddo |
---|
| 341 | enddo |
---|
| 342 | |
---|
| 343 | maxvalue=value(mpi_size-1) |
---|
| 344 | max_proc=index(mpi_size-1) |
---|
| 345 | |
---|
| 346 | do i=0,mpi_size-2 |
---|
| 347 | minvalue=value(i) |
---|
| 348 | min_proc=index(i) |
---|
| 349 | |
---|
| 350 | if (jj_nb_dissip(max_proc)>3) then |
---|
| 351 | if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) then |
---|
| 352 | jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 |
---|
| 353 | jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 |
---|
| 354 | exit |
---|
| 355 | else |
---|
| 356 | if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) & |
---|
| 357 | - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) then |
---|
| 358 | jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 |
---|
| 359 | jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 |
---|
| 360 | exit |
---|
| 361 | endif |
---|
| 362 | endif |
---|
| 363 | endif |
---|
| 364 | enddo |
---|
| 365 | |
---|
| 366 | deallocate(value) |
---|
| 367 | deallocate(index) |
---|
| 368 | |
---|
| 369 | CALL create_distrib(jj_nb_dissip,new_dist) |
---|
| 370 | |
---|
| 371 | end subroutine AdjustBands_dissip |
---|
| 372 | |
---|
| 373 | subroutine AdjustBands_physic |
---|
| 374 | use times |
---|
[1707] | 375 | #ifdef CPP_PHYS |
---|
| 376 | ! Ehouarn: what follows is only related to // physics |
---|
[1632] | 377 | USE mod_phys_lmdz_para, only : klon_mpi_para_nb |
---|
| 378 | #endif |
---|
[1864] | 379 | USE parallel_lmdz |
---|
[1632] | 380 | implicit none |
---|
| 381 | |
---|
| 382 | integer :: i,Index |
---|
| 383 | real,allocatable,dimension(:) :: value |
---|
| 384 | integer,allocatable,dimension(:) :: Inc |
---|
| 385 | real :: medium |
---|
| 386 | integer :: NbTot,sgn |
---|
| 387 | |
---|
| 388 | allocate(value(0:mpi_size-1)) |
---|
| 389 | allocate(Inc(0:mpi_size-1)) |
---|
| 390 | |
---|
| 391 | |
---|
| 392 | call allgather_timer_average |
---|
| 393 | |
---|
| 394 | medium=0 |
---|
| 395 | do i=0,mpi_size-1 |
---|
| 396 | value(i)=timer_average(jj_nb_physic(i),timer_physic,i) |
---|
| 397 | medium=medium+value(i) |
---|
| 398 | enddo |
---|
| 399 | |
---|
| 400 | medium=medium/mpi_size |
---|
| 401 | NbTot=0 |
---|
[1707] | 402 | #ifdef CPP_PHYS |
---|
[1632] | 403 | do i=0,mpi_size-1 |
---|
| 404 | Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i)) |
---|
| 405 | NbTot=NbTot+Inc(i) |
---|
| 406 | enddo |
---|
| 407 | |
---|
| 408 | if (NbTot>=0) then |
---|
| 409 | Sgn=1 |
---|
| 410 | else |
---|
| 411 | Sgn=-1 |
---|
| 412 | NbTot=-NbTot |
---|
| 413 | endif |
---|
| 414 | |
---|
| 415 | Index=0 |
---|
| 416 | do i=1,NbTot |
---|
| 417 | Inc(Index)=Inc(Index)-Sgn |
---|
| 418 | Index=Index+1 |
---|
| 419 | if (Index>mpi_size-1) Index=0 |
---|
| 420 | enddo |
---|
| 421 | |
---|
| 422 | do i=0,mpi_size-1 |
---|
| 423 | distrib_phys(i)=klon_mpi_para_nb(i)+inc(i) |
---|
| 424 | enddo |
---|
| 425 | #endif |
---|
| 426 | |
---|
| 427 | end subroutine AdjustBands_physic |
---|
| 428 | |
---|
| 429 | subroutine WriteBands |
---|
[1864] | 430 | USE parallel_lmdz |
---|
[1632] | 431 | implicit none |
---|
| 432 | include "dimensions.h" |
---|
| 433 | |
---|
| 434 | integer :: i,j |
---|
| 435 | character (len=4) :: siim,sjjm,sllm,sproc |
---|
| 436 | character (len=255) :: filename |
---|
| 437 | integer :: unit_number=10 |
---|
| 438 | integer :: ierr |
---|
| 439 | |
---|
| 440 | write(siim,'(i3)') iim |
---|
| 441 | write(sjjm,'(i3)') jjm |
---|
| 442 | write(sllm,'(i3)') llm |
---|
| 443 | write(sproc,'(i3)') mpi_size |
---|
| 444 | |
---|
| 445 | filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_' & |
---|
| 446 | //TRIM(ADJUSTL(sproc))//'prc.dat' |
---|
| 447 | |
---|
| 448 | OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr) |
---|
| 449 | |
---|
| 450 | if (ierr==0) then |
---|
| 451 | |
---|
| 452 | ! write (unit_number,*) '*** Bandes caldyn ***' |
---|
| 453 | do i=0,mpi_size-1 |
---|
| 454 | write (unit_number,*) i,jj_nb_caldyn(i) |
---|
| 455 | enddo |
---|
| 456 | |
---|
| 457 | ! write (unit_number,*) '*** Bandes vanleer ***' |
---|
| 458 | do i=0,mpi_size-1 |
---|
| 459 | write (unit_number,*) i,jj_nb_vanleer(i) |
---|
| 460 | enddo |
---|
| 461 | |
---|
| 462 | ! write (unit_number,*) '*** Bandes dissip ***' |
---|
| 463 | do i=0,mpi_size-1 |
---|
| 464 | write (unit_number,*) i,jj_nb_dissip(i) |
---|
| 465 | enddo |
---|
| 466 | |
---|
| 467 | do i=0,mpi_size-1 |
---|
| 468 | write (unit_number,*) i,distrib_phys(i) |
---|
| 469 | enddo |
---|
| 470 | |
---|
| 471 | CLOSE(unit_number) |
---|
| 472 | else |
---|
| 473 | print *,'probleme lors de l ecriture des bandes' |
---|
| 474 | endif |
---|
| 475 | |
---|
| 476 | end subroutine WriteBands |
---|
| 477 | |
---|
| 478 | end module Bands |
---|
| 479 | |
---|
| 480 | |
---|