Changeset 16b61a5 in flexpart.git for src/mpi_mod.f90


Ignore:
Timestamp:
Oct 14, 2016, 3:19:00 PM (8 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, univie
Children:
4c64400
Parents:
861805a
Message:

Reworked the domain-filling option (MPI). Fixed a slow loop which had errors in loop counter (MPI)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/mpi_mod.f90

    r861805a r16b61a5  
    122122  logical, parameter :: mp_dev_mode = .false.
    123123  logical, parameter :: mp_dbg_out = .false.
    124   logical, parameter :: mp_time_barrier=.false.
     124  logical, parameter :: mp_time_barrier=.true.
    125125  logical, parameter :: mp_measure_time=.false.
    126126  logical, parameter :: mp_exact_numpart=.true.
     
    251251        write(*,FMT='(80("#"))')
    252252      end if
    253       lmp_sync=.true. ! :DBG: eso fix this...
     253      lmp_sync=.true.
    254254    end if
    255255
     
    261261!**********************************************************************
    262262
     263!    id_read = min(mp_np-1, 1)
    263264    id_read = mp_np-1
    264265
     
    486487    integer :: i,jj, addone
    487488
     489! Exit if too many particles
     490    if (num_part.gt.maxpart_mpi) then
     491      write(*,*) '#####################################################'
     492      write(*,*) '#### ERROR - TOTAL NUMBER OF PARTICLES REQUIRED  ####'
     493      write(*,*) '#### EXCEEDS THE MAXIMUM ALLOWED NUMBER. REDUCE  ####'
     494      write(*,*) '#### EITHER NUMBER OF PARTICLES PER RELEASE POINT####'
     495      write(*,*) '#### OR INCREASE MAXPART.                        ####'
     496      write(*,*) '#####################################################'
     497!      call MPI_FINALIZE(mp_ierr)
     498      stop
     499    end if
     500
     501
    488502! Time for MPI communications
    489503!****************************
     
    527541
    528542    if (mp_measure_time) call mpif_mtime('commtime',1)
    529     write(*,*) "PID ", mp_partid, "ending MPI_Scatter operation"
    530543
    531544    goto 601
     
    535548
    536549! After the transfer of particles it it possible that the value of
    537 ! "numpart" is larger than the actual, so we reduce it if there are
     550! "numpart" is larger than the actual used, so we reduce it if there are
    538551! invalid particles at the end of the arrays
    539552601 do i=num_part, 1, -1
     
    628641        if ( numparticles_mpi(idx_arr(m)).gt.mp_min_redist.and.&
    629642             & real(num_trans)/real(numparticles_mpi(idx_arr(m))).gt.mp_redist_fract) then
     643! DBG
     644          ! write(*,*) 'mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, numparticles_mpi', &
     645          !      &mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, numparticles_mpi
     646! DBG
    630647          call mpif_redist_part(itime, idx_arr(m), idx_arr(i), num_trans/2)
    631648        end if
     
    661678    integer :: i, j, minpart, ipart, maxnumpart
    662679 
     680! Check for invalid input arguments
     681!**********************************
     682 if (src_proc.eq.dest_proc) then
     683   write(*,*) '<mpi_mod::mpif_redist_part>: Error: &
     684        &src_proc.eq.dest_proc'
     685   stop
     686 end if
     687
    663688! Measure time for MPI communications
    664689!************************************
     
    674699      ul=numpart
    675700
    676       ! if (mp_dev_mode) then
    677       !   write(*,FMT='(72("#"))')
    678       !   write(*,*) "Sending ", num_trans, "particles (from/to)", src_proc, dest_proc
    679       !   write(*,*) "numpart before: ", numpart
    680       ! end if
     701      if (mp_dev_mode) then
     702        write(*,FMT='(72("#"))')
     703        write(*,*) "Sending ", num_trans, "particles (from/to)", src_proc, dest_proc
     704        write(*,*) "numpart before: ", numpart
     705      end if
    681706
    682707      call MPI_SEND(nclass(ll:ul),num_trans,&
     
    717742      numpart = numpart-num_trans
    718743
    719       ! if (mp_dev_mode) then
    720       !   write(*,*) "numpart after: ", numpart
    721       !   write(*,FMT='(72("#"))')
    722       ! end if
     744      if (mp_dev_mode) then
     745        write(*,*) "numpart after: ", numpart
     746        write(*,FMT='(72("#"))')
     747      end if
    723748
    724749    else if (mp_partid.eq.dest_proc) then
     
    731756      ul=numpart+num_trans
    732757
    733       ! if (mp_dev_mode) then
    734       !   write(*,FMT='(72("#"))')
    735       !   write(*,*) "Receiving ", num_trans, "particles (from/to)", src_proc, dest_proc
    736       !   write(*,*) "numpart before: ", numpart
    737       ! end if
     758      if (mp_dev_mode) then
     759        write(*,FMT='(72("#"))')
     760        write(*,*) "Receiving ", num_trans, "particles (from/to)", src_proc, dest_proc
     761        write(*,*) "numpart before: ", numpart
     762      end if
    738763
    739764! We could receive the data directly at nclass(ll:ul) etc., but this leaves
     
    785810      do i=1, num_trans
    786811! Take into acount that we may have transferred invalid particles
    787         if (itra1_tmp(minpart).ne.itime) goto 200
     812        if (itra1_tmp(i).ne.itime) cycle
    788813        do ipart=minpart,maxnumpart
    789814          if (itra1(ipart).ne.itime) then
    790             itra1(ipart) = itra1_tmp(minpart)
    791             npoint(ipart) = npoint_tmp(minpart)
    792             nclass(ipart) = nclass_tmp(minpart)
    793             idt(ipart) = idt_tmp(minpart)
    794             itramem(ipart) = itramem_tmp(minpart)
    795             itrasplit(ipart) = itrasplit_tmp(minpart)
    796             xtra1(ipart) = xtra1_tmp(minpart)
    797             ytra1(ipart) = ytra1_tmp(minpart)
    798             ztra1(ipart) = ztra1_tmp(minpart)
    799             xmass1(ipart,:) = xmass1_tmp(minpart,:)
    800 ! Increase numpart, if necessary
    801             numpart=max(numpart,ipart)
     815            itra1(ipart) = itra1_tmp(i)
     816            npoint(ipart) = npoint_tmp(i)
     817            nclass(ipart) = nclass_tmp(i)
     818            idt(ipart) = idt_tmp(i)
     819            itramem(ipart) = itramem_tmp(i)
     820            itrasplit(ipart) = itrasplit_tmp(i)
     821            xtra1(ipart) = xtra1_tmp(i)
     822            ytra1(ipart) = ytra1_tmp(i)
     823            ztra1(ipart) = ztra1_tmp(i)
     824            xmass1(ipart,:) = xmass1_tmp(i,:)
    802825            goto 200 ! Storage space has been found, stop searching
    803826          end if
     827! :TODO: add check for if too many particles requiried
    804828        end do
    805 200     minpart=minpart+1
     829200     minpart=ipart+1
    806830      end do
     831! Increase numpart, if necessary
     832      numpart=max(numpart,ipart)
    807833
    808834      deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,itrasplit_tmp,&
    809835           & xtra1_tmp,ytra1_tmp,ztra1_tmp,xmass1_tmp)
    810836
    811       ! if (mp_dev_mode) then
    812       !   write(*,*) "numpart after: ", numpart
    813       !   write(*,FMT='(72("#"))')
    814       ! end if
     837      if (mp_dev_mode) then
     838        write(*,*) "numpart after: ", numpart
     839        write(*,FMT='(72("#"))')
     840      end if
    815841
    816842    else
     
    27272753! write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR VERTTRANSFORM:',&
    27282754!      & mp_vt_time_total
    2729 ! NB: the 'flush' function is possibly a gfortran-specific extension
    2730           call flush()
     2755! NB: the 'flush' function is possibly a gfortran-specific extension,
     2756! comment out if it gives problems
     2757!          call flush()
    27312758        end if
    27322759      end do
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG