Changeset d9f0585 in flexpart.git for src/FLEXPART_MPI.f90


Ignore:
Timestamp:
May 8, 2017, 8:34:40 AM (7 years ago)
Author:
Sabine <sabine.eckhardt@…>
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.
Message:

Merge branch 'dev' of git.nilu.no:flexpart/flexpart into dev

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/FLEXPART_MPI.f90

    r79abee9 rc8fc724  
    5555
    5656
    57 ! Initialize mpi
    58 !*********************
     57  ! Initialize mpi
     58  !*********************
    5959  call mpif_init
    6060
    6161  if (mp_measure_time) call mpif_mtime('flexpart',0)
    6262
    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
    6668
    6769  ! Generate a large number of random numbers
     
    7981  flexversion_major = '10' ! Major version number, also used for species file names
    8082!  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)'
    8284  verbosity=0
    8385
     
    306308  endif
    307309
    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
    311315
    312316  ! For continuation of previous run, read in particle positions
     
    318322    endif
    319323    ! readwind process skips this step
    320     if (lmp_use_reader.and..not.lmpreader) call readpartpositions
     324    if (.not.(lmpreader.and.lmp_use_reader)) call readpartpositions
    321325  else
    322326    if (verbosity.gt.0 .and. lroot) then
     
    425429  end do
    426430
     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
    427441
    428442! Calculate particle trajectories
     
    448462! NIK 16.02.2005
    449463  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, &
    451465         & 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, &
    453467         & mp_comm_used, mp_ierr)
    454468  else
    455469    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, &
    457471           & 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, &
    459473           & mp_comm_used, mp_ierr)
    460474    end if
     
    462476
    463477  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
    468487
    469488    write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG