Changeset d9f0585 in flexpart.git for src/FLEXPART_MPI.f90
- Timestamp:
- May 8, 2017, 8:34:40 AM (7 years ago)
- Branches:
- master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug, univie
- Children:
- 6985a98
- Parents:
- d404d98 (diff), c8fc724 (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the (diff) links above to see all the changes relative to each parent. - File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
src/FLEXPART_MPI.f90
r79abee9 rc8fc724 55 55 56 56 57 ! Initialize mpi58 !*********************57 ! Initialize mpi 58 !********************* 59 59 call mpif_init 60 60 61 61 if (mp_measure_time) call mpif_mtime('flexpart',0) 62 62 63 ! Initialize arrays in com_mod 64 !***************************** 65 call com_mod_allocate_part(maxpart_mpi) 63 ! Initialize arrays in com_mod 64 !***************************** 65 66 if(.not.(lmpreader.and.lmp_use_reader)) call com_mod_allocate_part(maxpart_mpi) 67 66 68 67 69 ! Generate a large number of random numbers … … 79 81 flexversion_major = '10' ! Major version number, also used for species file names 80 82 ! flexversion='Ver. 10 Beta MPI (2015-05-01)' 81 flexversion='Ver. '//trim(flexversion_major)//' Beta MPI (2015-05-01)'83 flexversion='Ver. '//trim(flexversion_major)//'.1beta MPI (2016-11-02)' 82 84 verbosity=0 83 85 … … 306 308 endif 307 309 308 do j=1, size(itra1) ! maxpart_mpi 309 itra1(j)=-999999999 310 end do 310 if (.not.(lmpreader.and.lmp_use_reader)) then 311 do j=1, size(itra1) ! maxpart_mpi 312 itra1(j)=-999999999 313 end do 314 end if 311 315 312 316 ! For continuation of previous run, read in particle positions … … 318 322 endif 319 323 ! readwind process skips this step 320 if ( lmp_use_reader.and..not.lmpreader) call readpartpositions324 if (.not.(lmpreader.and.lmp_use_reader)) call readpartpositions 321 325 else 322 326 if (verbosity.gt.0 .and. lroot) then … … 425 429 end do 426 430 431 ! Inform whether output kernel is used or not 432 !********************************************* 433 434 if (lroot) then 435 if (lnokernel) then 436 write(*,*) "Concentrations are calculated without using kernel" 437 else 438 write(*,*) "Concentrations are calculated using kernel" 439 end if 440 end if 427 441 428 442 ! Calculate particle trajectories … … 448 462 ! NIK 16.02.2005 449 463 if (lroot) then 450 call MPI_Reduce(MPI_IN_PLACE, tot_blc_count, 1, MPI_INTEGER8, MPI_SUM, id_root, &464 call MPI_Reduce(MPI_IN_PLACE, tot_blc_count, nspec, MPI_INTEGER8, MPI_SUM, id_root, & 451 465 & mp_comm_used, mp_ierr) 452 call MPI_Reduce(MPI_IN_PLACE, tot_inc_count, 1, MPI_INTEGER8, MPI_SUM, id_root, &466 call MPI_Reduce(MPI_IN_PLACE, tot_inc_count, nspec, MPI_INTEGER8, MPI_SUM, id_root, & 453 467 & mp_comm_used, mp_ierr) 454 468 else 455 469 if (mp_partgroup_pid.ge.0) then ! Skip for readwind process 456 call MPI_Reduce(tot_blc_count, 0, 1, MPI_INTEGER8, MPI_SUM, id_root, &470 call MPI_Reduce(tot_blc_count, 0, nspec, MPI_INTEGER8, MPI_SUM, id_root, & 457 471 & mp_comm_used, mp_ierr) 458 call MPI_Reduce(tot_inc_count, 0, 1, MPI_INTEGER8, MPI_SUM, id_root, &472 call MPI_Reduce(tot_inc_count, 0, nspec, MPI_INTEGER8, MPI_SUM, id_root, & 459 473 & mp_comm_used, mp_ierr) 460 474 end if … … 462 476 463 477 if (lroot) then 464 write(*,*) '**********************************************' 465 write(*,*) 'Total number of occurences of below-cloud scavenging', tot_blc_count 466 write(*,*) 'Total number of occurences of in-cloud scavenging', tot_inc_count 467 write(*,*) '**********************************************' 478 do i=1,nspec 479 write(*,*) '**********************************************' 480 write(*,*) 'Scavenging statistics for species ', species(i), ':' 481 write(*,*) 'Total number of occurences of below-cloud scavenging', & 482 & tot_blc_count(i) 483 write(*,*) 'Total number of occurences of in-cloud scavenging', & 484 & tot_inc_count(i) 485 write(*,*) '**********************************************' 486 end do 468 487 469 488 write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
Note: See TracChangeset
for help on using the changeset viewer.