Changeset 0a94e13 in flexpart.git for src/mpi_mod.f90


Ignore:
Timestamp:
May 6, 2019, 11:43:21 AM (5 years ago)
Author:
Espen Sollum ATMOS <eso@…>
Branches:
master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug
Children:
0c8c7f2
Parents:
328fdf9
Message:

Added ipout=3 option for time averaged particle output

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/mpi_mod.f90

    • Property mode changed from 100644 to 100755
    r328fdf9 r0a94e13  
    8888! Variables for MPI processes in the 'particle' group
    8989  integer, allocatable, dimension(:) :: mp_partgroup_rank
     90  integer, allocatable, dimension(:) :: npart_per_process
    9091  integer :: mp_partgroup_comm, mp_partgroup_pid, mp_partgroup_np
    9192
     
    125126! mp_time_barrier   Measure MPI barrier time
    126127! mp_exact_numpart  Use an extra MPI communication to give the exact number of particles
    127 !                   to standard output (this does *not* otherwise affect the simulation)
    128 ! mp_rebalance      Attempt to rebalance particle between processes
     128!                   to standard output (this does not otherwise affect the simulation)
    129129  logical, parameter :: mp_dbg_mode = .false.
    130130  logical, parameter :: mp_dev_mode = .false.
     
    133133  logical, parameter :: mp_measure_time=.false.
    134134  logical, parameter :: mp_exact_numpart=.true.
    135   logical, parameter :: mp_rebalance=.true.
    136135
    137136! for measuring CPU/Wall time
     
    146145  real(dp),private :: mp_getfields_wtime_beg, mp_getfields_wtime_end, mp_getfields_wtime_total=0.
    147146  real(sp),private :: mp_getfields_time_beg, mp_getfields_time_end, mp_getfields_time_total=0.
    148 !  real(dp),private :: mp_readwind_wtime_beg, mp_readwind_wtime_end, mp_readwind_wtime_total=0.
    149 !  real(sp),private :: mp_readwind_time_beg, mp_readwind_time_end, mp_readwind_time_total=0.
     147  real(dp),private :: mp_readwind_wtime_beg, mp_readwind_wtime_end, mp_readwind_wtime_total=0.
     148  real(sp),private :: mp_readwind_time_beg, mp_readwind_time_end, mp_readwind_time_total=0.
    150149  real(dp),private :: mp_io_wtime_beg, mp_io_wtime_end, mp_io_wtime_total=0.
    151150  real(sp),private :: mp_io_time_beg, mp_io_time_end, mp_io_time_total=0.
     
    367366    end if
    368367
     368! Allocate array for number of particles per process   
     369    allocate(npart_per_process(0:mp_partgroup_np-1))
     370
    369371! Write whether MPI_IN_PLACE is used or not
    370372#ifdef USE_MPIINPLACE
     
    607609    integer :: i,jj,nn, num_part=1,m,imin, num_trans
    608610    logical :: first_iter
    609     integer,allocatable,dimension(:) :: numparticles_mpi, idx_arr
     611    integer,allocatable,dimension(:) :: idx_arr
    610612    real,allocatable,dimension(:) :: sorted ! TODO: we don't really need this
    611613
     
    616618! All processes exchange information on number of particles
    617619!****************************************************************************
    618     allocate(numparticles_mpi(0:mp_partgroup_np-1), &
    619          &idx_arr(0:mp_partgroup_np-1), sorted(0:mp_partgroup_np-1))
    620 
    621     call MPI_Allgather(numpart, 1, MPI_INTEGER, numparticles_mpi, &
     620    allocate( idx_arr(0:mp_partgroup_np-1), sorted(0:mp_partgroup_np-1))
     621
     622    call MPI_Allgather(numpart, 1, MPI_INTEGER, npart_per_process, &
    622623         & 1, MPI_INTEGER, mp_comm_used, mp_ierr)
    623624
     
    625626! Sort from lowest to highest
    626627! Initial guess: correct order
    627     sorted(:) = numparticles_mpi(:)
     628    sorted(:) = npart_per_process(:)
    628629    do i=0, mp_partgroup_np-1
    629630      idx_arr(i) = i
    630631    end do
     632
     633! Do not rebalance particles for ipout=3   
     634    if (ipout.eq.3) return
    631635
    632636! For each successive element in index array, see if a lower value exists
     
    654658    m=mp_partgroup_np-1 ! index for last sorted process (most particles)
    655659    do i=0,mp_partgroup_np/2-1
    656       num_trans = numparticles_mpi(idx_arr(m)) - numparticles_mpi(idx_arr(i))
     660      num_trans = npart_per_process(idx_arr(m)) - npart_per_process(idx_arr(i))
    657661      if (mp_partid.eq.idx_arr(m).or.mp_partid.eq.idx_arr(i)) then
    658         if ( numparticles_mpi(idx_arr(m)).gt.mp_min_redist.and.&
    659              & real(num_trans)/real(numparticles_mpi(idx_arr(m))).gt.mp_redist_fract) then
     662        if ( npart_per_process(idx_arr(m)).gt.mp_min_redist.and.&
     663             & real(num_trans)/real(npart_per_process(idx_arr(m))).gt.mp_redist_fract) then
    660664! DBG
    661           ! write(*,*) 'mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, numparticles_mpi', &
    662           !      &mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, numparticles_mpi
     665          ! write(*,*) 'mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, npart_per_process', &
     666          !      &mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, npart_per_process
    663667! DBG
    664668          call mpif_redist_part(itime, idx_arr(m), idx_arr(i), num_trans/2)
     
    668672    end do
    669673
    670     deallocate(numparticles_mpi, idx_arr, sorted)
     674    deallocate(idx_arr, sorted)
    671675
    672676  end subroutine mpif_calculate_part_redist
     
    27162720      end if
    27172721
    2718 !    case ('readwind')
    2719 !      if (imode.eq.0) then
    2720 !        call cpu_time(mp_readwind_time_beg)
    2721 !        mp_readwind_wtime_beg = mpi_wtime()
    2722 !      else
    2723 !        call cpu_time(mp_readwind_time_end)
    2724 !        mp_readwind_wtime_end = mpi_wtime()
    2725 !
    2726 !        mp_readwind_time_total = mp_readwind_time_total + &
    2727 !             &(mp_readwind_time_end - mp_readwind_time_beg)
    2728 !        mp_readwind_wtime_total = mp_readwind_wtime_total + &
    2729 !             &(mp_readwind_wtime_end - mp_readwind_wtime_beg)
    2730 !      end if
     2722   case ('readwind')
     2723     if (imode.eq.0) then
     2724       call cpu_time(mp_readwind_time_beg)
     2725       mp_readwind_wtime_beg = mpi_wtime()
     2726     else
     2727       call cpu_time(mp_readwind_time_end)
     2728       mp_readwind_wtime_end = mpi_wtime()
     2729
     2730       mp_readwind_time_total = mp_readwind_time_total + &
     2731            &(mp_readwind_time_end - mp_readwind_time_beg)
     2732       mp_readwind_wtime_total = mp_readwind_wtime_total + &
     2733            &(mp_readwind_wtime_end - mp_readwind_wtime_beg)
     2734     end if
    27312735
    27322736    case ('commtime')
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG